Geant4-11
source
processes
electromagnetic
lowenergy
src
G4PhotoElectricAngularGeneratorSauterGavrila.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: G4PhotoElectricAngularGeneratorSauterGavrila
33
//
34
// Creation date: 10 May 2004
35
//
36
// Modifications:
37
// 10 May 2003 P. Rodrigues First implementation acording with new design
38
//
39
// Class Description:
40
//
41
// Concrete class for PhotoElectric Electron Angular Distribution Generation
42
// This model is a re-implementation of the Photolectric angular distribution
43
// developed my M. Maire for the Standard EM Physics G4PhotoElectricEffect
44
//
45
// Class Description: End
46
//
47
// -------------------------------------------------------------------
48
//
49
50
#include "
G4PhotoElectricAngularGeneratorSauterGavrila.hh
"
51
#include "
G4PhysicalConstants.hh
"
52
#include "
Randomize.hh
"
53
54
// -------------------------------------------------------------------
55
G4PhotoElectricAngularGeneratorSauterGavrila::G4PhotoElectricAngularGeneratorSauterGavrila
():
56
G4VEmAngularDistribution
(
"AngularGenSauterGavrilaLowE"
)
57
{}
58
59
// -------------------------------------------------------------------
60
61
G4PhotoElectricAngularGeneratorSauterGavrila::~G4PhotoElectricAngularGeneratorSauterGavrila
()
62
{}
63
64
// -------------------------------------------------------------------
65
66
G4ThreeVector
&
67
G4PhotoElectricAngularGeneratorSauterGavrila::SampleDirection
(
68
const
G4DynamicParticle
* dp,
69
G4double
,
G4int
,
const
G4Material
*)
70
{
71
72
// Compute Theta distribution of the emitted electron, with respect to the
73
// incident Gamma.
74
// The Sauter-Gavrila distribution for the K-shell is used.
75
G4double
costeta = 1.;
76
G4double
Phi =
twopi
*
G4UniformRand
();
77
G4double
cosphi = std::cos(Phi);
78
G4double
sinphi = std::sin(Phi);
79
G4double
sinteta = 0;
80
G4double
gamma = 1. + dp->
GetKineticEnergy
()/
electron_mass_c2
;
81
82
if
(gamma > 5.) {
83
fLocalDirection
= dp->
GetMomentumDirection
();
84
return
fLocalDirection
;
85
// Bugzilla 1120
86
// SI on 05/09/2010 as suggested by JG 04/09/10
87
}
88
89
G4double
beta
= std::sqrt((gamma - 1)*(gamma + 1))/gamma;
90
G4double
b = 0.5*gamma*(gamma - 1)*(gamma - 2);
91
92
G4double
rndm,term,greject,grejsup;
93
if
(gamma < 2.) grejsup = gamma*gamma*(1.+b-
beta
*b);
94
else
grejsup = gamma*gamma*(1.+b+
beta
*b);
95
96
do
{ rndm = 1.-2*
G4UniformRand
();
97
costeta = (rndm+
beta
)/(rndm*
beta
+1.);
98
term = 1.-
beta
*costeta;
99
greject = (1.-costeta*costeta)*(1.+b*term)/(term*term);
100
}
while
(greject <
G4UniformRand
()*grejsup);
101
102
sinteta = std::sqrt((1 - costeta)*(1 + costeta));
103
fLocalDirection
.
set
(sinteta*cosphi, sinteta*sinphi, costeta);
104
fLocalDirection
.
rotateUz
(dp->
GetMomentumDirection
());
105
return
fLocalDirection
;
106
}
107
108
// -------------------------------------------------------------------
109
110
void
G4PhotoElectricAngularGeneratorSauterGavrila::PrintGeneratorInformation
()
const
111
{
112
G4cout
<<
"\n"
<<
G4endl
;
113
G4cout
<<
""
<<
G4endl
;
114
G4cout
<<
"Re-implementation of the photolectric angular distribution"
<<
G4endl
;
115
G4cout
<<
"developed my M. Maire for the Standard EM Physics G4PhotoElectricEffect"
<<
G4endl
;
116
G4cout
<<
"It computes the theta distribution of the emitted electron, with respect to the"
<<
G4endl
;
117
G4cout
<<
"incident Gamma, using the Sauter-Gavrila distribution for the K-shell\n"
<<
G4endl
;
118
}
G4PhotoElectricAngularGeneratorSauterGavrila.hh
G4PhysicalConstants.hh
twopi
static constexpr double twopi
Definition:
G4SIunits.hh:56
G4double
double G4double
Definition:
G4Types.hh:83
G4int
int G4int
Definition:
G4Types.hh:85
G4endl
#define G4endl
Definition:
G4ios.hh:57
G4cout
G4GLOB_DLL std::ostream G4cout
Randomize.hh
G4UniformRand
#define G4UniformRand()
Definition:
Randomize.hh:52
CLHEP::Hep3Vector
Definition:
ThreeVector.h:36
CLHEP::Hep3Vector::set
void set(double x, double y, double z)
CLHEP::Hep3Vector::rotateUz
Hep3Vector & rotateUz(const Hep3Vector &)
Definition:
ThreeVector.cc:33
G4DynamicParticle
Definition:
G4DynamicParticle.hh:65
G4DynamicParticle::GetMomentumDirection
const G4ThreeVector & GetMomentumDirection() const
G4DynamicParticle::GetKineticEnergy
G4double GetKineticEnergy() const
G4Material
Definition:
G4Material.hh:118
G4PhotoElectricAngularGeneratorSauterGavrila::SampleDirection
G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double e=0.0, G4int shellId=0, const G4Material *mat=nullptr) override
Definition:
G4PhotoElectricAngularGeneratorSauterGavrila.cc:67
G4PhotoElectricAngularGeneratorSauterGavrila::~G4PhotoElectricAngularGeneratorSauterGavrila
~G4PhotoElectricAngularGeneratorSauterGavrila()
Definition:
G4PhotoElectricAngularGeneratorSauterGavrila.cc:61
G4PhotoElectricAngularGeneratorSauterGavrila::PrintGeneratorInformation
void PrintGeneratorInformation() const override
Definition:
G4PhotoElectricAngularGeneratorSauterGavrila.cc:110
G4PhotoElectricAngularGeneratorSauterGavrila::G4PhotoElectricAngularGeneratorSauterGavrila
G4PhotoElectricAngularGeneratorSauterGavrila()
Definition:
G4PhotoElectricAngularGeneratorSauterGavrila.cc:55
G4VEmAngularDistribution
Definition:
G4VEmAngularDistribution.hh:59
G4VEmAngularDistribution::fLocalDirection
G4ThreeVector fLocalDirection
Definition:
G4VEmAngularDistribution.hh:103
anonymous_namespace{G4PionRadiativeDecayChannel.cc}::beta
const G4double beta
Definition:
G4PionRadiativeDecayChannel.cc:44
source.hepunit.electron_mass_c2
float electron_mass_c2
Definition:
hepunit.py:273
Generated by
1.9.3