Geant4
10.01.p03
Main Page
Related Pages
Modules
Namespaces
Classes
Files
File List
File Members
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
// $Id: G4SauterGavrilaAngularDistribution.cc 78929 2014-02-04 09:00:46Z gcosmo $
27
//
28
// -------------------------------------------------------------------
29
//
30
// GEANT4 Class file
31
//
32
//
33
// File name: G4SauterGavrilaAngularDistribution
34
//
35
// Author: Vladimir Ivanchenko using Michel Maire algorithm
36
// developed for Geant3
37
//
38
// Creation date: 23 July 2012
39
//
40
//
41
// -------------------------------------------------------------------
42
//
43
44
#include "
G4SauterGavrilaAngularDistribution.hh
"
45
#include "
G4PhysicalConstants.hh
"
46
#include "
Randomize.hh
"
47
48
G4SauterGavrilaAngularDistribution::G4SauterGavrilaAngularDistribution
()
49
:
G4VEmAngularDistribution
(
"AngularGenSauterGavrila"
)
50
{}
51
52
G4SauterGavrilaAngularDistribution::~G4SauterGavrilaAngularDistribution
()
53
{}
54
55
G4ThreeVector
&
56
G4SauterGavrilaAngularDistribution::SampleDirection
(
57
const
G4DynamicParticle
* dp,
G4double
,
G4int
,
const
G4Material
*)
58
{
59
G4double
tau = dp->
GetKineticEnergy
()/electron_mass_c2;
60
static
const
G4double
taulimit = 50.0;
61
62
if
(tau > taulimit) {
63
fLocalDirection
= dp->
GetMomentumDirection
();
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
}
while
(g <
G4UniformRand
()*grej);
83
84
G4double
cost = 1 -
z
;
85
G4double
sint = std::sqrt(z*(2 - z));
86
G4double
phi = CLHEP::twopi*
G4UniformRand
();
87
88
fLocalDirection
.set(sint*std::cos(phi), sint*std::sin(phi), cost);
89
fLocalDirection
.rotateUz(dp->
GetMomentumDirection
());
90
}
91
return
fLocalDirection
;
92
}
93
94
void
G4SauterGavrilaAngularDistribution::PrintGeneratorInformation
()
const
95
{
96
G4cout
<<
"\n"
<<
G4endl
;
97
G4cout
<<
"Non-polarized photoelectric effect angular generator."
<<
G4endl
;
98
G4cout
<<
"The Sauter-Gavrila distribution for the K-shell is used."
<<
G4endl
;
99
G4cout
<<
"Originally developed by M.Maire for Geant3"
100
<<
G4endl
;
101
}
G4DynamicParticle::GetKineticEnergy
G4double GetKineticEnergy() const
G4ThreeVector
CLHEP::Hep3Vector G4ThreeVector
Definition:
G4ThreeVector.hh:42
z
G4double z
Definition:
TRTMaterials.hh:39
G4Material
Definition:
G4Material.hh:120
G4DynamicParticle
Definition:
G4DynamicParticle.hh:73
G4SauterGavrilaAngularDistribution::~G4SauterGavrilaAngularDistribution
virtual ~G4SauterGavrilaAngularDistribution()
Definition:
G4SauterGavrilaAngularDistribution.cc:52
G4SauterGavrilaAngularDistribution::SampleDirection
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double e=0.0, G4int shellId=0, const G4Material *mat=0)
Definition:
G4SauterGavrilaAngularDistribution.cc:56
G4int
int G4int
Definition:
G4Types.hh:78
G4SauterGavrilaAngularDistribution::PrintGeneratorInformation
void PrintGeneratorInformation() const
Definition:
G4SauterGavrilaAngularDistribution.cc:94
G4SauterGavrilaAngularDistribution.hh
G4UniformRand
#define G4UniformRand()
Definition:
Randomize.hh:93
G4cout
G4GLOB_DLL std::ostream G4cout
G4DynamicParticle::GetMomentumDirection
const G4ThreeVector & GetMomentumDirection() const
Randomize.hh
A
static const G4double A[nN]
Definition:
G4ElectroNuclearCrossSection.cc:115
G4PhysicalConstants.hh
G4SauterGavrilaAngularDistribution::G4SauterGavrilaAngularDistribution
G4SauterGavrilaAngularDistribution()
Definition:
G4SauterGavrilaAngularDistribution.cc:48
g
static const double g
Definition:
G4SIunits.hh:162
G4VEmAngularDistribution
Definition:
G4VEmAngularDistribution.hh:60
G4endl
#define G4endl
Definition:
G4ios.hh:61
G4VEmAngularDistribution::fLocalDirection
G4ThreeVector fLocalDirection
Definition:
G4VEmAngularDistribution.hh:90
G4double
double G4double
Definition:
G4Types.hh:76
geant4.10.01.p03
source
processes
electromagnetic
standard
src
G4SauterGavrilaAngularDistribution.cc
Generated on Fri Feb 19 2016 15:53:58 for Geant4 by
1.8.8