Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros 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//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class file
30//
31//
32// File name: G4ModifiedTsai
33//
34// Author: Andreia Trindade (andreia@lip.pt)
35// Pedro Rodrigues (psilva@lip.pt)
36// Luis Peralta (luis@lip.pt)
37//
38// Creation date: 21 March 2003
39//
40// Modifications:
41// 21 Mar 2003 A.Trindade First implementation acording with new design
42// 24 Mar 2003 A.Trindade Fix in Tsai generator in order to prevent theta
43// generation above pi
44// 13 Oct 2010 V.Ivanchenko Moved to standard and improved comment
45// 26.04.2011 V.Grichine Clean-up of PolarAngle method
46//
47// Class Description:
48//
49// Bremsstrahlung Angular Distribution Generation
50// suggested by L.Urban (Geant3 manual (1993) Phys211)
51// Derived from Tsai distribution (Rev Mod Phys 49,421(1977))
52//
53// Class Description: End
54//
55// -------------------------------------------------------------------
56//
57
58#include "G4ModifiedTsai.hh"
59#include "Randomize.hh"
60#include "G4Log.hh"
62
63//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
64
66 : G4VEmAngularDistribution("ModifiedTsai")
67{}
68
69//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
70
72
73//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
74
77 G4double, G4int, const G4Material*)
78{
79 // Sample gamma angle (Z - axis along the parent particle).
80 G4double cost = SampleCosTheta(dp->GetKineticEnergy());
81 G4double sint = std::sqrt((1 - cost)*(1 + cost));
82 G4double phi = CLHEP::twopi*G4UniformRand();
83
84 fLocalDirection.set(sint*std::cos(phi), sint*std::sin(phi), cost);
86
87 return fLocalDirection;
88}
89
90//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
91
92G4double G4ModifiedTsai::SampleCosTheta(G4double kinEnergy)
93{
94 // Universal distribution suggested by L. Urban (Geant3 manual (1993)
95 // Phys211) derived from Tsai distribution (Rev Mod Phys 49,421(1977))
96
97 G4double uMax = 2*(1. + kinEnergy/CLHEP::electron_mass_c2);
98
99 static const G4double a1 = 1.6;
100 static const G4double a2 = a1/3.;
101 static const G4double border = 0.25;
102 G4double u;
103 CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
104
105 do {
106 G4double uu = -G4Log(rndmEngine->flat()*rndmEngine->flat());
107 u = (border > rndmEngine->flat()) ? uu*a1 : uu*a2;
108
109 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
110 } while(u > uMax);
111
112 return 1.0 - 2.0*u*u/(uMax*uMax);
113}
114
115//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
116
118 G4double elecKinEnergy,
119 G4double posiKinEnergy,
120 G4ThreeVector& dirElectron,
121 G4ThreeVector& dirPositron,
122 G4int, const G4Material*)
123{
124 G4double phi = CLHEP::twopi * G4UniformRand();
125 G4double sinp = std::sin(phi);
126 G4double cosp = std::cos(phi);
127
128 G4double cost = SampleCosTheta(elecKinEnergy);
129 G4double sint = std::sqrt((1. - cost)*(1. + cost));
130
131 dirElectron.set(sint*cosp, sint*sinp, cost);
132 dirElectron.rotateUz(dp->GetMomentumDirection());
133
134 cost = SampleCosTheta(posiKinEnergy);
135 sint = std::sqrt((1. - cost)*(1. + cost));
136
137 dirPositron.set(-sint*cosp, -sint*sinp, cost);
138 dirPositron.rotateUz(dp->GetMomentumDirection());
139}
140
141//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
142
144{
145 G4cout << "\n" << G4endl;
146 G4cout << "Bremsstrahlung Angular Generator is Modified Tsai" << G4endl;
147 G4cout << "Distribution suggested by L.Urban (Geant3 manual (1993) Phys211)"
148 << G4endl;
149 G4cout << "Derived from Tsai distribution (Rev Mod Phys 49,421(1977)) \n"
150 << G4endl;
151}
152
153//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
G4double G4Log(G4double x)
Definition: G4Log.hh:227
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
virtual double flat()=0
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4ModifiedTsai(const G4String &name="")
G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double out_energy, G4int Z, const G4Material *mat=nullptr) final
void PrintGeneratorInformation() const final
~G4ModifiedTsai() override
void SamplePairDirections(const G4DynamicParticle *dp, G4double elecKinEnergy, G4double posiKinEnergy, G4ThreeVector &dirElectron, G4ThreeVector &dirPositron, G4int Z=0, const G4Material *mat=nullptr) final