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