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

#include <G4PhotoElectricAngularGeneratorSauterGavrila.hh>

+ Inheritance diagram for G4PhotoElectricAngularGeneratorSauterGavrila:

Public Member Functions

 G4PhotoElectricAngularGeneratorSauterGavrila ()
 
 ~G4PhotoElectricAngularGeneratorSauterGavrila ()
 
G4ThreeVectorSampleDirection (const G4DynamicParticle *dp, G4double e=0.0, G4int shellId=0, const G4Material *mat=nullptr) override
 
void PrintGeneratorInformation () const override
 
G4PhotoElectricAngularGeneratorSauterGavrilaoperator= (const G4PhotoElectricAngularGeneratorSauterGavrila &right)=delete
 
 G4PhotoElectricAngularGeneratorSauterGavrila (const G4PhotoElectricAngularGeneratorSauterGavrila &)=delete
 
- Public Member Functions inherited from G4VEmAngularDistribution
 G4VEmAngularDistribution (const G4String &name)
 
virtual ~G4VEmAngularDistribution ()
 
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

Constructor & Destructor Documentation

◆ G4PhotoElectricAngularGeneratorSauterGavrila() [1/2]

G4PhotoElectricAngularGeneratorSauterGavrila::G4PhotoElectricAngularGeneratorSauterGavrila ( )
explicit

Definition at line 55 of file G4PhotoElectricAngularGeneratorSauterGavrila.cc.

55 :
56 G4VEmAngularDistribution("AngularGenSauterGavrilaLowE")
57{}
G4VEmAngularDistribution(const G4String &name)

◆ ~G4PhotoElectricAngularGeneratorSauterGavrila()

G4PhotoElectricAngularGeneratorSauterGavrila::~G4PhotoElectricAngularGeneratorSauterGavrila ( )

Definition at line 61 of file G4PhotoElectricAngularGeneratorSauterGavrila.cc.

62{}

◆ G4PhotoElectricAngularGeneratorSauterGavrila() [2/2]

G4PhotoElectricAngularGeneratorSauterGavrila::G4PhotoElectricAngularGeneratorSauterGavrila ( const G4PhotoElectricAngularGeneratorSauterGavrila & )
delete

Member Function Documentation

◆ operator=()

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

◆ PrintGeneratorInformation()

void G4PhotoElectricAngularGeneratorSauterGavrila::PrintGeneratorInformation ( ) const
overridevirtual

Reimplemented from G4VEmAngularDistribution.

Definition at line 110 of file G4PhotoElectricAngularGeneratorSauterGavrila.cc.

111{
112 G4cout << "\n" << G4endl;
113 G4cout << "" << G4endl;
114 G4cout << "Re-implementation of the photolectric angular distribution" << G4endl;
115 G4cout << "developed my M. Maire for the Standard EM Physics G4PhotoElectricEffect" << G4endl;
116 G4cout << "It computes the theta distribution of the emitted electron, with respect to the" << G4endl;
117 G4cout << "incident Gamma, using the Sauter-Gavrila distribution for the K-shell\n" << G4endl;
118}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout

◆ SampleDirection()

G4ThreeVector & G4PhotoElectricAngularGeneratorSauterGavrila::SampleDirection ( const G4DynamicParticle * dp,
G4double e = 0.0,
G4int shellId = 0,
const G4Material * mat = nullptr )
overridevirtual

Implements G4VEmAngularDistribution.

Definition at line 67 of file G4PhotoElectricAngularGeneratorSauterGavrila.cc.

70{
71
72 // Compute Theta distribution of the emitted electron, with respect to the
73 // incident Gamma.
74 // The Sauter-Gavrila distribution for the K-shell is used.
75 G4double costeta = 1.;
76 G4double Phi = twopi * G4UniformRand();
77 G4double cosphi = std::cos(Phi);
78 G4double sinphi = std::sin(Phi);
79 G4double sinteta = 0;
80 G4double gamma = 1. + dp->GetKineticEnergy()/electron_mass_c2;
81
82 if (gamma > 5.) {
84 return fLocalDirection;
85 // Bugzilla 1120
86 // SI on 05/09/2010 as suggested by JG 04/09/10
87 }
88
89 G4double beta = std::sqrt((gamma - 1)*(gamma + 1))/gamma;
90 G4double b = 0.5*gamma*(gamma - 1)*(gamma - 2);
91
92 G4double rndm,term,greject,grejsup;
93 if (gamma < 2.) grejsup = gamma*gamma*(1.+b-beta*b);
94 else grejsup = gamma*gamma*(1.+b+beta*b);
95
96 do { rndm = 1.-2*G4UniformRand();
97 costeta = (rndm+beta)/(rndm*beta+1.);
98 term = 1.-beta*costeta;
99 greject = (1.-costeta*costeta)*(1.+b*term)/(term*term);
100 } while(greject < G4UniformRand()*grejsup);
101
102 sinteta = std::sqrt((1 - costeta)*(1 + costeta));
103 fLocalDirection.set(sinteta*cosphi, sinteta*sinphi, costeta);
105 return fLocalDirection;
106}
double G4double
Definition G4Types.hh:83
#define G4UniformRand()
Definition Randomize.hh:52
void set(double x, double y, double z)
Hep3Vector & rotateUz(const Hep3Vector &)
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const

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