ECCE @ EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4VEnergyLossProcess.hh
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file G4VEnergyLossProcess.hh
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 //
27 // -------------------------------------------------------------------
28 //
29 // GEANT4 Class header file
30 //
31 //
32 // File name: G4VEnergyLossProcess
33 //
34 // Author: Vladimir Ivanchenko on base of Laszlo Urban code
35 //
36 // Creation date: 03.01.2002
37 //
38 // Modifications: Vladimir Ivanchenko
39 //
40 // Class Description:
41 //
42 // It is the unified energy loss process it calculates the continuous
43 // energy loss for charged particles using a set of Energy Loss
44 // models valid for different energy regions. There are a possibility
45 // to create and access to dE/dx and range tables, or to calculate
46 // that information on fly.
47 
48 // -------------------------------------------------------------------
49 //
50 
51 #ifndef G4VEnergyLossProcess_h
52 #define G4VEnergyLossProcess_h 1
53 
55 #include "globals.hh"
56 #include "G4Material.hh"
57 #include "G4MaterialCutsCouple.hh"
58 #include "G4Track.hh"
59 #include "G4EmModelManager.hh"
60 #include "G4UnitsTable.hh"
62 #include "G4EmTableType.hh"
63 #include "G4PhysicsTable.hh"
64 #include "G4PhysicsVector.hh"
65 #include "G4EmParameters.hh"
66 
67 class G4Step;
69 class G4VEmModel;
71 class G4DataVector;
72 class G4Region;
73 class G4SafetyHelper;
75 class G4VSubCutProducer;
76 class G4EmBiasingManager;
77 class G4LossTableManager;
78 
79 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
80 
82 {
83 public:
84 
85  G4VEnergyLossProcess(const G4String& name = "EnergyLoss",
87 
88  virtual ~G4VEnergyLossProcess();
89 
90 private:
91  // clean vectors and arrays
92  void Clean();
93 
94  //------------------------------------------------------------------------
95  // Virtual methods to be implemented in concrete processes
96  //------------------------------------------------------------------------
97 
98 public:
99  virtual G4bool IsApplicable(const G4ParticleDefinition& p) override = 0;
100 
101  // obsolete to be removed
102  virtual void PrintInfo() {};
103 
104  virtual void ProcessDescription(std::ostream& outFile) const override;
105 
106 protected:
107 
108  virtual void StreamProcessInfo(std::ostream&) const {};
109 
111  const G4ParticleDefinition*) = 0;
112 
113  //------------------------------------------------------------------------
114  // Methods with standard implementation; may be overwritten if needed
115  //------------------------------------------------------------------------
116 
118  const G4Material*, G4double cut);
119 
120  //------------------------------------------------------------------------
121  // Virtual methods implementation common to all EM ContinuousDiscrete
122  // processes. Further inheritance is not assumed
123  //------------------------------------------------------------------------
124 
125 public:
126 
127  // prepare all tables
128  virtual void PreparePhysicsTable(const G4ParticleDefinition&) override;
129 
130  // build all tables
131  virtual void BuildPhysicsTable(const G4ParticleDefinition&) override;
132 
133  // build a table
135 
136  // build a table
138 
139  // Called before tracking of each new G4Track
140  virtual void StartTracking(G4Track*) override;
141 
142  // Step limit from AlongStep
144  const G4Track&,
145  G4double previousStepSize,
146  G4double currentMinimumStep,
147  G4double& currentSafety,
148  G4GPILSelection* selection) override;
149 
150  // Step limit from cross section
152  const G4Track& track,
153  G4double previousStepSize,
154  G4ForceCondition* condition) override;
155 
156  // AlongStep computations
157  virtual G4VParticleChange* AlongStepDoIt(const G4Track&,
158  const G4Step&) override;
159 
160  // Sampling of secondaries in vicinity of geometrical boundary
161  // Return sum of secodaries energy
162  G4double SampleSubCutSecondaries(std::vector<G4Track*>&, const G4Step&,
163  G4VEmModel* model, G4int matIdx);
164 
165  // PostStep sampling of secondaries
166  virtual G4VParticleChange* PostStepDoIt(const G4Track&,
167  const G4Step&) override;
168 
169  // Store all PhysicsTable in files.
170  // Return false in case of any fatal failure at I/O
172  const G4String& directory,
173  G4bool ascii = false) override;
174 
175  // Retrieve all Physics from a files.
176  // Return true if all the Physics Table are built.
177  // Return false if any fatal failure.
179  const G4String& directory,
180  G4bool ascii) override;
181 
182 private:
183 
184  // summary printout after initialisation
185  void StreamInfo(std::ostream& out, const G4ParticleDefinition& part,
186  G4bool rst=false) const;
187 
188  // store a table
190  G4PhysicsTable*, G4bool ascii,
191  const G4String& directory,
192  const G4String& tname);
193 
194  // retrieve a table
196  G4PhysicsTable*, G4bool ascii,
197  const G4String& directory,
198  const G4String& tname,
199  G4bool mandatory);
200 
201  //------------------------------------------------------------------------
202  // Public interface to cross section, mfp and sampling of fluctuations
203  // These methods are not used in run time
204  //------------------------------------------------------------------------
205 
206 public:
207 
208  // access to dispersion of restricted energy loss
210  const G4DynamicParticle* dp,
211  G4double length);
212 
213  // Access to cross section table
215  const G4MaterialCutsCouple* couple);
217  const G4MaterialCutsCouple* couple,
218  G4double logKineticEnergy);
219 
220  // access to cross section
222 
223  // access to step limit
225  G4double previousStepSize,
226  G4double currentMinimumStep,
227  G4double& currentSafety);
228 
229 protected:
230 
231  // implementation of the pure virtual method
232  virtual G4double GetMeanFreePath(const G4Track& track,
233  G4double previousStepSize,
234  G4ForceCondition* condition) override;
235 
236  // implementation of the pure virtual method
238  G4double previousStepSize,
239  G4double currentMinimumStep,
240  G4double& currentSafety) override;
241 
242  //------------------------------------------------------------------------
243  // Run time method which may be also used by derived processes
244  //------------------------------------------------------------------------
245 
246  // creeation of an empty vector for cross section
248  G4double cut);
249 
250  inline size_t CurrentMaterialCutsCoupleIndex() const;
251 
252  //------------------------------------------------------------------------
253  // Specific methods to set, access, modify models
254  //------------------------------------------------------------------------
255 
256  // Select model in run time
257  inline void SelectModel(G4double kinEnergy);
258 
259 public:
260  // Select model by energy and region index
261  inline G4VEmModel* SelectModelForMaterial(G4double kinEnergy,
262  size_t& idx) const;
263 
264  // Add EM model coupled with fluctuation model for region, smaller value
265  // of order defines which pair of models will be selected for a given
266  // energy interval
267  void AddEmModel(G4int, G4VEmModel*,
268  G4VEmFluctuationModel* fluc = 0,
269  const G4Region* region = nullptr);
270 
271  // Define new energy range for the model identified by the name
272  void UpdateEmModel(const G4String&, G4double, G4double);
273 
274  // Assign a model to a process local list, to enable the list in run time
275  // the derived process should execute AddEmModel(..) for all such models
276  void SetEmModel(G4VEmModel*, G4int index=0);
277 
278  // return a model from the local list
279  G4VEmModel* EmModel(size_t index=0) const;
280 
281  // Access to models
282  G4VEmModel* GetModelByIndex(G4int idx = 0, G4bool ver = false) const;
283 
284  G4int NumberOfModels() const;
285 
286  // Assign a fluctuation model to a process
287  inline void SetFluctModel(G4VEmFluctuationModel*);
288 
289  // return the assigned fluctuation model
291 
292  //------------------------------------------------------------------------
293  // Define and access particle type
294  //------------------------------------------------------------------------
295 
296 protected:
297  inline void SetParticle(const G4ParticleDefinition* p);
298  inline void SetSecondaryParticle(const G4ParticleDefinition* p);
299 
300 public:
301  inline void SetBaseParticle(const G4ParticleDefinition* p);
302  inline const G4ParticleDefinition* Particle() const;
303  inline const G4ParticleDefinition* BaseParticle() const;
304  inline const G4ParticleDefinition* SecondaryParticle() const;
305 
306  //------------------------------------------------------------------------
307  // Get/set parameters to configure the process at initialisation time
308  //------------------------------------------------------------------------
309 
310  // Add subcutoff option for the region
311  void ActivateSubCutoff(G4bool val, const G4Region* region = nullptr);
312 
313  // Activate biasing
314  void SetCrossSectionBiasingFactor(G4double f, G4bool flag = true);
315 
317  const G4String& region,
318  G4bool flag = true);
319 
320  void ActivateSecondaryBiasing(const G4String& region, G4double factor,
321  G4double energyLimit);
322 
323  // Add subcutoff process (bremsstrahlung) to sample secondary
324  // particle production in vicinity of the geometry boundary
326 
327  inline void SetLossFluctuations(G4bool val);
328 
329  inline void SetIntegral(G4bool val);
330  inline G4bool IsIntegral() const;
331 
332  // Set/Get flag "isIonisation"
333  void SetIonisation(G4bool val);
334  inline G4bool IsIonisationProcess() const;
335 
336  // Redefine parameteters for stepping control
337  void SetLinearLossLimit(G4double val);
338  void SetStepFunction(G4double v1, G4double v2, G4bool lock=true);
340 
341  inline G4int NumberOfSubCutoffRegions() const;
342 
343  //------------------------------------------------------------------------
344  // Specific methods to path Physics Tables to the process
345  //------------------------------------------------------------------------
346 
347  void SetDEDXTable(G4PhysicsTable* p, G4EmTableType tType);
348  void SetCSDARangeTable(G4PhysicsTable* pRange);
352 
355 
356  // Binning for dEdx, range, inverse range and labda tables
357  void SetDEDXBinning(G4int nbins);
358 
359  // Min kinetic energy for tables
360  void SetMinKinEnergy(G4double e);
361  inline G4double MinKinEnergy() const;
362 
363  // Max kinetic energy for tables
364  void SetMaxKinEnergy(G4double e);
365  inline G4double MaxKinEnergy() const;
366 
367  // Biasing parameters
368  inline G4double CrossSectionBiasingFactor() const;
369 
370  // Return values for given G4MaterialCutsCouple
371  inline G4double GetDEDX(G4double kineticEnergy, const G4MaterialCutsCouple*);
372  inline G4double GetDEDX(G4double kineticEnergy, const G4MaterialCutsCouple*,
373  G4double logKineticEnergy);
374  inline G4double GetDEDXForSubsec(G4double kineticEnergy,
375  const G4MaterialCutsCouple*);
376  inline G4double GetRange(G4double kineticEnergy, const G4MaterialCutsCouple*);
377  inline G4double GetRange(G4double kineticEnergy, const G4MaterialCutsCouple*,
378  G4double logKineticEnergy);
379  inline G4double GetCSDARange(G4double kineticEnergy,
380  const G4MaterialCutsCouple*);
381  inline G4double GetRangeForLoss(G4double kineticEnergy,
382  const G4MaterialCutsCouple*);
383  inline G4double GetRangeForLoss(G4double kineticEnergy,
384  const G4MaterialCutsCouple*,
385  G4double logKineticEnergy);
386  inline G4double GetKineticEnergy(G4double range,
387  const G4MaterialCutsCouple*);
388  inline G4double GetLambda(G4double kineticEnergy,const G4MaterialCutsCouple*);
389  inline G4double GetLambda(G4double kineticEnergy,const G4MaterialCutsCouple*,
390  G4double logKineticEnergy);
391 
392  inline G4bool TablesAreBuilt() const;
393 
394  // Access to specific tables
395  inline G4PhysicsTable* DEDXTable() const;
396  inline G4PhysicsTable* DEDXTableForSubsec() const;
397  inline G4PhysicsTable* DEDXunRestrictedTable() const;
398  inline G4PhysicsTable* IonisationTable() const;
400  inline G4PhysicsTable* CSDARangeTable() const;
401  inline G4PhysicsTable* SecondaryRangeTable() const;
402  inline G4PhysicsTable* RangeTableForLoss() const;
403  inline G4PhysicsTable* InverseRangeTable() const;
404  inline G4PhysicsTable* LambdaTable() const;
405  inline G4PhysicsTable* SubLambdaTable() const;
406 
407  //------------------------------------------------------------------------
408  // Run time method for simulation of ionisation
409  //------------------------------------------------------------------------
410 
411  // access atom on which interaction happens
412  const G4Element* GetCurrentElement() const;
413 
414  // Set scaling parameters for ions is needed to G4EmCalculator
415  inline void SetDynamicMassCharge(G4double massratio, G4double charge2ratio);
416 
417 private:
418 
420 
421  void PrintWarning(G4String, G4double val);
422 
423  // define material and indexes
424  inline void DefineMaterial(const G4MaterialCutsCouple* couple);
425 
426  //------------------------------------------------------------------------
427  // Compute values using scaling relation, mass and charge of based particle
428  //------------------------------------------------------------------------
429  inline G4double GetDEDXForScaledEnergy(G4double scaledKinEnergy);
430  inline G4double GetDEDXForScaledEnergy(G4double scaledKinEnergy,
431  G4double logScaledKinEnergy);
432  inline G4double GetSubDEDXForScaledEnergy(G4double scaledKinEnergy);
433  inline G4double GetIonisationForScaledEnergy(G4double scaledKinEnergy);
434  inline G4double GetSubIonisationForScaledEnergy(G4double scaledKinEnergy);
435  inline G4double GetScaledRangeForScaledEnergy(G4double scaledKinEnergy);
436  inline G4double GetScaledRangeForScaledEnergy(G4double scaledKinEnergy,
437  G4double logScaledKinEnergy);
438  inline G4double GetLimitScaledRangeForScaledEnergy(G4double scaledKinEnergy);
440  G4double logScaledKinEnergy);
442  inline G4double GetLambdaForScaledEnergy(G4double scaledKinEnergy);
443  inline G4double GetLambdaForScaledEnergy(G4double scaledKinEnergy,
444  G4double logScaledKinEnergy);
445  void ComputeLambdaForScaledEnergy(G4double scaledKinEnergy,
446  G4double logScaledKinEnergy);
447 
448  // hide assignment operator
451 
452  // ======== Parameters of the class fixed at construction =========
453 
459 
465 
466  // ======== Parameters of the class fixed at initialisation =======
467 
468  std::vector<G4VEmModel*> emModels;
472  std::vector<const G4Region*> scoffRegions;
475 
476  std::vector<G4VEnergyLossProcess*> scProcesses;
478 
479  // tables and vectors
491 
492  size_t idxDEDX;
493  size_t idxDEDXSub;
497  size_t idxRange;
498  size_t idxCSDA;
499  size_t idxSecRange;
501  size_t idxLambda;
502  size_t idxSubLambda;
503 
504  std::vector<G4double> theDEDXAtMaxEnergy;
505  std::vector<G4double> theRangeAtMaxEnergy;
506  std::vector<G4double> theEnergyOfCrossSectionMax;
507  std::vector<G4double> theCrossSectionMax;
508 
509  const std::vector<G4double>* theDensityFactor;
510  const std::vector<G4int>* theDensityIdx;
511 
514 
516 
519 
524 
531 
550 
551 protected:
552 
557 
567 
568  // ======== Cached values - may be state dependent ================
569 
570 private:
571 
572  std::vector<G4DynamicParticle*> secParticles;
573  std::vector<G4Track*> scTracks;
574 
576 
579  size_t lastIdx;
580 
586 
588 
592 };
593 
594 // ======== Run time inline methods ================
595 
597 {
598  return currentCoupleIndex;
599 }
600 
601 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
602 
604 {
607 }
608 
609 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
610 
612  G4double kinEnergy, size_t& idx) const
613 {
614  return modelManager->SelectModel(kinEnergy, idx);
615 }
616 
617 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
618 
619 inline void
621 {
622  if(couple != currentCouple) {
623  currentCouple = couple;
624  currentMaterial = couple->GetMaterial();
625  currentCoupleIndex = couple->GetIndex();
626  basedCoupleIndex = (*theDensityIdx)[currentCoupleIndex];
627  fFactor = chargeSqRatio*biasFactor*(*theDensityFactor)[currentCoupleIndex];
630  idxLambda = idxSubLambda = 0;
631  }
632 }
633 
634 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
635 
637  G4double charge2ratio)
638 {
639  massRatio = massratio;
641  fFactor = charge2ratio*biasFactor*(*theDensityFactor)[currentCoupleIndex];
642  chargeSqRatio = charge2ratio;
644 }
645 
646 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
647 
649 {
650  /*
651  G4cout << "G4VEnergyLossProcess::GetDEDX: Idx= "
652  << basedCoupleIndex << " E(MeV)= " << e
653  << " Emin= " << minKinEnergy << " Factor= " << fFactor
654  << " " << theDEDXTable << G4endl; */
655  G4double x = fFactor*(*theDEDXTable)[basedCoupleIndex]->Value(e, idxDEDX);
656  if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
657  return x;
658 }
659 
660 inline
662 {
663  /*
664  G4cout << "G4VEnergyLossProcess::GetDEDX: Idx= "
665  << basedCoupleIndex << " E(MeV)= " << e
666  << " Emin= " << minKinEnergy << " Factor= " << fFactor
667  << " " << theDEDXTable << G4endl; */
668  G4double x = fFactor*(*theDEDXTable)[basedCoupleIndex]->LogVectorValue(e,loge);
669  if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
670  return x;
671 }
672 
673 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
674 
676 {
677  G4double x =
678  fFactor*(*theDEDXSubTable)[basedCoupleIndex]->Value(e, idxDEDXSub);
679  if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
680  return x;
681 }
682 
683 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
684 
686 {
687  G4double x =
688  fFactor*(*theIonisationTable)[basedCoupleIndex]->Value(e, idxIonisation);
689  if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
690  return x;
691 }
692 
693 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
694 
695 inline
697 {
698  G4double x = fFactor*
699  (*theIonisationSubTable)[basedCoupleIndex]->Value(e, idxIonisationSub);
700  if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
701  return x;
702 }
703 
704 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
705 
707 {
708  //G4cout << "G4VEnergyLossProcess::GetScaledRange: Idx= "
709  // << basedCoupleIndex << " E(MeV)= " << e
710  // << " lastIdx= " << lastIdx << " " << theRangeTableForLoss << G4endl;
714  computedRange =
715  ((*theRangeTableForLoss)[basedCoupleIndex])->Value(e, idxRange);
716  if(e < minKinEnergy) { computedRange *= std::sqrt(e/minKinEnergy); }
717  }
718  //G4cout << "G4VEnergyLossProcess::GetScaledRange: Idx= "
719  // << basedCoupleIndex << " E(MeV)= " << e
720  // << " R= " << fRange << " " << theRangeTableForLoss << G4endl;
721 
722  return computedRange;
723 }
724 
725 inline G4double
727 {
728  //G4cout << "G4VEnergyLossProcess::GetScaledRange: Idx= "
729  // << basedCoupleIndex << " E(MeV)= " << e
730  // << " lastIdx= " << lastIdx << " " << theRangeTableForLoss << G4endl;
734  computedRange =
735  ((*theRangeTableForLoss)[basedCoupleIndex])->LogVectorValue(e, loge);
736  if(e < minKinEnergy) { computedRange *= std::sqrt(e/minKinEnergy); }
737  }
738  //G4cout << "G4VEnergyLossProcess::GetScaledRange: Idx= "
739  // << basedCoupleIndex << " E(MeV)= " << e
740  // << " R= " << fRange << " " << theRangeTableForLoss << G4endl;
741 
742  return computedRange;
743 }
744 
745 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
746 inline G4double
748 {
749  G4double x;
750  if (e < maxKinEnergyCSDA) {
751  x = ((*theCSDARangeTable)[basedCoupleIndex])->Value(e, idxCSDA);
752  if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
753  } else {
756  }
757  return x;
758 }
759 
760 inline G4double
762  G4double loge)
763 {
764  G4double x;
765  if (e < maxKinEnergyCSDA) {
766  x = ((*theCSDARangeTable)[basedCoupleIndex])->LogVectorValue(e, loge);
767  if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
768  } else {
771  }
772  return x;
773 }
774 
775 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
776 
778 {
779  //G4cout << "G4VEnergyLossProcess::GetEnergy: Idx= "
780  // << basedCoupleIndex << " R(mm)= " << r << " "
781  // << theInverseRangeTable << G4endl;
782  G4PhysicsVector* v = (*theInverseRangeTable)[basedCoupleIndex];
783  G4double rmin = v->Energy(0);
784  G4double e = 0.0;
785  if(r >= rmin) { e = v->Value(r, idxInverseRange); }
786  else if(r > 0.0) {
787  G4double x = r/rmin;
788  e = minKinEnergy*x*x;
789  }
790  return e;
791 }
792 
793 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
794 
796 {
797  return fFactor*((*theLambdaTable)[basedCoupleIndex])->Value(e, idxLambda);
798 }
799 
800 inline G4double
802 {
803  return fFactor*((*theLambdaTable)[basedCoupleIndex])->LogVectorValue(e,loge);
804 }
805 
806 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
807 
808 inline G4double
810  const G4MaterialCutsCouple* couple)
811 {
812  DefineMaterial(couple);
813  return GetDEDXForScaledEnergy(kinEnergy*massRatio);
814 }
815 
816 inline G4double
818  const G4MaterialCutsCouple* couple,
819  G4double logKinEnergy)
820 {
821  DefineMaterial(couple);
822  return GetDEDXForScaledEnergy(kinEnergy*massRatio, logKinEnergy+logMassRatio);
823 }
824 
825 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
826 
827 inline G4double
829  const G4MaterialCutsCouple* couple)
830 {
831  DefineMaterial(couple);
832  return GetSubDEDXForScaledEnergy(kineticEnergy*massRatio);
833 }
834 
835 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
836 
837 inline G4double
839  const G4MaterialCutsCouple* couple)
840 {
841  G4double x = fRange;
842  DefineMaterial(couple);
843  if(theCSDARangeTable) {
845  } else if(theRangeTableForLoss) {
847  }
848  return x;
849 }
850 
851 inline G4double
853  const G4MaterialCutsCouple* couple,
854  G4double logKinEnergy)
855 {
856  G4double x = fRange;
857  DefineMaterial(couple);
858  if(theCSDARangeTable) {
860  logKinEnergy+logMassRatio);
861  } else if(theRangeTableForLoss) {
863  logKinEnergy+logMassRatio);
864  }
865  return x;
866 }
867 
868 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
869 
870 inline G4double
872  const G4MaterialCutsCouple* couple)
873 {
874  DefineMaterial(couple);
875  return (theCSDARangeTable) ?
877  : DBL_MAX;
878 }
879 
880 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
881 
882 inline G4double
884  const G4MaterialCutsCouple* couple)
885 {
886  // G4cout << "GetRangeForLoss: Range from " << GetProcessName() << G4endl;
887  DefineMaterial(couple);
889 }
890 
891 inline G4double
893  const G4MaterialCutsCouple* couple,
894  G4double logKinEnergy)
895 {
896  // G4cout << "GetRangeForLoss: Range from " << GetProcessName() << G4endl;
897  DefineMaterial(couple);
899  logKinEnergy+logMassRatio);
900 }
901 
902 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
903 
904 inline G4double
906  const G4MaterialCutsCouple* couple)
907 {
908  DefineMaterial(couple);
910 }
911 
912 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
913 
914 inline G4double
916  const G4MaterialCutsCouple* couple)
917 {
918  DefineMaterial(couple);
919  return theLambdaTable ? GetLambdaForScaledEnergy(kinEnergy*massRatio) : 0.0;
920 }
921 
922 inline G4double
924  const G4MaterialCutsCouple* couple,
925  G4double logKinEnergy)
926 {
927  DefineMaterial(couple);
928  return theLambdaTable
929  ? GetLambdaForScaledEnergy(kinEnergy*massRatio, logKinEnergy+logMassRatio)
930  : 0.0;
931 }
932 
933 // ======== Get/Set inline methods used at initialisation ================
934 
936 {
937  fluctModel = p;
938 }
939 
940 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
941 
943 {
944  return fluctModel;
945 }
946 
947 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
948 
950 {
951  particle = p;
952 }
953 
954 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
955 
956 inline void
958 {
960 }
961 
962 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
963 
964 inline void
966 {
967  baseParticle = p;
968 }
969 
970 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
971 
973 {
974  return particle;
975 }
976 
977 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
978 
980 {
981  return baseParticle;
982 }
983 
984 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
985 
986 inline const G4ParticleDefinition*
988 {
989  return secondaryParticle;
990 }
991 
992 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
993 
995 {
996  lossFluctuationFlag = val;
997  actLossFluc = true;
998 }
999 
1000 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1001 
1003 {
1004  integral = val;
1005  actIntegral = true;
1006 }
1007 
1008 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1009 
1011 {
1012  return integral;
1013 }
1014 
1015 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1016 
1018 {
1019  return isIonisation;
1020 }
1021 
1022 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1023 
1025 {
1026  return nSCoffRegions;
1027 }
1028 
1029 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1030 
1032 {
1033  return minKinEnergy;
1034 }
1035 
1036 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1037 
1039 {
1040  return maxKinEnergy;
1041 }
1042 
1043 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1044 
1046 {
1047  return biasFactor;
1048 }
1049 
1050 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1051 
1053 {
1054  return tablesAreBuilt;
1055 }
1056 
1057 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1058 
1060 {
1061  return theDEDXTable;
1062 }
1063 
1064 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1065 
1067 {
1068  return theDEDXSubTable;
1069 }
1070 
1071 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1072 
1074 {
1075  return theDEDXunRestrictedTable;
1076 }
1077 
1078 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1079 
1081 {
1082  return theIonisationTable;
1083 }
1084 
1085 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1086 
1088 {
1089  return theIonisationSubTable;
1090 }
1091 
1092 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1093 
1095 {
1096  return theCSDARangeTable;
1097 }
1098 
1099 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1100 
1102 {
1103  return theSecondaryRangeTable;
1104 }
1105 
1106 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1107 
1109 {
1110  return theRangeTableForLoss;
1111 }
1112 
1113 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1114 
1116 {
1117  return theInverseRangeTable;
1118 }
1119 
1120 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1121 
1123 {
1124  return theLambdaTable;
1125 }
1126 
1127 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1128 
1130 {
1131  return theSubLambdaTable;
1132 }
1133 
1134 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1135 
1136 #endif