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

#include <G4DipBustGenerator.hh>

+ Inheritance diagram for G4DipBustGenerator:

Public Member Functions

 G4DipBustGenerator (const G4String &name="")
 
virtual ~G4DipBustGenerator ()
 
virtual G4ThreeVectorSampleDirection (const G4DynamicParticle *dp, G4double out_energy, G4int Z, const G4Material *mat=nullptr) final
 
virtual void SamplePairDirections (const G4DynamicParticle *dp, G4double elecKinEnergy, G4double posiKinEnergy, G4ThreeVector &dirElectron, G4ThreeVector &dirPositron, G4int Z=0, const G4Material *mat=nullptr) final
 
G4double PolarAngle (G4double initial_energy, G4double final_energy, G4int Z)
 
virtual void PrintGeneratorInformation () const final
 
- Public Member Functions inherited from G4VEmAngularDistribution
 G4VEmAngularDistribution (const G4String &name)
 
virtual ~G4VEmAngularDistribution ()
 
virtual G4ThreeVectorSampleDirection (const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
 
virtual G4ThreeVectorSampleDirectionForShell (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)
 
const G4StringGetName () const
 
G4VEmAngularDistributionoperator= (const G4VEmAngularDistribution &right)=delete
 
 G4VEmAngularDistribution (const G4VEmAngularDistribution &)=delete
 

Additional Inherited Members

- Protected Attributes inherited from G4VEmAngularDistribution
G4ThreeVector fLocalDirection
 
G4bool fPolarisation
 

Detailed Description

Definition at line 54 of file G4DipBustGenerator.hh.

Constructor & Destructor Documentation

◆ G4DipBustGenerator()

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

Definition at line 62 of file G4DipBustGenerator.cc.

◆ ~G4DipBustGenerator()

G4DipBustGenerator::~G4DipBustGenerator ( )
virtual

Definition at line 68 of file G4DipBustGenerator.cc.

69{}

Member Function Documentation

◆ PolarAngle()

G4double G4DipBustGenerator::PolarAngle ( G4double  initial_energy,
G4double  final_energy,
G4int  Z 
)

Definition at line 106 of file G4DipBustGenerator.cc.

109{
110 const G4double cosTheta = SampleCosTheta(eTkin);
111 G4double theta = std::acos(cosTheta);
112 theta = std::min(std::max(theta, 0.), CLHEP::pi);
113 return theta;
114}
double G4double
Definition: G4Types.hh:83

◆ PrintGeneratorInformation()

void G4DipBustGenerator::PrintGeneratorInformation ( ) const
finalvirtual

Definition at line 144 of file G4DipBustGenerator.cc.

145{
146 G4cout << "\n" << G4endl;
147 G4cout << "Angular Generator based on classical formula from" << G4endl;
148 G4cout << "J.D. Jackson, Classical Electrodynamics, Wiley, New York 1975"
149 << G4endl;
150}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout

◆ SampleDirection()

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

Implements G4VEmAngularDistribution.

Definition at line 90 of file G4DipBustGenerator.cc.

92{
93 const G4double cosTheta = SampleCosTheta(dp->GetKineticEnergy());
94
95 const G4double sinTheta = std::sqrt((1. - cosTheta)*(1. + cosTheta));
96 const G4double phi = CLHEP::twopi*G4UniformRand();
97
98 fLocalDirection.set(sinTheta*std::cos(phi), sinTheta*std::sin(phi),cosTheta);
100
101 return fLocalDirection;
102}
#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

◆ SamplePairDirections()

void G4DipBustGenerator::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 G4DipBustGenerator.cc.

124{
125 const G4double phi = CLHEP::twopi * G4UniformRand();
126 const G4double sinp = std::sin(phi);
127 const 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}

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