#include <G4ModifiedTsai.hh>
|
| G4ModifiedTsai (const G4String &name="") |
|
| ~G4ModifiedTsai () override |
|
G4ThreeVector & | SampleDirection (const G4DynamicParticle *dp, G4double out_energy, G4int Z, const G4Material *mat=nullptr) final |
|
void | SamplePairDirections (const G4DynamicParticle *dp, G4double elecKinEnergy, G4double posiKinEnergy, G4ThreeVector &dirElectron, G4ThreeVector &dirPositron, G4int Z=0, const G4Material *mat=nullptr) final |
|
void | PrintGeneratorInformation () const final |
|
G4ModifiedTsai & | operator= (const G4ModifiedTsai &right)=delete |
|
| G4ModifiedTsai (const G4ModifiedTsai &)=delete |
|
| G4VEmAngularDistribution (const G4String &name) |
|
virtual | ~G4VEmAngularDistribution () |
|
virtual G4ThreeVector & | SampleDirection (const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0 |
|
virtual G4ThreeVector & | SampleDirectionForShell (const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, G4int shellID, const G4Material *) |
|
virtual void | SamplePairDirections (const G4DynamicParticle *dp, G4double elecKinEnergy, G4double posiKinEnergy, G4ThreeVector &dirElectron, G4ThreeVector &dirPositron, G4int Z=0, const G4Material *mat=nullptr) |
|
virtual void | PrintGeneratorInformation () const |
|
const G4String & | GetName () const |
|
G4VEmAngularDistribution & | operator= (const G4VEmAngularDistribution &right)=delete |
|
| G4VEmAngularDistribution (const G4VEmAngularDistribution &)=delete |
|
Definition at line 62 of file G4ModifiedTsai.hh.
◆ G4ModifiedTsai() [1/2]
G4ModifiedTsai::G4ModifiedTsai |
( |
const G4String & |
name = "" | ) |
|
|
explicit |
◆ ~G4ModifiedTsai()
G4ModifiedTsai::~G4ModifiedTsai |
( |
| ) |
|
|
overridedefault |
◆ G4ModifiedTsai() [2/2]
◆ operator=()
◆ PrintGeneratorInformation()
void G4ModifiedTsai::PrintGeneratorInformation |
( |
| ) |
const |
|
finalvirtual |
Reimplemented from G4VEmAngularDistribution.
Definition at line 143 of file G4ModifiedTsai.cc.
144{
146 G4cout <<
"Bremsstrahlung Angular Generator is Modified Tsai" <<
G4endl;
147 G4cout <<
"Distribution suggested by L.Urban (Geant3 manual (1993) Phys211)"
149 G4cout <<
"Derived from Tsai distribution (Rev Mod Phys 49,421(1977)) \n"
151}
G4GLOB_DLL std::ostream G4cout
◆ SampleDirection()
Implements G4VEmAngularDistribution.
Definition at line 76 of file G4ModifiedTsai.cc.
78{
79
81 G4double sint = std::sqrt((1 - cost)*(1 + cost));
83
86
88}
void set(double x, double y, double z)
Hep3Vector & rotateUz(const Hep3Vector &)
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4ThreeVector fLocalDirection
◆ SamplePairDirections()
Reimplemented from G4VEmAngularDistribution.
Definition at line 117 of file G4ModifiedTsai.cc.
123{
127
128 G4double cost = SampleCosTheta(elecKinEnergy);
129 G4double sint = std::sqrt((1. - cost)*(1. + cost));
130
131 dirElectron.
set(sint*cosp, sint*sinp, cost);
133
134 cost = SampleCosTheta(posiKinEnergy);
135 sint = std::sqrt((1. - cost)*(1. + cost));
136
137 dirPositron.
set(-sint*cosp, -sint*sinp, cost);
139}
The documentation for this class was generated from the following files: