Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ModifiedTsai Class Reference

#include <G4ModifiedTsai.hh>

+ Inheritance diagram for G4ModifiedTsai:

Public Member Functions

 G4ModifiedTsai (const G4String &name="")
 
virtual ~G4ModifiedTsai ()
 
virtual G4ThreeVectorSampleDirection (const G4DynamicParticle *dp, G4double out_energy, G4int Z, const G4Material *mat=0)
 
void PrintGeneratorInformation () const
 
- Public Member Functions inherited from G4VEmAngularDistribution
 G4VEmAngularDistribution (const G4String &name)
 
virtual ~G4VEmAngularDistribution ()
 
virtual G4ThreeVectorSampleDirection (const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
 
const G4StringGetName () const
 

Additional Inherited Members

- Protected Attributes inherited from G4VEmAngularDistribution
G4ThreeVector fLocalDirection
 

Detailed Description

Definition at line 64 of file G4ModifiedTsai.hh.

Constructor & Destructor Documentation

◆ G4ModifiedTsai()

G4ModifiedTsai::G4ModifiedTsai ( const G4String name = "")

Definition at line 63 of file G4ModifiedTsai.cc.

◆ ~G4ModifiedTsai()

G4ModifiedTsai::~G4ModifiedTsai ( )
virtual

Definition at line 67 of file G4ModifiedTsai.cc.

68{}

Member Function Documentation

◆ PrintGeneratorInformation()

void G4ModifiedTsai::PrintGeneratorInformation ( ) const

Definition at line 103 of file G4ModifiedTsai.cc.

104{
105 G4cout << "\n" << G4endl;
106 G4cout << "Bremsstrahlung Angular Generator is Modified Tsai" << G4endl;
107 G4cout << "Distribution suggested by L.Urban (Geant3 manual (1993) Phys211)"
108 << G4endl;
109 G4cout << "Derived from Tsai distribution (Rev Mod Phys 49,421(1977)) \n"
110 << G4endl;
111}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout

◆ SampleDirection()

G4ThreeVector & G4ModifiedTsai::SampleDirection ( const G4DynamicParticle dp,
G4double  out_energy,
G4int  Z,
const G4Material mat = 0 
)
virtual

Implements G4VEmAngularDistribution.

Definition at line 71 of file G4ModifiedTsai.cc.

73{
74 // Sample gamma angle (Z - axis along the parent particle).
75 // Universal distribution suggested by L. Urban (Geant3 manual (1993)
76 // Phys211) derived from Tsai distribution (Rev Mod Phys 49,421(1977))
77
78 G4double uMax = 2*(1. + dp->GetKineticEnergy()/electron_mass_c2);
79
80 const G4double a1 = 0.625;
81 const G4double a2 = 1.875;
82 const G4double border = 0.25;
83 G4double u;
84
85 do {
86 u = - std::log(G4UniformRand()*G4UniformRand());
87
88 if ( border > G4UniformRand() ) { u /= a1; }
89 else { u /= a2; }
90
91 } while(u > uMax);
92
93 G4double cost = 1.0 - 2*u*u/(uMax*uMax);
94 G4double sint = std::sqrt((1 - cost)*(1 + cost));
95 G4double phi = CLHEP::twopi*G4UniformRand();
96
97 fLocalDirection.set(sint*std::cos(phi), sint*std::sin(phi), cost);
99
100 return fLocalDirection;
101}
double G4double
Definition: G4Types.hh:64
#define G4UniformRand()
Definition: Randomize.hh:53
void set(double x, double y, double z)
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const

The documentation for this class was generated from the following files: