Geant4  10.01.p03
G4IonParametrisedLossModel.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: G4IonParametrisedLossModel.cc 87443 2014-12-04 12:26:31Z gunter $
27 //
28 // ===========================================================================
29 // GEANT4 class source file
30 //
31 // Class: G4IonParametrisedLossModel
32 //
33 // Base class: G4VEmModel (utils)
34 //
35 // Author: Anton Lechner (Anton.Lechner@cern.ch)
36 //
37 // First implementation: 10. 11. 2008
38 //
39 // Modifications: 03. 02. 2009 - Bug fix iterators (AL)
40 // 11. 03. 2009 - Introduced new table handler(G4IonDEDXHandler)
41 // and modified method to add/remove tables
42 // (tables are now built in init. phase),
43 // Minor bug fix in ComputeDEDXPerVolume (AL)
44 // 11. 05. 2009 - Introduced scaling algorithm for heavier ions:
45 // G4IonDEDXScalingICRU73 (AL)
46 // 12. 11. 2009 - Moved from original ICRU 73 classes to new
47 // class (G4IonStoppingData), which is capable
48 // of reading stopping power data files stored
49 // in G4LEDATA (requires G4EMLOW6.8 or higher).
50 // Simultanesouly, the upper energy limit of
51 // ICRU 73 is increased to 1 GeV/nucleon.
52 // - Removed nuclear stopping from Corrections-
53 // AlongStep since dedicated process was created.
54 // - Added function for switching off scaling
55 // of heavy ions from ICRU 73 data
56 // - Minor fix in ComputeLossForStep function
57 // - Minor fix in ComputeDEDXPerVolume (AL)
58 // 23. 11. 2009 - Changed energy loss limit from 0.15 to 0.01
59 // to improve accuracy for large steps (AL)
60 // 24. 11. 2009 - Bug fix: Range calculation corrected if same
61 // materials appears with different cuts in diff.
62 // regions (added UpdateRangeCache function and
63 // modified BuildRangeVector, ComputeLossForStep
64 // functions accordingly, added new cache param.)
65 // - Removed GetRange function (AL)
66 // 04. 11. 2010 - Moved virtual methods to the source (VI)
67 //
68 //
69 // Class description:
70 // Model for computing the energy loss of ions by employing a
71 // parameterisation of dE/dx tables (by default ICRU 73 tables). For
72 // ion-material combinations and/or projectile energies not covered
73 // by this model, the G4BraggIonModel and G4BetheBloch models are
74 // employed.
75 //
76 // Comments:
77 //
78 // ===========================================================================
79 
80 
82 #include "G4PhysicalConstants.hh"
83 #include "G4SystemOfUnits.hh"
84 #include "G4LPhysicsFreeVector.hh"
85 #include "G4IonStoppingData.hh"
86 #include "G4VIonDEDXTable.hh"
89 #include "G4BraggIonModel.hh"
90 #include "G4BetheBlochModel.hh"
91 #include "G4ProductionCutsTable.hh"
93 #include "G4LossTableManager.hh"
94 #include "G4GenericIon.hh"
95 #include "G4Electron.hh"
96 #include "G4DeltaAngle.hh"
97 #include "Randomize.hh"
98 
99 //#define PRINT_TABLE_BUILT
100 
101 
102 // #########################################################################
103 
105  const G4ParticleDefinition*,
106  const G4String& nam)
107  : G4VEmModel(nam),
108  braggIonModel(0),
109  betheBlochModel(0),
110  nmbBins(90),
111  nmbSubBins(100),
112  particleChangeLoss(0),
113  corrFactor(1.0),
114  energyLossLimit(0.01),
115  cutEnergies(0)
116 {
118  genericIonPDGMass = genericIon -> GetPDGMass();
119  corrections = G4LossTableManager::Instance() -> EmCorrections();
120 
121  // The upper limit of the current model is set to 100 TeV
122  SetHighEnergyLimit(100.0 * TeV);
123 
124  // The Bragg ion and Bethe Bloch models are instantiated
127 
128  // By default ICRU 73 stopping power tables are loaded:
129  AddDEDXTable("ICRU73",
130  new G4IonStoppingData("ion_stopping_data/icru73"),
131  new G4IonDEDXScalingICRU73());
132 
133  // The boundaries for the range tables are set
134  lowerEnergyEdgeIntegr = 0.025 * MeV;
136 
137  // Cache parameters are set
138  cacheParticle = 0;
139  cacheMass = 0;
140  cacheElecMassRatio = 0;
141  cacheChargeSquare = 0;
142 
143  // Cache parameters are set
144  rangeCacheParticle = 0;
148 
149  // Cache parameters are set
150  dedxCacheParticle = 0;
151  dedxCacheMaterial = 0;
152  dedxCacheEnergyCut = 0;
157 
158  // default generator
160 }
161 
162 // #########################################################################
163 
165 
166  // dE/dx table objects are deleted and the container is cleared
167  LossTableList::iterator iterTables = lossTableList.begin();
168  LossTableList::iterator iterTables_end = lossTableList.end();
169 
170  for(;iterTables != iterTables_end; ++iterTables) { delete *iterTables; }
171  lossTableList.clear();
172 
173  // range table
174  RangeEnergyTable::iterator itr = r.begin();
175  RangeEnergyTable::iterator itr_end = r.end();
176  for(;itr != itr_end; ++itr) { delete itr->second; }
177  r.clear();
178 
179  // inverse range
180  EnergyRangeTable::iterator ite = E.begin();
181  EnergyRangeTable::iterator ite_end = E.end();
182  for(;ite != ite_end; ++ite) { delete ite->second; }
183  E.clear();
184 
185 }
186 
187 // #########################################################################
188 
190  const G4ParticleDefinition*,
191  const G4MaterialCutsCouple* couple) {
192 
193  return couple -> GetMaterial() -> GetIonisation() ->
194  GetMeanExcitationEnergy();
195 }
196 
197 // #########################################################################
198 
200  const G4ParticleDefinition* particle,
201  G4double kineticEnergy) {
202 
203  // ############## Maximum energy of secondaries ##########################
204  // Function computes maximum energy of secondary electrons which are
205  // released by an ion
206  //
207  // See Geant4 physics reference manual (version 9.1), section 9.1.1
208  //
209  // Ref.: W.M. Yao et al, Jour. of Phys. G 33 (2006) 1.
210  // C.Caso et al. (Part. Data Group), Europ. Phys. Jour. C 3 1 (1998).
211  // B. Rossi, High energy particles, New York, NY: Prentice-Hall (1952).
212  //
213  // (Implementation adapted from G4BraggIonModel)
214 
215  if(particle != cacheParticle) UpdateCache(particle);
216 
217  G4double tau = kineticEnergy/cacheMass;
218  G4double tmax = 2.0 * electron_mass_c2 * tau * (tau + 2.) /
219  (1. + 2.0 * (tau + 1.) * cacheElecMassRatio +
221 
222  return tmax;
223 }
224 
225 // #########################################################################
226 
228  const G4ParticleDefinition* particle,
229  const G4Material* material,
230  G4double kineticEnergy) { // Kinetic energy
231 
232  G4double chargeSquareRatio = corrections ->
233  EffectiveChargeSquareRatio(particle,
234  material,
235  kineticEnergy);
236  corrFactor = chargeSquareRatio *
237  corrections -> EffectiveChargeCorrection(particle,
238  material,
239  kineticEnergy);
240  return corrFactor;
241 }
242 
243 // #########################################################################
244 
246  const G4ParticleDefinition* particle,
247  const G4Material* material,
248  G4double kineticEnergy) { // Kinetic energy
249 
250  return corrections -> GetParticleCharge(particle, material, kineticEnergy);
251 }
252 
253 // #########################################################################
254 
256  const G4ParticleDefinition* particle,
257  const G4DataVector& cuts) {
258 
259  // Cached parameters are reset
260  cacheParticle = 0;
261  cacheMass = 0;
262  cacheElecMassRatio = 0;
263  cacheChargeSquare = 0;
264 
265  // Cached parameters are reset
266  rangeCacheParticle = 0;
270 
271  // Cached parameters are reset
272  dedxCacheParticle = 0;
273  dedxCacheMaterial = 0;
274  dedxCacheEnergyCut = 0;
279 
280  // The cache of loss tables is cleared
281  LossTableList::iterator iterTables = lossTableList.begin();
282  LossTableList::iterator iterTables_end = lossTableList.end();
283 
284  for(;iterTables != iterTables_end; iterTables++)
285  (*iterTables) -> ClearCache();
286 
287  // Range vs energy and energy vs range vectors from previous runs are
288  // cleared
289  RangeEnergyTable::iterator iterRange = r.begin();
290  RangeEnergyTable::iterator iterRange_end = r.end();
291 
292  for(;iterRange != iterRange_end; iterRange++) {
293  delete iterRange->second;
294  }
295  r.clear();
296 
297  EnergyRangeTable::iterator iterEnergy = E.begin();
298  EnergyRangeTable::iterator iterEnergy_end = E.end();
299 
300  for(;iterEnergy != iterEnergy_end; iterEnergy++) {
301  delete iterEnergy->second;
302  }
303  E.clear();
304 
305  // The cut energies are (re)loaded
306  size_t size = cuts.size();
307  cutEnergies.clear();
308  for(size_t i = 0; i < size; i++) cutEnergies.push_back(cuts[i]);
309 
310  // All dE/dx vectors are built
311  const G4ProductionCutsTable* coupleTable=
313  size_t nmbCouples = coupleTable -> GetTableSize();
314 
315 #ifdef PRINT_TABLE_BUILT
316  G4cout << "G4IonParametrisedLossModel::Initialise():"
317  << " Building dE/dx vectors:"
318  << G4endl;
319 #endif
320 
321  for (size_t i = 0; i < nmbCouples; i++) {
322 
323  const G4MaterialCutsCouple* couple =
324  coupleTable -> GetMaterialCutsCouple(i);
325 
326  const G4Material* material = couple -> GetMaterial();
327  // G4ProductionCuts* productionCuts = couple -> GetProductionCuts();
328 
329  for(G4int atomicNumberIon = 3; atomicNumberIon < 102; atomicNumberIon++) {
330 
331  LossTableList::iterator iter = lossTableList.begin();
332  LossTableList::iterator iter_end = lossTableList.end();
333 
334  for(;iter != iter_end; iter++) {
335 
336  if(*iter == 0) {
337  G4cout << "G4IonParametrisedLossModel::Initialise():"
338  << " Skipping illegal table."
339  << G4endl;
340  }
341 
342  G4bool isApplicable =
343  (*iter) -> BuildDEDXTable(atomicNumberIon, material);
344  if(isApplicable) {
345 
346 #ifdef PRINT_TABLE_BUILT
347  G4cout << " Atomic Number Ion = " << atomicNumberIon
348  << ", Material = " << material -> GetName()
349  << ", Table = " << (*iter) -> GetName()
350  << G4endl;
351 #endif
352  break;
353  }
354  }
355  }
356  }
357 
358  // The particle change object
359  if(! particleChangeLoss) {
363  }
364 
365  // The G4BraggIonModel and G4BetheBlochModel instances are initialised with
366  // the same settings as the current model:
367  braggIonModel -> Initialise(particle, cuts);
368  betheBlochModel -> Initialise(particle, cuts);
369 }
370 
371 // #########################################################################
372 
374  const G4ParticleDefinition* particle,
375  G4double kineticEnergy,
376  G4double atomicNumber,
377  G4double,
378  G4double cutEnergy,
379  G4double maxKinEnergy) {
380 
381  // ############## Cross section per atom ################################
382  // Function computes ionization cross section per atom
383  //
384  // See Geant4 physics reference manual (version 9.1), section 9.1.3
385  //
386  // Ref.: W.M. Yao et al, Jour. of Phys. G 33 (2006) 1.
387  // B. Rossi, High energy particles, New York, NY: Prentice-Hall (1952).
388  //
389  // (Implementation adapted from G4BraggIonModel)
390 
391  G4double crosssection = 0.0;
392  G4double tmax = MaxSecondaryEnergy(particle, kineticEnergy);
393  G4double maxEnergy = std::min(tmax, maxKinEnergy);
394 
395  if(cutEnergy < tmax) {
396 
397  G4double energy = kineticEnergy + cacheMass;
398  G4double betaSquared = kineticEnergy *
399  (energy + cacheMass) / (energy * energy);
400 
401  crosssection = 1.0 / cutEnergy - 1.0 / maxEnergy -
402  betaSquared * std::log(maxEnergy / cutEnergy) / tmax;
403 
404  crosssection *= twopi_mc2_rcl2 * cacheChargeSquare / betaSquared;
405  }
406 
407 #ifdef PRINT_DEBUG_CS
408  G4cout << "########################################################"
409  << G4endl
410  << "# G4IonParametrisedLossModel::ComputeCrossSectionPerAtom"
411  << G4endl
412  << "# particle =" << particle -> GetParticleName()
413  << G4endl
414  << "# cut(MeV) = " << cutEnergy/MeV
415  << G4endl;
416 
417  G4cout << "#"
418  << std::setw(13) << std::right << "E(MeV)"
419  << std::setw(14) << "CS(um)"
420  << std::setw(14) << "E_max_sec(MeV)"
421  << G4endl
422  << "# ------------------------------------------------------"
423  << G4endl;
424 
425  G4cout << std::setw(14) << std::right << kineticEnergy / MeV
426  << std::setw(14) << crosssection / (um * um)
427  << std::setw(14) << tmax / MeV
428  << G4endl;
429 #endif
430 
431  crosssection *= atomicNumber;
432 
433  return crosssection;
434 }
435 
436 // #########################################################################
437 
439  const G4Material* material,
440  const G4ParticleDefinition* particle,
441  G4double kineticEnergy,
442  G4double cutEnergy,
443  G4double maxEnergy) {
444 
445  G4double nbElecPerVolume = material -> GetTotNbOfElectPerVolume();
446  G4double cross = ComputeCrossSectionPerAtom(particle,
447  kineticEnergy,
448  nbElecPerVolume, 0,
449  cutEnergy,
450  maxEnergy);
451 
452  return cross;
453 }
454 
455 // #########################################################################
456 
458  const G4Material* material,
459  const G4ParticleDefinition* particle,
460  G4double kineticEnergy,
461  G4double cutEnergy) {
462 
463  // ############## dE/dx ##################################################
464  // Function computes dE/dx values, where following rules are adopted:
465  // A. If the ion-material pair is covered by any native ion data
466  // parameterisation, then:
467  // * This parameterization is used for energies below a given energy
468  // limit,
469  // * whereas above the limit the Bethe-Bloch model is applied, in
470  // combination with an effective charge estimate and high order
471  // correction terms.
472  // A smoothing procedure is applied to dE/dx values computed with
473  // the second approach. The smoothing factor is based on the dE/dx
474  // values of both approaches at the transition energy (high order
475  // correction terms are included in the calculation of the transition
476  // factor).
477  // B. If the particle is a generic ion, the BraggIon and Bethe-Bloch
478  // models are used and a smoothing procedure is applied to values
479  // obtained with the second approach.
480  // C. If the ion-material is not covered by any ion data parameterization
481  // then:
482  // * The BraggIon model is used for energies below a given energy
483  // limit,
484  // * whereas above the limit the Bethe-Bloch model is applied, in
485  // combination with an effective charge estimate and high order
486  // correction terms.
487  // Also in this case, a smoothing procedure is applied to dE/dx values
488  // computed with the second model.
489 
490  G4double dEdx = 0.0;
491  UpdateDEDXCache(particle, material, cutEnergy);
492 
493  LossTableList::iterator iter = dedxCacheIter;
494 
495  if(iter != lossTableList.end()) {
496 
497  G4double transitionEnergy = dedxCacheTransitionEnergy;
498 
499  if(transitionEnergy > kineticEnergy) {
500 
501  dEdx = (*iter) -> GetDEDX(particle, material, kineticEnergy);
502 
503  G4double dEdxDeltaRays = DeltaRayMeanEnergyTransferRate(material,
504  particle,
505  kineticEnergy,
506  cutEnergy);
507  dEdx -= dEdxDeltaRays;
508  }
509  else {
511 
512  G4double chargeSquare =
513  GetChargeSquareRatio(particle, material, kineticEnergy);
514 
515  G4double scaledKineticEnergy = kineticEnergy * massRatio;
516  G4double scaledTransitionEnergy = transitionEnergy * massRatio;
517 
518  G4double lowEnergyLimit = betheBlochModel -> LowEnergyLimit();
519 
520  if(scaledTransitionEnergy >= lowEnergyLimit) {
521 
523  material, genericIon,
524  scaledKineticEnergy, cutEnergy);
525 
526  dEdx *= chargeSquare;
527 
528  dEdx += corrections -> ComputeIonCorrections(particle,
529  material, kineticEnergy);
530 
532  kineticEnergy;
533 
534  dEdx *= factor;
535  }
536  }
537  }
538  else {
539  G4double massRatio = 1.0;
540  G4double chargeSquare = 1.0;
541 
542  if(particle != genericIon) {
543 
544  chargeSquare = GetChargeSquareRatio(particle, material, kineticEnergy);
545  massRatio = genericIonPDGMass / particle -> GetPDGMass();
546  }
547 
548  G4double scaledKineticEnergy = kineticEnergy * massRatio;
549 
550  G4double lowEnergyLimit = betheBlochModel -> LowEnergyLimit();
551  if(scaledKineticEnergy < lowEnergyLimit) {
553  material, genericIon,
554  scaledKineticEnergy, cutEnergy);
555 
556  dEdx *= chargeSquare;
557  }
558  else {
559  G4double dEdxLimitParam = braggIonModel -> ComputeDEDXPerVolume(
560  material, genericIon,
561  lowEnergyLimit, cutEnergy);
562 
563  G4double dEdxLimitBetheBloch = betheBlochModel -> ComputeDEDXPerVolume(
564  material, genericIon,
565  lowEnergyLimit, cutEnergy);
566 
567  if(particle != genericIon) {
568  G4double chargeSquareLowEnergyLimit =
569  GetChargeSquareRatio(particle, material,
570  lowEnergyLimit / massRatio);
571 
572  dEdxLimitParam *= chargeSquareLowEnergyLimit;
573  dEdxLimitBetheBloch *= chargeSquareLowEnergyLimit;
574 
575  dEdxLimitBetheBloch +=
576  corrections -> ComputeIonCorrections(particle,
577  material, lowEnergyLimit / massRatio);
578  }
579 
580  G4double factor = (1.0 + (dEdxLimitParam/dEdxLimitBetheBloch - 1.0)
581  * lowEnergyLimit / scaledKineticEnergy);
582 
584  material, genericIon,
585  scaledKineticEnergy, cutEnergy);
586 
587  dEdx *= chargeSquare;
588 
589  if(particle != genericIon) {
590  dEdx += corrections -> ComputeIonCorrections(particle,
591  material, kineticEnergy);
592  }
593 
594  dEdx *= factor;
595  }
596 
597  }
598 
599  if (dEdx < 0.0) dEdx = 0.0;
600 
601  return dEdx;
602 }
603 
604 // #########################################################################
605 
607  const G4ParticleDefinition* particle, // Projectile (ion)
608  const G4Material* material, // Absorber material
609  G4double lowerBoundary, // Minimum energy per nucleon
610  G4double upperBoundary, // Maximum energy per nucleon
611  G4int numBins, // Number of bins
612  G4bool logScaleEnergy) { // Logarithmic scaling of energy
613 
614  G4double atomicMassNumber = particle -> GetAtomicMass();
615  G4double materialDensity = material -> GetDensity();
616 
617  G4cout << "# dE/dx table for " << particle -> GetParticleName()
618  << " in material " << material -> GetName()
619  << " of density " << materialDensity / g * cm3
620  << " g/cm3"
621  << G4endl
622  << "# Projectile mass number A1 = " << atomicMassNumber
623  << G4endl
624  << "# ------------------------------------------------------"
625  << G4endl;
626  G4cout << "#"
627  << std::setw(13) << std::right << "E"
628  << std::setw(14) << "E/A1"
629  << std::setw(14) << "dE/dx"
630  << std::setw(14) << "1/rho*dE/dx"
631  << G4endl;
632  G4cout << "#"
633  << std::setw(13) << std::right << "(MeV)"
634  << std::setw(14) << "(MeV)"
635  << std::setw(14) << "(MeV/cm)"
636  << std::setw(14) << "(MeV*cm2/mg)"
637  << G4endl
638  << "# ------------------------------------------------------"
639  << G4endl;
640 
641  G4double energyLowerBoundary = lowerBoundary * atomicMassNumber;
642  G4double energyUpperBoundary = upperBoundary * atomicMassNumber;
643 
644  if(logScaleEnergy) {
645 
646  energyLowerBoundary = std::log(energyLowerBoundary);
647  energyUpperBoundary = std::log(energyUpperBoundary);
648  }
649 
650  G4double deltaEnergy = (energyUpperBoundary - energyLowerBoundary) /
651  G4double(nmbBins);
652 
653  for(int i = 0; i < numBins + 1; i++) {
654 
655  G4double energy = energyLowerBoundary + i * deltaEnergy;
656  if(logScaleEnergy) energy = std::exp(energy);
657 
658  G4double dedx = ComputeDEDXPerVolume(material, particle, energy, DBL_MAX);
659  G4cout.precision(6);
660  G4cout << std::setw(14) << std::right << energy / MeV
661  << std::setw(14) << energy / atomicMassNumber / MeV
662  << std::setw(14) << dedx / MeV * cm
663  << std::setw(14) << dedx / materialDensity / (MeV*cm2/(0.001*g))
664  << G4endl;
665  }
666 }
667 
668 // #########################################################################
669 
671  const G4ParticleDefinition* particle, // Projectile (ion)
672  const G4Material* material, // Absorber material
673  G4double lowerBoundary, // Minimum energy per nucleon
674  G4double upperBoundary, // Maximum energy per nucleon
675  G4int numBins, // Number of bins
676  G4bool logScaleEnergy) { // Logarithmic scaling of energy
677 
678  LossTableList::iterator iter = lossTableList.begin();
679  LossTableList::iterator iter_end = lossTableList.end();
680 
681  for(;iter != iter_end; iter++) {
682  G4bool isApplicable = (*iter) -> IsApplicable(particle, material);
683  if(isApplicable) {
684  (*iter) -> PrintDEDXTable(particle, material,
685  lowerBoundary, upperBoundary,
686  numBins,logScaleEnergy);
687  break;
688  }
689  }
690 }
691 
692 // #########################################################################
693 
695  std::vector<G4DynamicParticle*>* secondaries,
696  const G4MaterialCutsCouple* couple,
697  const G4DynamicParticle* particle,
698  G4double cutKinEnergySec,
699  G4double userMaxKinEnergySec) {
700 
701 
702  // ############## Sampling of secondaries #################################
703  // The probability density function (pdf) of the kinetic energy T of a
704  // secondary electron may be written as:
705  // pdf(T) = f(T) * g(T)
706  // where
707  // f(T) = (Tmax - Tcut) / (Tmax * Tcut) * (1 / T^2)
708  // g(T) = 1 - beta^2 * T / Tmax
709  // where Tmax is the maximum kinetic energy of the secondary, Tcut
710  // is the lower energy cut and beta is the kinetic energy of the
711  // projectile.
712  //
713  // Sampling of the kinetic energy of a secondary electron:
714  // 1) T0 is sampled from f(T) using the cumulated distribution function
715  // F(T) = int_Tcut^T f(T')dT'
716  // 2) T is accepted or rejected by evaluating the rejection function g(T)
717  // at the sampled energy T0 against a randomly sampled value
718  //
719  //
720  // See Geant4 physics reference manual (version 9.1), section 9.1.4
721  //
722  //
723  // Reference pdf: W.M. Yao et al, Jour. of Phys. G 33 (2006) 1.
724  //
725  // (Implementation adapted from G4BraggIonModel)
726 
727  G4double rossiMaxKinEnergySec = MaxSecondaryKinEnergy(particle);
728  G4double maxKinEnergySec =
729  std::min(rossiMaxKinEnergySec, userMaxKinEnergySec);
730 
731  if(cutKinEnergySec >= maxKinEnergySec) return;
732 
733  G4double kineticEnergy = particle -> GetKineticEnergy();
734 
735  G4double energy = kineticEnergy + cacheMass;
736  G4double betaSquared = kineticEnergy * (energy + cacheMass)
737  / (energy * energy);
738 
739  G4double kinEnergySec;
740  G4double grej;
741 
742  do {
743 
744  // Sampling kinetic energy from f(T) (using F(T)):
745  G4double xi = G4UniformRand();
746  kinEnergySec = cutKinEnergySec * maxKinEnergySec /
747  (maxKinEnergySec * (1.0 - xi) + cutKinEnergySec * xi);
748 
749  // Deriving the value of the rejection function at the obtained kinetic
750  // energy:
751  grej = 1.0 - betaSquared * kinEnergySec / rossiMaxKinEnergySec;
752 
753  if(grej > 1.0) {
754  G4cout << "G4IonParametrisedLossModel::SampleSecondary Warning: "
755  << "Majorant 1.0 < "
756  << grej << " for e= " << kinEnergySec
757  << G4endl;
758  }
759 
760  } while( G4UniformRand() >= grej );
761 
762  const G4Material* mat = couple->GetMaterial();
763  G4int Z = SelectRandomAtomNumber(mat);
764 
766 
767  G4DynamicParticle* delta = new G4DynamicParticle(electron,
768  GetAngularDistribution()->SampleDirection(particle, kinEnergySec,
769  Z, mat),
770  kinEnergySec);
771 
772 
773  secondaries->push_back(delta);
774 
775  // Change kinematics of primary particle
776  G4ThreeVector direction = particle ->GetMomentumDirection();
777  G4double totalMomentum = std::sqrt(kineticEnergy*(energy + cacheMass));
778 
779  G4ThreeVector finalP = totalMomentum*direction - delta->GetMomentum();
780  finalP = finalP.unit();
781 
782  kineticEnergy -= kinEnergySec;
783 
786 }
787 
788 // #########################################################################
789 
791  const G4ParticleDefinition* particle,
792  const G4MaterialCutsCouple* matCutsCouple) {
793 
794  // ############## Caching ##################################################
795  // If the ion-material-cut combination is covered by any native ion data
796  // parameterisation (for low energies), range vectors are computed
797 
798  if(particle == rangeCacheParticle &&
799  matCutsCouple == rangeCacheMatCutsCouple) {
800  }
801  else{
802  rangeCacheParticle = particle;
803  rangeCacheMatCutsCouple = matCutsCouple;
804 
805  const G4Material* material = matCutsCouple -> GetMaterial();
806  LossTableList::iterator iter = IsApplicable(particle, material);
807 
808  // If any table is applicable, the transition factor is computed:
809  if(iter != lossTableList.end()) {
810 
811  // Build range-energy and energy-range vectors if they don't exist
812  IonMatCouple ionMatCouple = std::make_pair(particle, matCutsCouple);
813  RangeEnergyTable::iterator iterRange = r.find(ionMatCouple);
814 
815  if(iterRange == r.end()) BuildRangeVector(particle, matCutsCouple);
816 
817  rangeCacheEnergyRange = E[ionMatCouple];
818  rangeCacheRangeEnergy = r[ionMatCouple];
819  }
820  else {
823  }
824  }
825 }
826 
827 // #########################################################################
828 
830  const G4ParticleDefinition* particle,
831  const G4Material* material,
832  G4double cutEnergy) {
833 
834  // ############## Caching ##################################################
835  // If the ion-material combination is covered by any native ion data
836  // parameterisation (for low energies), a transition factor is computed
837  // which is applied to Bethe-Bloch results at higher energies to
838  // guarantee a smooth transition.
839  // This factor only needs to be calculated for the first step an ion
840  // performs inside a certain material.
841 
842  if(particle == dedxCacheParticle &&
843  material == dedxCacheMaterial &&
844  cutEnergy == dedxCacheEnergyCut) {
845  }
846  else {
847 
848  dedxCacheParticle = particle;
849  dedxCacheMaterial = material;
850  dedxCacheEnergyCut = cutEnergy;
851 
852  G4double massRatio = genericIonPDGMass / particle -> GetPDGMass();
853  dedxCacheGenIonMassRatio = massRatio;
854 
855  LossTableList::iterator iter = IsApplicable(particle, material);
856  dedxCacheIter = iter;
857 
858  // If any table is applicable, the transition factor is computed:
859  if(iter != lossTableList.end()) {
860 
861  // Retrieving the transition energy from the parameterisation table
862  G4double transitionEnergy =
863  (*iter) -> GetUpperEnergyEdge(particle, material);
864  dedxCacheTransitionEnergy = transitionEnergy;
865 
866  // Computing dE/dx from low-energy parameterisation at
867  // transition energy
868  G4double dEdxParam = (*iter) -> GetDEDX(particle, material,
869  transitionEnergy);
870 
871  G4double dEdxDeltaRays = DeltaRayMeanEnergyTransferRate(material,
872  particle,
873  transitionEnergy,
874  cutEnergy);
875  dEdxParam -= dEdxDeltaRays;
876 
877  // Computing dE/dx from Bethe-Bloch formula at transition
878  // energy
879  G4double transitionChargeSquare =
880  GetChargeSquareRatio(particle, material, transitionEnergy);
881 
882  G4double scaledTransitionEnergy = transitionEnergy * massRatio;
883 
884  G4double dEdxBetheBloch =
886  material, genericIon,
887  scaledTransitionEnergy, cutEnergy);
888  dEdxBetheBloch *= transitionChargeSquare;
889 
890  // Additionally, high order corrections are added
891  dEdxBetheBloch +=
892  corrections -> ComputeIonCorrections(particle,
893  material, transitionEnergy);
894 
895  // Computing transition factor from both dE/dx values
897  (dEdxParam - dEdxBetheBloch)/dEdxBetheBloch
898  * transitionEnergy;
899  }
900  else {
901 
902  dedxCacheParticle = particle;
903  dedxCacheMaterial = material;
904  dedxCacheEnergyCut = cutEnergy;
905 
907  genericIonPDGMass / particle -> GetPDGMass();
908 
911  }
912  }
913 }
914 
915 // #########################################################################
916 
918  const G4MaterialCutsCouple* couple,
919  const G4DynamicParticle* dynamicParticle,
920  G4double& eloss,
921  G4double&,
922  G4double length) {
923 
924  // ############## Corrections for along step energy loss calculation ######
925  // The computed energy loss (due to electronic stopping) is overwritten
926  // by this function if an ion data parameterization is available for the
927  // current ion-material pair.
928  // No action on the energy loss (due to electronic stopping) is performed
929  // if no parameterization is available. In this case the original
930  // generic ion tables (in combination with the effective charge) are used
931  // in the along step DoIt function.
932  //
933  //
934  // (Implementation partly adapted from G4BraggIonModel/G4BetheBlochModel)
935 
936  const G4ParticleDefinition* particle = dynamicParticle -> GetDefinition();
937  const G4Material* material = couple -> GetMaterial();
938 
939  G4double kineticEnergy = dynamicParticle -> GetKineticEnergy();
940 
941  if(kineticEnergy == eloss) { return; }
942 
943  G4double cutEnergy = DBL_MAX;
944  size_t cutIndex = couple -> GetIndex();
945  cutEnergy = cutEnergies[cutIndex];
946 
947  UpdateDEDXCache(particle, material, cutEnergy);
948 
949  LossTableList::iterator iter = dedxCacheIter;
950 
951  // If parameterization for ions is available the electronic energy loss
952  // is overwritten
953  if(iter != lossTableList.end()) {
954 
955  // The energy loss is calculated using the ComputeDEDXPerVolume function
956  // and the step length (it is assumed that dE/dx does not change
957  // considerably along the step)
958  eloss =
959  length * ComputeDEDXPerVolume(material, particle,
960  kineticEnergy, cutEnergy);
961 
962 #ifdef PRINT_DEBUG
963  G4cout.precision(6);
964  G4cout << "########################################################"
965  << G4endl
966  << "# G4IonParametrisedLossModel::CorrectionsAlongStep"
967  << G4endl
968  << "# cut(MeV) = " << cutEnergy/MeV
969  << G4endl;
970 
971  G4cout << "#"
972  << std::setw(13) << std::right << "E(MeV)"
973  << std::setw(14) << "l(um)"
974  << std::setw(14) << "l*dE/dx(MeV)"
975  << std::setw(14) << "(l*dE/dx)/E"
976  << G4endl
977  << "# ------------------------------------------------------"
978  << G4endl;
979 
980  G4cout << std::setw(14) << std::right << kineticEnergy / MeV
981  << std::setw(14) << length / um
982  << std::setw(14) << eloss / MeV
983  << std::setw(14) << eloss / kineticEnergy * 100.0
984  << G4endl;
985 #endif
986 
987  // If the energy loss exceeds a certain fraction of the kinetic energy
988  // (the fraction is indicated by the parameter "energyLossLimit") then
989  // the range tables are used to derive a more accurate value of the
990  // energy loss
991  if(eloss > energyLossLimit * kineticEnergy) {
992 
993  eloss = ComputeLossForStep(couple, particle,
994  kineticEnergy,length);
995 
996 #ifdef PRINT_DEBUG
997  G4cout << "# Correction applied:"
998  << G4endl;
999 
1000  G4cout << std::setw(14) << std::right << kineticEnergy / MeV
1001  << std::setw(14) << length / um
1002  << std::setw(14) << eloss / MeV
1003  << std::setw(14) << eloss / kineticEnergy * 100.0
1004  << G4endl;
1005 #endif
1006 
1007  }
1008  }
1009 
1010  // For all corrections below a kinetic energy between the Pre- and
1011  // Post-step energy values is used
1012  G4double energy = kineticEnergy - eloss * 0.5;
1013  if(energy < 0.0) energy = kineticEnergy * 0.5;
1014 
1015  G4double chargeSquareRatio = corrections ->
1016  EffectiveChargeSquareRatio(particle,
1017  material,
1018  energy);
1019  GetModelOfFluctuations() -> SetParticleAndCharge(particle,
1020  chargeSquareRatio);
1021 
1022  // A correction is applied considering the change of the effective charge
1023  // along the step (the parameter "corrFactor" refers to the effective
1024  // charge at the beginning of the step). Note: the correction is not
1025  // applied for energy loss values deriving directly from parameterized
1026  // ion stopping power tables
1027  G4double transitionEnergy = dedxCacheTransitionEnergy;
1028 
1029  if(iter != lossTableList.end() && transitionEnergy < kineticEnergy) {
1030  chargeSquareRatio *= corrections -> EffectiveChargeCorrection(particle,
1031  material,
1032  energy);
1033 
1034  G4double chargeSquareRatioCorr = chargeSquareRatio/corrFactor;
1035  eloss *= chargeSquareRatioCorr;
1036  }
1037  else if (iter == lossTableList.end()) {
1038 
1039  chargeSquareRatio *= corrections -> EffectiveChargeCorrection(particle,
1040  material,
1041  energy);
1042 
1043  G4double chargeSquareRatioCorr = chargeSquareRatio/corrFactor;
1044  eloss *= chargeSquareRatioCorr;
1045  }
1046 
1047  // Ion high order corrections are applied if the current model does not
1048  // overwrite the energy loss (i.e. when the effective charge approach is
1049  // used)
1050  if(iter == lossTableList.end()) {
1051 
1052  G4double scaledKineticEnergy = kineticEnergy * dedxCacheGenIonMassRatio;
1053  G4double lowEnergyLimit = betheBlochModel -> LowEnergyLimit();
1054 
1055  // Corrections are only applied in the Bethe-Bloch energy region
1056  if(scaledKineticEnergy > lowEnergyLimit)
1057  eloss += length *
1058  corrections -> IonHighOrderCorrections(particle, couple, energy);
1059  }
1060 }
1061 
1062 // #########################################################################
1063 
1065  const G4ParticleDefinition* particle,
1066  const G4MaterialCutsCouple* matCutsCouple) {
1067 
1068  G4double cutEnergy = DBL_MAX;
1069  size_t cutIndex = matCutsCouple -> GetIndex();
1070  cutEnergy = cutEnergies[cutIndex];
1071 
1072  const G4Material* material = matCutsCouple -> GetMaterial();
1073 
1074  G4double massRatio = genericIonPDGMass / particle -> GetPDGMass();
1075 
1076  G4double lowerEnergy = lowerEnergyEdgeIntegr / massRatio;
1077  G4double upperEnergy = upperEnergyEdgeIntegr / massRatio;
1078 
1079  G4double logLowerEnergyEdge = std::log(lowerEnergy);
1080  G4double logUpperEnergyEdge = std::log(upperEnergy);
1081 
1082  G4double logDeltaEnergy = (logUpperEnergyEdge - logLowerEnergyEdge) /
1083  G4double(nmbBins);
1084 
1085  G4double logDeltaIntegr = logDeltaEnergy / G4double(nmbSubBins);
1086 
1087  G4LPhysicsFreeVector* energyRangeVector =
1089  lowerEnergy,
1090  upperEnergy);
1091 
1092  G4double dedxLow = ComputeDEDXPerVolume(material,
1093  particle,
1094  lowerEnergy,
1095  cutEnergy);
1096 
1097  G4double range = 2.0 * lowerEnergy / dedxLow;
1098 
1099  energyRangeVector -> PutValues(0, lowerEnergy, range);
1100 
1101  G4double logEnergy = std::log(lowerEnergy);
1102  for(size_t i = 1; i < nmbBins+1; i++) {
1103 
1104  G4double logEnergyIntegr = logEnergy;
1105 
1106  for(size_t j = 0; j < nmbSubBins; j++) {
1107 
1108  G4double binLowerBoundary = std::exp(logEnergyIntegr);
1109  logEnergyIntegr += logDeltaIntegr;
1110 
1111  G4double binUpperBoundary = std::exp(logEnergyIntegr);
1112  G4double deltaIntegr = binUpperBoundary - binLowerBoundary;
1113 
1114  G4double energyIntegr = binLowerBoundary + 0.5 * deltaIntegr;
1115 
1116  G4double dedxValue = ComputeDEDXPerVolume(material,
1117  particle,
1118  energyIntegr,
1119  cutEnergy);
1120 
1121  if(dedxValue > 0.0) range += deltaIntegr / dedxValue;
1122 
1123 #ifdef PRINT_DEBUG_DETAILS
1124  G4cout << " E = "<< energyIntegr/MeV
1125  << " MeV -> dE = " << deltaIntegr/MeV
1126  << " MeV -> dE/dx = " << dedxValue/MeV*mm
1127  << " MeV/mm -> dE/(dE/dx) = " << deltaIntegr /
1128  dedxValue / mm
1129  << " mm -> range = " << range / mm
1130  << " mm " << G4endl;
1131 #endif
1132  }
1133 
1134  logEnergy += logDeltaEnergy;
1135 
1136  G4double energy = std::exp(logEnergy);
1137 
1138  energyRangeVector -> PutValues(i, energy, range);
1139 
1140 #ifdef PRINT_DEBUG_DETAILS
1141  G4cout << "G4IonParametrisedLossModel::BuildRangeVector() bin = "
1142  << i <<", E = "
1143  << energy / MeV << " MeV, R = "
1144  << range / mm << " mm"
1145  << G4endl;
1146 #endif
1147 
1148  }
1149  energyRangeVector -> SetSpline(true);
1150 
1151  G4double lowerRangeEdge =
1152  energyRangeVector -> Value(lowerEnergy);
1153  G4double upperRangeEdge =
1154  energyRangeVector -> Value(upperEnergy);
1155 
1156  G4LPhysicsFreeVector* rangeEnergyVector
1157  = new G4LPhysicsFreeVector(nmbBins+1,
1158  lowerRangeEdge,
1159  upperRangeEdge);
1160 
1161  for(size_t i = 0; i < nmbBins+1; i++) {
1162  G4double energy = energyRangeVector -> Energy(i);
1163  rangeEnergyVector ->
1164  PutValues(i, energyRangeVector -> Value(energy), energy);
1165  }
1166 
1167  rangeEnergyVector -> SetSpline(true);
1168 
1169 #ifdef PRINT_DEBUG_TABLES
1170  G4cout << *energyLossVector
1171  << *energyRangeVector
1172  << *rangeEnergyVector << G4endl;
1173 #endif
1174 
1175  IonMatCouple ionMatCouple = std::make_pair(particle, matCutsCouple);
1176 
1177  E[ionMatCouple] = energyRangeVector;
1178  r[ionMatCouple] = rangeEnergyVector;
1179 }
1180 
1181 // #########################################################################
1182 
1184  const G4MaterialCutsCouple* matCutsCouple,
1185  const G4ParticleDefinition* particle,
1186  G4double kineticEnergy,
1187  G4double stepLength) {
1188 
1189  G4double loss = 0.0;
1190 
1191  UpdateRangeCache(particle, matCutsCouple);
1192 
1193  G4PhysicsVector* energyRange = rangeCacheEnergyRange;
1194  G4PhysicsVector* rangeEnergy = rangeCacheRangeEnergy;
1195 
1196  if(energyRange != 0 && rangeEnergy != 0) {
1197 
1198  G4double lowerEnEdge = energyRange -> Energy( 0 );
1199  G4double lowerRangeEdge = rangeEnergy -> Energy( 0 );
1200 
1201  // Computing range for pre-step kinetic energy:
1202  G4double range = energyRange -> Value(kineticEnergy);
1203 
1204  // Energy below vector boundary:
1205  if(kineticEnergy < lowerEnEdge) {
1206 
1207  range = energyRange -> Value(lowerEnEdge);
1208  range *= std::sqrt(kineticEnergy / lowerEnEdge);
1209  }
1210 
1211 #ifdef PRINT_DEBUG
1212  G4cout << "G4IonParametrisedLossModel::ComputeLossForStep() range = "
1213  << range / mm << " mm, step = " << stepLength / mm << " mm"
1214  << G4endl;
1215 #endif
1216 
1217  // Remaining range:
1218  G4double remRange = range - stepLength;
1219 
1220  // If range is smaller than step length, the loss is set to kinetic
1221  // energy
1222  if(remRange < 0.0) loss = kineticEnergy;
1223  else if(remRange < lowerRangeEdge) {
1224 
1225  G4double ratio = remRange / lowerRangeEdge;
1226  loss = kineticEnergy - ratio * ratio * lowerEnEdge;
1227  }
1228  else {
1229 
1230  G4double energy = rangeEnergy -> Value(range - stepLength);
1231  loss = kineticEnergy - energy;
1232  }
1233  }
1234 
1235  if(loss < 0.0) loss = 0.0;
1236 
1237  return loss;
1238 }
1239 
1240 // #########################################################################
1241 
1243  const G4String& nam,
1244  G4VIonDEDXTable* table,
1245  G4VIonDEDXScalingAlgorithm* algorithm) {
1246 
1247  if(table == 0) {
1248  G4cerr << "G4IonParametrisedLossModel::AddDEDXTable() Cannot "
1249  << " add table: Invalid pointer."
1250  << G4endl;
1251 
1252  return false;
1253  }
1254 
1255  // Checking uniqueness of name
1256  LossTableList::iterator iter = lossTableList.begin();
1257  LossTableList::iterator iter_end = lossTableList.end();
1258 
1259  for(;iter != iter_end; iter++) {
1260  G4String tableName = (*iter) -> GetName();
1261 
1262  if(tableName == nam) {
1263  G4cerr << "G4IonParametrisedLossModel::AddDEDXTable() Cannot "
1264  << " add table: Name already exists."
1265  << G4endl;
1266 
1267  return false;
1268  }
1269  }
1270 
1271  G4VIonDEDXScalingAlgorithm* scalingAlgorithm = algorithm;
1272  if(scalingAlgorithm == 0)
1273  scalingAlgorithm = new G4VIonDEDXScalingAlgorithm;
1274 
1275  G4IonDEDXHandler* handler =
1276  new G4IonDEDXHandler(table, scalingAlgorithm, nam);
1277 
1278  lossTableList.push_front(handler);
1279 
1280  return true;
1281 }
1282 
1283 // #########################################################################
1284 
1286  const G4String& nam) {
1287 
1288  LossTableList::iterator iter = lossTableList.begin();
1289  LossTableList::iterator iter_end = lossTableList.end();
1290 
1291  for(;iter != iter_end; iter++) {
1292  G4String tableName = (*iter) -> GetName();
1293 
1294  if(tableName == nam) {
1295  delete (*iter);
1296 
1297  // Remove from table list
1298  lossTableList.erase(iter);
1299 
1300  // Range vs energy and energy vs range vectors are cleared
1301  RangeEnergyTable::iterator iterRange = r.begin();
1302  RangeEnergyTable::iterator iterRange_end = r.end();
1303 
1304  for(;iterRange != iterRange_end; iterRange++)
1305  delete iterRange -> second;
1306  r.clear();
1307 
1308  EnergyRangeTable::iterator iterEnergy = E.begin();
1309  EnergyRangeTable::iterator iterEnergy_end = E.end();
1310 
1311  for(;iterEnergy != iterEnergy_end; iterEnergy++)
1312  delete iterEnergy -> second;
1313  E.clear();
1314 
1315  return true;
1316  }
1317  }
1318 
1319  return false;
1320 }
1321 
1322 // #########################################################################
1323 
1325 
1326  RemoveDEDXTable("ICRU73");
1327  AddDEDXTable("ICRU73", new G4IonStoppingData("ion_stopping_data/icru73"));
1328 }
1329 
1330 // #########################################################################
virtual G4double GetChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double)
static const double cm
Definition: G4SIunits.hh:106
G4double MaxSecondaryKinEnergy(const G4DynamicParticle *dynParticle)
Definition: G4VEmModel.hh:464
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:623
static const double MeV
Definition: G4SIunits.hh:193
static G4LossTableManager * Instance()
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:120
static const double cm2
Definition: G4SIunits.hh:107
void BuildRangeVector(const G4ParticleDefinition *, const G4MaterialCutsCouple *)
CLHEP::Hep3Vector G4ThreeVector
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:616
static G4GenericIon * Definition()
Definition: G4GenericIon.cc:49
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:599
G4bool AddDEDXTable(const G4String &name, G4VIonDEDXTable *table, G4VIonDEDXScalingAlgorithm *algorithm=0)
G4VEmFluctuationModel * GetModelOfFluctuations()
Definition: G4VEmModel.hh:592
int G4int
Definition: G4Types.hh:78
const G4ParticleDefinition * cacheParticle
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double, G4double, G4double, G4double, G4double)
void UpdateRangeCache(const G4ParticleDefinition *, const G4MaterialCutsCouple *)
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:707
const G4MaterialCutsCouple * rangeCacheMatCutsCouple
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double, G4double)
virtual G4double GetParticleCharge(const G4ParticleDefinition *, const G4Material *, G4double)
#define G4UniformRand()
Definition: Randomize.hh:93
G4GLOB_DLL std::ostream G4cout
G4ParticleChangeForLoss * particleChangeLoss
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=0)
Definition: G4VEmModel.cc:402
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double, G4double)
bool G4bool
Definition: G4Types.hh:79
const G4ThreeVector & GetMomentumDirection() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
static const double cm3
Definition: G4SIunits.hh:108
static const double second
Definition: G4SIunits.hh:138
void SetProposedMomentumDirection(const G4ThreeVector &dir)
G4double DeltaRayMeanEnergyTransferRate(const G4Material *, const G4ParticleDefinition *, G4double, G4double)
void PrintDEDXTableHandlers(const G4ParticleDefinition *, const G4Material *, G4double, G4double, G4int, G4bool)
static G4ProductionCutsTable * GetProductionCutsTable()
G4bool RemoveDEDXTable(const G4String &name)
void UpdateCache(const G4ParticleDefinition *)
static const G4double factor
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double)
G4double energy(const ThreeVector &p, const G4double m)
std::pair< const G4ParticleDefinition *, const G4MaterialCutsCouple * > IonMatCouple
static const double g
Definition: G4SIunits.hh:162
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:606
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double, G4double, G4double)
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
const G4String & GetName() const
Definition: G4VEmModel.hh:777
static const double TeV
Definition: G4SIunits.hh:197
G4IonParametrisedLossModel(const G4ParticleDefinition *particle=0, const G4String &name="ParamICRU73")
LossTableList::iterator IsApplicable(const G4ParticleDefinition *, const G4Material *)
void UpdateDEDXCache(const G4ParticleDefinition *, const G4Material *, G4double cutEnergy)
double G4double
Definition: G4Types.hh:76
const G4ParticleDefinition * dedxCacheParticle
void PrintDEDXTable(const G4ParticleDefinition *, const G4Material *, G4double, G4double, G4int, G4bool)
G4double ComputeLossForStep(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double, G4double)
#define DBL_MAX
Definition: templates.hh:83
const G4ParticleDefinition * rangeCacheParticle
static const double mm
Definition: G4SIunits.hh:102
G4ThreeVector GetMomentum() const
virtual G4double Value(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:361
const G4Material * GetMaterial() const
G4GLOB_DLL std::ostream G4cerr
G4int SelectRandomAtomNumber(const G4Material *)
Definition: G4VEmModel.hh:546
virtual void CorrectionsAlongStep(const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &, G4double &, G4double)