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

#include <G4SauterGavrilaAngularDistribution.hh>

Inheritance diagram for G4SauterGavrilaAngularDistribution:
G4VEmAngularDistribution

Public Member Functions

 G4SauterGavrilaAngularDistribution ()
 
 G4SauterGavrilaAngularDistribution (const G4SauterGavrilaAngularDistribution &)=delete
 
const G4StringGetName () const
 
G4SauterGavrilaAngularDistributionoperator= (const G4SauterGavrilaAngularDistribution &right)=delete
 
void PrintGeneratorInformation () const override
 
G4ThreeVectorSampleDirection (const G4DynamicParticle *dp, G4double e=0.0, G4int shellId=0, const G4Material *mat=0) final
 
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)
 
 ~G4SauterGavrilaAngularDistribution () override
 

Protected Attributes

G4ThreeVector fLocalDirection
 
G4bool fPolarisation
 

Private Attributes

G4String fName
 

Detailed Description

Definition at line 52 of file G4SauterGavrilaAngularDistribution.hh.

Constructor & Destructor Documentation

◆ G4SauterGavrilaAngularDistribution() [1/2]

G4SauterGavrilaAngularDistribution::G4SauterGavrilaAngularDistribution ( )
explicit

Definition at line 48 of file G4SauterGavrilaAngularDistribution.cc.

49 : G4VEmAngularDistribution("SauterGavrila")
50{}
G4VEmAngularDistribution(const G4String &name)

◆ ~G4SauterGavrilaAngularDistribution()

G4SauterGavrilaAngularDistribution::~G4SauterGavrilaAngularDistribution ( )
override

Definition at line 52 of file G4SauterGavrilaAngularDistribution.cc.

53{}

◆ G4SauterGavrilaAngularDistribution() [2/2]

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

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

◆ PrintGeneratorInformation()

void G4SauterGavrilaAngularDistribution::PrintGeneratorInformation ( ) const
overridevirtual

Reimplemented from G4VEmAngularDistribution.

Definition at line 110 of file G4SauterGavrilaAngularDistribution.cc.

111{
112 G4cout << "\n" << G4endl;
113 G4cout << "Non-polarized photoelectric effect angular generator." << G4endl;
114 G4cout << "The Sauter-Gavrila distribution for the K-shell is used."<<G4endl;
115 G4cout << "Originally developed by M.Maire for Geant3"
116 << G4endl;
117}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout

References G4cout, and G4endl.

◆ SampleDirection()

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

Implements G4VEmAngularDistribution.

Definition at line 55 of file G4SauterGavrilaAngularDistribution.cc.

57{
58 static const G4double emin = 1*CLHEP::eV;
59 static const G4double emax = 100*CLHEP::MeV;
60
62 if (energy > emax) {
64 } else {
65 // Initial algorithm according Penelope 2008 manual and
66 // F.Sauter Ann. Physik 9, 217(1931); 11, 454(1931).
67 // Modified according Penelope 2014 manual
68 G4double costheta = 0.0;
69
70 // 1) initialize energy-dependent variables
71 // Variable naming according to Eq. (2.24) of Penelope Manual
72 // (pag. 44)
74 G4double gamma = 1.0 + tau;
75 G4double beta = std::sqrt(tau*(tau + 2.0))/gamma;
76
77 // ac corresponds to "A" of Eq. (2.31)
78 //
79 G4double ac = (1.0 - beta)/beta;
80 G4double a1 = 0.5*beta*gamma*tau*(gamma-2.0);
81 G4double a2 = ac + 2.0;
82 // gtmax = maximum of the rejection function according to Eq. (2.28),
83 // obtained for tsam=0
84 G4double gtmax = 2.0*(a1 + 1.0/ac);
85
86 G4double tsam = 0.0;
87 G4double gtr = 0.0;
88
89 //2) sampling. Eq. (2.31) of Penelope Manual
90 // tsam = 1-std::cos(theta)
91 // gtr = rejection function according to Eq. (2.28)
92 do{
93 G4double rand = G4UniformRand();
94 tsam = 2.0*ac * (2.0*rand + a2*std::sqrt(rand)) / (a2*a2 - 4.0*rand);
95 gtr = (2.0 - tsam) * (a1 + 1.0/(ac+tsam));
96 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
97 } while(G4UniformRand()*gtmax > gtr);
98
99 costheta = 1.0 - tsam;
100
101 G4double sint = std::sqrt(tsam*(2.0 - tsam));
103
104 fLocalDirection.set(sint*std::cos(phi), sint*std::sin(phi), costheta);
106 }
107 return fLocalDirection;
108}
static const G4double emax
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
G4double GetKineticEnergy() const
static constexpr double twopi
Definition: SystemOfUnits.h:56
static constexpr double MeV
static constexpr double eV
G4double energy(const ThreeVector &p, const G4double m)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
float electron_mass_c2
Definition: hepunit.py:273

References anonymous_namespace{G4PionRadiativeDecayChannel.cc}::beta, source.hepunit::electron_mass_c2, emax, G4INCL::KinematicsUtils::energy(), CLHEP::eV, G4VEmAngularDistribution::fLocalDirection, G4UniformRand, G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMomentumDirection(), G4INCL::Math::max(), CLHEP::MeV, CLHEP::Hep3Vector::rotateUz(), 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 G4VEmAngularDistribution::SamplePairDirections ( const G4DynamicParticle dp,
G4double  elecKinEnergy,
G4double  posiKinEnergy,
G4ThreeVector dirElectron,
G4ThreeVector dirPositron,
G4int  Z = 0,
const G4Material mat = nullptr 
)
virtualinherited

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: