Geant4
10.03.p03
Main Page
Related Pages
Modules
Namespaces
Classes
Files
File List
File Members
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Groups
Pages
G4ModifiedTsai.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: G4ModifiedTsai.cc 91726 2015-08-03 15:41:36Z gcosmo $
27
//
28
// -------------------------------------------------------------------
29
//
30
// GEANT4 Class file
31
//
32
//
33
// File name: G4ModifiedTsai
34
//
35
// Author: Andreia Trindade (andreia@lip.pt)
36
// Pedro Rodrigues (psilva@lip.pt)
37
// Luis Peralta (luis@lip.pt)
38
//
39
// Creation date: 21 March 2003
40
//
41
// Modifications:
42
// 21 Mar 2003 A.Trindade First implementation acording with new design
43
// 24 Mar 2003 A.Trindade Fix in Tsai generator in order to prevent theta
44
// generation above pi
45
// 13 Oct 2010 V.Ivanchenko Moved to standard and improved comment
46
// 26.04.2011 V.Grichine Clean-up of PolarAngle method
47
//
48
// Class Description:
49
//
50
// Bremsstrahlung Angular Distribution Generation
51
// suggested by L.Urban (Geant3 manual (1993) Phys211)
52
// Derived from Tsai distribution (Rev Mod Phys 49,421(1977))
53
//
54
// Class Description: End
55
//
56
// -------------------------------------------------------------------
57
//
58
59
#include "
G4ModifiedTsai.hh
"
60
#include "
G4PhysicalConstants.hh
"
61
#include "
Randomize.hh
"
62
#include "
G4Log.hh
"
63
64
G4ModifiedTsai::G4ModifiedTsai
(
const
G4String
&)
65
:
G4VEmAngularDistribution
(
"AngularGenUrban"
)
66
{}
67
68
G4ModifiedTsai::~G4ModifiedTsai
()
69
{}
70
71
G4ThreeVector
&
72
G4ModifiedTsai::SampleDirection
(
const
G4DynamicParticle
* dp,
73
G4double
,
G4int
,
const
G4Material
*)
74
{
75
// Sample gamma angle (Z - axis along the parent particle).
76
// Universal distribution suggested by L. Urban (Geant3 manual (1993)
77
// Phys211) derived from Tsai distribution (Rev Mod Phys 49,421(1977))
78
79
G4double
uMax = 2*(1. + dp->
GetKineticEnergy
()/
electron_mass_c2
);
80
81
static
const
G4double
a1 = 0.625;
82
static
const
G4double
a2 = 1.875;
83
static
const
G4double
border
= 0.25;
84
G4double
u;
85
86
do
{
87
u = -
G4Log
(
G4UniformRand
()*
G4UniformRand
());
88
89
if
( border >
G4UniformRand
() ) { u /= a1; }
90
else
{ u /= a2; }
91
92
// Loop checking, 03-Aug-2015, Vladimir Ivanchenko
93
}
while
(u > uMax);
94
95
G4double
cost = 1.0 - 2*u*u/(uMax*uMax);
96
G4double
sint = std::sqrt((1 - cost)*(1 + cost));
97
G4double
phi =
CLHEP::twopi
*
G4UniformRand
();
98
99
fLocalDirection
.
set
(sint*std::cos(phi), sint*std::sin(phi), cost);
100
fLocalDirection
.
rotateUz
(dp->
GetMomentumDirection
());
101
102
return
fLocalDirection
;
103
}
104
105
void
G4ModifiedTsai::PrintGeneratorInformation
()
const
106
{
107
G4cout
<<
"\n"
<<
G4endl
;
108
G4cout
<<
"Bremsstrahlung Angular Generator is Modified Tsai"
<<
G4endl
;
109
G4cout
<<
"Distribution suggested by L.Urban (Geant3 manual (1993) Phys211)"
110
<<
G4endl
;
111
G4cout
<<
"Derived from Tsai distribution (Rev Mod Phys 49,421(1977)) \n"
112
<<
G4endl
;
113
}
CLHEP::Hep3Vector::set
void set(double x, double y, double z)
CLHEP::Hep3Vector
Definition:
ThreeVector.h:41
border
static unsigned border[]
Definition:
csz_inflate.cc:298
G4DynamicParticle::GetKineticEnergy
G4double GetKineticEnergy() const
G4ModifiedTsai::SampleDirection
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double out_energy, G4int Z, const G4Material *mat=nullptr) final
Definition:
G4ModifiedTsai.cc:72
G4Material
Definition:
G4Material.hh:120
G4DynamicParticle
Definition:
G4DynamicParticle.hh:73
G4int
int G4int
Definition:
G4Types.hh:78
G4UniformRand
#define G4UniformRand()
Definition:
Randomize.hh:97
G4cout
G4GLOB_DLL std::ostream G4cout
G4ModifiedTsai::PrintGeneratorInformation
virtual void PrintGeneratorInformation() const final
Definition:
G4ModifiedTsai.cc:105
G4DynamicParticle::GetMomentumDirection
const G4ThreeVector & GetMomentumDirection() const
CLHEP::Hep3Vector::rotateUz
Hep3Vector & rotateUz(const Hep3Vector &)
Definition:
ThreeVector.cc:38
python.hepunit.electron_mass_c2
float electron_mass_c2
Definition:
hepunit.py:274
Randomize.hh
G4Log.hh
G4PhysicalConstants.hh
G4Log
G4double G4Log(G4double x)
Definition:
G4Log.hh:230
G4ModifiedTsai::G4ModifiedTsai
G4ModifiedTsai(const G4String &name="")
Definition:
G4ModifiedTsai.cc:64
G4ModifiedTsai::~G4ModifiedTsai
virtual ~G4ModifiedTsai()
Definition:
G4ModifiedTsai.cc:68
G4VEmAngularDistribution
Definition:
G4VEmAngularDistribution.hh:60
G4endl
#define G4endl
Definition:
G4ios.hh:61
G4ModifiedTsai.hh
G4VEmAngularDistribution::fLocalDirection
G4ThreeVector fLocalDirection
Definition:
G4VEmAngularDistribution.hh:90
G4double
double G4double
Definition:
G4Types.hh:76
CLHEP::twopi
static constexpr double twopi
Definition:
SystemOfUnits.h:55
G4String
Definition:
G4String.hh:45
source
geant4.10.03.p03
source
processes
electromagnetic
standard
src
G4ModifiedTsai.cc
Generated on Tue Nov 28 2017 21:44:07 for Geant4 by
1.8.5