Geant4  10.02
G4EmCalculator.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: G4EmCalculator.cc 93291 2015-10-15 10:04:27Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4EmCalculator
34 //
35 // Author: Vladimir Ivanchenko
36 //
37 // Creation date: 28.06.2004
38 //
39 // Modifications:
40 // 12.09.2004 Add verbosity (V.Ivanchenko)
41 // 17.11.2004 Change signature of methods, add new methods (V.Ivanchenko)
42 // 08.04.2005 Major optimisation of internal interfaces (V.Ivantchenko)
43 // 08.05.2005 Use updated interfaces (V.Ivantchenko)
44 // 23.10.2005 Fix computations for ions (V.Ivantchenko)
45 // 11.01.2006 Add GetCSDARange (V.Ivantchenko)
46 // 26.01.2006 Rename GetRange -> GetRangeFromRestricteDEDX (V.Ivanchenko)
47 // 14.03.2006 correction in GetCrossSectionPerVolume (mma)
48 // suppress GetCrossSectionPerAtom
49 // elm->GetA() in ComputeCrossSectionPerAtom
50 // 22.03.2006 Add ComputeElectronicDEDX and ComputeTotalDEDX (V.Ivanchenko)
51 // 13.05.2006 Add Corrections for ion stopping (V.Ivanchenko)
52 // 29.09.2006 Uncomment computation of smoothing factor (V.Ivanchenko)
53 // 27.10.2006 Change test energy to access lowEnergy model from
54 // 10 keV to 1 keV (V. Ivanchenko)
55 // 15.03.2007 Add ComputeEnergyCutFromRangeCut methods (V.Ivanchenko)
56 // 21.04.2008 Updated computations for ions (V.Ivanchenko)
57 //
58 // Class Description:
59 //
60 // -------------------------------------------------------------------
61 //
62 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
63 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
64 
65 #include "G4EmCalculator.hh"
66 #include "G4SystemOfUnits.hh"
67 #include "G4LossTableManager.hh"
68 #include "G4NistManager.hh"
69 #include "G4VEmProcess.hh"
70 #include "G4VEnergyLossProcess.hh"
71 #include "G4VMultipleScattering.hh"
72 #include "G4Material.hh"
73 #include "G4MaterialCutsCouple.hh"
74 #include "G4ParticleDefinition.hh"
75 #include "G4ParticleTable.hh"
76 #include "G4IonTable.hh"
77 #include "G4PhysicsTable.hh"
78 #include "G4ProductionCutsTable.hh"
79 #include "G4ProcessManager.hh"
80 #include "G4ionEffectiveCharge.hh"
81 #include "G4RegionStore.hh"
82 #include "G4Element.hh"
83 #include "G4EmCorrections.hh"
84 #include "G4GenericIon.hh"
85 #include "G4ProcessVector.hh"
86 #include "G4Gamma.hh"
87 #include "G4Electron.hh"
88 #include "G4Positron.hh"
89 
90 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
91 
93 {
97  nLocalMaterials = 0;
98  verbose = 0;
100  currentCouple = nullptr;
101  currentMaterial = cutMaterial = nullptr;
102  currentParticle = nullptr;
103  lambdaParticle = nullptr;
104  baseParticle = nullptr;
105  currentLambda = nullptr;
106  currentModel = nullptr;
107  currentProcess = nullptr;
108  loweModel = nullptr;
109  chargeSquare = 1.0;
110  massRatio = 1.0;
111  mass = 0.0;
112  currentCut = 0.0;
113  cutenergy[0] = cutenergy[1] = cutenergy[2] = DBL_MAX;
116  currentName = "";
117  lambdaName = "";
121  isIon = false;
122  isApplicable = false;
123 }
124 
125 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
126 
128 {
129  delete ionEffCharge;
130  for (G4int i=0; i<nLocalMaterials; ++i) {
131  delete localCouples[i];
132  }
133 }
134 
135 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
136 
138  const G4ParticleDefinition* p,
139  const G4Material* mat,
140  const G4Region* region)
141 {
142  G4double res = 0.0;
143  const G4MaterialCutsCouple* couple = FindCouple(mat, region);
144  if(couple && UpdateParticle(p, kinEnergy) ) {
145  res = manager->GetDEDX(p, kinEnergy, couple);
146 
147  if(isIon) {
148  if(FindEmModel(p, currentProcessName, kinEnergy)) {
149  G4double length = CLHEP::nm;
150  G4double eloss = res*length;
151  //G4cout << "### GetDEDX: E= " << kinEnergy << " dedx0= " << res
152  // << " de= " << eloss << G4endl;;
153  G4double niel = 0.0;
154  dynParticle.SetKineticEnergy(kinEnergy);
155  currentModel->GetChargeSquareRatio(p, mat, kinEnergy);
156  currentModel->CorrectionsAlongStep(couple,&dynParticle,eloss,niel,length);
157  res = eloss/length;
158  //G4cout << " de1= " << eloss << " res1= " << res
159  // << " " << p->GetParticleName() <<G4endl;;
160  }
161  }
162 
163  if(verbose>0) {
164  G4cout << "G4EmCalculator::GetDEDX: E(MeV)= " << kinEnergy/MeV
165  << " DEDX(MeV/mm)= " << res*mm/MeV
166  << " DEDX(MeV*cm^2/g)= " << res*gram/(MeV*cm2*mat->GetDensity())
167  << " " << p->GetParticleName()
168  << " in " << mat->GetName()
169  << " isIon= " << isIon
170  << G4endl;
171  }
172  }
173  return res;
174 }
175 
176 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
177 
179  const G4ParticleDefinition* p,
180  const G4Material* mat,
181  const G4Region* region)
182 {
183  G4double res = 0.0;
184  const G4MaterialCutsCouple* couple = FindCouple(mat,region);
185  if(couple && UpdateParticle(p, kinEnergy)) {
186  res = manager->GetRangeFromRestricteDEDX(p, kinEnergy, couple);
187  if(verbose>0) {
188  G4cout << "G4EmCalculator::GetRange: E(MeV)= " << kinEnergy/MeV
189  << " range(mm)= " << res/mm
190  << " " << p->GetParticleName()
191  << " in " << mat->GetName()
192  << G4endl;
193  }
194  }
195  return res;
196 }
197 
198 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
199 
201  const G4ParticleDefinition* p,
202  const G4Material* mat,
203  const G4Region* region)
204 {
205  G4double res = 0.0;
206  if(!G4LossTableManager::Instance()->BuildCSDARange()) {
208  ed << "G4EmCalculator::GetCSDARange: CSDA table is not built; "
209  << " use UI command: /process/eLoss/CSDARange true";
210  G4Exception("G4EmCalculator::GetCSDARange", "em0077",
211  JustWarning, ed);
212  return res;
213  }
214 
215  const G4MaterialCutsCouple* couple = FindCouple(mat,region);
216  if(couple && UpdateParticle(p, kinEnergy)) {
217  res = manager->GetCSDARange(p, kinEnergy, couple);
218  if(verbose>0) {
219  G4cout << "G4EmCalculator::GetRange: E(MeV)= " << kinEnergy/MeV
220  << " range(mm)= " << res/mm
221  << " " << p->GetParticleName()
222  << " in " << mat->GetName()
223  << G4endl;
224  }
225  }
226  return res;
227 }
228 
229 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
230 
232  const G4ParticleDefinition* p,
233  const G4Material* mat,
234  const G4Region* region)
235 {
236  G4double res = 0.0;
237  const G4MaterialCutsCouple* couple = FindCouple(mat,region);
238  if(couple && UpdateParticle(p, kinEnergy)) {
239  res = manager->GetRange(p, kinEnergy, couple);
240  if(verbose>0) {
241  G4cout << "G4EmCalculator::GetRange: E(MeV)= " << kinEnergy/MeV
242  << " range(mm)= " << res/mm
243  << " " << p->GetParticleName()
244  << " in " << mat->GetName()
245  << G4endl;
246  }
247  }
248  return res;
249 }
250 
251 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
252 
254  const G4ParticleDefinition* p,
255  const G4Material* mat,
256  const G4Region* region)
257 {
258  G4double res = 0.0;
259  const G4MaterialCutsCouple* couple = FindCouple(mat,region);
260  if(couple && UpdateParticle(p, 1.0*GeV)) {
261  res = manager->GetEnergy(p, range, couple);
262  if(verbose>0) {
263  G4cout << "G4EmCalculator::GetKinEnergy: Range(mm)= " << range/mm
264  << " KinE(MeV)= " << res/MeV
265  << " " << p->GetParticleName()
266  << " in " << mat->GetName()
267  << G4endl;
268  }
269  }
270  return res;
271 }
272 
273 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
274 
276  const G4ParticleDefinition* p,
277  const G4String& processName,
278  const G4Material* mat,
279  const G4Region* region)
280 {
281  G4double res = 0.0;
282  const G4MaterialCutsCouple* couple = FindCouple(mat,region);
283 
284  if(couple && UpdateParticle(p, kinEnergy)) {
285  G4int idx = couple->GetIndex();
286  FindLambdaTable(p, processName, kinEnergy);
287 
288  if(currentLambda) {
289  G4double e = kinEnergy*massRatio;
290  res = (((*currentLambda)[idx])->Value(e))*chargeSquare;
291  } else {
292  res = ComputeCrossSectionPerVolume(kinEnergy, p, processName, mat,
293  kinEnergy);
294  }
295  if(verbose>0) {
296  G4cout << "G4EmCalculator::GetXSPerVolume: E(MeV)= " << kinEnergy/MeV
297  << " cross(cm-1)= " << res*cm
298  << " " << p->GetParticleName()
299  << " in " << mat->GetName();
300  if(verbose>1)
301  G4cout << " idx= " << idx << " Escaled((MeV)= "
302  << kinEnergy*massRatio
303  << " q2= " << chargeSquare;
304  G4cout << G4endl;
305  }
306  }
307  return res;
308 }
309 
310 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
311 
313  const G4String& particle,
314  G4int Z,
316  G4double kinEnergy)
317 {
318  G4double res = 0.0;
319  const G4ParticleDefinition* p = FindParticle(particle);
321  if(p && ad) {
322  res = ad->GetShellIonisationCrossSectionPerAtom(p, Z, shell, kinEnergy);
323  }
324  return res;
325 }
326 
327 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
328 
330  const G4ParticleDefinition* p,
331  const G4String& processName,
332  const G4Material* mat,
333  const G4Region* region)
334 {
335  G4double res = DBL_MAX;
336  G4double x = GetCrossSectionPerVolume(kinEnergy,p, processName, mat,region);
337  if(x > 0.0) { res = 1.0/x; }
338  if(verbose>1) {
339  G4cout << "G4EmCalculator::GetMeanFreePath: E(MeV)= " << kinEnergy/MeV
340  << " MFP(mm)= " << res/mm
341  << " " << p->GetParticleName()
342  << " in " << mat->GetName()
343  << G4endl;
344  }
345  return res;
346 }
347 
348 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
349 
351 {
353  G4cout << "##### DEDX Table for " << p->GetParticleName() << G4endl;
354  if(elp) G4cout << *(elp->DEDXTable()) << G4endl;
355 }
356 
357 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
358 
360 {
362  G4cout << "##### Range Table for " << p->GetParticleName() << G4endl;
363  if(elp) G4cout << *(elp->RangeTableForLoss()) << G4endl;
364 }
365 
366 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
367 
369 {
371  G4cout << "### G4EmCalculator: Inverse Range Table for "
372  << p->GetParticleName() << G4endl;
373  if(elp) G4cout << *(elp->InverseRangeTable()) << G4endl;
374 }
375 
376 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
377 
379  const G4ParticleDefinition* p,
380  const G4String& processName,
381  const G4Material* mat,
382  G4double cut)
383 {
384  SetupMaterial(mat);
385  G4double res = 0.0;
386  if(verbose > 1) {
387  G4cout << "### G4EmCalculator::ComputeDEDX: " << p->GetParticleName()
388  << " in " << currentMaterialName
389  << " e(MeV)= " << kinEnergy/MeV << " cut(MeV)= " << cut/MeV
390  << G4endl;
391  }
392  if(UpdateParticle(p, kinEnergy)) {
393  if(FindEmModel(p, processName, kinEnergy)) {
394 
395  // Special case of ICRU'73 model
396  if(currentModel->GetName() == "ParamICRU73") {
397  res = currentModel->ComputeDEDXPerVolume(mat, p, kinEnergy, cut);
398  if(verbose > 1) {
399  G4cout << " ICRU73 ion E(MeV)= " << kinEnergy << " ";
400  G4cout << currentModel->GetName() << ": DEDX(MeV/mm)= " << res*mm/MeV
401  << " DEDX(MeV*cm^2/g)= "
402  << res*gram/(MeV*cm2*mat->GetDensity())
403  << G4endl;
404  }
405  } else {
406 
407  G4double escaled = kinEnergy*massRatio;
408  if(baseParticle) {
410  mat, baseParticle, escaled, cut) * chargeSquare;
411  if(verbose > 1) {
413  << " Escaled(MeV)= " << escaled;
414  }
415  } else {
416  res = currentModel->ComputeDEDXPerVolume(mat, p, kinEnergy, cut);
417  if(verbose > 1) { G4cout << " no basePart E(MeV)= " << kinEnergy << " "; }
418  }
419  if(verbose > 1) {
420  G4cout << currentModel->GetName() << ": DEDX(MeV/mm)= " << res*mm/MeV
421  << " DEDX(MeV*cm^2/g)= "
422  << res*gram/(MeV*cm2*mat->GetDensity())
423  << G4endl;
424  }
425 
426  // emulate smoothing procedure
428  // G4cout << "massRatio= " << massRatio << " eth= " << eth << G4endl;
429  if(loweModel) {
430  G4double res0 = 0.0;
431  G4double res1 = 0.0;
432  if(baseParticle) {
433  res1 = currentModel->ComputeDEDXPerVolume(mat, baseParticle, eth, cut)
434  * chargeSquare;
435  res0 = loweModel->ComputeDEDXPerVolume(mat, baseParticle, eth, cut)
436  * chargeSquare;
437  } else {
438  res1 = currentModel->ComputeDEDXPerVolume(mat, p, eth, cut);
439  res0 = loweModel->ComputeDEDXPerVolume(mat, p, eth, cut);
440  }
441  if(verbose > 1) {
442  G4cout << "At boundary energy(MeV)= " << eth/MeV
443  << " DEDX(MeV/mm)= " << res1*mm/MeV
444  << G4endl;
445  }
446 
447  //G4cout << "eth= " << eth << " escaled= " << escaled
448  // << " res0= " << res0 << " res1= "
449  // << res1 << " q2= " << chargeSquare << G4endl;
450 
451  if(res1 > 0.0 && escaled > 0.0) {
452  res *= (1.0 + (res0/res1 - 1.0)*eth/escaled);
453  }
454  }
455 
456  // low energy correction for ions
457  if(isIon) {
458  G4double length = CLHEP::nm;
459  const G4Region* r = 0;
460  const G4MaterialCutsCouple* couple = FindCouple(mat, r);
461  G4double eloss = res*length;
462  G4double niel = 0.0;
463  dynParticle.SetKineticEnergy(kinEnergy);
464  currentModel->GetChargeSquareRatio(p, mat, kinEnergy);
465  currentModel->CorrectionsAlongStep(couple,&dynParticle,eloss,niel,length);
466  res = eloss/length;
467 
468  if(verbose > 1) {
469  G4cout << "After Corrections: DEDX(MeV/mm)= " << res*mm/MeV
470  << " DEDX(MeV*cm^2/g)= " << res*gram/(MeV*cm2*mat->GetDensity())
471  << G4endl;
472  }
473  }
474  }
475  }
476  if(verbose > 0) {
477  G4cout << "Sum: E(MeV)= " << kinEnergy/MeV
478  << " DEDX(MeV/mm)= " << res*mm/MeV
479  << " DEDX(MeV*cm^2/g)= " << res*gram/(MeV*cm2*mat->GetDensity())
480  << " cut(MeV)= " << cut/MeV
481  << " " << p->GetParticleName()
482  << " in " << currentMaterialName
483  << " Zi^2= " << chargeSquare
484  << " isIon=" << isIon
485  << G4endl;
486  }
487  }
488  return res;
489 }
490 
491 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
492 
494  const G4ParticleDefinition* part,
495  const G4Material* mat,
496  G4double cut)
497 {
498  SetupMaterial(mat);
499  G4double dedx = 0.0;
500  if(UpdateParticle(part, kinEnergy)) {
501 
503  const std::vector<G4VEnergyLossProcess*> vel =
504  lManager->GetEnergyLossProcessVector();
505  G4int n = vel.size();
506 
507  //G4cout << "ComputeElectronicDEDX for " << part->GetParticleName()
508  // << " n= " << n << G4endl;
509 
510  for(G4int i=0; i<n; ++i) {
511  if(vel[i]) {
512  G4VProcess* p = reinterpret_cast<G4VProcess*>(vel[i]);
513  if(ActiveForParticle(part, p)) {
514  //G4cout << "idx= " << i << " " << (vel[i])->GetProcessName()
515  // << " " << (vel[i])->Particle()->GetParticleName() << G4endl;
516  dedx += ComputeDEDX(kinEnergy,part,(vel[i])->GetProcessName(),mat,cut);
517  }
518  }
519  }
520  }
521  return dedx;
522 }
523 
524 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
525 
527  const G4ParticleDefinition* part,
528  const G4Material* mat,
529  G4double rangecut)
530 {
531  SetupMaterial(mat);
532  G4double dedx = 0.0;
533  if(UpdateParticle(part, kinEnergy)) {
534 
536  const std::vector<G4VEnergyLossProcess*> vel =
537  lManager->GetEnergyLossProcessVector();
538  G4int n = vel.size();
539 
540  if(mat != cutMaterial) {
541  cutMaterial = mat;
545  }
546 
547  //G4cout << "ComputeElectronicDEDX for " << part->GetParticleName()
548  // << " n= " << n << G4endl;
549 
550  for(G4int i=0; i<n; ++i) {
551  if(vel[i]) {
552  G4VProcess* p = reinterpret_cast<G4VProcess*>(vel[i]);
553  if(ActiveForParticle(part, p)) {
554  //G4cout << "idx= " << i << " " << (vel[i])->GetProcessName()
555  // << " " << (vel[i])->Particle()->GetParticleName() << G4endl;
556  const G4ParticleDefinition* sec = (vel[i])->SecondaryParticle();
557  G4int idx = 0;
558  if(sec == G4Electron::Electron()) { idx = 1; }
559  else if(sec == G4Positron::Positron()) { idx = 2; }
560 
561  dedx += ComputeDEDX(kinEnergy,part,(vel[i])->GetProcessName(),
562  mat,cutenergy[idx]);
563  }
564  }
565  }
566  }
567  return dedx;
568 }
569 
570 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
571 
573  const G4ParticleDefinition* part,
574  const G4Material* mat,
575  G4double cut)
576 {
577  G4double dedx = ComputeElectronicDEDX(kinEnergy,part,mat,cut);
578  if(mass > 700.*MeV) { dedx += ComputeNuclearDEDX(kinEnergy,part,mat); }
579  return dedx;
580 }
581 
582 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
583 
585  const G4ParticleDefinition* p,
586  const G4Material* mat)
587 {
588 
589  G4double res = corr->NuclearDEDX(p, mat, kinEnergy, false);
590 
591  if(verbose > 1) {
592  G4cout << p->GetParticleName() << " E(MeV)= " << kinEnergy/MeV
593  << " NuclearDEDX(MeV/mm)= " << res*mm/MeV
594  << " NuclearDEDX(MeV*cm^2/g)= "
595  << res*gram/(MeV*cm2*mat->GetDensity())
596  << G4endl;
597  }
598  return res;
599 }
600 
601 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
602 
604  G4double kinEnergy,
605  const G4ParticleDefinition* p,
606  const G4String& processName,
607  const G4Material* mat,
608  G4double cut)
609 {
610  SetupMaterial(mat);
611  G4double res = 0.0;
612  if(UpdateParticle(p, kinEnergy)) {
613  if(FindEmModel(p, processName, kinEnergy)) {
614  G4double e = kinEnergy;
615  if(baseParticle) {
616  e *= kinEnergy*massRatio;
618  mat, baseParticle, e, cut, e) * chargeSquare;
619  } else {
620  res = currentModel->CrossSectionPerVolume(mat, p, e, cut, e);
621  }
622  if(verbose>0) {
623  G4cout << "G4EmCalculator::ComputeXSPerVolume: E(MeV)= " << kinEnergy/MeV
624  << " cross(cm-1)= " << res*cm
625  << " cut(keV)= " << cut/keV
626  << " " << p->GetParticleName()
627  << " in " << mat->GetName()
628  << G4endl;
629  }
630  }
631  }
632  return res;
633 }
634 
635 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
636 
638  G4double kinEnergy,
639  const G4ParticleDefinition* p,
640  const G4String& processName,
641  G4double Z, G4double A,
642  G4double cut)
643 {
644  G4double res = 0.0;
645  if(UpdateParticle(p, kinEnergy)) {
647  if(FindEmModel(p, processName, kinEnergy)) {
648  G4double e = kinEnergy;
649  if(baseParticle) {
650  e *= kinEnergy*massRatio;
653  baseParticle, e, Z, A, cut) * chargeSquare;
654  } else {
656  res = currentModel->ComputeCrossSectionPerAtom(p, e, Z, A, cut);
657  }
658  if(verbose>0) {
659  G4cout << "E(MeV)= " << kinEnergy/MeV
660  << " cross(barn)= " << res/barn
661  << " " << p->GetParticleName()
662  << " Z= " << Z << " A= " << A/(g/mole) << " g/mole"
663  << G4endl;
664  }
665  }
666  }
667  return res;
668 }
669 
670 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
671 
673  const G4ParticleDefinition* p,
674  const G4String& processName,
675  G4int Z, G4int shellIdx,
676  G4double cut)
677 {
678  G4double res = 0.0;
679  if(UpdateParticle(p, kinEnergy)) {
680  CheckMaterial(Z);
681  if(FindEmModel(p, processName, kinEnergy)) {
682  G4double e = kinEnergy;
683  if(baseParticle) {
684  e *= kinEnergy*massRatio;
687  e, cut) * chargeSquare;
688  } else {
690  res = currentModel->ComputeCrossSectionPerAtom(p, Z, shellIdx, e, cut);
691  }
692  if(verbose>0) {
693  G4cout << "E(MeV)= " << kinEnergy/MeV
694  << " cross(barn)= " << res/barn
695  << " " << p->GetParticleName()
696  << " Z= " << Z << " shellIdx= " << shellIdx << G4endl;
697  }
698  }
699  }
700  return res;
701 }
702 
703 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
704 
705 G4double
707  const G4Material* mat)
708 {
709  G4double res = 0.0;
710  const G4ParticleDefinition* gamma = G4Gamma::Gamma();
711  res += ComputeCrossSectionPerVolume(kinEnergy, gamma, "conv", mat, 0.0);
712  res += ComputeCrossSectionPerVolume(kinEnergy, gamma, "compt", mat, 0.0);
713  res += ComputeCrossSectionPerVolume(kinEnergy, gamma, "phot", mat, 0.0);
714  res += ComputeCrossSectionPerVolume(kinEnergy, gamma, "Rayl", mat, 0.0);
715  if(res > 0.0) { res = 1.0/res; }
716  return res;
717 }
718 
719 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
720 
722  const G4String& particle,
723  G4int Z,
725  G4double kinEnergy,
726  const G4Material* mat)
727 {
728  G4double res = 0.0;
729  const G4ParticleDefinition* p = FindParticle(particle);
731  if(p && ad) {
732  res = ad->ComputeShellIonisationCrossSectionPerAtom(p, Z, shell,
733  kinEnergy, mat);
734  }
735  return res;
736 }
737 
738 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
739 
741  const G4ParticleDefinition* p,
742  const G4String& processName,
743  const G4Material* mat,
744  G4double cut)
745 {
746  G4double mfp = DBL_MAX;
747  G4double x = ComputeCrossSectionPerVolume(kinEnergy, p, processName, mat, cut);
748  if(x > 0.0) { mfp = 1.0/x; }
749  if(verbose>1) {
750  G4cout << "E(MeV)= " << kinEnergy/MeV
751  << " MFP(mm)= " << mfp/mm
752  << " " << p->GetParticleName()
753  << " in " << mat->GetName()
754  << G4endl;
755  }
756  return mfp;
757 }
758 
759 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
760 
762  G4double range,
763  const G4ParticleDefinition* part,
764  const G4Material* mat)
765 {
767  ConvertRangeToEnergy(part, mat, range);
768 }
769 
770 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
771 
773  G4double kinEnergy)
774 {
775  if(p != currentParticle) {
776 
777  // new particle
778  currentParticle = p;
779  dynParticle.SetDefinition(const_cast<G4ParticleDefinition*>(p));
780  dynParticle.SetKineticEnergy(kinEnergy);
781  baseParticle = 0;
783  massRatio = 1.0;
784  mass = p->GetPDGMass();
785  chargeSquare = 1.0;
787  currentProcessName = "";
788  isIon = false;
789 
790  // ionisation process exist
791  if(currentProcess) {
794 
795  // base particle is used
796  if(baseParticle) {
799  chargeSquare = q*q;
800  }
801 
802  if(p->GetParticleType() == "nucleus"
803  && currentParticleName != "deuteron"
804  && currentParticleName != "triton"
805  && currentParticleName != "alpha+"
806  && currentParticleName != "helium"
807  && currentParticleName != "hydrogen"
808  ) {
809  isIon = true;
812  // G4cout << p->GetParticleName()
813  // << " in " << currentMaterial->GetName()
814  // << " e= " << kinEnergy << G4endl;
815  }
816  }
817  }
818 
819  // Effective charge for ions
820  if(isIon) {
821  chargeSquare =
824  if(currentProcess) {
826  //G4cout <<"NewP: massR= "<< massRatio << " q2= "
827  // << chargeSquare << G4endl;
828  }
829  }
830  return true;
831 }
832 
833 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
834 
836 {
837  const G4ParticleDefinition* p = 0;
838  if(name != currentParticleName) {
840  if(!p) {
841  G4cout << "### WARNING: G4EmCalculator::FindParticle fails to find "
842  << name << G4endl;
843  }
844  } else {
845  p = currentParticle;
846  }
847  return p;
848 }
849 
850 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
851 
853 {
854  const G4ParticleDefinition* p = ionTable->GetIon(Z,A,0);
855  return p;
856 }
857 
858 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
859 
861 {
862  if(name != currentMaterialName) {
864  if(!currentMaterial) {
865  G4cout << "### WARNING: G4EmCalculator::FindMaterial fails to find "
866  << name << G4endl;
867  }
868  }
869  return currentMaterial;
870 }
871 
872 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
873 
875 {
876  const G4Region* r = 0;
877  if(reg != "" && reg != "world") {
879  } else {
880  r = G4RegionStore::GetInstance()->GetRegion("DefaultRegionForTheWorld");
881  }
882  return r;
883 }
884 
885 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
886 
888  const G4Material* material,
889  const G4Region* region)
890 {
891  const G4MaterialCutsCouple* couple = 0;
892  SetupMaterial(material);
893  if(currentMaterial) {
894  // Access to materials
895  const G4ProductionCutsTable* theCoupleTable=
897  const G4Region* r = region;
898  if(r) {
899  couple = theCoupleTable->GetMaterialCutsCouple(material,
900  r->GetProductionCuts());
901  } else {
903  size_t nr = store->size();
904  if(0 < nr) {
905  for(size_t i=0; i<nr; ++i) {
906  couple = theCoupleTable->GetMaterialCutsCouple(
907  material, ((*store)[i])->GetProductionCuts());
908  if(couple) { break; }
909  }
910  }
911  }
912  }
913  if(!couple) {
915  ed << "G4EmCalculator::FindCouple: fail for material <"
916  << currentMaterialName << ">";
917  if(region) { ed << " and region " << region->GetName(); }
918  G4Exception("G4EmCalculator::FindCouple", "em0078",
919  FatalException, ed);
920  }
921  return couple;
922 }
923 
924 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
925 
927 {
928  SetupMaterial(material);
929  if(!currentMaterial) { return false; }
930  for (G4int i=0; i<nLocalMaterials; ++i) {
931  if(material == localMaterials[i] && cut == localCuts[i]) {
934  currentCut = cut;
935  return true;
936  }
937  }
938  const G4MaterialCutsCouple* cc = new G4MaterialCutsCouple(material);
939  localMaterials.push_back(material);
940  localCouples.push_back(cc);
941  localCuts.push_back(cut);
942  nLocalMaterials++;
943  currentCouple = cc;
945  currentCut = cut;
946  return true;
947 }
948 
949 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
950 
952  const G4String& processName,
953  G4double kinEnergy)
954 {
955  // Search for the process
956  if (!currentLambda || p != lambdaParticle || processName != lambdaName) {
957  lambdaName = processName;
958  currentLambda = 0;
959  lambdaParticle = p;
960 
961  const G4ParticleDefinition* part = p;
962  if(isIon) { part = theGenericIon; }
963 
964  // Search for energy loss process
965  currentName = processName;
966  currentModel = 0;
967  loweModel = 0;
968 
969  G4VEnergyLossProcess* elproc = FindEnLossProcess(part, processName);
970  if(elproc) {
971  currentLambda = elproc->LambdaTable();
972  if(currentLambda) {
973  isApplicable = true;
974  if(verbose>1) {
975  G4cout << "G4VEnergyLossProcess is found out: " << currentName
976  << G4endl;
977  }
978  }
979  return;
980  }
981 
982  // Search for discrete process
983  G4VEmProcess* proc = FindDiscreteProcess(part, processName);
984  if(proc) {
985  currentLambda = proc->LambdaTable();
986  if(currentLambda) {
987  isApplicable = true;
988  if(verbose>1) {
989  G4cout << "G4VEmProcess is found out: " << currentName << G4endl;
990  }
991  }
992  return;
993  }
994 
995  // Search for msc process
996  G4VMultipleScattering* msc = FindMscProcess(part, processName);
997  if(msc) {
998  currentModel = msc->SelectModel(kinEnergy,0);
999  if(currentModel) {
1001  if(currentLambda) {
1002  isApplicable = true;
1003  if(verbose>1) {
1004  G4cout << "G4VMultipleScattering is found out: " << currentName
1005  << G4endl;
1006  }
1007  }
1008  }
1009  }
1010  }
1011 }
1012 
1013 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1014 
1016  const G4String& processName,
1017  G4double kinEnergy)
1018 {
1019  isApplicable = false;
1020  if(!p || !currentMaterial) {
1021  G4cout << "G4EmCalculator::FindEmModel WARNING: no particle"
1022  << " or materail defined; particle: " << p << G4endl;
1023  return isApplicable;
1024  }
1025  G4String partname = p->GetParticleName();
1026  const G4ParticleDefinition* part = p;
1027  G4double scaledEnergy = kinEnergy*massRatio;
1028  if(isIon) { part = theGenericIon; }
1029 
1030  if(verbose > 1) {
1031  G4cout << "## G4EmCalculator::FindEmModel for " << partname
1032  << " (type= " << p->GetParticleType()
1033  << ") and " << processName << " at E(MeV)= " << scaledEnergy
1034  << G4endl;
1035  if(p != part) { G4cout << " GenericIon is the base particle" << G4endl; }
1036  }
1037 
1038  // Search for energy loss process
1039  currentName = processName;
1040  currentModel = 0;
1041  loweModel = 0;
1042  size_t idx = 0;
1043 
1044  G4VEnergyLossProcess* elproc = FindEnLossProcess(part, processName);
1045  if(elproc) {
1046  currentModel = elproc->SelectModelForMaterial(scaledEnergy, idx);
1048  currentModel->SetupForMaterial(part, currentMaterial, scaledEnergy);
1050  if(eth > 0.0) {
1051  loweModel = elproc->SelectModelForMaterial(eth - CLHEP::eV, idx);
1052  if(loweModel == currentModel) { loweModel = nullptr; }
1053  else {
1056  }
1057  }
1058  }
1059 
1060  // Search for discrete process
1061  if(!currentModel) {
1062  G4VEmProcess* proc = FindDiscreteProcess(part, processName);
1063  if(proc) {
1064  currentModel = proc->SelectModelForMaterial(kinEnergy, idx);
1066  currentModel->SetupForMaterial(part, currentMaterial, kinEnergy);
1068  if(eth > 0.0) {
1069  loweModel = proc->SelectModelForMaterial(eth - CLHEP::eV, idx);
1070  if(loweModel == currentModel) { loweModel = 0; }
1071  else {
1074  }
1075  }
1076  }
1077  }
1078 
1079  // Search for msc process
1080  if(!currentModel) {
1081  G4VMultipleScattering* proc = FindMscProcess(part, processName);
1082  if(proc) {
1083  currentModel = proc->SelectModel(kinEnergy, idx);
1084  loweModel = 0;
1085  }
1086  }
1087  if(currentModel) {
1088  if(loweModel == currentModel) { loweModel = 0; }
1089  isApplicable = true;
1091  if(loweModel) {
1093  }
1094  if(verbose > 1) {
1095  G4cout << " Model <" << currentModel->GetName()
1096  << "> Emin(MeV)= " << currentModel->LowEnergyLimit()/MeV
1097  << " for " << part->GetParticleName();
1098  if(elproc) {
1099  G4cout << " and " << elproc->GetProcessName() << " " << elproc
1100  << G4endl;
1101  }
1102  if(loweModel) {
1103  G4cout << " LowEnergy model <" << loweModel->GetName() << ">";
1104  }
1105  G4cout << G4endl;
1106  }
1107  }
1108  return isApplicable;
1109 }
1110 
1111 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1112 
1114  const G4ParticleDefinition* p)
1115 {
1116  G4VEnergyLossProcess* elp = 0;
1117  G4String partname = p->GetParticleName();
1118  const G4ParticleDefinition* part = p;
1119 
1120  if(p->GetParticleType() == "nucleus"
1121  && currentParticleName != "deuteron"
1122  && currentParticleName != "triton"
1123  && currentParticleName != "alpha+"
1124  && currentParticleName != "helium"
1125  && currentParticleName != "hydrogen"
1126  ) { part = theGenericIon; }
1127 
1129  return elp;
1130 }
1131 
1132 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1133 
1136  const G4String& processName)
1137 {
1138  G4VEnergyLossProcess* proc = 0;
1139  const std::vector<G4VEnergyLossProcess*> v =
1141  G4int n = v.size();
1142  for(G4int i=0; i<n; ++i) {
1143  if((v[i])->GetProcessName() == processName) {
1144  G4VProcess* p = reinterpret_cast<G4VProcess*>(v[i]);
1145  if(ActiveForParticle(part, p)) {
1146  proc = v[i];
1147  break;
1148  }
1149  }
1150  }
1151  return proc;
1152 }
1153 
1154 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1155 
1156 G4VEmProcess*
1158  const G4String& processName)
1159 {
1160  G4VEmProcess* proc = 0;
1161  const std::vector<G4VEmProcess*> v =
1163  G4int n = v.size();
1164  for(G4int i=0; i<n; ++i) {
1165  if((v[i])->GetProcessName() == processName) {
1166  G4VProcess* p = reinterpret_cast<G4VProcess*>(v[i]);
1167  if(ActiveForParticle(part, p)) {
1168  proc = v[i];
1169  break;
1170  }
1171  }
1172  }
1173  return proc;
1174 }
1175 
1176 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1177 
1180  const G4String& processName)
1181 {
1182  G4VMultipleScattering* proc = 0;
1183  const std::vector<G4VMultipleScattering*> v =
1185  G4int n = v.size();
1186  for(G4int i=0; i<n; ++i) {
1187  if((v[i])->GetProcessName() == processName) {
1188  G4VProcess* p = reinterpret_cast<G4VProcess*>(v[i]);
1189  if(ActiveForParticle(part, p)) {
1190  proc = v[i];
1191  break;
1192  }
1193  }
1194  }
1195  return proc;
1196 }
1197 
1198 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1199 
1201  const G4String& processName)
1202 {
1203  G4VProcess* proc = 0;
1204  const G4ProcessManager* procman = part->GetProcessManager();
1205  G4ProcessVector* pv = procman->GetProcessList();
1206  G4int nproc = pv->size();
1207  for(G4int i=0; i<nproc; ++i) {
1208  if(processName == (*pv)[i]->GetProcessName()) {
1209  proc = (*pv)[i];
1210  break;
1211  }
1212  }
1213  return proc;
1214 }
1215 
1216 
1217 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1218 
1220  G4VProcess* proc)
1221 {
1222  G4ProcessManager* pm = part->GetProcessManager();
1223  G4ProcessVector* pv = pm->GetProcessList();
1224  G4int n = pv->size();
1225  G4bool res = false;
1226  for(G4int i=0; i<n; ++i) {
1227  if((*pv)[i] == proc) {
1228  if(pm->GetProcessActivation(i)) { res = true; }
1229  break;
1230  }
1231  }
1232  return res;
1233 }
1234 
1235 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1236 
1238 {
1239  if(mat) {
1240  currentMaterial = mat;
1241  currentMaterialName = mat->GetName();
1242  } else {
1243  currentMaterial = 0;
1244  currentMaterialName = "";
1245  }
1246 }
1247 
1248 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1249 
1251 {
1253 }
1254 
1255 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1256 
1258 {
1259  G4bool isFound = false;
1260  if(currentMaterial) {
1262  for(size_t i=0; i<nn; ++i) {
1263  if(Z == G4lrint(currentMaterial->GetElement(i)->GetZ())) {
1264  isFound = true;
1265  break;
1266  }
1267  }
1268  }
1269  if(!isFound) {
1271  }
1272 }
1273 
1274 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1275 
1277 {
1278  verbose = verb;
1279 }
1280 
1281 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1282 
static const double cm
Definition: G4SIunits.hh:118
G4PhysicsTable * LambdaTable() const
static const double gram
Definition: G4SIunits.hh:175
G4VEmModel * loweModel
G4VEmProcess * FindDiscreteProcess(const G4ParticleDefinition *, const G4String &processName)
G4ProductionCuts * GetProductionCuts() const
const G4PhysicsTable * currentLambda
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
void FindLambdaTable(const G4ParticleDefinition *, const G4String &processName, G4double kinEnergy)
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:258
G4Material * FindOrBuildMaterial(const G4String &name, G4bool isotopes=true, G4bool warning=false)
const G4ParticleDefinition * theGenericIon
static const double MeV
Definition: G4SIunits.hh:211
G4VEmModel * SelectModel(G4double kinEnergy, size_t idx)
G4Region * GetRegion(const G4String &name, G4bool verbose=true) const
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)
Definition: G4VEmModel.cc:244
const G4String & GetName() const
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4LossTableManager * Instance()
static const double cm2
Definition: G4SIunits.hh:119
virtual void CorrectionsAlongStep(const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
Definition: G4VEmModel.cc:362
void SetDynamicMassCharge(G4double massratio, G4double charge2ratio)
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4double GetDEDX(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
const G4Material * FindMaterial(const G4String &)
G4Material * FindOrBuildSimpleMaterial(G4int Z, G4bool warning=false)
const G4Material * cutMaterial
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:604
std::vector< const G4MaterialCutsCouple * > localCouples
const G4ParticleDefinition * currentParticle
G4PhysicsTable * RangeTableForLoss() const
G4double EffectiveChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4String name
Definition: TRTMaterials.hh:40
const G4String & GetName() const
Definition: G4Material.hh:178
void PrintRangeTable(const G4ParticleDefinition *)
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:491
virtual G4double GetShellIonisationCrossSectionPerAtom(const G4ParticleDefinition *, G4int Z, G4AtomicShellEnumerator shell, G4double kinE, const G4Material *mat=0)=0
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:403
G4double GetZ() const
Definition: G4Element.hh:131
G4bool UpdateCouple(const G4Material *, G4double cut)
G4double GetDensity() const
Definition: G4Material.hh:180
G4double ComputeTotalDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, G4double cut=DBL_MAX)
G4double ComputeElectronicDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *mat, G4double cut=DBL_MAX)
G4PhysicsTable * GetCrossSectionTable()
Definition: G4VEmModel.hh:826
G4double EffectiveChargeCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4EmCorrections * corr
G4ionEffectiveCharge * ionEffCharge
const G4MaterialCutsCouple * currentCouple
G4double GetCrossSectionPerVolume(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, const G4Region *r=0)
G4double GetMeanFreePath(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, const G4Region *r=0)
void SetupMaterial(const G4Material *)
void PrintInverseRangeTable(const G4ParticleDefinition *)
G4double GetRange(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=0)
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:202
G4ProcessManager * GetProcessManager() const
const std::vector< G4VEmProcess * > & GetEmProcessVector()
int G4int
Definition: G4Types.hh:78
G4DynamicParticle dynParticle
static G4NistManager * Instance()
const G4String & GetParticleName() const
G4double ComputeShellIonisationCrossSectionPerAtom(const G4String &part, G4int Z, G4AtomicShellEnumerator shell, G4double kinEnergy, const G4Material *mat=0)
G4double GetCSDARange(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4double GetDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=0)
static const G4double reg
static G4RegionStore * GetInstance()
G4double GetCSDARange(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=0)
G4double cutenergy[3]
G4double GetShellIonisationCrossSectionPerAtom(const G4String &part, G4int Z, G4AtomicShellEnumerator shell, G4double kinEnergy)
G4PhysicsTable * LambdaTable() const
G4VEnergyLossProcess * FindEnergyLossProcess(const G4ParticleDefinition *)
G4IonTable * GetIonTable() const
virtual void InitialiseForMaterial(const G4ParticleDefinition *, const G4Material *)
Definition: G4VEmModel.cc:225
G4GLOB_DLL std::ostream G4cout
double A(double temperature)
G4VEnergyLossProcess * currentProcess
const G4MaterialCutsCouple * FindCouple(const G4Material *, const G4Region *r=0)
G4PhysicsTable * DEDXTable() const
G4double ComputeDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, G4double cut=DBL_MAX)
G4String currentName
G4double NuclearDEDX(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy, G4bool fluct=true)
const G4ParticleDefinition * baseParticle
bool G4bool
Definition: G4Types.hh:79
G4EmCorrections * EmCorrections()
G4bool FindEmModel(const G4ParticleDefinition *, const G4String &processName, G4double kinEnergy)
static const double nm
Definition: G4SIunits.hh:111
G4double ComputeNuclearDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *)
G4double GetRangeFromRestricteDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=0)
G4DataVector localCuts
static const double GeV
Definition: G4SIunits.hh:214
const G4String & GetParticleType() const
const G4ParticleDefinition * BaseParticle() const
const G4int n
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
G4double GetEnergy(const G4ParticleDefinition *aParticle, G4double range, const G4MaterialCutsCouple *couple)
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
void CheckMaterial(G4int Z)
G4bool ActiveForParticle(const G4ParticleDefinition *part, G4VProcess *proc)
const std::vector< G4VEnergyLossProcess * > & GetEnergyLossProcessVector()
G4VMultipleScattering * FindMscProcess(const G4ParticleDefinition *, const G4String &processName)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:313
void SetKineticEnergy(G4double aEnergy)
virtual G4double ComputeShellIonisationCrossSectionPerAtom(const G4ParticleDefinition *, G4int Z, G4AtomicShellEnumerator shell, G4double kinE, const G4Material *mat=0)=0
G4double GetRangeFromRestricteDEDX(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
std::vector< const G4Material * > localMaterials
G4VEnergyLossProcess * FindEnLossProcess(const G4ParticleDefinition *, const G4String &processName)
G4int size() const
virtual G4double GetChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:345
G4VEnergyLossProcess * GetEnergyLossProcess(const G4ParticleDefinition *)
static G4ProductionCutsTable * GetProductionCutsTable()
static const double eV
Definition: G4SIunits.hh:212
static G4Positron * Positron()
Definition: G4Positron.cc:94
static G4GenericIon * GenericIon()
Definition: G4GenericIon.cc:93
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
Definition: G4VEmModel.cc:249
G4VEmModel * SelectModelForMaterial(G4double kinEnergy, size_t &idxRegion) const
const G4ParticleDefinition * FindParticle(const G4String &)
G4PhysicsTable * InverseRangeTable() const
G4double GetPDGMass() const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
G4double chargeSquare
static G4ParticleTable * GetParticleTable()
int G4lrint(double ad)
Definition: templates.hh:163
G4double ComputeMeanFreePath(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, G4double cut=0.0)
G4double ComputeCrossSectionPerShell(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, G4int Z, G4int shellIdx, G4double cut=0.0)
G4double GetKinEnergy(G4double range, const G4ParticleDefinition *, const G4Material *, const G4Region *r=0)
static const double g
Definition: G4SIunits.hh:180
const G4Region * FindRegion(const G4String &)
G4IonTable * ionTable
G4String currentMaterialName
const G4ParticleDefinition * FindIon(G4int Z, G4int A)
G4String currentProcessName
G4String currentParticleName
const G4double x[NPOINTSGL]
G4double ComputeCrossSectionPerVolume(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, G4double cut=0.0)
virtual G4double ComputeCrossSectionPerShell(const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:323
G4double ComputeGammaAttenuationLength(G4double kinEnergy, const G4Material *)
static const double mole
Definition: G4SIunits.hh:283
G4bool GetProcessActivation(G4VProcess *aProcess) const
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
G4double ComputeCrossSectionPerAtom(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, G4double Z, G4double A, G4double cut=0.0)
const G4String & GetName() const
Definition: G4VEmModel.hh:795
G4VProcess * FindProcess(const G4ParticleDefinition *part, const G4String &processName)
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
void SetVerbose(G4int val)
void PrintDEDXTable(const G4ParticleDefinition *)
G4VAtomDeexcitation * AtomDeexcitation()
static const double barn
Definition: G4SIunits.hh:104
static const double keV
Definition: G4SIunits.hh:213
const G4ParticleDefinition * lambdaParticle
double G4double
Definition: G4Types.hh:76
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
G4double ComputeEnergyCutFromRangeCut(G4double range, const G4ParticleDefinition *, const G4Material *)
const G4Material * currentMaterial
G4VEmModel * currentModel
G4double GetPDGCharge() const
G4bool UpdateParticle(const G4ParticleDefinition *, G4double kinEnergy)
#define DBL_MAX
Definition: templates.hh:83
G4double GetRange(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
static const double mm
Definition: G4SIunits.hh:114
G4LossTableManager * manager
G4double ComputeDEDXForCutInRange(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *mat, G4double rangecut=DBL_MAX)
G4AtomicShellEnumerator
const std::vector< G4VMultipleScattering * > & GetMultipleScatteringVector()
G4VEmModel * SelectModelForMaterial(G4double kinEnergy, size_t &idx) const
G4ProcessVector * GetProcessList() const
G4NistManager * nist