Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4LivermoreNuclearGammaConversionModel.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$
27//
28// Authors: G.Depaola & F.Longo
29//
30
33#include "G4SystemOfUnits.hh"
34
35//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
36
37using namespace std;
38
39//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
40
42 const G4String& nam)
43 :G4VEmModel(nam),fParticleChange(0),smallEnergy(2.*MeV),
44 isInitialised(false),
45 crossSectionHandler(0),meanFreePathTable(0)
46{
47 lowEnergyLimit = 2.0*electron_mass_c2;
48 highEnergyLimit = 100 * GeV;
49 SetHighEnergyLimit(highEnergyLimit);
50
51 verboseLevel= 0;
52 // Verbosity scale:
53 // 0 = nothing
54 // 1 = warning for energy non-conservation
55 // 2 = details of energy budget
56 // 3 = calculation of cross sections, file openings, sampling of atoms
57 // 4 = entering in methods
58
59 if(verboseLevel > 0) {
60 G4cout << "Livermore Nuclear Gamma conversion is constructed " << G4endl
61 << "Energy range: "
62 << lowEnergyLimit / MeV << " MeV - "
63 << highEnergyLimit / GeV << " GeV"
64 << G4endl;
65 }
66}
67
68//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
69
71{
72 if (crossSectionHandler) delete crossSectionHandler;
73}
74
75//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
76
77void
79 const G4DataVector&)
80{
81 if (verboseLevel > 3)
82 G4cout << "Calling G4LivermoreNuclearGammaConversionModel::Initialise()" << G4endl;
83
84 if (crossSectionHandler)
85 {
86 crossSectionHandler->Clear();
87 delete crossSectionHandler;
88 }
89
90 // Read data tables for all materials
91
92 crossSectionHandler = new G4CrossSectionHandler();
93 crossSectionHandler->Initialise(0,lowEnergyLimit,100.*GeV,400);
94 G4String crossSectionFile = "pairdata/pp-pair-cs-"; // here only pair in nuclear field cs should be used
95 crossSectionHandler->LoadData(crossSectionFile);
96
97 //
98
99 if (verboseLevel > 0) {
100 G4cout << "Loaded cross section files for Livermore GammaConversion" << G4endl;
101 G4cout << "To obtain the total cross section this should be used only " << G4endl
102 << "in connection with G4ElectronGammaConversion " << G4endl;
103 }
104
105 if (verboseLevel > 0) {
106 G4cout << "Livermore Nuclear Gamma Conversion model is initialized " << G4endl
107 << "Energy range: "
108 << LowEnergyLimit() / MeV << " MeV - "
109 << HighEnergyLimit() / GeV << " GeV"
110 << G4endl;
111 }
112
113 if(isInitialised) return;
115 isInitialised = true;
116}
117
118//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
119
122 G4double GammaEnergy,
125{
126 if (verboseLevel > 3) {
127 G4cout << "Calling ComputeCrossSectionPerAtom() of G4LivermoreNuclearGammaConversionModel"
128 << G4endl;
129 }
130 if (GammaEnergy < lowEnergyLimit || GammaEnergy > highEnergyLimit) return 0;
131
132 G4double cs = crossSectionHandler->FindValue(G4int(Z), GammaEnergy);
133 return cs;
134}
135
136//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
137
138void G4LivermoreNuclearGammaConversionModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
139 const G4MaterialCutsCouple* couple,
140 const G4DynamicParticle* aDynamicGamma,
141 G4double,
142 G4double)
143{
144
145// The energies of the e+ e- secondaries are sampled using the Bethe - Heitler
146// cross sections with Coulomb correction. A modified version of the random
147// number techniques of Butcher & Messel is used (Nuc Phys 20(1960),15).
148
149// Note 1 : Effects due to the breakdown of the Born approximation at low
150// energy are ignored.
151// Note 2 : The differential cross section implicitly takes account of
152// pair creation in both nuclear and atomic electron fields. However triplet
153// prodution is not generated.
154
155 if (verboseLevel > 3)
156 G4cout << "Calling SampleSecondaries() of G4LivermoreNuclearGammaConversionModel" << G4endl;
157
158 G4double photonEnergy = aDynamicGamma->GetKineticEnergy();
159 G4ParticleMomentum photonDirection = aDynamicGamma->GetMomentumDirection();
160
161 G4double epsilon ;
162 G4double epsilon0Local = electron_mass_c2 / photonEnergy ;
163
164 // Do it fast if photon energy < 2. MeV
165 if (photonEnergy < smallEnergy )
166 {
167 epsilon = epsilon0Local + (0.5 - epsilon0Local) * G4UniformRand();
168 }
169 else
170 {
171 // Select randomly one element in the current material
172 //const G4Element* element = crossSectionHandler->SelectRandomElement(couple,photonEnergy);
173 const G4ParticleDefinition* particle = aDynamicGamma->GetDefinition();
174 const G4Element* element = SelectRandomAtom(couple,particle,photonEnergy);
175
176 if (element == 0)
177 {
178 G4cout << "G4LivermoreNuclearGammaConversionModel::SampleSecondaries - element = 0"
179 << G4endl;
180 return;
181 }
182 G4IonisParamElm* ionisation = element->GetIonisation();
183 if (ionisation == 0)
184 {
185 G4cout << "G4LivermoreNuclearGammaConversionModel::SampleSecondaries - ionisation = 0"
186 << G4endl;
187 return;
188 }
189
190 // Extract Coulomb factor for this Element
191 G4double fZ = 8. * (ionisation->GetlogZ3());
192 if (photonEnergy > 50. * MeV) fZ += 8. * (element->GetfCoulomb());
193
194 // Limits of the screening variable
195 G4double screenFactor = 136. * epsilon0Local / (element->GetIonisation()->GetZ3()) ;
196 G4double screenMax = std::exp ((42.24 - fZ)/8.368) - 0.952 ;
197 G4double screenMin = std::min(4.*screenFactor,screenMax) ;
198
199 // Limits of the energy sampling
200 G4double epsilon1 = 0.5 - 0.5 * std::sqrt(1. - screenMin / screenMax) ;
201 G4double epsilonMin = std::max(epsilon0Local,epsilon1);
202 G4double epsilonRange = 0.5 - epsilonMin ;
203
204 // Sample the energy rate of the created electron (or positron)
205 G4double screen;
206 G4double gReject ;
207
208 G4double f10 = ScreenFunction1(screenMin) - fZ;
209 G4double f20 = ScreenFunction2(screenMin) - fZ;
210 G4double normF1 = std::max(f10 * epsilonRange * epsilonRange,0.);
211 G4double normF2 = std::max(1.5 * f20,0.);
212
213 do {
214 if (normF1 / (normF1 + normF2) > G4UniformRand() )
215 {
216 epsilon = 0.5 - epsilonRange * std::pow(G4UniformRand(), 0.3333) ;
217 screen = screenFactor / (epsilon * (1. - epsilon));
218 gReject = (ScreenFunction1(screen) - fZ) / f10 ;
219 }
220 else
221 {
222 epsilon = epsilonMin + epsilonRange * G4UniformRand();
223 screen = screenFactor / (epsilon * (1 - epsilon));
224 gReject = (ScreenFunction2(screen) - fZ) / f20 ;
225 }
226 } while ( gReject < G4UniformRand() );
227
228 } // End of epsilon sampling
229
230 // Fix charges randomly
231
232 G4double electronTotEnergy;
233 G4double positronTotEnergy;
234
235 if (G4int(2*G4UniformRand()))
236 {
237 electronTotEnergy = (1. - epsilon) * photonEnergy;
238 positronTotEnergy = epsilon * photonEnergy;
239 }
240 else
241 {
242 positronTotEnergy = (1. - epsilon) * photonEnergy;
243 electronTotEnergy = epsilon * photonEnergy;
244 }
245
246 // Scattered electron (positron) angles. ( Z - axis along the parent photon)
247 // Universal distribution suggested by L. Urban (Geant3 manual (1993) Phys211),
248 // derived from Tsai distribution (Rev. Mod. Phys. 49, 421 (1977)
249
250 G4double u;
251 const G4double a1 = 0.625;
252 G4double a2 = 3. * a1;
253 // G4double d = 27. ;
254
255 // if (9. / (9. + d) > G4UniformRand())
256 if (0.25 > G4UniformRand())
257 {
258 u = - std::log(G4UniformRand() * G4UniformRand()) / a1 ;
259 }
260 else
261 {
262 u = - std::log(G4UniformRand() * G4UniformRand()) / a2 ;
263 }
264
265 G4double thetaEle = u*electron_mass_c2/electronTotEnergy;
266 G4double thetaPos = u*electron_mass_c2/positronTotEnergy;
267 G4double phi = twopi * G4UniformRand();
268
269 G4double dxEle= std::sin(thetaEle)*std::cos(phi),dyEle= std::sin(thetaEle)*std::sin(phi),dzEle=std::cos(thetaEle);
270 G4double dxPos=-std::sin(thetaPos)*std::cos(phi),dyPos=-std::sin(thetaPos)*std::sin(phi),dzPos=std::cos(thetaPos);
271
272
273 // Kinematics of the created pair:
274 // the electron and positron are assumed to have a symetric angular
275 // distribution with respect to the Z axis along the parent photon
276
277 G4double electronKineEnergy = std::max(0.,electronTotEnergy - electron_mass_c2) ;
278
279 // SI - The range test has been removed wrt original G4LowEnergyGammaconversion class
280
281 G4ThreeVector electronDirection (dxEle, dyEle, dzEle);
282 electronDirection.rotateUz(photonDirection);
283
285 electronDirection,
286 electronKineEnergy);
287
288 // The e+ is always created (even with kinetic energy = 0) for further annihilation
289 G4double positronKineEnergy = std::max(0.,positronTotEnergy - electron_mass_c2) ;
290
291 // SI - The range test has been removed wrt original G4LowEnergyGammaconversion class
292
293 G4ThreeVector positronDirection (dxPos, dyPos, dzPos);
294 positronDirection.rotateUz(photonDirection);
295
296 // Create G4DynamicParticle object for the particle2
298 positronDirection, positronKineEnergy);
299 // Fill output vector
300
301 fvect->push_back(particle1);
302 fvect->push_back(particle2);
303
304 // kill incident photon
307
308}
309
310//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
311
312G4double G4LivermoreNuclearGammaConversionModel::ScreenFunction1(G4double screenVariable)
313{
314 // Compute the value of the screening function 3*phi1 - phi2
315
316 G4double value;
317
318 if (screenVariable > 1.)
319 value = 42.24 - 8.368 * std::log(screenVariable + 0.952);
320 else
321 value = 42.392 - screenVariable * (7.796 - 1.961 * screenVariable);
322
323 return value;
324}
325
326//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
327
328G4double G4LivermoreNuclearGammaConversionModel::ScreenFunction2(G4double screenVariable)
329{
330 // Compute the value of the screening function 1.5*phi1 - 0.5*phi2
331
332 G4double value;
333
334 if (screenVariable > 1.)
335 value = 42.24 - 8.368 * std::log(screenVariable + 0.952);
336 else
337 value = 41.405 - screenVariable * (5.828 - 0.8945 * screenVariable);
338
339 return value;
340}
341
@ fStopAndKill
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4double GetfCoulomb() const
Definition: G4Element.hh:201
G4IonisParamElm * GetIonisation() const
Definition: G4Element.hh:209
G4double GetlogZ3() const
G4double GetZ3() const
G4LivermoreNuclearGammaConversionModel(const G4ParticleDefinition *p=0, const G4String &nam="LivermoreNuclearGammaConversion")
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
void SetProposedKineticEnergy(G4double proposedKinEnergy)
static G4Positron * Positron()
Definition: G4Positron.cc:94
G4double FindValue(G4int Z, G4double e) const
void LoadData(const G4String &dataFile)
void Initialise(G4VDataSetAlgorithm *interpolation=0, G4double minE=250 *CLHEP::eV, G4double maxE=100 *CLHEP::GeV, G4int numberOfBins=200, G4double unitE=CLHEP::MeV, G4double unitData=CLHEP::barn, G4int minZ=1, G4int maxZ=99)
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:585
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:109
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:529
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:522
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:459
void ProposeTrackStatus(G4TrackStatus status)