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