Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4SauterGavrilaAngularDistribution.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class file
30//
31//
32// File name: G4SauterGavrilaAngularDistribution
33//
34// Author: Vladimir Ivanchenko using Michel Maire algorithm
35// developed for Geant3
36//
37// Creation date: 23 July 2012
38//
39// Modified:
40// 04.05.2017 Marilena Bandieramonte implemented Penelope 2014 algorithm
41// -------------------------------------------------------------------
42//
43
46#include "Randomize.hh"
47
49 : G4VEmAngularDistribution("SauterGavrila")
50{}
51
53 = default;
54
56 const G4DynamicParticle* dp, G4double, G4int, const G4Material*)
57{
58 static const G4double emin = 1*CLHEP::eV;
59 static const G4double emax = 100*CLHEP::MeV;
60
61 G4double energy = std::max(dp->GetKineticEnergy(), emin);
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)
73 G4double tau = energy/electron_mass_c2;
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));
102 G4double phi = CLHEP::twopi*G4UniformRand();
103
104 fLocalDirection.set(sint*std::cos(phi), sint*std::sin(phi), costheta);
106 }
107 return fLocalDirection;
108}
109
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}
118
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#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
G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double e=0.0, G4int shellId=0, const G4Material *mat=nullptr) final