ECCE @ EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4eBremParametrizedModel.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file G4eBremParametrizedModel.cc
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: G4eBremParametrizedModel
33 //
34 // Author: Andreas Schaelicke
35 //
36 // Creation date: 06.04.2011
37 //
38 // Modifications:
39 //
40 // Main References:
41 // - based on G4eBremsstrahlungModel and G4eBremsstrahlungRelModel
42 // -------------------------------------------------------------------
43 //
44 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
45 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
46 
48 #include "G4PhysicalConstants.hh"
49 #include "G4SystemOfUnits.hh"
50 #include "G4Electron.hh"
51 #include "G4Positron.hh"
52 #include "G4Gamma.hh"
53 #include "Randomize.hh"
54 #include "G4Material.hh"
55 #include "G4Element.hh"
56 #include "G4ElementVector.hh"
57 #include "G4ProductionCutsTable.hh"
59 #include "G4LossTableManager.hh"
60 #include "G4ModifiedTsai.hh"
61 #include "G4Exp.hh"
62 #include "G4Log.hh"
63 
64 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
65 
66 const G4double G4eBremParametrizedModel::xgi[]={ 0.0199, 0.1017, 0.2372, 0.4083,
67  0.5917, 0.7628, 0.8983, 0.9801 };
68 const G4double G4eBremParametrizedModel::wgi[]={ 0.0506, 0.1112, 0.1569, 0.1813,
69  0.1813, 0.1569, 0.1112, 0.0506 };
70 
71 static const G4double tlow = 1.*CLHEP::MeV;
72 
73 //
74 // GEANT4 internal units.
75 //
76 static const G4double
77  ah10 = 4.67733E+00, ah11 =-6.19012E-01, ah12 = 2.02225E-02,
78  ah20 =-7.34101E+00, ah21 = 1.00462E+00, ah22 =-3.20985E-02,
79  ah30 = 2.93119E+00, ah31 =-4.03761E-01, ah32 = 1.25153E-02;
80 
81 static const G4double
82  bh10 = 4.23071E+00, bh11 =-6.10995E-01, bh12 = 1.95531E-02,
83  bh20 =-7.12527E+00, bh21 = 9.69160E-01, bh22 =-2.74255E-02,
84  bh30 = 2.69925E+00, bh31 =-3.63283E-01, bh32 = 9.55316E-03;
85 
86 static const G4double
87  al00 =-2.05398E+00, al01 = 2.38815E-02, al02 = 5.25483E-04,
88  al10 =-7.69748E-02, al11 =-6.91499E-02, al12 = 2.22453E-03,
89  al20 = 4.06463E-02, al21 =-1.01281E-02, al22 = 3.40919E-04;
90 
91 static const G4double
92  bl00 = 1.04133E+00, bl01 =-9.43291E-03, bl02 =-4.54758E-04,
93  bl10 = 1.19253E-01, bl11 = 4.07467E-02, bl12 =-1.30718E-03,
94  bl20 =-1.59391E-02, bl21 = 7.27752E-03, bl22 =-1.94405E-04;
95 
96 using namespace std;
97 
99  const G4String& nam)
100  : G4VEmModel(nam),
101  particle(0),
102  isElectron(true),
105  isInitialised(false)
106 {
108 
109  minThreshold = 0.1*keV;
110  lowKinEnergy = 10.*MeV;
112 
114 
116 
118  = densityFactor = densityCorr = fMax = fCoulomb = 0.;
119 
121  if(p) { SetParticle(p); }
122 }
123 
124 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
125 
127 {
128  facFel = G4Log(184.15);
129  facFinel = G4Log(1194.);
130 }
131 
132 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
133 
135 {
136 }
137 
138 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
139 
141 {
142  particle = p;
143  particleMass = p->GetPDGMass();
144  if(p == G4Electron::Electron()) { isElectron = true; }
145  else { isElectron = false;}
146 }
147 
148 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
149 
151  const G4MaterialCutsCouple*)
152 {
153  return minThreshold;
154 }
155 
156 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
157 
159  const G4Material* mat,
160  G4double kineticEnergy)
161 {
163 
164  // calculate threshold for density effect
165  kinEnergy = kineticEnergy;
166  totalEnergy = kineticEnergy + particleMass;
168 }
169 
170 
171 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
172 
174  const G4DataVector& cuts)
175 {
176  if(p) { SetParticle(p); }
177 
179 
180  currentZ = 0.;
181 
182  if(IsMaster()) { InitialiseElementSelectors(p, cuts); }
183 
184  if(isInitialised) { return; }
186  isInitialised = true;
187 }
188 
189 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
190 
192  G4VEmModel* masterModel)
193 {
195 }
196 
197 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
198 
200  const G4Material* material,
201  const G4ParticleDefinition* p,
202  G4double kineticEnergy,
203  G4double cutEnergy)
204 {
205  if(!particle) { SetParticle(p); }
206  if(kineticEnergy < lowKinEnergy) { return 0.0; }
207  G4double cut = std::min(cutEnergy, kineticEnergy);
208  if(cut == 0.0) { return 0.0; }
209 
210  SetupForMaterial(particle, material,kineticEnergy);
211 
212  const G4ElementVector* theElementVector = material->GetElementVector();
213  const G4double* theAtomicNumDensityVector = material->GetAtomicNumDensityVector();
214 
215  G4double dedx = 0.0;
216 
217  // loop for elements in the material
218  for (size_t i=0; i<material->GetNumberOfElements(); i++) {
219 
220  G4VEmModel::SetCurrentElement((*theElementVector)[i]);
221  SetCurrentElement((*theElementVector)[i]->GetZ());
222 
223  dedx += theAtomicNumDensityVector[i]*currentZ*currentZ*ComputeBremLoss(cut);
224  }
225  dedx *= bremFactor;
226 
227  return dedx;
228 }
229 
230 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
231 
233 {
234  G4double loss = 0.0;
235 
236  // number of intervals and integration step
237  G4double vcut = cut/totalEnergy;
238  G4int n = (G4int)(20*vcut) + 3;
239  G4double delta = vcut/G4double(n);
240 
241  G4double e0 = 0.0;
242  G4double xs;
243 
244  // integration
245  for(G4int l=0; l<n; l++) {
246 
247  for(G4int i=0; i<8; i++) {
248 
249  G4double eg = (e0 + xgi[i]*delta)*totalEnergy;
250 
251  xs = ComputeDXSectionPerAtom(eg);
252 
253  loss += wgi[i]*xs/(1.0 + densityCorr/(eg*eg));
254  }
255  e0 += delta;
256  }
257 
258  loss *= delta*totalEnergy;
259 
260  return loss;
261 }
262 
263 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
264 
266  const G4ParticleDefinition* p,
267  G4double kineticEnergy,
269  G4double cutEnergy,
270  G4double maxEnergy)
271 {
272  if(!particle) { SetParticle(p); }
273  if(kineticEnergy < lowKinEnergy) { return 0.0; }
274  G4double cut = std::min(cutEnergy, kineticEnergy);
275  G4double tmax = std::min(maxEnergy, kineticEnergy);
276 
277  if(cut >= tmax) { return 0.0; }
278 
280 
282 
283  // allow partial integration
284  if(tmax < kinEnergy) { cross -= ComputeXSectionPerAtom(tmax); }
285 
286  cross *= Z*Z*bremFactor;
287 
288  return cross;
289 }
290 
291 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
292 
294 {
295  G4double cross = 0.0;
296 
297  // number of intervals and integration step
298  G4double vcut = G4Log(cut/totalEnergy);
300  G4int n = (G4int)(0.45*(vmax - vcut)) + 4;
301  // n=1; // integration test
302  G4double delta = (vmax - vcut)/G4double(n);
303 
304  G4double e0 = vcut;
305  G4double xs;
306 
307  // integration
308  for(G4int l=0; l<n; l++) {
309 
310  for(G4int i=0; i<8; i++) {
311 
312  G4double eg = G4Exp(e0 + xgi[i]*delta)*totalEnergy;
313 
314  xs = ComputeDXSectionPerAtom(eg);
315 
316  cross += wgi[i]*xs/(1.0 + densityCorr/(eg*eg));
317  }
318  e0 += delta;
319  }
320 
321  cross *= delta;
322 
323  return cross;
324 }
325 
326 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
327 
328 // compute the value of the screening function 3*PHI1 - PHI2
329 
331 {
332  G4double screenVal;
333 
334  if (ScreenVariable > 1.)
335  screenVal = 42.24 - 8.368*G4Log(ScreenVariable+0.952);
336  else
337  screenVal = 42.392 - ScreenVariable* (7.796 - 1.961*ScreenVariable);
338 
339  return screenVal;
340 }
341 
342 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
343 
344 // compute the value of the screening function 1.5*PHI1 - 0.5*PHI2
345 
347 {
348  G4double screenVal;
349 
350  if (ScreenVariable > 1.)
351  screenVal = 42.24 - 8.368*G4Log(ScreenVariable+0.952);
352  else
353  screenVal = 41.734 - ScreenVariable* (6.484 - 1.250*ScreenVariable);
354 
355  return screenVal;
356 }
357 
358 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
359 
360 // Parametrized cross section
362  G4double kineticEnergy,
363  G4double gammaEnergy, G4double Z)
364 {
366  G4double FZ = lnZ* (4.- 0.55*lnZ);
367  G4double Z3 = z13;
368  G4double ZZ = z13*nist->GetZ13(G4lrint(Z)+1);
369 
370  totalEnergy = kineticEnergy + electron_mass_c2;
371 
372  // G4double x, epsil, greject, migdal, grejmax, q;
373  G4double epsil, greject;
374  G4double U = G4Log(kineticEnergy/electron_mass_c2);
375  G4double U2 = U*U;
376 
377  // precalculated parameters
378  G4double ah, bh;
379 
380  if (kineticEnergy > tlow) {
381 
382  G4double ah1 = ah10 + ZZ* (ah11 + ZZ* ah12);
383  G4double ah2 = ah20 + ZZ* (ah21 + ZZ* ah22);
384  G4double ah3 = ah30 + ZZ* (ah31 + ZZ* ah32);
385 
386  G4double bh1 = bh10 + ZZ* (bh11 + ZZ* bh12);
387  G4double bh2 = bh20 + ZZ* (bh21 + ZZ* bh22);
388  G4double bh3 = bh30 + ZZ* (bh31 + ZZ* bh32);
389 
390  ah = 1. + (ah1*U2 + ah2*U + ah3) / (U2*U);
391  bh = 0.75 + (bh1*U2 + bh2*U + bh3) / (U2*U);
392 
393  // limit of the screening variable
394  G4double screenfac =
395  136.*electron_mass_c2/(Z3*totalEnergy);
396 
397  epsil = gammaEnergy/totalEnergy; // epsil = x*kineticEnergy/totalEnergy;
398  G4double screenvar = screenfac*epsil/(1.0-epsil);
399  G4double F1 = max(ScreenFunction1(screenvar) - FZ ,0.);
400  G4double F2 = max(ScreenFunction2(screenvar) - FZ ,0.);
401 
402 
403  greject = (F1 - epsil* (ah*F1 - bh*epsil*F2))/8.; // 1./(42.392 - FZ);
404 
405  std::cout << " yy = "<<epsil<<std::endl;
406  std::cout << " F1/(...) "<<F1/(42.392 - FZ)<<std::endl;
407  std::cout << " F2/(...) "<<F2/(42.392 - FZ)<<std::endl;
408  std::cout << " (42.392 - FZ) " << (42.392 - FZ) <<std::endl;
409 
410  } else {
411 
412  G4double al0 = al00 + ZZ* (al01 + ZZ* al02);
413  G4double al1 = al10 + ZZ* (al11 + ZZ* al12);
414  G4double al2 = al20 + ZZ* (al21 + ZZ* al22);
415 
416  G4double bl0 = bl00 + ZZ* (bl01 + ZZ* bl02);
417  G4double bl1 = bl10 + ZZ* (bl11 + ZZ* bl12);
418  G4double bl2 = bl20 + ZZ* (bl21 + ZZ* bl22);
419 
420  ah = al0 + al1*U + al2*U2;
421  bh = bl0 + bl1*U + bl2*U2;
422 
423  G4double x=gammaEnergy/kineticEnergy;
424  greject=(1. + x* (ah + bh*x));
425 
426  /*
427  // Compute the maximum of the rejection function
428  grejmax = max(1. + xmin* (ah + bh*xmin), 1.+ah+bh);
429  G4double xm = -ah/(2.*bh);
430  if ( xmin < xm && xm < xmax) grejmax = max(grejmax, 1.+ xm* (ah + bh*xm));
431  */
432  }
433 
434  return greject;
435 }
436 
437 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
438 
440 {
441 
442  if(gammaEnergy < 0.0) { return 0.0; }
443 
444  G4double y = gammaEnergy/totalEnergy;
445 
446  G4double main=0.;
447  //secondTerm=0.;
448 
449  // ** form factors complete screening case **
450  // only valid for high energies (and if LPM suppression does not play a role)
451  main = (3./4.*y*y - y + 1.) * ( (Fel-fCoulomb) + Finel/currentZ );
452  // secondTerm = (1.-y)/12.*(1.+1./currentZ);
453 
454  std::cout<<" F1(0) "<<ScreenFunction1(0.) <<std::endl;
455  std::cout<<" F1(0) "<<ScreenFunction2(0.) <<std::endl;
456  std::cout<<"Ekin = "<<kinEnergy<<std::endl;
457  std::cout<<"Z = "<<currentZ<<std::endl;
458  std::cout<<"main = "<<main<<std::endl;
459  std::cout<<" y = "<<y<<std::endl;
460  std::cout<<" Fel-fCoulomb "<< (Fel-fCoulomb) <<std::endl;
461 
463  std::cout<<"main2 = "<<main2<<std::endl;
464  std::cout<<"main2tot = "<<main2 * ( (Fel-fCoulomb) + Finel/currentZ )/(Fel-fCoulomb);
465 
466  G4double cross = main2; //main+secondTerm;
467  return cross;
468 }
469 
470 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
471 
473  std::vector<G4DynamicParticle*>* vdp,
474  const G4MaterialCutsCouple* couple,
475  const G4DynamicParticle* dp,
476  G4double cutEnergy,
477  G4double maxEnergy)
478 {
479  G4double kineticEnergy = dp->GetKineticEnergy();
480  if(kineticEnergy < lowKinEnergy) { return; }
481  G4double cut = std::min(cutEnergy, kineticEnergy);
482  G4double emax = std::min(maxEnergy, kineticEnergy);
483  if(cut >= emax) { return; }
484 
485  SetupForMaterial(particle, couple->GetMaterial(),kineticEnergy);
486 
487  const G4Element* elm = SelectTargetAtom(couple,particle,kineticEnergy,
488  dp->GetLogKineticEnergy(),cut,emax);
489  SetCurrentElement(elm->GetZ());
490 
491  kinEnergy = kineticEnergy;
492  totalEnergy = kineticEnergy + particleMass;
494 
495  G4double xmin = G4Log(cut*cut + densityCorr);
496  G4double xmax = G4Log(emax*emax + densityCorr);
497  G4double gammaEnergy, f, x;
498 
499  CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
500 
501  do {
502  x = G4Exp(xmin + rndmEngine->flat()*(xmax - xmin)) - densityCorr;
503  if(x < 0.0) x = 0.0;
504  gammaEnergy = sqrt(x);
505  f = ComputeDXSectionPerAtom(gammaEnergy);
506 
507  if ( f > fMax ) {
508  G4cout << "### G4eBremParametrizedModel Warning: Majoranta exceeded! "
509  << f << " > " << fMax
510  << " Egamma(MeV)= " << gammaEnergy
511  << " E(mEV)= " << kineticEnergy
512  << G4endl;
513  }
514 
515  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
516  } while (f < fMax*rndmEngine->flat());
517 
518  //
519  // angles of the emitted gamma. ( Z - axis along the parent particle)
520  // use general interface
521  //
522  G4ThreeVector gammaDirection =
523  GetAngularDistribution()->SampleDirection(dp, totalEnergy-gammaEnergy,
524  G4lrint(currentZ),
525  couple->GetMaterial());
526 
527  // create G4DynamicParticle object for the Gamma
528  G4DynamicParticle* gamma = new G4DynamicParticle(theGamma,gammaDirection,
529  gammaEnergy);
530  vdp->push_back(gamma);
531 
532  G4double totMomentum = sqrt(kineticEnergy*(totalEnergy + electron_mass_c2));
533  G4ThreeVector direction = (totMomentum*dp->GetMomentumDirection()
534  - gammaEnergy*gammaDirection).unit();
535 
536  // energy of primary
537  G4double finalE = kineticEnergy - gammaEnergy;
538 
539  // stop tracking and create new secondary instead of primary
540  if(gammaEnergy > SecondaryThreshold()) {
543  G4DynamicParticle* el =
544  new G4DynamicParticle(const_cast<G4ParticleDefinition*>(particle),
545  direction, finalE);
546  vdp->push_back(el);
547 
548  // continue tracking
549  } else {
552  }
553 }
554 
555 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
556 
557