Geant4  10.00.p01
G4VEnergyLossProcess.hh
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: G4VEnergyLossProcess.hh 76333 2013-11-08 14:31:50Z gcosmo $
27 // GEANT4 tag $Name:
28 //
29 // -------------------------------------------------------------------
30 //
31 // GEANT4 Class header file
32 //
33 //
34 // File name: G4VEnergyLossProcess
35 //
36 // Author: Vladimir Ivanchenko on base of Laszlo Urban code
37 //
38 // Creation date: 03.01.2002
39 //
40 // Modifications:
41 //
42 // 26-12-02 Secondary production moved to derived classes (V.Ivanchenko)
43 // 20-01-03 Migrade to cut per region (V.Ivanchenko)
44 // 24-01-03 Make models region aware (V.Ivanchenko)
45 // 05-02-03 Fix compilation warnings (V.Ivanchenko)
46 // 13-02-03 SubCutoffProcessors defined for regions (V.Ivanchenko)
47 // 17-02-03 Fix problem of store/restore tables (V.Ivanchenko)
48 // 26-02-03 Region dependent step limit (V.Ivanchenko)
49 // 26-03-03 Add GetDEDXDispersion (V.Ivanchenko)
50 // 09-04-03 Fix problem of negative range limit for non integral (V.Ivanchenko)
51 // 13-05-03 Add calculation of precise range (V.Ivanchenko)
52 // 21-07-03 Add UpdateEmModel method (V.Ivanchenko)
53 // 12-11-03 G4EnergyLossSTD -> G4EnergyLossProcess (V.Ivanchenko)
54 // 14-01-04 Activate precise range calculation (V.Ivanchenko)
55 // 10-03-04 Fix problem of step limit calculation (V.Ivanchenko)
56 // 30-06-04 make destructor virtual (V.Ivanchenko)
57 // 05-07-04 fix problem of GenericIons seen at small cuts (V.Ivanchenko)
58 // 03-08-04 Add DEDX table to all processes for control on integral range(VI)
59 // 06-08-04 Clear up names of member functions (V.Ivanchenko)
60 // 27-08-04 Add NeedBuildTables method (V.Ivanchneko)
61 // 09-09-04 Bug fix for the integral mode with 2 peaks (V.Ivanchneko)
62 // 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivanchenko)
63 // 08-04-05 Major optimisation of internal interfaces (V.Ivanchenko)
64 // 11-04-05 Use MaxSecondaryEnergy from a model (V.Ivanchenko)
65 // 10-01-05 Remove SetStepLimits (V.Ivanchenko)
66 // 10-01-06 PreciseRange -> CSDARange (V.Ivantchenko)
67 // 13-01-06 Remove AddSubCutSecondaries and cleanup (V.Ivantchenko)
68 // 20-01-06 Introduce G4EmTableType and reducing number of methods (VI)
69 // 26-01-06 Add public method GetCSDARange (V.Ivanchenko)
70 // 22-03-06 Add SetDynamicMassCharge (V.Ivanchenko)
71 // 23-03-06 Use isIonisation flag (V.Ivanchenko)
72 // 13-05-06 Add method to access model by index (V.Ivanchenko)
73 // 14-01-07 add SetEmModel(index) and SetFluctModel() (mma)
74 // 15-01-07 Add separate ionisation tables and reorganise get/set methods for
75 // dedx tables (V.Ivanchenko)
76 // 13-03-07 use SafetyHelper instead of navigator (V.Ivanchenko)
77 // 27-07-07 use stl vector for emModels instead of C-array (V.Ivanchenko)
78 // 25-09-07 More accurate handling zero xsect in
79 // PostStepGetPhysicalInteractionLength (V.Ivanchenko)
80 // 27-10-07 Virtual functions moved to source (V.Ivanchenko)
81 // 15-07-08 Reorder class members for further multi-thread development (VI)
82 //
83 // Class Description:
84 //
85 // It is the unified energy loss process it calculates the continuous
86 // energy loss for charged particles using a set of Energy Loss
87 // models valid for different energy regions. There are a possibility
88 // to create and access to dE/dx and range tables, or to calculate
89 // that information on fly.
90 
91 // -------------------------------------------------------------------
92 //
93 
94 #ifndef G4VEnergyLossProcess_h
95 #define G4VEnergyLossProcess_h 1
96 
98 #include "globals.hh"
99 #include "G4Material.hh"
100 #include "G4MaterialCutsCouple.hh"
101 #include "G4Track.hh"
102 #include "G4EmModelManager.hh"
103 #include "G4UnitsTable.hh"
105 #include "G4EmTableType.hh"
106 #include "G4PhysicsTable.hh"
107 #include "G4PhysicsVector.hh"
108 
109 class G4Step;
111 class G4VEmModel;
113 class G4DataVector;
114 class G4Region;
115 class G4SafetyHelper;
116 class G4VAtomDeexcitation;
117 class G4EmBiasingManager;
118 class G4LossTableManager;
119 
120 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
121 
123 {
124 public:
125 
126  G4VEnergyLossProcess(const G4String& name = "EnergyLoss",
128 
129  virtual ~G4VEnergyLossProcess();
130 
131 private:
132  // clean vectors and arrays
133  void Clean();
134 
135  //------------------------------------------------------------------------
136  // Virtual methods to be implemented in concrete processes
137  //------------------------------------------------------------------------
138 
139 public:
140  virtual G4bool IsApplicable(const G4ParticleDefinition& p) = 0;
141 
142  virtual void PrintInfo() = 0;
143 
144 protected:
145 
147  const G4ParticleDefinition*) = 0;
148 
149  //------------------------------------------------------------------------
150  // Methods with standard implementation; may be overwritten if needed
151  //------------------------------------------------------------------------
152 
154  const G4Material*, G4double cut);
155 
156  //------------------------------------------------------------------------
157  // Virtual methods implementation common to all EM ContinuousDiscrete
158  // processes. Further inheritance is not assumed
159  //------------------------------------------------------------------------
160 
161 public:
162 
163  // prepare all tables
165 
166  // build all tables
168 
169  // build a table
171 
172  // build a table
174 
175  // summary printout after initialisation
176  void PrintInfoDefinition(const G4ParticleDefinition& part);
177 
178  // Called before tracking of each new G4Track
179  void StartTracking(G4Track*);
180 
181  // Step limit from AlongStep
183  G4double previousStepSize,
184  G4double currentMinimumStep,
185  G4double& currentSafety,
186  G4GPILSelection* selection);
187 
188  // Step limit from cross section
190  G4double previousStepSize,
192 
193  // AlongStep computations
195 
196  // Sampling of secondaries in vicinity of geometrical boundary
197  // Return sum of secodaries energy
198  G4double SampleSubCutSecondaries(std::vector<G4Track*>&, const G4Step&,
199  G4VEmModel* model, G4int matIdx);
200 
201  // PostStep sampling of secondaries
202  G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step&);
203 
204  // Store all PhysicsTable in files.
205  // Return false in case of any fatal failure at I/O
207  const G4String& directory,
208  G4bool ascii = false);
209 
210  // Retrieve all Physics from a files.
211  // Return true if all the Physics Table are built.
212  // Return false if any fatal failure.
214  const G4String& directory,
215  G4bool ascii);
216 
217 private:
218  // store a table
220  G4PhysicsTable*, G4bool ascii,
221  const G4String& directory,
222  const G4String& tname);
223 
224  // retrieve a table
226  G4PhysicsTable*, G4bool ascii,
227  const G4String& directory,
228  const G4String& tname,
229  G4bool mandatory);
230 
231  //------------------------------------------------------------------------
232  // Public interface to cross section, mfp and sampling of fluctuations
233  // These methods are not used in run time
234  //------------------------------------------------------------------------
235 
236 public:
237 
238  // access to dispersion of restricted energy loss
240  const G4DynamicParticle* dp,
241  G4double length);
242 
243  // Access to cross section table
245  const G4MaterialCutsCouple* couple);
246 
247  // access to cross section
248  G4double MeanFreePath(const G4Track& track);
249 
250  // access to step limit
251  G4double ContinuousStepLimit(const G4Track& track,
252  G4double previousStepSize,
253  G4double currentMinimumStep,
254  G4double& currentSafety);
255 
256 protected:
257 
258  // implementation of the pure virtual method
259  G4double GetMeanFreePath(const G4Track& track,
260  G4double previousStepSize,
261  G4ForceCondition* condition);
262 
263  // implementation of the pure virtual method
265  G4double previousStepSize,
266  G4double currentMinimumStep,
267  G4double& currentSafety);
268 
269  //------------------------------------------------------------------------
270  // Run time method which may be also used by derived processes
271  //------------------------------------------------------------------------
272 
273  // creeation of an empty vector for cross section
275  G4double cut);
276 
277  inline size_t CurrentMaterialCutsCoupleIndex() const;
278 
279  //------------------------------------------------------------------------
280  // Specific methods to set, access, modify models
281  //------------------------------------------------------------------------
282 
283  // Select model in run time
284  inline void SelectModel(G4double kinEnergy);
285 
286 public:
287  // Select model by energy and region index
288  inline G4VEmModel* SelectModelForMaterial(G4double kinEnergy,
289  size_t& idx) const;
290 
291  // Add EM model coupled with fluctuation model for region, smaller value
292  // of order defines which pair of models will be selected for a given
293  // energy interval
294  void AddEmModel(G4int, G4VEmModel*,
295  G4VEmFluctuationModel* fluc = 0,
296  const G4Region* region = 0);
297 
298  // Define new energy range for the model identified by the name
299  void UpdateEmModel(const G4String&, G4double, G4double);
300 
301  // Assign a model to a process
302  void SetEmModel(G4VEmModel*, G4int index=1);
303 
304  // return the assigned model
305  G4VEmModel* EmModel(G4int index=1) const;
306 
307  // Access to models
308  G4VEmModel* GetModelByIndex(G4int idx = 0, G4bool ver = false) const;
309 
310  G4int NumberOfModels() const;
311 
312  // Assign a fluctuation model to a process
314 
315  // return the assigned fluctuation model
317 
318  //------------------------------------------------------------------------
319  // Define and access particle type
320  //------------------------------------------------------------------------
321 
322 protected:
323  inline void SetParticle(const G4ParticleDefinition* p);
324  inline void SetSecondaryParticle(const G4ParticleDefinition* p);
325 
326 public:
327  inline void SetBaseParticle(const G4ParticleDefinition* p);
328  inline const G4ParticleDefinition* Particle() const;
329  inline const G4ParticleDefinition* BaseParticle() const;
330  inline const G4ParticleDefinition* SecondaryParticle() const;
331 
332  //------------------------------------------------------------------------
333  // Get/set parameters to configure the process at initialisation time
334  //------------------------------------------------------------------------
335 
336  // Add subcutoff option for the region
337  void ActivateSubCutoff(G4bool val, const G4Region* region = 0);
338 
339  // Activate biasing
340  void SetCrossSectionBiasingFactor(G4double f, G4bool flag = true);
341 
342  void ActivateForcedInteraction(G4double length = 0.0,
343  const G4String& region = "",
344  G4bool flag = true);
345 
346  void ActivateSecondaryBiasing(const G4String& region, G4double factor,
347  G4double energyLimit);
348 
349  // Add subcutoff process (bremsstrahlung) to sample secondary
350  // particle production in vicinity of the geometry boundary
352 
353  inline void SetLossFluctuations(G4bool val);
354  inline void SetRandomStep(G4bool val);
355 
356  inline void SetIntegral(G4bool val);
357  inline G4bool IsIntegral() const;
358 
359  // Set/Get flag "isIonisation"
360  inline void SetIonisation(G4bool val);
361  inline G4bool IsIonisationProcess() const;
362 
363  // Redefine parameteters for stepping control
364  inline void SetLinearLossLimit(G4double val);
365  inline void SetMinSubRange(G4double val);
366  inline void SetLambdaFactor(G4double val);
367  inline void SetStepFunction(G4double v1, G4double v2);
368  inline void SetLowestEnergyLimit(G4double);
369 
370  inline G4int NumberOfSubCutoffRegions() const;
371 
372  //------------------------------------------------------------------------
373  // Specific methods to path Physics Tables to the process
374  //------------------------------------------------------------------------
375 
376  void SetDEDXTable(G4PhysicsTable* p, G4EmTableType tType);
377  void SetCSDARangeTable(G4PhysicsTable* pRange);
381 
384 
385  // Binning for dEdx, range, inverse range and labda tables
386  inline void SetDEDXBinning(G4int nbins);
387  inline void SetLambdaBinning(G4int nbins);
388 
389  // Binning for dEdx, range, and inverse range tables
390  inline void SetDEDXBinningForCSDARange(G4int nbins);
391 
392  // Min kinetic energy for tables
393  inline void SetMinKinEnergy(G4double e);
394  inline G4double MinKinEnergy() const;
395 
396  // Max kinetic energy for tables
397  inline void SetMaxKinEnergy(G4double e);
398  inline G4double MaxKinEnergy() const;
399 
400  // Max kinetic energy for tables
401  inline void SetMaxKinEnergyForCSDARange(G4double e);
402 
403  // Biasing parameters
404  inline G4double CrossSectionBiasingFactor() const;
405 
406  // Return values for given G4MaterialCutsCouple
407  inline G4double GetDEDX(G4double& kineticEnergy, const G4MaterialCutsCouple*);
408  inline G4double GetDEDXForSubsec(G4double& kineticEnergy,
409  const G4MaterialCutsCouple*);
410  inline G4double GetRange(G4double& kineticEnergy, const G4MaterialCutsCouple*);
411  inline G4double GetCSDARange(G4double& kineticEnergy, const G4MaterialCutsCouple*);
412  inline G4double GetRangeForLoss(G4double& kineticEnergy, const G4MaterialCutsCouple*);
414  inline G4double GetLambda(G4double& kineticEnergy, const G4MaterialCutsCouple*);
415 
416  inline G4bool TablesAreBuilt() const;
417 
418  // Access to specific tables
419  inline G4PhysicsTable* DEDXTable() const;
420  inline G4PhysicsTable* DEDXTableForSubsec() const;
421  inline G4PhysicsTable* DEDXunRestrictedTable() const;
422  inline G4PhysicsTable* IonisationTable() const;
424  inline G4PhysicsTable* CSDARangeTable() const;
425  inline G4PhysicsTable* SecondaryRangeTable() const;
426  inline G4PhysicsTable* RangeTableForLoss() const;
427  inline G4PhysicsTable* InverseRangeTable() const;
428  inline G4PhysicsTable* LambdaTable() const;
429  inline G4PhysicsTable* SubLambdaTable() const;
430 
431  //------------------------------------------------------------------------
432  // Run time method for simulation of ionisation
433  //------------------------------------------------------------------------
434 
435  // access atom on which interaction happens
436  const G4Element* GetCurrentElement() const;
437 
438  // sample range at the end of a step
439  // inline G4double SampleRange();
440 
441  // Set scaling parameters for ions is needed to G4EmCalculator
442  inline void SetDynamicMassCharge(G4double massratio, G4double charge2ratio);
443 
444 private:
445 
446  void FillSecondariesAlongStep(G4double& eloss, G4double& weight);
447 
448  // define material and indexes
449  inline void DefineMaterial(const G4MaterialCutsCouple* couple);
450 
451  //------------------------------------------------------------------------
452  // Compute values using scaling relation, mass and charge of based particle
453  //------------------------------------------------------------------------
454 
455  inline G4double GetDEDXForScaledEnergy(G4double scaledKinEnergy);
456  inline G4double GetSubDEDXForScaledEnergy(G4double scaledKinEnergy);
457  inline G4double GetIonisationForScaledEnergy(G4double scaledKinEnergy);
458  inline G4double GetSubIonisationForScaledEnergy(G4double scaledKinEnergy);
459  inline G4double GetScaledRangeForScaledEnergy(G4double scaledKinEnergy);
460  inline G4double GetLimitScaledRangeForScaledEnergy(G4double scaledKinEnergy);
462  inline G4double GetLambdaForScaledEnergy(G4double scaledKinEnergy);
463  inline void ComputeLambdaForScaledEnergy(G4double scaledKinEnergy);
464 
465  // hide assignment operator
468 
469  // ======== Parameters of the class fixed at construction =========
470 
475 
481 
482  // G4PhysicsVector* vstrag;
483 
484  // ======== Parameters of the class fixed at initialisation =======
485 
486  std::vector<G4VEmModel*> emModels;
489  std::vector<const G4Region*> scoffRegions;
492 
493  std::vector<G4VEnergyLossProcess*> scProcesses;
495 
496  // tables and vectors
508 
509  size_t idxDEDX;
510  size_t idxDEDXSub;
514  size_t idxRange;
515  size_t idxCSDA;
516  size_t idxSecRange;
518  size_t idxLambda;
519  size_t idxSubLambda;
520 
521  std::vector<G4double> theDEDXAtMaxEnergy;
522  std::vector<G4double> theRangeAtMaxEnergy;
523  std::vector<G4double> theEnergyOfCrossSectionMax;
524  std::vector<G4double> theCrossSectionMax;
525 
526  const std::vector<G4double>* theDensityFactor;
527  const std::vector<G4int>* theDensityIdx;
528 
531 
533 
536 
541 
548 
560 
561 protected:
562 
564 
565  // ======== Cached values - may be state dependent ================
566 
567 private:
568 
569  std::vector<G4DynamicParticle*> secParticles;
570  std::vector<G4Track*> scTracks;
571 
573 
579  size_t lastIdx;
580 
582 
587 
595 
597 
601 };
602 
603 // ======== Run time inline methods ================
604 
606 {
607  return currentCoupleIndex;
608 }
609 
610 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
611 
613 {
616 }
617 
618 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
619 
621  G4double kinEnergy, size_t& idx) const
622 {
623  return modelManager->SelectModel(kinEnergy, idx);
624 }
625 
626 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
627 
628 inline void
630 {
631  if(couple != currentCouple) {
632  currentCouple = couple;
633  currentMaterial = couple->GetMaterial();
634  currentCoupleIndex = couple->GetIndex();
635  basedCoupleIndex = (*theDensityIdx)[currentCoupleIndex];
636  fFactor = chargeSqRatio*biasFactor*(*theDensityFactor)[currentCoupleIndex];
639  idxLambda = idxSubLambda = 0;
640  }
641 }
642 
643 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
644 
646  G4double charge2ratio)
647 {
648  massRatio = massratio;
649  fFactor = charge2ratio*biasFactor*(*theDensityFactor)[currentCoupleIndex];
650  chargeSqRatio = charge2ratio;
652 }
653 
654 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
655 
657 {
658  /*
659  G4cout << "G4VEnergyLossProcess::GetDEDX: Idx= "
660  << basedCoupleIndex << " E(MeV)= " << e
661  << " Emin= " << minKinEnergy << " Factor= " << fFactor
662  << " " << theDEDXTable << G4endl; */
663  G4double x = fFactor*(*theDEDXTable)[basedCoupleIndex]->Value(e, idxDEDX);
664  if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
665  return x;
666 }
667 
668 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
669 
671 {
672  G4double x =
673  fFactor*(*theDEDXSubTable)[basedCoupleIndex]->Value(e, idxDEDXSub);
674  if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
675  return x;
676 }
677 
678 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
679 
681 {
682  G4double x =
683  fFactor*(*theIonisationTable)[basedCoupleIndex]->Value(e, idxIonisation);
684  if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
685  return x;
686 }
687 
688 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
689 
690 inline
692 {
693  G4double x = fFactor*
694  (*theIonisationSubTable)[basedCoupleIndex]->Value(e, idxIonisationSub);
695  if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
696  return x;
697 }
698 
699 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
700 
702 {
703  //G4cout << "G4VEnergyLossProcess::GetScaledRange: Idx= "
704  // << basedCoupleIndex << " E(MeV)= " << e
705  // << " lastIdx= " << lastIdx << " " << theRangeTableForLoss << G4endl;
708  preStepRangeEnergy = e;
709  computedRange =
710  ((*theRangeTableForLoss)[basedCoupleIndex])->Value(e, idxRange);
711  if(e < minKinEnergy) { computedRange *= std::sqrt(e/minKinEnergy); }
712  }
713  //G4cout << "G4VEnergyLossProcess::GetScaledRange: Idx= "
714  // << basedCoupleIndex << " E(MeV)= " << e
715  // << " R= " << fRange << " " << theRangeTableForLoss << G4endl;
716 
717  return computedRange;
718 }
719 
720 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
721 
722 inline G4double
724 {
725  G4double x;
726  if (e < maxKinEnergyCSDA) {
727  x = ((*theCSDARangeTable)[basedCoupleIndex])->Value(e, idxCSDA);
728  if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
729  } else {
732  }
733  return x;
734 }
735 
736 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
737 
739 {
740  //G4cout << "G4VEnergyLossProcess::GetEnergy: Idx= "
741  // << basedCoupleIndex << " R(mm)= " << r << " "
742  // << theInverseRangeTable << G4endl;
743  G4PhysicsVector* v = (*theInverseRangeTable)[basedCoupleIndex];
744  G4double rmin = v->Energy(0);
745  G4double e = 0.0;
746  if(r >= rmin) { e = v->Value(r, idxInverseRange); }
747  else if(r > 0.0) {
748  G4double x = r/rmin;
749  e = minKinEnergy*x*x;
750  }
751  return e;
752 }
753 
754 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
755 
757 {
758  return fFactor*((*theLambdaTable)[basedCoupleIndex])->Value(e, idxLambda);
759 }
760 
761 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
762 
763 inline G4double
765  const G4MaterialCutsCouple* couple)
766 {
767  DefineMaterial(couple);
768  return GetDEDXForScaledEnergy(kineticEnergy*massRatio);
769 }
770 
771 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
772 
773 inline G4double
775  const G4MaterialCutsCouple* couple)
776 {
777  DefineMaterial(couple);
778  return GetSubDEDXForScaledEnergy(kineticEnergy*massRatio);
779 }
780 
781 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
782 
783 inline G4double
785  const G4MaterialCutsCouple* couple)
786 {
787  G4double x = fRange;
788  DefineMaterial(couple);
789  if(theCSDARangeTable) {
791  * reduceFactor;
792  } else if(theRangeTableForLoss) {
794  }
795  return x;
796 }
797 
798 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
799 
800 inline G4double
802  const G4MaterialCutsCouple* couple)
803 {
804  DefineMaterial(couple);
805  G4double x = DBL_MAX;
806  if(theCSDARangeTable) {
808  }
809  return x;
810 }
811 
812 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
813 
814 inline G4double
816  const G4MaterialCutsCouple* couple)
817 {
818  // G4cout << "GetRangeForLoss: Range from " << GetProcessName() << G4endl;
819  DefineMaterial(couple);
820  G4double x =
822  //G4cout << "GetRangeForLoss: Range from " << GetProcessName()
823  // << " e= " << kineticEnergy << " r= " << x << G4endl;
824  return x;
825 }
826 
827 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
828 
829 inline G4double
831  const G4MaterialCutsCouple* couple)
832 {
833  DefineMaterial(couple);
835 }
836 
837 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
838 
839 inline G4double
841  const G4MaterialCutsCouple* couple)
842 {
843  DefineMaterial(couple);
844  G4double x = 0.0;
845  if(theLambdaTable) { x = GetLambdaForScaledEnergy(kineticEnergy*massRatio); }
846  return x;
847 }
848 
849 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
850 
852 {
854  if (e <= mfpKinEnergy) {
856 
857  } else {
859  if(e1 > mfpKinEnergy) {
861  G4double preStepLambda1 = GetLambdaForScaledEnergy(e1);
862  if(preStepLambda1 > preStepLambda) {
863  mfpKinEnergy = e1;
864  preStepLambda = preStepLambda1;
865  }
866  } else {
868  }
869  }
870 }
871 
872 // ======== Get/Set inline methods used at initialisation ================
873 
875 {
876  fluctModel = p;
877 }
878 
879 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
880 
882 {
883  return fluctModel;
884 }
885 
886 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
887 
889 {
890  particle = p;
891 }
892 
893 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
894 
896 {
897  secondaryParticle = p;
898 }
899 
900 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
901 
903 {
904  baseParticle = p;
905 }
906 
907 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
908 
910 {
911  return particle;
912 }
913 
914 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
915 
917 {
918  return baseParticle;
919 }
920 
921 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
922 
924 {
925  return secondaryParticle;
926 }
927 
928 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
929 
931 {
932  lossFluctuationFlag = val;
933 }
934 
935 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
936 
938 {
939  rndmStepFlag = val;
940 }
941 
942 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
943 
945 {
946  integral = val;
947 }
948 
949 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
950 
952 {
953  return integral;
954 }
955 
956 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
957 
959 {
960  isIonisation = val;
963 }
964 
965 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
966 
968 {
969  return isIonisation;
970 }
971 
972 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
973 
975 {
976  linLossLimit = val;
977 }
978 
979 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
980 
982 {
983  minSubRange = val;
984 }
985 
986 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
987 
989 {
990  if(val > 0.0 && val <= 1.0) { lambdaFactor = val; }
991 }
992 
993 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
994 
996 {
997  dRoverRange = v1;
998  finalRange = v2;
999  if (dRoverRange > 0.999) { dRoverRange = 1.0; }
1000 }
1001 
1002 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1003 
1005 {
1006  lowestKinEnergy = val;
1007 }
1008 
1009 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1010 
1012 {
1013  return nSCoffRegions;
1014 }
1015 
1016 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1017 
1019 {
1020  nBins = nbins;
1021 }
1022 
1023 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1024 
1026 {
1027  nBins = nbins;
1028 }
1029 
1030 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1031 
1033 {
1034  nBinsCSDA = nbins;
1035 }
1036 
1037 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1038 
1040 {
1041  minKinEnergy = e;
1042 }
1043 
1044 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1045 
1047 {
1048  return minKinEnergy;
1049 }
1050 
1051 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1052 
1054 {
1055  maxKinEnergy = e;
1056  if(e < maxKinEnergyCSDA) { maxKinEnergyCSDA = e; }
1057 }
1058 
1059 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1060 
1062 {
1063  return maxKinEnergy;
1064 }
1065 
1066 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1067 
1069 {
1070  maxKinEnergyCSDA = e;
1071 }
1072 
1073 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1074 
1076 {
1077  return biasFactor;
1078 }
1079 
1080 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1081 
1083 {
1084  return tablesAreBuilt;
1085 }
1086 
1087 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1088 
1090 {
1091  return theDEDXTable;
1092 }
1093 
1094 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1095 
1097 {
1098  return theDEDXSubTable;
1099 }
1100 
1101 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1102 
1104 {
1105  return theDEDXunRestrictedTable;
1106 }
1107 
1108 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1109 
1111 {
1112  return theIonisationTable;
1113 }
1114 
1115 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1116 
1118 {
1119  return theIonisationSubTable;
1120 }
1121 
1122 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1123 
1125 {
1126  return theCSDARangeTable;
1127 }
1128 
1129 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1130 
1132 {
1133  return theSecondaryRangeTable;
1134 }
1135 
1136 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1137 
1139 {
1140  return theRangeTableForLoss;
1141 }
1142 
1143 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1144 
1146 {
1147  return theInverseRangeTable;
1148 }
1149 
1150 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1151 
1153 {
1154  return theLambdaTable;
1155 }
1156 
1157 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1158 
1160 {
1161  return theSubLambdaTable;
1162 }
1163 
1164 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1165 
1166 #endif
G4double GetKineticEnergy(G4double &range, const G4MaterialCutsCouple *)
G4VEnergyLossProcess(const G4String &name="EnergyLoss", G4ProcessType type=fElectromagnetic)
G4double condition(const G4ErrorSymMatrix &m)
G4LossTableManager * lManager
std::vector< G4double > theEnergyOfCrossSectionMax
void SetDEDXTable(G4PhysicsTable *p, G4EmTableType tType)
void SetIntegral(G4bool val)
void SetMaxKinEnergyForCSDARange(G4double e)
G4GPILSelection aGPILSelection
void SetDynamicMassCharge(G4double massratio, G4double charge2ratio)
void PrintInfoDefinition(const G4ParticleDefinition &part)
void SetIonisation(G4bool val)
void BuildPhysicsTable(const G4ParticleDefinition &)
G4double GetDEDX(G4double &kineticEnergy, const G4MaterialCutsCouple *)
G4PhysicsTable * SubLambdaTable() const
G4VEnergyLossProcess & operator=(const G4VEnergyLossProcess &right)
G4PhysicsTable * RangeTableForLoss() const
G4double GetContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)
G4String name
Definition: TRTMaterials.hh:40
const std::vector< G4double > * theDensityFactor
G4PhysicsTable * SecondaryRangeTable() const
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
G4double GetRangeForLoss(G4double &kineticEnergy, const G4MaterialCutsCouple *)
void AddCollaborativeProcess(G4VEnergyLossProcess *)
const std::vector< G4int > * theDensityIdx
void ActivateForcedInteraction(G4double length=0.0, const G4String &region="", G4bool flag=true)
void SetLinearLossLimit(G4double val)
const G4ParticleDefinition * secondaryParticle
G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
void SetStepFunction(G4double v1, G4double v2)
G4PhysicsTable * CSDARangeTable() const
G4PhysicsTable * theSecondaryRangeTable
G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false)
G4VEmFluctuationModel * fluctModel
G4ParticleChangeForLoss fParticleChange
G4double GetCSDARange(G4double &kineticEnergy, const G4MaterialCutsCouple *)
void SelectModel(G4double kinEnergy)
void SetLambdaBinning(G4int nbins)
G4PhysicsTable * theIonisationTable
G4VEmModel * GetModelByIndex(G4int idx=0, G4bool ver=false) const
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
G4PhysicsTable * IonisationTableForSubsec() const
G4double GetDEDXForSubsec(G4double &kineticEnergy, const G4MaterialCutsCouple *)
void ActivateSecondaryBiasing(const G4String &region, G4double factor, G4double energyLimit)
void SetFluctModel(G4VEmFluctuationModel *)
G4PhysicsVector * LambdaPhysicsVector(const G4MaterialCutsCouple *, G4double cut)
G4PhysicsTable * IonisationTable() const
G4VAtomDeexcitation * atomDeexcitation
int G4int
Definition: G4Types.hh:78
G4PhysicsTable * BuildLambdaTable(G4EmTableType tType=fRestricted)
const G4ParticleDefinition * theGamma
const G4MaterialCutsCouple * currentCouple
const G4DataVector * theCuts
G4PhysicsTable * theSubLambdaTable
G4PhysicsTable * LambdaTable() const
void SetInverseRangeTable(G4PhysicsTable *p)
const G4ParticleDefinition * SecondaryParticle() const
void SetCrossSectionBiasingFactor(G4double f, G4bool flag=true)
G4VEmFluctuationModel * FluctModel()
G4double SampleSubCutSecondaries(std::vector< G4Track * > &, const G4Step &, G4VEmModel *model, G4int matIdx)
const G4ParticleDefinition * baseParticle
void FillSecondariesAlongStep(G4double &eloss, G4double &weight)
const G4DataVector * theSubCuts
G4PhysicsTable * DEDXTable() const
G4double ScaledKinEnergyForLoss(G4double range)
G4double GetDEDXForScaledEnergy(G4double scaledKinEnergy)
const G4Element * GetCurrentElement() const
const G4ParticleDefinition * theGenericIon
const G4ParticleDefinition * thePositron
void ComputeLambdaForScaledEnergy(G4double scaledKinEnergy)
bool G4bool
Definition: G4Types.hh:79
std::vector< G4Track * > scTracks
void SetLossFluctuations(G4bool val)
void SetDEDXBinningForCSDARange(G4int nbins)
G4double GetIonisationForScaledEnergy(G4double scaledKinEnergy)
G4double GetSubDEDXForScaledEnergy(G4double scaledKinEnergy)
G4PhysicsTable * theCSDARangeTable
G4int NumberOfSubCutoffRegions() const
void SetMaxKinEnergy(G4double e)
G4double GetLambda(G4double &kineticEnergy, const G4MaterialCutsCouple *)
std::vector< G4double > theDEDXAtMaxEnergy
G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &)
const G4ParticleDefinition * BaseParticle() const
Definition: G4Step.hh:76
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=0)
const G4ParticleDefinition * particle
G4PhysicsTable * DEDXTableForSubsec() const
G4double Energy(size_t index) const
G4double Value(G4double theEnergy, size_t &lastidx) const
G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
void SetSecondaryParticle(const G4ParticleDefinition *p)
static const G4double e1
std::vector< const G4Region * > scoffRegions
G4PhysicsTable * theLambdaTable
G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii)
void PreparePhysicsTable(const G4ParticleDefinition &)
const G4ParticleDefinition * theElectron
const G4Material * currentMaterial
void SetLambdaTable(G4PhysicsTable *p)
G4PhysicsTable * theInverseRangeTable
G4SafetyHelper * safetyHelper
std::vector< G4double > theCrossSectionMax
virtual void PrintInfo()=0
G4PhysicsTable * InverseRangeTable() const
G4double CrossSectionBiasingFactor() const
const G4ParticleDefinition * Particle() const
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:419
G4double MeanFreePath(const G4Track &track)
G4PhysicsTable * theDEDXunRestrictedTable
void UpdateEmModel(const G4String &, G4double, G4double)
G4bool StoreTable(const G4ParticleDefinition *p, G4PhysicsTable *, G4bool ascii, const G4String &directory, const G4String &tname)
virtual G4bool IsApplicable(const G4ParticleDefinition &p)=0
G4PhysicsTable * theRangeTableForLoss
void SetEmModel(G4VEmModel *, G4int index=1)
G4double ContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)
G4double MinKinEnergy() const
std::vector< G4VEmModel * > emModels
G4double GetLimitScaledRangeForScaledEnergy(G4double scaledKinEnergy)
G4double MaxKinEnergy() const
G4VEmModel * SelectModel(G4double &energy, size_t &index)
void SetCSDARangeTable(G4PhysicsTable *pRange)
G4double GetDEDXDispersion(const G4MaterialCutsCouple *couple, const G4DynamicParticle *dp, G4double length)
G4EmTableType
G4double GetSubIonisationForScaledEnergy(G4double scaledKinEnergy)
G4double GetLambdaForScaledEnergy(G4double scaledKinEnergy)
void SetMinSubRange(G4double val)
void SetSecondaryRangeTable(G4PhysicsTable *p)
G4PhysicsTable * DEDXunRestrictedTable() const
std::vector< G4double > theRangeAtMaxEnergy
G4bool RetrieveTable(const G4ParticleDefinition *p, G4PhysicsTable *, G4bool ascii, const G4String &directory, const G4String &tname, G4bool mandatory)
G4PhysicsTable * theDEDXTable
void SetLambdaFactor(G4double val)
double G4double
Definition: G4Types.hh:76
void SetLowestEnergyLimit(G4double)
void SetSubLambdaTable(G4PhysicsTable *p)
virtual G4double MinPrimaryEnergy(const G4ParticleDefinition *, const G4Material *, G4double cut)
G4ForceCondition
G4EmBiasingManager * biasManager
void ActivateSubCutoff(G4bool val, const G4Region *region=0)
void SetRangeTableForLoss(G4PhysicsTable *p)
G4PhysicsTable * theDEDXSubTable
#define DBL_MAX
Definition: templates.hh:83
G4double GetScaledRangeForScaledEnergy(G4double scaledKinEnergy)
void SetParticle(const G4ParticleDefinition *p)
void SetBaseParticle(const G4ParticleDefinition *p)
std::vector< G4VEnergyLossProcess * > scProcesses
G4PhysicsTable * BuildDEDXTable(G4EmTableType tType=fRestricted)
G4VEmModel * EmModel(G4int index=1) const
void DefineMaterial(const G4MaterialCutsCouple *couple)
G4VEmModel * SelectModelForMaterial(G4double kinEnergy, size_t &idx) const
void SetDEDXBinning(G4int nbins)
G4double GetRange(G4double &kineticEnergy, const G4MaterialCutsCouple *)
virtual void InitialiseEnergyLossProcess(const G4ParticleDefinition *, const G4ParticleDefinition *)=0
G4GPILSelection
const G4Material * GetMaterial() const
void SetRandomStep(G4bool val)
G4PhysicsTable * theIonisationSubTable
size_t CurrentMaterialCutsCoupleIndex() const
G4double CrossSectionPerVolume(G4double kineticEnergy, const G4MaterialCutsCouple *couple)
std::vector< G4DynamicParticle * > secParticles
G4ProcessType
G4bool IsIonisationProcess() const
void SetMinKinEnergy(G4double e)
G4EmModelManager * modelManager