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

#include <G4ModifiedMephi.hh>

Inheritance diagram for G4ModifiedMephi:
G4VEmAngularDistribution

Public Member Functions

 G4ModifiedMephi (const G4ModifiedMephi &)=delete
 
 G4ModifiedMephi (const G4String &name="")
 
const G4StringGetName () const
 
G4ModifiedMephioperator= (const G4ModifiedMephi &right)=delete
 
void PrintGeneratorInformation () const override
 
G4ThreeVectorSampleDirection (const G4DynamicParticle *dp, G4double gEnergy, 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
 
 ~G4ModifiedMephi () override
 

Protected Attributes

G4ThreeVector fLocalDirection
 
G4bool fPolarisation
 

Private Member Functions

G4double SampleCosTheta (G4double primKinEnergy, G4double gEnergy, G4double mass)
 

Private Attributes

G4String fName
 

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.

52 : G4VEmAngularDistribution("ModifiedMephi")
53{}
G4VEmAngularDistribution(const G4String &name)

◆ ~G4ModifiedMephi()

G4ModifiedMephi::~G4ModifiedMephi ( )
override

Definition at line 57 of file G4ModifiedMephi.cc.

58{}

◆ G4ModifiedMephi() [2/2]

G4ModifiedMephi::G4ModifiedMephi ( const G4ModifiedMephi )
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=()

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

◆ PrintGeneratorInformation()

void G4ModifiedMephi::PrintGeneratorInformation ( ) const
overridevirtual

Reimplemented from G4VEmAngularDistribution.

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

References G4cout, and G4endl.

◆ SampleCosTheta()

G4double G4ModifiedMephi::SampleCosTheta ( G4double  primKinEnergy,
G4double  gEnergy,
G4double  mass 
)
private

Definition at line 81 of file G4ModifiedMephi.cc.

84{
85 G4double gam = 1.0 + primKinEnergy/mass;
86 G4double rmax = gam*CLHEP::halfpi*std::min(1.0, gam*mass/gEnergy - 1.0);
87 G4double rmax2= rmax*rmax;
88 G4double x = G4UniformRand()*rmax2/(1.0 + rmax2);
89
90 return std::cos(std::sqrt(x/(1.0 - x))/gam);
91}
double G4double
Definition: G4Types.hh:83
#define G4UniformRand()
Definition: Randomize.hh:52
static constexpr double halfpi
Definition: SystemOfUnits.h:57
T min(const T t1, const T t2)
brief Return the smallest of the two arguments

References G4UniformRand, G4InuclParticleNames::gam, CLHEP::halfpi, and G4INCL::Math::min().

Referenced by SampleDirection(), and SamplePairDirections().

◆ 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));
72
73 fLocalDirection.set(sint*std::cos(phi), sint*std::sin(phi), cost);
75
76 return fLocalDirection;
77}
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
G4double SampleCosTheta(G4double primKinEnergy, G4double gEnergy, G4double mass)
static constexpr double twopi
Definition: SystemOfUnits.h:56

References G4VEmAngularDistribution::fLocalDirection, G4UniformRand, G4DynamicParticle::GetDefinition(), G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMomentumDirection(), G4ParticleDefinition::GetPDGMass(), 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 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{
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}

References G4UniformRand, G4DynamicParticle::GetDefinition(), G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMomentumDirection(), G4ParticleDefinition::GetPDGMass(), 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: