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
G4LivermoreGammaConversionModel.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// Author: Sebastien Incerti
27// 22 January 2012
28// on base of G4LivermoreGammaConversionModel
29
32#include "G4SystemOfUnits.hh"
33
34//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
35
36using namespace std;
37
38//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
39
41 const G4String& nam)
42:G4VEmModel(nam),smallEnergy(2.*MeV),isInitialised(false),maxZ(99)
43{
44 fParticleChange = 0;
45
46 lowEnergyLimit = 2.0*electron_mass_c2;
47 data.resize(maxZ+1,0);
48
49 verboseLevel= 0;
50 // Verbosity scale for debugging purposes:
51 // 0 = nothing
52 // 1 = calculation of cross sections, file openings...
53 // 2 = entering in methods
54
55 if(verboseLevel > 0)
56 {
57 G4cout << "G4LivermoreGammaConversionModel is constructed " << G4endl;
58 }
59}
60
61//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
62
64{
65 for(G4int i=0; i<=maxZ; ++i) { delete data[i]; }
66}
67
68//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
69
70void
72 const G4DataVector& cuts)
73{
74 if (verboseLevel > 1)
75 {
76 G4cout << "Calling Initialise() of G4LivermoreGammaConversionModel." << G4endl
77 << "Energy range: "
78 << LowEnergyLimit() / MeV << " MeV - "
79 << HighEnergyLimit() / GeV << " GeV"
80 << G4endl;
81 }
82
83 // Initialise element selector
84
85 InitialiseElementSelectors(particle, cuts);
86
87 // Access to elements
88
89 char* path = getenv("G4LEDATA");
90
91 G4ProductionCutsTable* theCoupleTable =
93 G4int numOfCouples = theCoupleTable->GetTableSize();
94
95 for(G4int i=0; i<numOfCouples; ++i)
96 {
97 const G4Material* material =
98 theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
99 const G4ElementVector* theElementVector = material->GetElementVector();
100 G4int nelm = material->GetNumberOfElements();
101
102 for (G4int j=0; j<nelm; ++j)
103 {
104
105 G4int Z = (G4int)(*theElementVector)[j]->GetZ();
106 if(Z < 1) { Z = 1; }
107 else if(Z > maxZ) { Z = maxZ; }
108 if(!data[Z]) { ReadData(Z, path); }
109 }
110 }
111 //
112
113 if(isInitialised) { return; }
114 fParticleChange = GetParticleChangeForGamma();
115 isInitialised = true;
116}
117
118//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
119
120void G4LivermoreGammaConversionModel::ReadData(size_t Z, const char* path)
121{
122 if (verboseLevel > 1)
123 {
124 G4cout << "Calling ReadData() of G4LivermoreGammaConversionModel"
125 << G4endl;
126 }
127
128 if(data[Z]) { return; }
129
130 const char* datadir = path;
131
132 if(!datadir)
133 {
134 datadir = getenv("G4LEDATA");
135 if(!datadir)
136 {
137 G4Exception("G4LivermoreGammaConversionModel::ReadData()",
138 "em0006",FatalException,
139 "Environment variable G4LEDATA not defined");
140 return;
141 }
142 }
143
144 //
145
146 data[Z] = new G4LPhysicsFreeVector();
147
148 // Activation of spline interpolation
149 data[Z] ->SetSpline(true);
150 //
151
152 std::ostringstream ost;
153 ost << datadir << "/livermore/pair/pp-cs-" << Z <<".dat";
154 std::ifstream fin(ost.str().c_str());
155
156 if( !fin.is_open())
157 {
159 ed << "G4LivermoreGammaConversionModel data file <" << ost.str().c_str()
160 << "> is not opened!" << G4endl;
161 G4Exception("G4LivermoreGammaConversionModel::ReadData()",
162 "em0003",FatalException,
163 ed,"G4LEDATA version should be G4EMLOW6.27 or later.");
164 return;
165 }
166
167 else
168 {
169
170 if(verboseLevel > 3) { G4cout << "File " << ost.str()
171 << " is opened by G4LivermoreGammaConversionModel" << G4endl;}
172
173 data[Z]->Retrieve(fin, true);
174 }
175
176
177}
178
179//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
180
183 G4double GammaEnergy,
186{
187 if (verboseLevel > 1)
188 {
189 G4cout << "Calling ComputeCrossSectionPerAtom() of G4LivermoreGammaConversionModel"
190 << G4endl;
191 }
192
193 if (GammaEnergy < lowEnergyLimit) { return 0.0; }
194
195 G4double xs = 0.0;
196
197 G4int intZ=G4int(Z);
198
199 if(intZ < 1 || intZ > maxZ) { return xs; }
200
201 G4LPhysicsFreeVector* pv = data[intZ];
202
203 // element was not initialised
204 if(!pv)
205 {
206 char* path = getenv("G4LEDATA");
207 ReadData(intZ, path);
208 pv = data[intZ];
209 if(!pv) { return xs; }
210 }
211 // x-section is taken from the table
212 xs = pv->Value(GammaEnergy);
213
214 if(verboseLevel > 0)
215 {
216 G4int n = pv->GetVectorLength() - 1;
217 G4cout << "****** DEBUG: tcs value for Z=" << Z << " at energy (MeV)=" << GammaEnergy/MeV << G4endl;
218 G4cout << " cs (Geant4 internal unit)=" << xs << G4endl;
219 G4cout << " -> first cs value in EADL data file (iu) =" << (*pv)[0] << G4endl;
220 G4cout << " -> last cs value in EADL data file (iu) =" << (*pv)[n] << G4endl;
221 G4cout << "*********************************************************" << G4endl;
222 }
223
224 return xs;
225
226}
227
228//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
229
231 std::vector<G4DynamicParticle*>* fvect,
232 const G4MaterialCutsCouple* couple,
233 const G4DynamicParticle* aDynamicGamma,
235{
236
237// The energies of the e+ e- secondaries are sampled using the Bethe - Heitler
238// cross sections with Coulomb correction. A modified version of the random
239// number techniques of Butcher & Messel is used (Nuc Phys 20(1960),15).
240
241// Note 1 : Effects due to the breakdown of the Born approximation at low
242// energy are ignored.
243// Note 2 : The differential cross section implicitly takes account of
244// pair creation in both nuclear and atomic electron fields. However triplet
245// prodution is not generated.
246
247 if (verboseLevel > 1)
248 G4cout << "Calling SampleSecondaries() of G4LivermoreGammaConversionModel" << G4endl;
249
250 G4double photonEnergy = aDynamicGamma->GetKineticEnergy();
251 G4ParticleMomentum photonDirection = aDynamicGamma->GetMomentumDirection();
252
253 G4double epsilon ;
254 G4double epsilon0Local = electron_mass_c2 / photonEnergy ;
255
256 // Do it fast if photon energy < 2. MeV
257 if (photonEnergy < smallEnergy )
258 {
259 epsilon = epsilon0Local + (0.5 - epsilon0Local) * G4UniformRand();
260 }
261 else
262 {
263 // Select randomly one element in the current material
264
265 const G4ParticleDefinition* particle = aDynamicGamma->GetDefinition();
266 const G4Element* element = SelectRandomAtom(couple,particle,photonEnergy);
267
268 if (element == 0)
269 {
270 G4cout << "G4LivermoreGammaConversionModel::SampleSecondaries - element = 0"
271 << G4endl;
272 return;
273 }
274 G4IonisParamElm* ionisation = element->GetIonisation();
275 if (ionisation == 0)
276 {
277 G4cout << "G4LivermoreGammaConversionModel::SampleSecondaries - ionisation = 0"
278 << G4endl;
279 return;
280 }
281
282 // Extract Coulomb factor for this Element
283 G4double fZ = 8. * (ionisation->GetlogZ3());
284 if (photonEnergy > 50. * MeV) fZ += 8. * (element->GetfCoulomb());
285
286 // Limits of the screening variable
287 G4double screenFactor = 136. * epsilon0Local / (element->GetIonisation()->GetZ3()) ;
288 G4double screenMax = std::exp ((42.24 - fZ)/8.368) - 0.952 ;
289 G4double screenMin = std::min(4.*screenFactor,screenMax) ;
290
291 // Limits of the energy sampling
292 G4double epsilon1 = 0.5 - 0.5 * std::sqrt(1. - screenMin / screenMax) ;
293 G4double epsilonMin = std::max(epsilon0Local,epsilon1);
294 G4double epsilonRange = 0.5 - epsilonMin ;
295
296 // Sample the energy rate of the created electron (or positron)
297 G4double screen;
298 G4double gReject ;
299
300 G4double f10 = ScreenFunction1(screenMin) - fZ;
301 G4double f20 = ScreenFunction2(screenMin) - fZ;
302 G4double normF1 = std::max(f10 * epsilonRange * epsilonRange,0.);
303 G4double normF2 = std::max(1.5 * f20,0.);
304
305 do
306 {
307 if (normF1 / (normF1 + normF2) > G4UniformRand() )
308 {
309 epsilon = 0.5 - epsilonRange * std::pow(G4UniformRand(), 0.333333) ;
310 screen = screenFactor / (epsilon * (1. - epsilon));
311 gReject = (ScreenFunction1(screen) - fZ) / f10 ;
312 }
313 else
314 {
315 epsilon = epsilonMin + epsilonRange * G4UniformRand();
316 screen = screenFactor / (epsilon * (1 - epsilon));
317 gReject = (ScreenFunction2(screen) - fZ) / f20 ;
318 }
319 } while ( gReject < G4UniformRand() );
320
321 } // End of epsilon sampling
322
323 // Fix charges randomly
324
325 G4double electronTotEnergy;
326 G4double positronTotEnergy;
327
328 if (G4UniformRand() > 0.5)
329 {
330 electronTotEnergy = (1. - epsilon) * photonEnergy;
331 positronTotEnergy = epsilon * photonEnergy;
332 }
333 else
334 {
335 positronTotEnergy = (1. - epsilon) * photonEnergy;
336 electronTotEnergy = epsilon * photonEnergy;
337 }
338
339 // Scattered electron (positron) angles. ( Z - axis along the parent photon)
340 // Universal distribution suggested by L. Urban (Geant3 manual (1993) Phys211),
341 // derived from Tsai distribution (Rev. Mod. Phys. 49, 421 (1977)
342
343 G4double u;
344 const G4double a1 = 0.625;
345 G4double a2 = 3. * a1;
346 // G4double d = 27. ;
347
348 // if (9. / (9. + d) > G4UniformRand())
349 if (0.25 > G4UniformRand())
350 {
351 u = - std::log(G4UniformRand() * G4UniformRand()) / a1 ;
352 }
353 else
354 {
355 u = - std::log(G4UniformRand() * G4UniformRand()) / a2 ;
356 }
357
358 G4double thetaEle = u*electron_mass_c2/electronTotEnergy;
359 G4double thetaPos = u*electron_mass_c2/positronTotEnergy;
360 G4double phi = twopi * G4UniformRand();
361
362 G4double dxEle= std::sin(thetaEle)*std::cos(phi),dyEle= std::sin(thetaEle)*std::sin(phi),dzEle=std::cos(thetaEle);
363 G4double dxPos=-std::sin(thetaPos)*std::cos(phi),dyPos=-std::sin(thetaPos)*std::sin(phi),dzPos=std::cos(thetaPos);
364
365
366 // Kinematics of the created pair:
367 // the electron and positron are assumed to have a symetric angular
368 // distribution with respect to the Z axis along the parent photon
369
370 G4double electronKineEnergy = std::max(0.,electronTotEnergy - electron_mass_c2) ;
371
372 G4ThreeVector electronDirection (dxEle, dyEle, dzEle);
373 electronDirection.rotateUz(photonDirection);
374
376 electronDirection,
377 electronKineEnergy);
378
379 // The e+ is always created
380 G4double positronKineEnergy = std::max(0.,positronTotEnergy - electron_mass_c2) ;
381
382 G4ThreeVector positronDirection (dxPos, dyPos, dzPos);
383 positronDirection.rotateUz(photonDirection);
384
385 // Create G4DynamicParticle object for the particle2
387 positronDirection,
388 positronKineEnergy);
389 // Fill output vector
390 fvect->push_back(particle1);
391 fvect->push_back(particle2);
392
393 // kill incident photon
394 fParticleChange->SetProposedKineticEnergy(0.);
395 fParticleChange->ProposeTrackStatus(fStopAndKill);
396
397}
398
399//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
400
402G4LivermoreGammaConversionModel::ScreenFunction1(G4double screenVariable)
403{
404 // Compute the value of the screening function 3*phi1 - phi2
405
406 G4double value;
407
408 if (screenVariable > 1.)
409 value = 42.24 - 8.368 * std::log(screenVariable + 0.952);
410 else
411 value = 42.392 - screenVariable * (7.796 - 1.961 * screenVariable);
412
413 return value;
414}
415
416//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
417
419G4LivermoreGammaConversionModel::ScreenFunction2(G4double screenVariable)
420{
421 // Compute the value of the screening function 1.5*phi1 - 0.5*phi2
422
423 G4double value;
424
425 if (screenVariable > 1.)
426 value = 42.24 - 8.368 * std::log(screenVariable + 0.952);
427 else
428 value = 41.405 - screenVariable * (5.828 - 0.8945 * screenVariable);
429
430 return value;
431}
432
std::vector< G4Element * > G4ElementVector
@ FatalException
@ 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
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)
G4LivermoreGammaConversionModel(const G4ParticleDefinition *p=0, const G4String &nam="LivermoreConversion")
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
const G4Material * GetMaterial() const
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:189
size_t GetNumberOfElements() const
Definition: G4Material.hh:185
void SetProposedKineticEnergy(G4double proposedKinEnergy)
G4double Value(G4double theEnergy)
size_t GetVectorLength() const
static G4Positron * Positron()
Definition: G4Positron.cc:94
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
static G4ProductionCutsTable * GetProductionCutsTable()
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 InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:123
void ProposeTrackStatus(G4TrackStatus status)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76