Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 79268 2014-02-20 16:46:31Z 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 "G4VEmProcess.hh"
69 #include "G4VEnergyLossProcess.hh"
70 #include "G4VMultipleScattering.hh"
71 #include "G4Material.hh"
72 #include "G4MaterialCutsCouple.hh"
73 #include "G4ParticleDefinition.hh"
74 #include "G4ParticleTable.hh"
75 #include "G4IonTable.hh"
76 #include "G4PhysicsTable.hh"
77 #include "G4ProductionCutsTable.hh"
78 #include "G4ProcessManager.hh"
79 #include "G4ionEffectiveCharge.hh"
80 #include "G4RegionStore.hh"
81 #include "G4Element.hh"
82 #include "G4EmCorrections.hh"
83 #include "G4GenericIon.hh"
84 #include "G4ProcessVector.hh"
85 #include "G4Gamma.hh"
86 
87 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
88 
90 {
91  manager = G4LossTableManager::Instance();
92  corr = manager->EmCorrections();
93  nLocalMaterials = 0;
94  verbose = 0;
95  currentCoupleIndex = 0;
96  currentCouple = 0;
97  currentMaterial = 0;
98  currentParticle = 0;
99  lambdaParticle = 0;
100  baseParticle = 0;
101  currentLambda = 0;
102  currentModel = 0;
103  currentProcess = 0;
104  loweModel = 0;
105  chargeSquare = 1.0;
106  massRatio = 1.0;
107  mass = 0.0;
108  currentCut = 0.0;
109  currentParticleName= "";
110  currentMaterialName= "";
111  currentName = "";
112  lambdaName = "";
113  theGenericIon = G4GenericIon::GenericIon();
114  ionEffCharge = new G4ionEffectiveCharge();
116  isIon = false;
117  isApplicable = false;
118 }
119 
120 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
121 
123 {
124  delete ionEffCharge;
125  for (G4int i=0; i<nLocalMaterials; ++i) {
126  delete localCouples[i];
127  }
128 }
129 
130 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
131 
133  const G4ParticleDefinition* p,
134  const G4Material* mat,
135  const G4Region* region)
136 {
137  G4double res = 0.0;
138  const G4MaterialCutsCouple* couple = FindCouple(mat, region);
139  if(couple && UpdateParticle(p, kinEnergy) ) {
140  res = manager->GetDEDX(p, kinEnergy, couple);
141 
142  if(isIon) {
143  if(FindEmModel(p, currentProcessName, kinEnergy)) {
144  G4double length = CLHEP::nm;
145  G4double eloss = res*length;
146  //G4cout << "### GetDEDX: E= " << kinEnergy << " dedx0= " << res
147  // << " de= " << eloss << G4endl;;
148  G4double niel = 0.0;
149  dynParticle.SetKineticEnergy(kinEnergy);
150  currentModel->GetChargeSquareRatio(p, mat, kinEnergy);
151  currentModel->CorrectionsAlongStep(couple,&dynParticle,eloss,niel,length);
152  res = eloss/length;
153  //G4cout << " de1= " << eloss << " res1= " << res
154  // << " " << p->GetParticleName() <<G4endl;;
155  }
156  }
157 
158  if(verbose>0) {
159  G4cout << "G4EmCalculator::GetDEDX: E(MeV)= " << kinEnergy/MeV
160  << " DEDX(MeV/mm)= " << res*mm/MeV
161  << " DEDX(MeV*cm^2/g)= " << res*gram/(MeV*cm2*mat->GetDensity())
162  << " " << p->GetParticleName()
163  << " in " << mat->GetName()
164  << " isIon= " << isIon
165  << G4endl;
166  }
167  }
168  return res;
169 }
170 
171 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
172 
174  const G4String& material, const G4String& reg)
175 {
176  return GetDEDX(kinEnergy,FindParticle(particle),
177  FindMaterial(material),FindRegion(reg));
178 }
179 
180 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
181 
183  const G4ParticleDefinition* p,
184  const G4Material* mat,
185  const G4Region* region)
186 {
187  G4double res = 0.0;
188  const G4MaterialCutsCouple* couple = FindCouple(mat,region);
189  if(couple && UpdateParticle(p, kinEnergy)) {
190  res = manager->GetRangeFromRestricteDEDX(p, kinEnergy, couple);
191  if(verbose>0) {
192  G4cout << "G4EmCalculator::GetRange: E(MeV)= " << kinEnergy/MeV
193  << " range(mm)= " << res/mm
194  << " " << p->GetParticleName()
195  << " in " << mat->GetName()
196  << G4endl;
197  }
198  }
199  return res;
200 }
201 
202 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
203 
205  const G4ParticleDefinition* p,
206  const G4Material* mat,
207  const G4Region* region)
208 {
209  G4double res = 0.0;
210  if(!G4LossTableManager::Instance()->BuildCSDARange()) {
212  ed << "G4EmCalculator::GetCSDARange: CSDA table is not built; "
213  << " use UI command: /process/eLoss/CSDARange true";
214  G4Exception("G4EmCalculator::GetCSDARange", "em0077",
215  JustWarning, ed);
216  return res;
217  }
218 
219  const G4MaterialCutsCouple* couple = FindCouple(mat,region);
220  if(couple && UpdateParticle(p, kinEnergy)) {
221  res = manager->GetCSDARange(p, kinEnergy, couple);
222  if(verbose>0) {
223  G4cout << "G4EmCalculator::GetRange: E(MeV)= " << kinEnergy/MeV
224  << " range(mm)= " << res/mm
225  << " " << p->GetParticleName()
226  << " in " << mat->GetName()
227  << G4endl;
228  }
229  }
230  return res;
231 }
232 
233 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
234 
236  const G4ParticleDefinition* p,
237  const G4Material* mat,
238  const G4Region* region)
239 {
240  G4double res = 0.0;
241  const G4MaterialCutsCouple* couple = FindCouple(mat,region);
242  if(couple && UpdateParticle(p, kinEnergy)) {
243  res = manager->GetRange(p, kinEnergy, couple);
244  if(verbose>0) {
245  G4cout << "G4EmCalculator::GetRange: E(MeV)= " << kinEnergy/MeV
246  << " range(mm)= " << res/mm
247  << " " << p->GetParticleName()
248  << " in " << mat->GetName()
249  << G4endl;
250  }
251  }
252  return res;
253 }
254 
255 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
256 
258  const G4String& particle,
259  const G4String& material,
260  const G4String& reg)
261 {
262  return GetRangeFromRestricteDEDX(kinEnergy,FindParticle(particle),
263  FindMaterial(material),FindRegion(reg));
264 }
265 
266 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
267 
269  const G4String& particle,
270  const G4String& material,
271  const G4String& reg)
272 {
273  return GetCSDARange(kinEnergy,FindParticle(particle),
274  FindMaterial(material),FindRegion(reg));
275 }
276 
277 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
278 
280  const G4String& particle,
281  const G4String& material,
282  const G4String& reg)
283 {
284  return GetRange(kinEnergy,FindParticle(particle),
285  FindMaterial(material),FindRegion(reg));
286 }
287 
288 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
289 
291  const G4ParticleDefinition* p,
292  const G4Material* mat,
293  const G4Region* region)
294 {
295  G4double res = 0.0;
296  const G4MaterialCutsCouple* couple = FindCouple(mat,region);
297  if(couple && UpdateParticle(p, 1.0*GeV)) {
298  res = manager->GetEnergy(p, range, couple);
299  if(verbose>0) {
300  G4cout << "G4EmCalculator::GetKinEnergy: Range(mm)= " << range/mm
301  << " KinE(MeV)= " << res/MeV
302  << " " << p->GetParticleName()
303  << " in " << mat->GetName()
304  << G4endl;
305  }
306  }
307  return res;
308 }
309 
310 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
311 
313  const G4String& material, const G4String& reg)
314 {
315  return GetKinEnergy(range,FindParticle(particle),
316  FindMaterial(material),FindRegion(reg));
317 }
318 
319 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
320 
322  const G4ParticleDefinition* p,
323  const G4String& processName,
324  const G4Material* mat,
325  const G4Region* region)
326 {
327  G4double res = 0.0;
328  const G4MaterialCutsCouple* couple = FindCouple(mat,region);
329 
330  if(couple && UpdateParticle(p, kinEnergy)) {
331  G4int idx = couple->GetIndex();
332  FindLambdaTable(p, processName, kinEnergy);
333 
334  if(currentLambda) {
335  G4double e = kinEnergy*massRatio;
336  res = (((*currentLambda)[idx])->Value(e))*chargeSquare;
337  } else {
338  res = ComputeCrossSectionPerVolume(kinEnergy, p, processName, mat,
339  kinEnergy);
340  }
341  if(verbose>0) {
342  G4cout << "G4EmCalculator::GetXSPerVolume: E(MeV)= " << kinEnergy/MeV
343  << " cross(cm-1)= " << res*cm
344  << " " << p->GetParticleName()
345  << " in " << mat->GetName();
346  if(verbose>1)
347  G4cout << " idx= " << idx << " Escaled((MeV)= "
348  << kinEnergy*massRatio
349  << " q2= " << chargeSquare;
350  G4cout << G4endl;
351  }
352  }
353  return res;
354 }
355 
356 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
357 
359  const G4String& particle,
360  const G4String& processName,
361  const G4String& material,
362  const G4String& reg)
363 {
364  return GetCrossSectionPerVolume(kinEnergy,FindParticle(particle),processName,
365  FindMaterial(material),FindRegion(reg));
366 }
367 
368 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
369 
371  const G4String& particle,
372  G4int Z,
374  G4double kinEnergy)
375 {
376  G4double res = 0.0;
377  const G4ParticleDefinition* p = FindParticle(particle);
378  G4VAtomDeexcitation* ad = manager->AtomDeexcitation();
379  if(p && ad) {
380  res = ad->GetShellIonisationCrossSectionPerAtom(p, Z, shell, kinEnergy);
381  }
382  return res;
383 }
384 
385 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
386 
388  const G4ParticleDefinition* p,
389  const G4String& processName,
390  const G4Material* mat,
391  const G4Region* region)
392 {
393  G4double res = DBL_MAX;
394  G4double x = GetCrossSectionPerVolume(kinEnergy,p, processName, mat,region);
395  if(x > 0.0) { res = 1.0/x; }
396  if(verbose>1) {
397  G4cout << "G4EmCalculator::GetMeanFreePath: E(MeV)= " << kinEnergy/MeV
398  << " MFP(mm)= " << res/mm
399  << " " << p->GetParticleName()
400  << " in " << mat->GetName()
401  << G4endl;
402  }
403  return res;
404 }
405 
406 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
407 
409  const G4String& particle,
410  const G4String& processName,
411  const G4String& material,
412  const G4String& reg)
413 {
414  return GetMeanFreePath(kinEnergy,FindParticle(particle),processName,
415  FindMaterial(material),FindRegion(reg));
416 }
417 
418 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
419 
421 {
422  const G4VEnergyLossProcess* elp = FindEnergyLossProcess(p);
423  G4cout << "##### DEDX Table for " << p->GetParticleName() << G4endl;
424  if(elp) G4cout << *(elp->DEDXTable()) << G4endl;
425 }
426 
427 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
428 
430 {
431  const G4VEnergyLossProcess* elp = FindEnergyLossProcess(p);
432  G4cout << "##### Range Table for " << p->GetParticleName() << G4endl;
433  if(elp) G4cout << *(elp->RangeTableForLoss()) << G4endl;
434 }
435 
436 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
437 
439 {
440  const G4VEnergyLossProcess* elp = FindEnergyLossProcess(p);
441  G4cout << "### G4EmCalculator: Inverse Range Table for "
442  << p->GetParticleName() << G4endl;
443  if(elp) G4cout << *(elp->InverseRangeTable()) << G4endl;
444 }
445 
446 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
447 
449  const G4ParticleDefinition* p,
450  const G4String& processName,
451  const G4Material* mat,
452  G4double cut)
453 {
454  currentMaterial = mat;
455  currentMaterialName = mat->GetName();
456  G4double res = 0.0;
457  if(verbose > 1) {
458  G4cout << "### G4EmCalculator::ComputeDEDX: " << p->GetParticleName()
459  << " in " << currentMaterialName
460  << " e(MeV)= " << kinEnergy/MeV << " cut(MeV)= " << cut/MeV
461  << G4endl;
462  }
463  if(UpdateParticle(p, kinEnergy)) {
464  if(FindEmModel(p, processName, kinEnergy)) {
465  G4double escaled = kinEnergy*massRatio;
466  if(baseParticle) {
467  res = currentModel->ComputeDEDXPerVolume(
468  mat, baseParticle, escaled, cut) * chargeSquare;
469  if(verbose > 1) {
470  G4cout << baseParticle->GetParticleName()
471  << " Escaled(MeV)= " << escaled;
472  }
473  } else {
474  res = currentModel->ComputeDEDXPerVolume(mat, p, kinEnergy, cut);
475  if(verbose > 1) { G4cout << " no basePart E(MeV)= " << kinEnergy << " "; }
476  }
477  if(verbose > 1) {
478  G4cout << currentModel->GetName() << ": DEDX(MeV/mm)= " << res*mm/MeV
479  << " DEDX(MeV*cm^2/g)= "
480  << res*gram/(MeV*cm2*mat->GetDensity())
481  << G4endl;
482  }
483 
484  // emulate smoothing procedure
485  G4double eth = currentModel->LowEnergyLimit();
486  // G4cout << "massRatio= " << massRatio << " eth= " << eth << G4endl;
487  if(loweModel) {
488  G4double res0 = 0.0;
489  G4double res1 = 0.0;
490  if(baseParticle) {
491  res1 = currentModel->ComputeDEDXPerVolume(mat, baseParticle, eth, cut)
492  * chargeSquare;
493  res0 = loweModel->ComputeDEDXPerVolume(mat, baseParticle, eth, cut)
494  * chargeSquare;
495  } else {
496  res1 = currentModel->ComputeDEDXPerVolume(mat, p, eth, cut);
497  res0 = loweModel->ComputeDEDXPerVolume(mat, p, eth, cut);
498  }
499  if(verbose > 1) {
500  G4cout << "At boundary energy(MeV)= " << eth/MeV
501  << " DEDX(MeV/mm)= " << res1*mm/MeV
502  << G4endl;
503  }
504 
505  //G4cout << "eth= " << eth << " escaled= " << escaled
506  // << " res0= " << res0 << " res1= "
507  // << res1 << " q2= " << chargeSquare << G4endl;
508 
509  if(res1 > 0.0 && escaled > 0.0) {
510  res *= (1.0 + (res0/res1 - 1.0)*eth/escaled);
511  }
512  }
513 
514  // low energy correction for ions
515  if(isIon) {
516  G4double length = CLHEP::nm;
517  const G4Region* r = 0;
518  const G4MaterialCutsCouple* couple = FindCouple(mat, r);
519  G4double eloss = res*length;
520  G4double niel = 0.0;
521  dynParticle.SetKineticEnergy(kinEnergy);
522  currentModel->GetChargeSquareRatio(p, mat, kinEnergy);
523  currentModel->CorrectionsAlongStep(couple,&dynParticle,eloss,niel,length);
524  res = eloss/length;
525 
526  if(verbose > 1) {
527  G4cout << "After Corrections: DEDX(MeV/mm)= " << res*mm/MeV
528  << " DEDX(MeV*cm^2/g)= " << res*gram/(MeV*cm2*mat->GetDensity())
529  << G4endl;
530  }
531  }
532  }
533 
534  if(verbose > 0) {
535  G4cout << "Sum: E(MeV)= " << kinEnergy/MeV
536  << " DEDX(MeV/mm)= " << res*mm/MeV
537  << " DEDX(MeV*cm^2/g)= " << res*gram/(MeV*cm2*mat->GetDensity())
538  << " cut(MeV)= " << cut/MeV
539  << " " << p->GetParticleName()
540  << " in " << currentMaterialName
541  << " Zi^2= " << chargeSquare
542  << " isIon=" << isIon
543  << G4endl;
544  }
545  }
546  return res;
547 }
548 
549 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
550 
552  const G4ParticleDefinition* part,
553  const G4Material* mat,
554  G4double cut)
555 {
556  currentMaterial = mat;
557  currentMaterialName = mat->GetName();
558  G4double dedx = 0.0;
559  if(UpdateParticle(part, kinEnergy)) {
560 
562  const std::vector<G4VEnergyLossProcess*> vel =
563  lManager->GetEnergyLossProcessVector();
564  G4int n = vel.size();
565 
566  //G4cout << "ComputeElectronicDEDX for " << part->GetParticleName()
567  // << " n= " << n << G4endl;
568 
569  for(G4int i=0; i<n; ++i) {
570  if(vel[i]) {
571  G4VProcess* p = reinterpret_cast<G4VProcess*>(vel[i]);
572  if(ActiveForParticle(part, p)) {
573  //G4cout << "idx= " << i << " " << (vel[i])->GetProcessName()
574  // << " " << (vel[i])->Particle()->GetParticleName() << G4endl;
575  dedx += ComputeDEDX(kinEnergy,part,(vel[i])->GetProcessName(),mat,cut);
576  }
577  }
578  }
579  }
580  return dedx;
581 }
582 
583 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
584 
586  const G4String& mat, G4double cut)
587 {
588  return ComputeElectronicDEDX(kinEnergy,FindParticle(part),FindMaterial(mat),cut);
589 }
590 
591 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
592 
594  const G4ParticleDefinition* part,
595  const G4Material* mat,
596  G4double cut)
597 {
598  G4double dedx = ComputeElectronicDEDX(kinEnergy,part,mat,cut);
599  if(mass > 700.*MeV) { dedx += ComputeNuclearDEDX(kinEnergy,part,mat); }
600  return dedx;
601 }
602 
603 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
604 
606  const G4String& part,
607  const G4String& mat,
608  G4double cut)
609 {
610  return ComputeTotalDEDX(kinEnergy,FindParticle(part),FindMaterial(mat),cut);
611 }
612 
613 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
614 
616  const G4String& particle,
617  const G4String& processName,
618  const G4String& material,
619  G4double cut)
620 {
621  return ComputeDEDX(kinEnergy,FindParticle(particle),processName,
622  FindMaterial(material),cut);
623 }
624 
625 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
626 
628  const G4ParticleDefinition* p,
629  const G4Material* mat)
630 {
631 
632  G4double res = corr->NuclearDEDX(p, mat, kinEnergy, false);
633 
634  if(verbose > 1) {
635  G4cout << p->GetParticleName() << " E(MeV)= " << kinEnergy/MeV
636  << " NuclearDEDX(MeV/mm)= " << res*mm/MeV
637  << " NuclearDEDX(MeV*cm^2/g)= "
638  << res*gram/(MeV*cm2*mat->GetDensity())
639  << G4endl;
640  }
641  return res;
642 }
643 
644 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
645 
647  const G4String& particle,
648  const G4String& material)
649 {
650  return ComputeNuclearDEDX(kinEnergy,FindParticle(particle),
651  FindMaterial(material));
652 }
653 
654 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
655 
657  G4double kinEnergy,
658  const G4ParticleDefinition* p,
659  const G4String& processName,
660  const G4Material* mat,
661  G4double cut)
662 {
663  currentMaterial = mat;
664  currentMaterialName = mat->GetName();
665  G4double res = 0.0;
666  if(UpdateParticle(p, kinEnergy)) {
667  if(FindEmModel(p, processName, kinEnergy)) {
668  G4double e = kinEnergy;
669  if(baseParticle) {
670  e *= kinEnergy*massRatio;
671  res = currentModel->CrossSectionPerVolume(
672  mat, baseParticle, e, cut, e) * chargeSquare;
673  } else {
674  res = currentModel->CrossSectionPerVolume(mat, p, e, cut, e);
675  }
676  if(verbose>0) {
677  G4cout << "G4EmCalculator::ComputeXSPerVolume: E(MeV)= " << kinEnergy/MeV
678  << " cross(cm-1)= " << res*cm
679  << " cut(keV)= " << cut/keV
680  << " " << p->GetParticleName()
681  << " in " << mat->GetName()
682  << G4endl;
683  }
684  }
685  }
686  return res;
687 }
688 
689 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
690 
692  G4double kinEnergy,
693  const G4String& particle,
694  const G4String& processName,
695  const G4String& material,
696  G4double cut)
697 {
698  return ComputeCrossSectionPerVolume(kinEnergy,FindParticle(particle),
699  processName,
700  FindMaterial(material),cut);
701 }
702 
703 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
704 
706  G4double kinEnergy,
707  const G4ParticleDefinition* p,
708  const G4String& processName,
709  G4double Z, G4double A,
710  G4double cut)
711 {
712  G4double res = 0.0;
713  if(UpdateParticle(p, kinEnergy)) {
714  if(FindEmModel(p, processName, kinEnergy)) {
715  G4double e = kinEnergy;
716  if(baseParticle) {
717  e *= kinEnergy*massRatio;
718  currentModel->InitialiseForElement(baseParticle, G4lrint(Z));
719  res = currentModel->ComputeCrossSectionPerAtom(
720  baseParticle, e, Z, A, cut) * chargeSquare;
721  } else {
722  currentModel->InitialiseForElement(p, G4lrint(Z));
723  res = currentModel->ComputeCrossSectionPerAtom(p, e, Z, A, cut);
724  }
725  if(verbose>0) {
726  G4cout << "E(MeV)= " << kinEnergy/MeV
727  << " cross(barn)= " << res/barn
728  << " " << p->GetParticleName()
729  << " Z= " << Z << " A= " << A/(g/mole) << " g/mole"
730  << G4endl;
731  }
732  }
733  }
734  return res;
735 }
736 
737 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
738 
740  const G4String& particle,
741  const G4String& processName,
742  const G4Element* elm,
743  G4double cut)
744 {
745  return ComputeCrossSectionPerAtom(kinEnergy,FindParticle(particle),
746  processName,
747  elm->GetZ(),elm->GetN(),cut);
748 }
749 
750 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
751 
752 G4double
754  const G4Material* mat)
755 {
756  G4double res = 0.0;
757  const G4ParticleDefinition* gamma = G4Gamma::Gamma();
758  res += ComputeCrossSectionPerVolume(kinEnergy, gamma, "conv", mat, 0.0);
759  res += ComputeCrossSectionPerVolume(kinEnergy, gamma, "compt", mat, 0.0);
760  res += ComputeCrossSectionPerVolume(kinEnergy, gamma, "phot", mat, 0.0);
761  res += ComputeCrossSectionPerVolume(kinEnergy, gamma, "Rayl", mat, 0.0);
762  if(res > 0.0) { res = 1.0/res; }
763  return res;
764 }
765 
766 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
767 
769  const G4String& particle,
770  G4int Z,
772  G4double kinEnergy,
773  const G4Material* mat)
774 {
775  G4double res = 0.0;
776  const G4ParticleDefinition* p = FindParticle(particle);
777  G4VAtomDeexcitation* ad = manager->AtomDeexcitation();
778  if(p && ad) {
779  res = ad->ComputeShellIonisationCrossSectionPerAtom(p, Z, shell,
780  kinEnergy, mat);
781  }
782  return res;
783 }
784 
785 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
786 
788  const G4ParticleDefinition* p,
789  const G4String& processName,
790  const G4Material* mat,
791  G4double cut)
792 {
793  G4double mfp = DBL_MAX;
794  G4double x = ComputeCrossSectionPerVolume(kinEnergy, p, processName, mat, cut);
795  if(x > 0.0) { mfp = 1.0/x; }
796  if(verbose>1) {
797  G4cout << "E(MeV)= " << kinEnergy/MeV
798  << " MFP(mm)= " << mfp/mm
799  << " " << p->GetParticleName()
800  << " in " << mat->GetName()
801  << G4endl;
802  }
803  return mfp;
804 }
805 
806 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
807 
809  const G4String& particle,
810  const G4String& processName,
811  const G4String& material,
812  G4double cut)
813 {
814  return ComputeMeanFreePath(kinEnergy,FindParticle(particle),processName,
815  FindMaterial(material),cut);
816 }
817 
818 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
819 
821  G4double range,
822  const G4ParticleDefinition* part,
823  const G4Material* mat)
824 {
826  ConvertRangeToEnergy(part, mat, range);
827 }
828 
829 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
830 
832  G4double range,
833  const G4String& particle,
834  const G4String& material)
835 {
836  return ComputeEnergyCutFromRangeCut(range,FindParticle(particle),
837  FindMaterial(material));
838 }
839 
840 
841 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
842 
843 G4bool G4EmCalculator::UpdateParticle(const G4ParticleDefinition* p,
844  G4double kinEnergy)
845 {
846  if(p != currentParticle) {
847 
848  // new particle
849  currentParticle = p;
850  dynParticle.SetDefinition(const_cast<G4ParticleDefinition*>(p));
851  dynParticle.SetKineticEnergy(kinEnergy);
852  baseParticle = 0;
853  currentParticleName = p->GetParticleName();
854  massRatio = 1.0;
855  mass = p->GetPDGMass();
856  chargeSquare = 1.0;
857  currentProcess = FindEnergyLossProcess(p);
858  currentProcessName = "";
859  isIon = false;
860 
861  // ionisation process exist
862  if(currentProcess) {
863  currentProcessName = currentProcess->GetProcessName();
864  baseParticle = currentProcess->BaseParticle();
865 
866  // base particle is used
867  if(baseParticle) {
868  massRatio = baseParticle->GetPDGMass()/p->GetPDGMass();
869  G4double q = p->GetPDGCharge()/baseParticle->GetPDGCharge();
870  chargeSquare = q*q;
871  }
872 
873  if(p->GetParticleType() == "nucleus"
874  && currentParticleName != "deuteron"
875  && currentParticleName != "triton"
876  && currentParticleName != "alpha+"
877  && currentParticleName != "helium"
878  && currentParticleName != "hydrogen"
879  ) {
880  isIon = true;
881  massRatio = theGenericIon->GetPDGMass()/p->GetPDGMass();
882  baseParticle = theGenericIon;
883  // G4cout << p->GetParticleName()
884  // << " in " << currentMaterial->GetName()
885  // << " e= " << kinEnergy << G4endl;
886  }
887  }
888  }
889 
890  // Effective charge for ions
891  if(isIon) {
892  chargeSquare =
893  corr->EffectiveChargeSquareRatio(p, currentMaterial, kinEnergy)
894  * corr->EffectiveChargeCorrection(p,currentMaterial,kinEnergy);
895  if(currentProcess) {
896  currentProcess->SetDynamicMassCharge(massRatio,chargeSquare);
897  //G4cout << "NewP: massR= " << massRatio << " q2= " << chargeSquare << G4endl;
898  }
899  }
900  return true;
901 }
902 
903 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
904 
906 {
907  const G4ParticleDefinition* p = 0;
908  if(name != currentParticleName) {
910  if(!p) {
911  G4cout << "### WARNING: G4EmCalculator::FindParticle fails to find "
912  << name << G4endl;
913  }
914  } else {
915  p = currentParticle;
916  }
917  return p;
918 }
919 
920 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
921 
923 {
924  const G4ParticleDefinition* p = ionTable->GetIon(Z,A,0);
925  return p;
926 }
927 
928 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
929 
931 {
932  if(name != currentMaterialName) {
933  currentMaterial = G4Material::GetMaterial(name, false);
934  currentMaterialName = name;
935  if(!currentMaterial) {
936  G4cout << "### WARNING: G4EmCalculator::FindMaterial fails to find "
937  << name << G4endl;
938  }
939  }
940  return currentMaterial;
941 }
942 
943 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
944 
946 {
947  const G4Region* r = 0;
948  if(reg != "" && reg != "world") {
950  } else {
951  r = G4RegionStore::GetInstance()->GetRegion("DefaultRegionForTheWorld");
952  }
953  return r;
954 }
955 
956 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
957 
959  const G4Material* material,
960  const G4Region* region)
961 {
962  const G4MaterialCutsCouple* couple = 0;
963  if(material) {
964  currentMaterial = material;
965  currentMaterialName = material->GetName();
966  // Access to materials
967  const G4ProductionCutsTable* theCoupleTable=
969  const G4Region* r = region;
970  if(r) {
971  couple = theCoupleTable->GetMaterialCutsCouple(material,
972  r->GetProductionCuts());
973  } else {
975  size_t nr = store->size();
976  if(0 < nr) {
977  for(size_t i=0; i<nr; ++i) {
978  couple = theCoupleTable->GetMaterialCutsCouple(
979  material, ((*store)[i])->GetProductionCuts());
980  if(couple) { break; }
981  }
982  }
983  }
984  }
985  if(!couple) {
987  ed << "G4EmCalculator::FindCouple: fail for material " << material
988  << " <" << currentMaterialName << " > and region " << region;
989  G4Exception("G4EmCalculator::FindCouple", "em0078",
990  FatalException, ed);
991  }
992  return couple;
993 }
994 
995 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
996 
997 G4bool G4EmCalculator::UpdateCouple(const G4Material* material, G4double cut)
998 {
999  if(!material) return false;
1000  currentMaterial = material;
1001  currentMaterialName = material->GetName();
1002  for (G4int i=0; i<nLocalMaterials; ++i) {
1003  if(material == localMaterials[i] && cut == localCuts[i]) {
1004  currentCouple = localCouples[i];
1005  currentCoupleIndex = currentCouple->GetIndex();
1006  currentCut = cut;
1007  return true;
1008  }
1009  }
1010  const G4MaterialCutsCouple* cc = new G4MaterialCutsCouple(material);
1011  localMaterials.push_back(material);
1012  localCouples.push_back(cc);
1013  localCuts.push_back(cut);
1014  nLocalMaterials++;
1015  currentCouple = cc;
1016  currentCoupleIndex = currentCouple->GetIndex();
1017  currentCut = cut;
1018  return true;
1019 }
1020 
1021 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1022 
1023 void G4EmCalculator::FindLambdaTable(const G4ParticleDefinition* p,
1024  const G4String& processName,
1025  G4double kinEnergy)
1026 {
1027  // Search for the process
1028  if (!currentLambda || p != lambdaParticle || processName != lambdaName) {
1029  lambdaName = processName;
1030  currentLambda = 0;
1031  lambdaParticle = p;
1032 
1033  const G4ParticleDefinition* part = p;
1034  if(isIon) { part = theGenericIon; }
1035 
1036  // Search for energy loss process
1037  currentName = processName;
1038  currentModel = 0;
1039  loweModel = 0;
1040 
1041  G4VEnergyLossProcess* elproc = FindEnLossProcess(part, processName);
1042  if(elproc) {
1043  currentLambda = elproc->LambdaTable();
1044  if(currentLambda) {
1045  isApplicable = true;
1046  if(verbose>1) {
1047  G4cout << "G4VEnergyLossProcess is found out: " << currentName
1048  << G4endl;
1049  }
1050  }
1051  return;
1052  }
1053 
1054  // Search for discrete process
1055  G4VEmProcess* proc = FindDiscreteProcess(part, processName);
1056  if(proc) {
1057  currentLambda = proc->LambdaTable();
1058  if(currentLambda) {
1059  isApplicable = true;
1060  if(verbose>1) {
1061  G4cout << "G4VEmProcess is found out: " << currentName << G4endl;
1062  }
1063  }
1064  return;
1065  }
1066 
1067  // Search for msc process
1068  G4VMultipleScattering* msc = FindMscProcess(part, processName);
1069  if(msc) {
1070  currentModel = msc->SelectModel(kinEnergy,0);
1071  /*
1072  if(currentModel) {
1073  currentLambda = currentModel->GetCrossSectionTable();
1074  if(currentLambda) {
1075  isApplicable = true;
1076  if(verbose>1) {
1077  G4cout << "G4VMultipleScattering is found out: " << currentName
1078  << G4endl;
1079  }
1080  }
1081  }
1082  */
1083  }
1084  }
1085 }
1086 
1087 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1088 
1089 G4bool G4EmCalculator::FindEmModel(const G4ParticleDefinition* p,
1090  const G4String& processName,
1091  G4double kinEnergy)
1092 {
1093  isApplicable = false;
1094  if(!p) {
1095  G4cout << "G4EmCalculator::FindEmModel WARNING: no particle defined"
1096  << G4endl;
1097  return isApplicable;
1098  }
1099  G4String partname = p->GetParticleName();
1100  const G4ParticleDefinition* part = p;
1101  G4double scaledEnergy = kinEnergy*massRatio;
1102  if(isIon) { part = theGenericIon; }
1103 
1104  if(verbose > 1) {
1105  G4cout << "## G4EmCalculator::FindEmModel for " << partname
1106  << " (type= " << p->GetParticleType()
1107  << ") and " << processName << " at E(MeV)= " << scaledEnergy
1108  << G4endl;
1109  if(p != part) { G4cout << " GenericIon is the base particle" << G4endl; }
1110  }
1111 
1112  // Search for energy loss process
1113  currentName = processName;
1114  currentModel = 0;
1115  loweModel = 0;
1116  size_t idx = 0;
1117 
1118  G4VEnergyLossProcess* elproc = FindEnLossProcess(part, processName);
1119  if(elproc) {
1120  currentModel = elproc->SelectModelForMaterial(scaledEnergy, idx);
1121  G4double eth = currentModel->LowEnergyLimit();
1122  if(eth > 0.0) {
1123  loweModel = elproc->SelectModelForMaterial(eth - CLHEP::eV, idx);
1124  if(loweModel == currentModel) { loweModel = 0; }
1125  }
1126  }
1127 
1128  // Search for discrete process
1129  if(!currentModel) {
1130  G4VEmProcess* proc = FindDiscreteProcess(part, processName);
1131  if(proc) {
1132  currentModel = proc->SelectModelForMaterial(kinEnergy, idx);
1133  G4double eth = currentModel->LowEnergyLimit();
1134  if(eth > 0.0) {
1135  loweModel = proc->SelectModelForMaterial(eth - CLHEP::eV, idx);
1136  if(loweModel == currentModel) { loweModel = 0; }
1137  }
1138  }
1139  }
1140 
1141  // Search for msc process
1142  if(!currentModel) {
1143  G4VMultipleScattering* proc = FindMscProcess(part, processName);
1144  if(proc) {
1145  currentModel = proc->SelectModel(kinEnergy, idx);
1146  loweModel = 0;
1147  }
1148  }
1149  if(currentModel) {
1150  if(loweModel == currentModel) { loweModel = 0; }
1151  isApplicable = true;
1152  currentModel->InitialiseForMaterial(part, currentMaterial);
1153  if(loweModel) {
1154  loweModel->InitialiseForMaterial(part, currentMaterial);
1155  }
1156  if(verbose > 1) {
1157  G4cout << " Model <" << currentModel->GetName()
1158  << "> Emin(MeV)= " << currentModel->LowEnergyLimit()/MeV
1159  << " for " << part->GetParticleName();
1160  if(elproc) {
1161  G4cout << " and " << elproc->GetProcessName() << " " << elproc
1162  << G4endl;
1163  }
1164  if(loweModel) {
1165  G4cout << " LowEnergy model <" << loweModel->GetName() << ">";
1166  }
1167  G4cout << G4endl;
1168  }
1169  }
1170  return isApplicable;
1171 }
1172 
1173 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1174 
1175 G4VEnergyLossProcess* G4EmCalculator::FindEnergyLossProcess(
1176  const G4ParticleDefinition* p)
1177 {
1178  G4VEnergyLossProcess* elp = 0;
1179  G4String partname = p->GetParticleName();
1180  const G4ParticleDefinition* part = p;
1181 
1182  if(p->GetParticleType() == "nucleus"
1183  && currentParticleName != "deuteron"
1184  && currentParticleName != "triton"
1185  && currentParticleName != "alpha+"
1186  && currentParticleName != "helium"
1187  && currentParticleName != "hydrogen"
1188  ) { part = theGenericIon; }
1189 
1191  return elp;
1192 }
1193 
1194 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1195 
1197 G4EmCalculator::FindEnLossProcess(const G4ParticleDefinition* part,
1198  const G4String& processName)
1199 {
1200  G4VEnergyLossProcess* proc = 0;
1201  const std::vector<G4VEnergyLossProcess*> v =
1203  G4int n = v.size();
1204  for(G4int i=0; i<n; ++i) {
1205  if((v[i])->GetProcessName() == processName) {
1206  G4VProcess* p = reinterpret_cast<G4VProcess*>(v[i]);
1207  if(ActiveForParticle(part, p)) {
1208  proc = v[i];
1209  break;
1210  }
1211  }
1212  }
1213  return proc;
1214 }
1215 
1216 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1217 
1218 G4VEmProcess*
1219 G4EmCalculator::FindDiscreteProcess(const G4ParticleDefinition* part,
1220  const G4String& processName)
1221 {
1222  G4VEmProcess* proc = 0;
1223  const std::vector<G4VEmProcess*> v =
1225  G4int n = v.size();
1226  for(G4int i=0; i<n; ++i) {
1227  if((v[i])->GetProcessName() == processName) {
1228  G4VProcess* p = reinterpret_cast<G4VProcess*>(v[i]);
1229  if(ActiveForParticle(part, p)) {
1230  proc = v[i];
1231  break;
1232  }
1233  }
1234  }
1235  return proc;
1236 }
1237 
1238 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1239 
1241 G4EmCalculator::FindMscProcess(const G4ParticleDefinition* part,
1242  const G4String& processName)
1243 {
1244  G4VMultipleScattering* proc = 0;
1245  const std::vector<G4VMultipleScattering*> v =
1247  G4int n = v.size();
1248  for(G4int i=0; i<n; ++i) {
1249  if((v[i])->GetProcessName() == processName) {
1250  G4VProcess* p = reinterpret_cast<G4VProcess*>(v[i]);
1251  if(ActiveForParticle(part, p)) {
1252  proc = v[i];
1253  break;
1254  }
1255  }
1256  }
1257  return proc;
1258 }
1259 
1260 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1261 
1262 G4bool G4EmCalculator::ActiveForParticle(const G4ParticleDefinition* part,
1263  G4VProcess* proc)
1264 {
1265  G4ProcessManager* pm = part->GetProcessManager();
1266  G4ProcessVector* pv = pm->GetProcessList();
1267  G4int n = pv->size();
1268  G4bool res = false;
1269  for(G4int i=0; i<n; ++i) {
1270  if((*pv)[i] == proc) {
1271  if(pm->GetProcessActivation(i)) { res = true; }
1272  break;
1273  }
1274  }
1275  return res;
1276 }
1277 
1278 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1279 
1281 {
1282  verbose = verb;
1283 }
1284 
1285 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1286 
G4PhysicsTable * LambdaTable() const
G4ProductionCuts * GetProductionCuts() const
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:599
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:245
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:231
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4LossTableManager * Instance()
virtual void CorrectionsAlongStep(const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
Definition: G4VEmModel.cc:339
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 &)
G4double GetN() const
Definition: G4Element.hh:134
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:578
G4PhysicsTable * RangeTableForLoss() const
G4double EffectiveChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
const char * p
Definition: xmltok.h:285
const G4String & GetName() const
Definition: G4Material.hh:176
void PrintRangeTable(const G4ParticleDefinition *)
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:449
virtual G4double GetShellIonisationCrossSectionPerAtom(const G4ParticleDefinition *, G4int Z, G4AtomicShellEnumerator shell, G4double kinE, const G4Material *mat=0)=0
G4double GetZ() const
Definition: G4Element.hh:131
G4double GetDensity() const
Definition: G4Material.hh:178
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)
const XML_Char * name
G4double EffectiveChargeCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
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 PrintInverseRangeTable(const G4ParticleDefinition *)
G4double GetRange(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=0)
G4ProcessManager * GetProcessManager() const
const std::vector< G4VEmProcess * > & GetEmProcessVector()
int G4int
Definition: G4Types.hh:78
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 G4RegionStore * GetInstance()
string material
Definition: eplot.py:19
G4double GetCSDARange(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=0)
G4double GetShellIonisationCrossSectionPerAtom(const G4String &part, G4int Z, G4AtomicShellEnumerator shell, G4double kinEnergy)
G4PhysicsTable * LambdaTable() const
function g(Y1, Y2, PT2)
Definition: hijing1.383.f:5205
G4IonTable * GetIonTable() const
virtual void InitialiseForMaterial(const G4ParticleDefinition *, const G4Material *)
Definition: G4VEmModel.cc:212
G4GLOB_DLL std::ostream G4cout
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)
G4double NuclearDEDX(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy, G4bool fluct=true)
bool G4bool
Definition: G4Types.hh:79
G4EmCorrections * EmCorrections()
G4double ComputeNuclearDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *)
G4double GetRangeFromRestricteDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=0)
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
const std::vector< G4VEnergyLossProcess * > & GetEnergyLossProcessVector()
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:300
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
G4int size() const
virtual G4double GetChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:322
G4VEnergyLossProcess * GetEnergyLossProcess(const G4ParticleDefinition *)
static G4ProductionCutsTable * GetProductionCutsTable()
static G4GenericIon * GenericIon()
Definition: G4GenericIon.cc:93
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
Definition: G4VEmModel.cc:236
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
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 GetKinEnergy(G4double range, const G4ParticleDefinition *, const G4Material *, const G4Region *r=0)
const G4Region * FindRegion(const G4String &)
const G4ParticleDefinition * FindIon(G4int Z, G4int A)
G4double ComputeCrossSectionPerVolume(G4double kinEnergy, const G4ParticleDefinition *, const G4String &processName, const G4Material *, G4double cut=0.0)
G4double ComputeGammaAttenuationLength(G4double kinEnergy, const G4Material *)
G4bool GetProcessActivation(G4VProcess *aProcess) const
#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:753
void SetVerbose(G4int val)
void PrintDEDXTable(const G4ParticleDefinition *)
G4VAtomDeexcitation * AtomDeexcitation()
double G4double
Definition: G4Types.hh:76
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
G4double ComputeEnergyCutFromRangeCut(G4double range, const G4ParticleDefinition *, const G4Material *)
G4double GetPDGCharge() const
#define DBL_MAX
Definition: templates.hh:83
G4double GetRange(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4AtomicShellEnumerator
const std::vector< G4VMultipleScattering * > & GetMultipleScatteringVector()
G4VEmModel * SelectModelForMaterial(G4double kinEnergy, size_t &idx) const
G4ProcessVector * GetProcessList() const