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

#include <G4ModifiedMephi.hh>

+ Inheritance diagram for G4ModifiedMephi:

Public Member Functions

 G4ModifiedMephi (const G4String &name="")
 
virtual ~G4ModifiedMephi ()
 
virtual G4ThreeVectorSampleDirection (const G4DynamicParticle *dp, G4double gEnergy, 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
 
virtual void PrintGeneratorInformation () const final
 
G4ModifiedMephioperator= (const G4ModifiedMephi &right)=delete
 
 G4ModifiedMephi (const G4ModifiedMephi &)=delete
 
- 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 55 of file G4ModifiedMephi.hh.

Constructor & Destructor Documentation

◆ G4ModifiedMephi() [1/2]

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

Definition at line 51 of file G4ModifiedMephi.cc.

◆ ~G4ModifiedMephi()

G4ModifiedMephi::~G4ModifiedMephi ( )
virtual

Definition at line 57 of file G4ModifiedMephi.cc.

58{}

◆ G4ModifiedMephi() [2/2]

G4ModifiedMephi::G4ModifiedMephi ( const G4ModifiedMephi )
delete

Member Function Documentation

◆ operator=()

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

◆ PrintGeneratorInformation()

void G4ModifiedMephi::PrintGeneratorInformation ( ) const
finalvirtual

Definition at line 125 of file G4ModifiedMephi.cc.

126{
127 G4cout << "\n" << G4endl;
128 G4cout << "Angular Generator is Modified Mephi" << G4endl;
129}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout

◆ SampleDirection()

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

Implements G4VEmAngularDistribution.

Definition at line 63 of file G4ModifiedMephi.cc.

66{
67 // Sample gamma angle (Z - axis along the parent particle).
68 G4double cost = SampleCosTheta(dp->GetKineticEnergy(), gEnergy,
69 dp->GetDefinition()->GetPDGMass());
70 G4double sint = std::sqrt((1.0 - cost)*(1.0 + cost));
71 G4double phi = CLHEP::twopi*G4UniformRand();
72
73 fLocalDirection.set(sint*std::cos(phi), sint*std::sin(phi), cost);
75
76 return fLocalDirection;
77}
double G4double
Definition: G4Types.hh:83
#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
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const

◆ SamplePairDirections()

void G4ModifiedMephi::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 95 of file G4ModifiedMephi.cc.

101{
102 G4double phi = CLHEP::twopi * G4UniformRand();
103 G4double sinp = std::sin(phi);
104 G4double cosp = std::cos(phi);
105
106 G4double ekin = dp->GetKineticEnergy();
107 G4double etwo = elecKinEnergy + posiKinEnergy;
108 G4double mass = dp->GetDefinition()->GetPDGMass();
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);
114 dirElectron.rotateUz(dp->GetMomentumDirection());
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);
120 dirPositron.rotateUz(dp->GetMomentumDirection());
121}

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