Geant4_10
G4PenelopeGammaConversionModel.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: G4PenelopeGammaConversionModel.cc 76220 2013-11-08 10:15:00Z gcosmo $
27 //
28 // Author: Luciano Pandola
29 //
30 // History:
31 // --------
32 // 13 Jan 2010 L Pandola First implementation (updated to Penelope08)
33 // 24 May 2011 L Pandola Renamed (make v2008 as default Penelope)
34 // 18 Sep 2013 L Pandola Migration to MT paradigm. Only master model deals with
35 // data and creates shared tables
36 //
37 
39 #include "G4PhysicalConstants.hh"
40 #include "G4SystemOfUnits.hh"
41 #include "G4ParticleDefinition.hh"
42 #include "G4MaterialCutsCouple.hh"
43 #include "G4ProductionCutsTable.hh"
44 #include "G4DynamicParticle.hh"
45 #include "G4Element.hh"
46 #include "G4Gamma.hh"
47 #include "G4Electron.hh"
48 #include "G4Positron.hh"
49 #include "G4PhysicsFreeVector.hh"
50 #include "G4MaterialCutsCouple.hh"
51 #include "G4AutoLock.hh"
52 
53 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
54 
56  const G4String& nam)
57  :G4VEmModel(nam),fParticleChange(0),fParticle(0),
58  logAtomicCrossSection(0),
59  fEffectiveCharge(0),fMaterialInvScreeningRadius(0),
60  fScreeningFunction(0),isInitialised(false),fLocalTable(false)
61 
62 {
63  fIntrinsicLowEnergyLimit = 2.0*electron_mass_c2;
64  fIntrinsicHighEnergyLimit = 100.0*GeV;
65  fSmallEnergy = 1.1*MeV;
66 
67  InitializeScreeningRadii();
68 
69  if (part)
70  SetParticle(part);
71 
72  // SetLowEnergyLimit(fIntrinsicLowEnergyLimit);
73  SetHighEnergyLimit(fIntrinsicHighEnergyLimit);
74  //
75  verboseLevel= 0;
76  // Verbosity scale:
77  // 0 = nothing
78  // 1 = warning for energy non-conservation
79  // 2 = details of energy budget
80  // 3 = calculation of cross sections, file openings, sampling of atoms
81  // 4 = entering in methods
82 }
83 
84 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
85 
87 {
88  //Delete shared tables, they exist only in the master model
89  if (IsMaster() || fLocalTable)
90  {
91  //std::map <G4int,G4PhysicsFreeVector*>::iterator i;
92  if (logAtomicCrossSection)
93  {
94  /*
95  for (i=logAtomicCrossSection->begin();i != logAtomicCrossSection->end();i++)
96  if (i->second) delete i->second;
97  */
98  delete logAtomicCrossSection;
99  }
100  if (fEffectiveCharge)
101  delete fEffectiveCharge;
102  if (fMaterialInvScreeningRadius)
103  delete fMaterialInvScreeningRadius;
104  if (fScreeningFunction)
105  delete fScreeningFunction;
106  }
107 }
108 
109 
110 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
111 
113  const G4DataVector&)
114 {
115  if (verboseLevel > 3)
116  G4cout << "Calling G4PenelopeGammaConversionModel::Initialise()" << G4endl;
117 
118  SetParticle(part);
119 
120  //Only the master model creates/fills/destroys the tables
121  if (IsMaster() && part == fParticle)
122  {
123  // logAtomicCrossSection is created only once, since it is never cleared
124  if (!logAtomicCrossSection)
125  logAtomicCrossSection = new std::map<G4int,G4PhysicsFreeVector*>;
126 
127  //delete old material data...
128  if (fEffectiveCharge)
129  {
130  delete fEffectiveCharge;
131  fEffectiveCharge = 0;
132  }
133  if (fMaterialInvScreeningRadius)
134  {
135  delete fMaterialInvScreeningRadius;
136  fMaterialInvScreeningRadius = 0;
137  }
138  if (fScreeningFunction)
139  {
140  delete fScreeningFunction;
141  fScreeningFunction = 0;
142  }
143  //and create new ones
144  fEffectiveCharge = new std::map<const G4Material*,G4double>;
145  fMaterialInvScreeningRadius = new std::map<const G4Material*,G4double>;
146  fScreeningFunction = new std::map<const G4Material*,std::pair<G4double,G4double> >;
147 
148  G4ProductionCutsTable* theCoupleTable =
150 
151  for (size_t i=0;i<theCoupleTable->GetTableSize();i++)
152  {
153  const G4Material* material =
154  theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
155  const G4ElementVector* theElementVector = material->GetElementVector();
156 
157  for (size_t j=0;j<material->GetNumberOfElements();j++)
158  {
159  G4int iZ = (G4int) theElementVector->at(j)->GetZ();
160  //read data files only in the master
161  if (!logAtomicCrossSection->count(iZ))
162  ReadDataFile(iZ);
163  }
164 
165  //check if material data are available
166  if (!fEffectiveCharge->count(material))
167  InitializeScreeningFunctions(material);
168  }
169 
170 
171  if (verboseLevel > 0) {
172  G4cout << "Penelope Gamma Conversion model v2008 is initialized " << G4endl
173  << "Energy range: "
174  << LowEnergyLimit() / MeV << " MeV - "
175  << HighEnergyLimit() / GeV << " GeV"
176  << G4endl;
177  }
178 
179  }
180 
181 
182  if(isInitialised) return;
184  isInitialised = true;
185 }
186 
187 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
188 
190  G4VEmModel *masterModel)
191 {
192  if (verboseLevel > 3)
193  G4cout << "Calling G4PenelopeGammaConversionModel::InitialiseLocal()" << G4endl;
194 
195  //
196  //Check that particle matches: one might have multiple master models (e.g.
197  //for e+ and e-).
198  //
199  if (part == fParticle)
200  {
201  //Get the const table pointers from the master to the workers
202  const G4PenelopeGammaConversionModel* theModel =
203  static_cast<G4PenelopeGammaConversionModel*> (masterModel);
204 
205  //Copy pointers to the data tables
206  fEffectiveCharge = theModel->fEffectiveCharge;
207  fMaterialInvScreeningRadius = theModel->fMaterialInvScreeningRadius;
208  fScreeningFunction = theModel->fScreeningFunction;
209  logAtomicCrossSection = theModel->logAtomicCrossSection;
210 
211  //Same verbosity for all workers, as the master
212  verboseLevel = theModel->verboseLevel;
213  }
214 
215  return;
216 }
217 
218 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
219 namespace { G4Mutex PenelopeGammaConversionModelMutex = G4MUTEX_INITIALIZER; }
220 
222  const G4ParticleDefinition*,
226 {
227  //
228  // Penelope model v2008.
229  // Cross section (including triplet production) read from database and managed
230  // through the G4CrossSectionHandler utility. Cross section data are from
231  // M.J. Berger and J.H. Hubbel (XCOM), Report NBSIR 887-3598
232  //
233 
234  if (energy < fIntrinsicLowEnergyLimit)
235  return 0;
236 
237  G4int iZ = (G4int) Z;
238 
239  //Either Initialize() was not called, or we are in a slave and InitializeLocal() was
240  //not invoked
241  if (!logAtomicCrossSection)
242  {
243  //create a **thread-local** version of the table. Used only for G4EmCalculator and
244  //Unit Tests
245  fLocalTable = true;
246  logAtomicCrossSection = new std::map<G4int,G4PhysicsFreeVector*>;
247  }
248  //now it should be ok
249  if (!logAtomicCrossSection->count(iZ))
250  {
251  //If we are here, it means that Initialize() was inkoved, but the MaterialTable was
252  //not filled up. This can happen in a UnitTest or via G4EmCalculator
254  ed << "Unable to retrieve the cross section table for Z=" << iZ << G4endl;
255  ed << "This can happen only in Unit Tests or via G4EmCalculator" << G4endl;
256  G4Exception("G4PenelopeGammaConversionModel::ComputeCrossSectionPerAtom()",
257  "em2018",JustWarning,ed);
258  //protect file reading via autolock
259  G4AutoLock lock(&PenelopeGammaConversionModelMutex);
260  ReadDataFile(iZ);
261  lock.unlock();
262  }
263 
264  G4double cs = 0;
265  G4double logene = std::log(energy);
266 
267  G4PhysicsFreeVector* theVec = logAtomicCrossSection->find(iZ)->second;
268 
269  G4double logXS = theVec->Value(logene);
270  cs = std::exp(logXS);
271 
272  if (verboseLevel > 2)
273  G4cout << "Gamma conversion cross section at " << energy/MeV << " MeV for Z=" << Z <<
274  " = " << cs/barn << " barn" << G4endl;
275  return cs;
276 }
277 
278 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
279 
280 void
281 G4PenelopeGammaConversionModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
282  const G4MaterialCutsCouple* couple,
283  const G4DynamicParticle* aDynamicGamma,
284  G4double,
285  G4double)
286 {
287  //
288  // Penelope model v2008.
289  // Final state is sampled according to the Bethe-Heitler model with Coulomb
290  // corrections, according to the semi-empirical model of
291  // J. Baro' et al., Radiat. Phys. Chem. 44 (1994) 531.
292  //
293  // The model uses the high energy Coulomb correction from
294  // H. Davies et al., Phys. Rev. 93 (1954) 788
295  // and atomic screening radii tabulated from
296  // J.H. Hubbel et al., J. Phys. Chem. Ref. Data 9 (1980) 1023
297  // for Z= 1 to 92.
298  //
299  if (verboseLevel > 3)
300  G4cout << "Calling SamplingSecondaries() of G4PenelopeGammaConversionModel" << G4endl;
301 
302  G4double photonEnergy = aDynamicGamma->GetKineticEnergy();
303 
304  // Always kill primary
307 
308  if (photonEnergy <= fIntrinsicLowEnergyLimit)
309  {
311  return ;
312  }
313 
314  G4ParticleMomentum photonDirection = aDynamicGamma->GetMomentumDirection();
315  const G4Material* mat = couple->GetMaterial();
316 
317  if (!fEffectiveCharge->count(mat))
318  {
320  ed << "Unable to allocate the EffectiveCharge data for " <<
321  mat->GetName() << G4endl;
322  G4Exception("G4PenelopeGammaConversion::SampleSecondaries()",
323  "em2019",FatalException,ed);
324  }
325 
326  // eps is the fraction of the photon energy assigned to e- (including rest mass)
327  G4double eps = 0;
328  G4double eki = electron_mass_c2/photonEnergy;
329 
330  //Do it fast for photon energy < 1.1 MeV (close to threshold)
331  if (photonEnergy < fSmallEnergy)
332  eps = eki + (1.0-2.0*eki)*G4UniformRand();
333  else
334  {
335  //Complete calculation
336  G4double effC = fEffectiveCharge->find(mat)->second;
337  G4double alz = effC*fine_structure_const;
338  G4double T = std::sqrt(2.0*eki);
339  G4double F00=(-1.774-1.210e1*alz+1.118e1*alz*alz)*T
340  +(8.523+7.326e1*alz-4.441e1*alz*alz)*T*T
341  -(1.352e1+1.211e2*alz-9.641e1*alz*alz)*T*T*T
342  +(8.946+6.205e1*alz-6.341e1*alz*alz)*T*T*T*T;
343 
344  G4double F0b = fScreeningFunction->find(mat)->second.second;
345  G4double g0 = F0b + F00;
346  G4double invRad = fMaterialInvScreeningRadius->find(mat)->second;
347  G4double bmin = 4.0*eki/invRad;
348  std::pair<G4double,G4double> scree = GetScreeningFunctions(bmin);
349  G4double g1 = scree.first;
350  G4double g2 = scree.second;
351  G4double g1min = g1+g0;
352  G4double g2min = g2+g0;
353  G4double xr = 0.5-eki;
354  G4double a1 = 2.*g1min*xr*xr/3.;
355  G4double p1 = a1/(a1+g2min);
356 
357  G4bool loopAgain = false;
358  //Random sampling of eps
359  do{
360  loopAgain = false;
361  if (G4UniformRand() <= p1)
362  {
363  G4double ru2m1 = 2.0*G4UniformRand()-1.0;
364  if (ru2m1 < 0)
365  eps = 0.5-xr*std::pow(std::abs(ru2m1),1./3.);
366  else
367  eps = 0.5+xr*std::pow(ru2m1,1./3.);
368  G4double B = eki/(invRad*eps*(1.0-eps));
369  scree = GetScreeningFunctions(B);
370  g1 = scree.first;
371  g1 = std::max(g1+g0,0.);
372  if (G4UniformRand()*g1min > g1)
373  loopAgain = true;
374  }
375  else
376  {
377  eps = eki+2.0*xr*G4UniformRand();
378  G4double B = eki/(invRad*eps*(1.0-eps));
379  scree = GetScreeningFunctions(B);
380  g2 = scree.second;
381  g2 = std::max(g2+g0,0.);
382  if (G4UniformRand()*g2min > g2)
383  loopAgain = true;
384  }
385  }while(loopAgain);
386 
387  }
388  if (verboseLevel > 4)
389  G4cout << "Sampled eps = " << eps << G4endl;
390 
391  G4double electronTotEnergy = eps*photonEnergy;
392  G4double positronTotEnergy = (1.0-eps)*photonEnergy;
393 
394  // Scattered electron (positron) angles. ( Z - axis along the parent photon)
395 
396  //electron kinematics
397  G4double electronKineEnergy = std::max(0.,electronTotEnergy - electron_mass_c2) ;
398  G4double costheta_el = G4UniformRand()*2.0-1.0;
399  G4double kk = std::sqrt(electronKineEnergy*(electronKineEnergy+2.*electron_mass_c2));
400  costheta_el = (costheta_el*electronTotEnergy+kk)/(electronTotEnergy+costheta_el*kk);
401  G4double phi_el = twopi * G4UniformRand() ;
402  G4double dirX_el = std::sqrt(1.-costheta_el*costheta_el) * std::cos(phi_el);
403  G4double dirY_el = std::sqrt(1.-costheta_el*costheta_el) * std::sin(phi_el);
404  G4double dirZ_el = costheta_el;
405 
406  //positron kinematics
407  G4double positronKineEnergy = std::max(0.,positronTotEnergy - electron_mass_c2) ;
408  G4double costheta_po = G4UniformRand()*2.0-1.0;
409  kk = std::sqrt(positronKineEnergy*(positronKineEnergy+2.*electron_mass_c2));
410  costheta_po = (costheta_po*positronTotEnergy+kk)/(positronTotEnergy+costheta_po*kk);
411  G4double phi_po = twopi * G4UniformRand() ;
412  G4double dirX_po = std::sqrt(1.-costheta_po*costheta_po) * std::cos(phi_po);
413  G4double dirY_po = std::sqrt(1.-costheta_po*costheta_po) * std::sin(phi_po);
414  G4double dirZ_po = costheta_po;
415 
416  // Kinematics of the created pair:
417  // the electron and positron are assumed to have a symetric angular
418  // distribution with respect to the Z axis along the parent photon
419  G4double localEnergyDeposit = 0. ;
420 
421  if (electronKineEnergy > 0.0)
422  {
423  G4ThreeVector electronDirection ( dirX_el, dirY_el, dirZ_el);
424  electronDirection.rotateUz(photonDirection);
426  electronDirection,
427  electronKineEnergy);
428  fvect->push_back(electron);
429  }
430  else
431  {
432  localEnergyDeposit += electronKineEnergy;
433  electronKineEnergy = 0;
434  }
435 
436  //Generate the positron. Real particle in any case, because it will annihilate. If below
437  //threshold, produce it at rest
438  if (positronKineEnergy < 0.0)
439  {
440  localEnergyDeposit += positronKineEnergy;
441  positronKineEnergy = 0; //produce it at rest
442  }
443  G4ThreeVector positronDirection(dirX_po,dirY_po,dirZ_po);
444  positronDirection.rotateUz(photonDirection);
446  positronDirection, positronKineEnergy);
447  fvect->push_back(positron);
448 
449  //Add rest of energy to the local energy deposit
450  fParticleChange->ProposeLocalEnergyDeposit(localEnergyDeposit);
451 
452  if (verboseLevel > 1)
453  {
454  G4cout << "-----------------------------------------------------------" << G4endl;
455  G4cout << "Energy balance from G4PenelopeGammaConversion" << G4endl;
456  G4cout << "Incoming photon energy: " << photonEnergy/keV << " keV" << G4endl;
457  G4cout << "-----------------------------------------------------------" << G4endl;
458  if (electronKineEnergy)
459  G4cout << "Electron (explicitely produced) " << electronKineEnergy/keV << " keV"
460  << G4endl;
461  if (positronKineEnergy)
462  G4cout << "Positron (not at rest) " << positronKineEnergy/keV << " keV" << G4endl;
463  G4cout << "Rest masses of e+/- " << 2.0*electron_mass_c2/keV << " keV" << G4endl;
464  if (localEnergyDeposit)
465  G4cout << "Local energy deposit " << localEnergyDeposit/keV << " keV" << G4endl;
466  G4cout << "Total final state: " << (electronKineEnergy+positronKineEnergy+
467  localEnergyDeposit+2.0*electron_mass_c2)/keV <<
468  " keV" << G4endl;
469  G4cout << "-----------------------------------------------------------" << G4endl;
470  }
471  if (verboseLevel > 0)
472  {
473  G4double energyDiff = std::fabs(electronKineEnergy+positronKineEnergy+
474  localEnergyDeposit+2.0*electron_mass_c2-photonEnergy);
475  if (energyDiff > 0.05*keV)
476  G4cout << "Warning from G4PenelopeGammaConversion: problem with energy conservation: "
477  << (electronKineEnergy+positronKineEnergy+
478  localEnergyDeposit+2.0*electron_mass_c2)/keV
479  << " keV (final) vs. " << photonEnergy/keV << " keV (initial)" << G4endl;
480  }
481 }
482 
483 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
484 
485 void G4PenelopeGammaConversionModel::ReadDataFile(const G4int Z)
486 {
487  if (!IsMaster())
488  //Should not be here!
489  G4Exception("G4PenelopeGammaConversionModel::ReadDataFile()",
490  "em0100",FatalException,"Worker thread in this method");
491 
492  if (verboseLevel > 2)
493  {
494  G4cout << "G4PenelopeGammaConversionModel::ReadDataFile()" << G4endl;
495  G4cout << "Going to read Gamma Conversion data files for Z=" << Z << G4endl;
496  }
497 
498  char* path = getenv("G4LEDATA");
499  if (!path)
500  {
501  G4String excep =
502  "G4PenelopeGammaConversionModel - G4LEDATA environment variable not set!";
503  G4Exception("G4PenelopeGammaConversionModel::ReadDataFile()",
504  "em0006",FatalException,excep);
505  return;
506  }
507 
508  /*
509  Read the cross section file
510  */
511  std::ostringstream ost;
512  if (Z>9)
513  ost << path << "/penelope/pairproduction/pdgpp" << Z << ".p08";
514  else
515  ost << path << "/penelope/pairproduction/pdgpp0" << Z << ".p08";
516  std::ifstream file(ost.str().c_str());
517  if (!file.is_open())
518  {
519  G4String excep = "G4PenelopeGammaConversionModel - data file " +
520  G4String(ost.str()) + " not found!";
521  G4Exception("G4PenelopeGammaConversionModel::ReadDataFile()",
522  "em0003",FatalException,excep);
523  }
524 
525  //I have to know in advance how many points are in the data list
526  //to initialize the G4PhysicsFreeVector()
527  size_t ndata=0;
528  G4String line;
529  while( getline(file, line) )
530  ndata++;
531  ndata -= 1; //remove one header line
532  //G4cout << "Found: " << ndata << " lines" << G4endl;
533 
534  file.clear();
535  file.close();
536  file.open(ost.str().c_str());
537  G4int readZ =0;
538  file >> readZ;
539 
540  if (verboseLevel > 3)
541  G4cout << "Element Z=" << Z << G4endl;
542 
543  //check the right file is opened.
544  if (readZ != Z)
545  {
547  ed << "Corrupted data file for Z=" << Z << G4endl;
548  G4Exception("G4PenelopeGammaConversionModel::ReadDataFile()",
549  "em0005",FatalException,ed);
550  }
551 
552  G4PhysicsFreeVector* theVec = new G4PhysicsFreeVector(ndata);
553  G4double ene=0,xs=0;
554  for (size_t i=0;i<ndata;i++)
555  {
556  file >> ene >> xs;
557  //dimensional quantities
558  ene *= eV;
559  xs *= barn;
560  if (xs < 1e-40*cm2) //protection against log(0)
561  xs = 1e-40*cm2;
562  theVec->PutValue(i,std::log(ene),std::log(xs));
563  }
564  file.close();
565 
566  if (!logAtomicCrossSection)
567  {
569  ed << "Problem with allocation of logAtomicCrossSection data table " << G4endl;
570  G4Exception("G4PenelopeGammaConversionModel::ReadDataFile()",
571  "em2020",FatalException,ed);
572  delete theVec;
573  return;
574  }
575  logAtomicCrossSection->insert(std::make_pair(Z,theVec));
576 
577  return;
578 
579 }
580 
581 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
582 
583 void G4PenelopeGammaConversionModel::InitializeScreeningRadii()
584 {
585  G4double temp[99] = {1.2281e+02,7.3167e+01,6.9228e+01,6.7301e+01,6.4696e+01,
586  6.1228e+01,5.7524e+01,5.4033e+01,5.0787e+01,4.7851e+01,4.6373e+01,
587  4.5401e+01,4.4503e+01,4.3815e+01,4.3074e+01,4.2321e+01,4.1586e+01,
588  4.0953e+01,4.0524e+01,4.0256e+01,3.9756e+01,3.9144e+01,3.8462e+01,
589  3.7778e+01,3.7174e+01,3.6663e+01,3.5986e+01,3.5317e+01,3.4688e+01,
590  3.4197e+01,3.3786e+01,3.3422e+01,3.3068e+01,3.2740e+01,3.2438e+01,
591  3.2143e+01,3.1884e+01,3.1622e+01,3.1438e+01,3.1142e+01,3.0950e+01,
592  3.0758e+01,3.0561e+01,3.0285e+01,3.0097e+01,2.9832e+01,2.9581e+01,
593  2.9411e+01,2.9247e+01,2.9085e+01,2.8930e+01,2.8721e+01,2.8580e+01,
594  2.8442e+01,2.8312e+01,2.8139e+01,2.7973e+01,2.7819e+01,2.7675e+01,
595  2.7496e+01,2.7285e+01,2.7093e+01,2.6911e+01,2.6705e+01,2.6516e+01,
596  2.6304e+01,2.6108e+01,2.5929e+01,2.5730e+01,2.5577e+01,2.5403e+01,
597  2.5245e+01,2.5100e+01,2.4941e+01,2.4790e+01,2.4655e+01,2.4506e+01,
598  2.4391e+01,2.4262e+01,2.4145e+01,2.4039e+01,2.3922e+01,2.3813e+01,
599  2.3712e+01,2.3621e+01,2.3523e+01,2.3430e+01,2.3331e+01,2.3238e+01,
600  2.3139e+01,2.3048e+01,2.2967e+01,2.2833e+01,2.2694e+01,2.2624e+01,
601  2.2545e+01,2.2446e+01,2.2358e+01,2.2264e+01};
602 
603  //copy temporary vector in class data member
604  for (G4int i=0;i<99;i++)
605  fAtomicScreeningRadius[i] = temp[i];
606 }
607 
608 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
609 
610 void G4PenelopeGammaConversionModel::InitializeScreeningFunctions(const G4Material* material)
611 {
612  if (!IsMaster())
613  //Should not be here!
614  G4Exception("G4PenelopeGammaConversionModel::InitializeScreeningFunctions()",
615  "em01001",FatalException,"Worker thread in this method");
616 
617  // This is subroutine GPPa0 of Penelope
618  //
619  // 1) calculate the effective Z for the purpose
620  //
621  G4double zeff = 0;
622  G4int intZ = 0;
623  G4int nElements = material->GetNumberOfElements();
624  const G4ElementVector* elementVector = material->GetElementVector();
625 
626  //avoid calculations if only one building element!
627  if (nElements == 1)
628  {
629  zeff = (*elementVector)[0]->GetZ();
630  intZ = (G4int) zeff;
631  }
632  else // many elements...let's do the calculation
633  {
634  const G4double* fractionVector = material->GetVecNbOfAtomsPerVolume();
635 
636  G4double atot = 0;
637  for (G4int i=0;i<nElements;i++)
638  {
639  G4double Zelement = (*elementVector)[i]->GetZ();
640  G4double Aelement = (*elementVector)[i]->GetAtomicMassAmu();
641  atot += Aelement*fractionVector[i];
642  zeff += Zelement*Aelement*fractionVector[i]; //average with the number of nuclei
643  }
644  atot /= material->GetTotNbOfAtomsPerVolume();
645  zeff /= (material->GetTotNbOfAtomsPerVolume()*atot);
646 
647  intZ = (G4int) (zeff+0.25);
648  if (intZ <= 0)
649  intZ = 1;
650  if (intZ > 99)
651  intZ = 99;
652  }
653 
654  if (fEffectiveCharge)
655  fEffectiveCharge->insert(std::make_pair(material,zeff));
656 
657  //
658  // 2) Calculate Coulomb Correction
659  //
660  G4double alz = fine_structure_const*zeff;
661  G4double alzSquared = alz*alz;
662  G4double fc = alzSquared*(0.202059-alzSquared*
663  (0.03693-alzSquared*
664  (0.00835-alzSquared*(0.00201-alzSquared*
665  (0.00049-alzSquared*
666  (0.00012-alzSquared*0.00003)))))
667  +1.0/(alzSquared+1.0));
668  //
669  // 3) Screening functions and low-energy corrections
670  //
671  G4double matRadius = 2.0/ fAtomicScreeningRadius[intZ-1];
672  if (fMaterialInvScreeningRadius)
673  fMaterialInvScreeningRadius->insert(std::make_pair(material,matRadius));
674 
675  std::pair<G4double,G4double> myPair(0,0);
676  G4double f0a = 4.0*std::log(fAtomicScreeningRadius[intZ-1]);
677  G4double f0b = f0a - 4.0*fc;
678  myPair.first = f0a;
679  myPair.second = f0b;
680 
681  if (fScreeningFunction)
682  fScreeningFunction->insert(std::make_pair(material,myPair));
683 
684  if (verboseLevel > 2)
685  {
686  G4cout << "Average Z for material " << material->GetName() << " = " <<
687  zeff << G4endl;
688  G4cout << "Effective radius for material " << material->GetName() << " = " <<
689  fAtomicScreeningRadius[intZ-1] << " m_e*c/hbar --> BCB = " <<
690  matRadius << G4endl;
691  G4cout << "Screening parameters F0 for material " << material->GetName() << " = " <<
692  f0a << "," << f0b << G4endl;
693  }
694  return;
695 }
696 
697 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
698 
699 std::pair<G4double,G4double>
700 G4PenelopeGammaConversionModel::GetScreeningFunctions(G4double B)
701 {
702  // This is subroutine SCHIFF of Penelope
703  //
704  // Screening Functions F1(B) and F2(B) in the Bethe-Heitler differential cross
705  // section for pair production
706  //
707  std::pair<G4double,G4double> result(0.,0.);
708  G4double BSquared = B*B;
709  G4double f1 = 2.0-2.0*std::log(1.0+BSquared);
710  G4double f2 = f1 - 6.66666666e-1; // (-2/3)
711  if (B < 1.0e-10)
712  f1 = f1-twopi*B;
713  else
714  {
715  G4double a0 = 4.0*B*std::atan(1./B);
716  f1 = f1 - a0;
717  f2 += 2.0*BSquared*(4.0-a0-3.0*std::log((1.0+BSquared)/BSquared));
718  }
719  G4double g1 = 0.5*(3.0*f1-f2);
720  G4double g2 = 0.25*(3.0*f1+f2);
721 
722  result.first = g1;
723  result.second = g2;
724 
725  return result;
726 }
727 
728 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...
729 
730 void G4PenelopeGammaConversionModel::SetParticle(const G4ParticleDefinition* p)
731 {
732  if(!fParticle) {
733  fParticle = p;
734  }
735 }
736 
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:599
std::vector< G4Element * > G4ElementVector
void PutValue(size_t binNumber, G4double binValue, G4double dataValue)
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4double GetKineticEnergy() const
#define F00
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:592
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX)
const char * p
Definition: xmltok.h:285
const G4String & GetName() const
Definition: G4Material.hh:176
G4PenelopeGammaConversionModel(const G4ParticleDefinition *p=0, const G4String &processName="PenConversion")
G4double G4NeutronHPJENDLHEData::G4double result
G4bool IsMaster() const
Definition: G4VEmModel.hh:676
TFile * file
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
int G4int
Definition: G4Types.hh:78
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:158
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:683
Float_t mat
Definition: plot.C:40
string material
Definition: eplot.py:19
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:204
Float_t f1
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
double energy
Definition: plottest35.C:25
Float_t f2
Float_t Z
Definition: plot.C:39
bool G4bool
Definition: G4Types.hh:79
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
const G4ThreeVector & GetMomentumDirection() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
TString part[npart]
Definition: Style.C:32
float electron_mass_c2
Definition: hepunit.py:274
G4double Value(G4double theEnergy, size_t &lastidx) const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:207
static G4ProductionCutsTable * GetProductionCutsTable()
G4int G4Mutex
Definition: G4Threading.hh:156
static G4Positron * Positron()
Definition: G4Positron.cc:94
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
static G4Electron * Electron()
Definition: G4Electron.cc:94
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:121
const G4Material * GetMaterial() const