#include <G4DipBustGenerator.hh>
|
| G4DipBustGenerator (const G4String &name="") |
|
virtual | ~G4DipBustGenerator () |
|
virtual G4ThreeVector & | SampleDirection (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 |
|
| 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) |
|
const G4String & | GetName () const |
|
G4VEmAngularDistribution & | operator= (const G4VEmAngularDistribution &right)=delete |
|
| G4VEmAngularDistribution (const G4VEmAngularDistribution &)=delete |
|
Definition at line 54 of file G4DipBustGenerator.hh.
◆ G4DipBustGenerator()
G4DipBustGenerator::G4DipBustGenerator |
( |
const G4String & |
name = "" | ) |
|
|
explicit |
◆ ~G4DipBustGenerator()
G4DipBustGenerator::~G4DipBustGenerator |
( |
| ) |
|
|
virtual |
◆ PolarAngle()
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}
◆ PrintGeneratorInformation()
void G4DipBustGenerator::PrintGeneratorInformation |
( |
| ) |
const |
|
finalvirtual |
Definition at line 144 of file G4DipBustGenerator.cc.
145{
147 G4cout <<
"Angular Generator based on classical formula from" <<
G4endl;
148 G4cout <<
"J.D. Jackson, Classical Electrodynamics, Wiley, New York 1975"
150}
G4GLOB_DLL std::ostream G4cout
◆ SampleDirection()
Implements G4VEmAngularDistribution.
Definition at line 90 of file G4DipBustGenerator.cc.
92{
94
95 const G4double sinTheta = std::sqrt((1. - cosTheta)*(1. + cosTheta));
97
100
102}
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 118 of file G4DipBustGenerator.cc.
124{
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);
134
135 cost = SampleCosTheta(posiKinEnergy);
136 sint = std::sqrt((1. - cost)*(1. + cost));
137
138 dirPositron.
set(-sint*cosp, -sint*sinp, cost);
140}
The documentation for this class was generated from the following files: