Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4SauterGavrilaAngularDistribution Class Reference

#include <G4SauterGavrilaAngularDistribution.hh>

Inheritance diagram for G4SauterGavrilaAngularDistribution:
Collaboration 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) final
 
virtual void PrintGeneratorInformation () const final
 
- 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 *)
 
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 ( )
explicit

Definition at line 48 of file G4SauterGavrilaAngularDistribution.cc.

49  : G4VEmAngularDistribution("AngularGenSauterGavrila")
50 {}
G4VEmAngularDistribution(const G4String &name)
G4SauterGavrilaAngularDistribution::~G4SauterGavrilaAngularDistribution ( )
virtual

Definition at line 52 of file G4SauterGavrilaAngularDistribution.cc.

53 {}

Member Function Documentation

void G4SauterGavrilaAngularDistribution::PrintGeneratorInformation ( ) const
finalvirtual

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 }
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
G4ThreeVector & G4SauterGavrilaAngularDistribution::SampleDirection ( const G4DynamicParticle dp,
G4double  e = 0.0,
G4int  shellId = 0,
const G4Material mat = 0 
)
finalvirtual

Implements G4VEmAngularDistribution.

Definition at line 56 of file G4SauterGavrilaAngularDistribution.cc.

58 {
60  static const G4double taulimit = 50.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  // algorithm according Penelope 2008 manual and
68  // F.Sauter Ann. Physik 9, 217(1931); 11, 454(1931).
69 
70  G4double gamma = tau + 1;
71  G4double beta = std::sqrt(tau*(tau + 2))/gamma;
72  G4double A = (1 - beta)/beta;
73  G4double Ap2 = A + 2;
74  G4double B = 0.5*beta*gamma*(gamma - 1)*(gamma - 2);
75  G4double grej = 2*(1 + A*B)/A;
76  G4double z, g;
77  do {
78  G4double q = G4UniformRand();
79  z = 2*A*(2*q + Ap2*std::sqrt(q))/(Ap2*Ap2 - 4*q);
80  g = (2 - z)*(1.0/(A + z) + B);
81 
82  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
83  } while(g < G4UniformRand()*grej);
84 
85  G4double cost = 1 - z;
86  G4double sint = std::sqrt(z*(2 - z));
88 
89  fLocalDirection.set(sint*std::cos(phi), sint*std::sin(phi), cost);
91  }
92  return fLocalDirection;
93 }
void set(double x, double y, double z)
G4double GetKineticEnergy() const
double B(double temperature)
static constexpr double g
Definition: G4SIunits.hh:183
static constexpr double electron_mass_c2
#define G4UniformRand()
Definition: Randomize.hh:97
double A(double temperature)
const G4ThreeVector & GetMomentumDirection() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:38
double G4double
Definition: G4Types.hh:76
static constexpr double twopi
Definition: SystemOfUnits.h:55

Here is the call graph for this function:


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