Geant4-11
Public Member Functions | Protected Attributes | Private Member Functions | Private Attributes
G4ModifiedTsai Class Reference

#include <G4ModifiedTsai.hh>

Inheritance diagram for G4ModifiedTsai:
G4VEmAngularDistribution

Public Member Functions

 G4ModifiedTsai (const G4ModifiedTsai &)=delete
 
 G4ModifiedTsai (const G4String &name="")
 
const G4StringGetName () const
 
G4ModifiedTsaioperator= (const G4ModifiedTsai &right)=delete
 
void PrintGeneratorInformation () const final
 
G4ThreeVectorSampleDirection (const G4DynamicParticle *dp, G4double out_energy, G4int Z, const G4Material *mat=nullptr) final
 
virtual G4ThreeVectorSampleDirectionForShell (const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, G4int shellID, const G4Material *)
 
void SamplePairDirections (const G4DynamicParticle *dp, G4double elecKinEnergy, G4double posiKinEnergy, G4ThreeVector &dirElectron, G4ThreeVector &dirPositron, G4int Z=0, const G4Material *mat=nullptr) final
 
 ~G4ModifiedTsai () override
 

Protected Attributes

G4ThreeVector fLocalDirection
 
G4bool fPolarisation
 

Private Member Functions

G4double SampleCosTheta (G4double kinEnergy)
 

Private Attributes

G4String fName
 

Detailed Description

Definition at line 62 of file G4ModifiedTsai.hh.

Constructor & Destructor Documentation

◆ G4ModifiedTsai() [1/2]

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

Definition at line 65 of file G4ModifiedTsai.cc.

66 : G4VEmAngularDistribution("ModifiedTsai")
67{}
G4VEmAngularDistribution(const G4String &name)

◆ ~G4ModifiedTsai()

G4ModifiedTsai::~G4ModifiedTsai ( )
override

Definition at line 71 of file G4ModifiedTsai.cc.

72{}

◆ G4ModifiedTsai() [2/2]

G4ModifiedTsai::G4ModifiedTsai ( const G4ModifiedTsai )
delete

Member Function Documentation

◆ GetName()

const G4String & G4VEmAngularDistribution::GetName ( ) const
inlineinherited

Definition at line 111 of file G4VEmAngularDistribution.hh.

112{
113 return fName;
114}

References G4VEmAngularDistribution::fName.

◆ operator=()

G4ModifiedTsai & G4ModifiedTsai::operator= ( const G4ModifiedTsai right)
delete

◆ PrintGeneratorInformation()

void G4ModifiedTsai::PrintGeneratorInformation ( ) const
finalvirtual

Reimplemented from G4VEmAngularDistribution.

Definition at line 144 of file G4ModifiedTsai.cc.

145{
146 G4cout << "\n" << G4endl;
147 G4cout << "Bremsstrahlung Angular Generator is Modified Tsai" << G4endl;
148 G4cout << "Distribution suggested by L.Urban (Geant3 manual (1993) Phys211)"
149 << G4endl;
150 G4cout << "Derived from Tsai distribution (Rev Mod Phys 49,421(1977)) \n"
151 << G4endl;
152}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout

References G4cout, and G4endl.

◆ SampleCosTheta()

G4double G4ModifiedTsai::SampleCosTheta ( G4double  kinEnergy)
private

Definition at line 93 of file G4ModifiedTsai.cc.

94{
95 // Universal distribution suggested by L. Urban (Geant3 manual (1993)
96 // Phys211) derived from Tsai distribution (Rev Mod Phys 49,421(1977))
97
98 G4double uMax = 2*(1. + kinEnergy/CLHEP::electron_mass_c2);
99
100 static const G4double a1 = 1.6;
101 static const G4double a2 = a1/3.;
102 static const G4double border = 0.25;
103 G4double u;
104 CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
105
106 do {
107 G4double uu = -G4Log(rndmEngine->flat()*rndmEngine->flat());
108 u = (border > rndmEngine->flat()) ? uu*a1 : uu*a2;
109
110 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
111 } while(u > uMax);
112
113 return 1.0 - 2.0*u*u/(uMax*uMax);
114}
G4double G4Log(G4double x)
Definition: G4Log.hh:226
double G4double
Definition: G4Types.hh:83
virtual double flat()=0
static constexpr double electron_mass_c2

References CLHEP::electron_mass_c2, CLHEP::HepRandomEngine::flat(), and G4Log().

Referenced by SampleDirection(), and SamplePairDirections().

◆ SampleDirection()

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

Implements G4VEmAngularDistribution.

Definition at line 77 of file G4ModifiedTsai.cc.

79{
80 // Sample gamma angle (Z - axis along the parent particle).
82 G4double sint = std::sqrt((1 - cost)*(1 + cost));
84
85 fLocalDirection.set(sint*std::cos(phi), sint*std::sin(phi), cost);
87
88 return fLocalDirection;
89}
#define G4UniformRand()
Definition: Randomize.hh:52
void set(double x, double y, double z)
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4double SampleCosTheta(G4double kinEnergy)
static constexpr double twopi
Definition: SystemOfUnits.h:56

References G4VEmAngularDistribution::fLocalDirection, G4UniformRand, G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMomentumDirection(), CLHEP::Hep3Vector::rotateUz(), SampleCosTheta(), CLHEP::Hep3Vector::set(), and CLHEP::twopi.

◆ SampleDirectionForShell()

G4ThreeVector & G4VEmAngularDistribution::SampleDirectionForShell ( const G4DynamicParticle dp,
G4double  finalTotalEnergy,
G4int  Z,
G4int  shellID,
const G4Material mat 
)
virtualinherited

◆ SamplePairDirections()

void G4ModifiedTsai::SamplePairDirections ( const G4DynamicParticle dp,
G4double  elecKinEnergy,
G4double  posiKinEnergy,
G4ThreeVector dirElectron,
G4ThreeVector dirPositron,
G4int  Z = 0,
const G4Material mat = nullptr 
)
finalvirtual

Reimplemented from G4VEmAngularDistribution.

Definition at line 118 of file G4ModifiedTsai.cc.

124{
126 G4double sinp = std::sin(phi);
127 G4double cosp = std::cos(phi);
128
129 G4double cost = SampleCosTheta(elecKinEnergy);
130 G4double sint = std::sqrt((1. - cost)*(1. + cost));
131
132 dirElectron.set(sint*cosp, sint*sinp, cost);
133 dirElectron.rotateUz(dp->GetMomentumDirection());
134
135 cost = SampleCosTheta(posiKinEnergy);
136 sint = std::sqrt((1. - cost)*(1. + cost));
137
138 dirPositron.set(-sint*cosp, -sint*sinp, cost);
139 dirPositron.rotateUz(dp->GetMomentumDirection());
140}

References G4UniformRand, G4DynamicParticle::GetMomentumDirection(), CLHEP::Hep3Vector::rotateUz(), SampleCosTheta(), CLHEP::Hep3Vector::set(), and CLHEP::twopi.

Field Documentation

◆ fLocalDirection

G4ThreeVector G4VEmAngularDistribution::fLocalDirection
protectedinherited

◆ fName

G4String G4VEmAngularDistribution::fName
privateinherited

Definition at line 108 of file G4VEmAngularDistribution.hh.

Referenced by G4VEmAngularDistribution::GetName().

◆ fPolarisation

G4bool G4VEmAngularDistribution::fPolarisation
protectedinherited

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