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

#include <G4SauterGavrilaAngularDistribution.hh>

+ Inheritance diagram for G4SauterGavrilaAngularDistribution:

Public Member Functions

 G4SauterGavrilaAngularDistribution ()
 
virtual ~G4SauterGavrilaAngularDistribution ()
 
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

Definition at line 55 of file G4SauterGavrilaAngularDistribution.hh.

Constructor & Destructor Documentation

◆ G4SauterGavrilaAngularDistribution()

G4SauterGavrilaAngularDistribution::G4SauterGavrilaAngularDistribution ( )

Definition at line 48 of file G4SauterGavrilaAngularDistribution.cc.

49 : G4VEmAngularDistribution("AngularGenSauterGavrila")
50{}

◆ ~G4SauterGavrilaAngularDistribution()

G4SauterGavrilaAngularDistribution::~G4SauterGavrilaAngularDistribution ( )
virtual

Definition at line 52 of file G4SauterGavrilaAngularDistribution.cc.

53{}

Member Function Documentation

◆ PrintGeneratorInformation()

void G4SauterGavrilaAngularDistribution::PrintGeneratorInformation ( ) const

Definition at line 95 of file G4SauterGavrilaAngularDistribution.cc.

96{
97 G4cout << "\n" << G4endl;
98 G4cout << "Non-polarized photoelectric effect angular generator." << G4endl;
99 G4cout << "The Sauter-Gavrila distribution for the K-shell is used."<<G4endl;
100 G4cout << "Originally developed by M.Maire for Geant3"
101 << G4endl;
102}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout

◆ SampleDirection()

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

Implements G4VEmAngularDistribution.

Definition at line 56 of file G4SauterGavrilaAngularDistribution.cc.

58{
59 G4double tau = dp->GetKineticEnergy()/electron_mass_c2;
60 const G4double taulimit = 30.0;
61
62 if (tau > taulimit) {
64 // Bugzilla 1120
65 // SI on 05/09/2010 as suggested by JG 04/09/10
66 } else {
67
68 G4double invgamma = 1.0/(tau + 1.0);
69 G4double beta = std::sqrt(tau*(tau + 2.0))*invgamma;
70 G4double b = 0.5*tau*(tau*tau - 1.0);
71 G4double invgamma2 = invgamma*invgamma;
72
73 G4double rndm,term,greject,grejsup,costeta,sint2;
74 if (tau < 1.) { grejsup = (1.+b-beta*b)/invgamma2; }
75 else { grejsup = (1.+b+beta*b)/invgamma2; }
76
77 do {
78 rndm = 1 - 2*G4UniformRand();
79 costeta = (rndm + beta)/(rndm*beta + 1);
80 term = invgamma2/(1 + beta*rndm);
81 sint2 = (1 - costeta)*(1 + costeta);
82 greject = sint2*(1 + b*term)/(term*term);
83
84 } while(greject < G4UniformRand()*grejsup);
85
86 G4double sinteta = std::sqrt(sint2);
87 G4double phi = CLHEP::twopi*G4UniformRand();
88
89 fLocalDirection.set(sinteta*std::cos(phi), sinteta*std::sin(phi), costeta);
91 }
92 return fLocalDirection;
93}
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: