Geant4 9.6.0
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 ()
 
virtual G4ThreeVectorSampleDirection (const G4DynamicParticle *dp, G4double e=0.0, G4int shellId=0, const G4Material *mat=0)
 
void PrintGeneratorInformation () const
 
- Public Member Functions inherited from G4VEmAngularDistribution
 G4VEmAngularDistribution (const G4String &name)
 
virtual ~G4VEmAngularDistribution ()
 
virtual G4ThreeVectorSampleDirection (const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
 
const G4StringGetName () const
 

Additional Inherited Members

- Protected Attributes inherited from G4VEmAngularDistribution
G4ThreeVector fLocalDirection
 

Detailed Description

Constructor & Destructor Documentation

◆ G4PhotoElectricAngularGeneratorSauterGavrila()

G4PhotoElectricAngularGeneratorSauterGavrila::G4PhotoElectricAngularGeneratorSauterGavrila ( )

Definition at line 54 of file G4PhotoElectricAngularGeneratorSauterGavrila.cc.

54 :
55 G4VEmAngularDistribution("AngularGenSauterGavrilaLowE")
56{}

◆ ~G4PhotoElectricAngularGeneratorSauterGavrila()

G4PhotoElectricAngularGeneratorSauterGavrila::~G4PhotoElectricAngularGeneratorSauterGavrila ( )

Definition at line 58 of file G4PhotoElectricAngularGeneratorSauterGavrila.cc.

59{}

Member Function Documentation

◆ PrintGeneratorInformation()

void G4PhotoElectricAngularGeneratorSauterGavrila::PrintGeneratorInformation ( ) const

Definition at line 104 of file G4PhotoElectricAngularGeneratorSauterGavrila.cc.

105{
106 G4cout << "\n" << G4endl;
107 G4cout << "" << G4endl;
108 G4cout << "Re-implementation of the photolectric angular distribution" << G4endl;
109 G4cout << "developed my M. Maire for the Standard EM Physics G4PhotoElectricEffect" << G4endl;
110 G4cout << "It computes the theta distribution of the emitted electron, with respect to the" << G4endl;
111 G4cout << "incident Gamma, using the Sauter-Gavrila distribution for the K-shell\n" << G4endl;
112}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout

◆ SampleDirection()

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

Implements G4VEmAngularDistribution.

Definition at line 62 of file G4PhotoElectricAngularGeneratorSauterGavrila.cc.

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

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