#include <G4ModifiedMephi.hh>
|
| G4ModifiedMephi (const G4String &name="") |
|
| ~G4ModifiedMephi () override |
|
G4ThreeVector & | SampleDirection (const G4DynamicParticle *dp, G4double gEnergy, 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 override |
|
G4ModifiedMephi & | operator= (const G4ModifiedMephi &right)=delete |
|
| G4ModifiedMephi (const G4ModifiedMephi &)=delete |
|
| G4VEmAngularDistribution (const G4String &name) |
|
virtual | ~G4VEmAngularDistribution () |
|
virtual G4ThreeVector & | SampleDirectionForShell (const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, G4int shellID, const G4Material *) |
|
const G4String & | GetName () const |
|
G4VEmAngularDistribution & | operator= (const G4VEmAngularDistribution &right)=delete |
|
| G4VEmAngularDistribution (const G4VEmAngularDistribution &)=delete |
|
Definition at line 55 of file G4ModifiedMephi.hh.
◆ G4ModifiedMephi() [1/2]
G4ModifiedMephi::G4ModifiedMephi |
( |
const G4String & | name = "" | ) |
|
|
explicit |
Definition at line 51 of file G4ModifiedMephi.cc.
53{}
G4VEmAngularDistribution(const G4String &name)
◆ ~G4ModifiedMephi()
G4ModifiedMephi::~G4ModifiedMephi |
( |
| ) |
|
|
override |
◆ G4ModifiedMephi() [2/2]
◆ operator=()
◆ PrintGeneratorInformation()
void G4ModifiedMephi::PrintGeneratorInformation |
( |
| ) |
const |
|
overridevirtual |
◆ SampleDirection()
Implements G4VEmAngularDistribution.
Definition at line 63 of file G4ModifiedMephi.cc.
66{
67
70 G4double sint = std::sqrt((1.0 - cost)*(1.0 + cost));
72
75
77}
void set(double x, double y, double z)
Hep3Vector & rotateUz(const Hep3Vector &)
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double GetPDGMass() const
G4ThreeVector fLocalDirection
◆ SamplePairDirections()
Reimplemented from G4VEmAngularDistribution.
Definition at line 95 of file G4ModifiedMephi.cc.
101{
105
107 G4double etwo = elecKinEnergy + posiKinEnergy;
109
110 G4double cost = SampleCosTheta(ekin, etwo, mass);
111 G4double sint = std::sqrt((1.0 - cost)*(1.0 + cost));
112
113 dirElectron.
set(sint*cosp, sint*sinp, cost);
115
116 cost = SampleCosTheta(ekin, etwo, mass);
117 sint = std::sqrt((1.0 - cost)*(1.0 + cost));
118
119 dirPositron.
set(-sint*cosp, -sint*sinp, cost);
121}
The documentation for this class was generated from the following files: