Geant4  10.02.p02
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 93264 2015-10-14 09:30:04Z 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 #include "G4EmParameters.hh"
109 
110 class G4Step;
112 class G4VEmModel;
114 class G4DataVector;
115 class G4Region;
116 class G4SafetyHelper;
117 class G4VAtomDeexcitation;
118 class G4VSubCutProducer;
119 class G4EmBiasingManager;
120 class G4LossTableManager;
121 
122 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
123 
125 {
126 public:
127 
128  G4VEnergyLossProcess(const G4String& name = "EnergyLoss",
130 
131  virtual ~G4VEnergyLossProcess();
132 
133 private:
134  // clean vectors and arrays
135  void Clean();
136 
137  //------------------------------------------------------------------------
138  // Virtual methods to be implemented in concrete processes
139  //------------------------------------------------------------------------
140 
141 public:
142  virtual G4bool IsApplicable(const G4ParticleDefinition& p) = 0;
143 
144  virtual void PrintInfo() = 0;
145 
146  virtual void ProcessDescription(std::ostream& outFile) const; // = 0;
147 
148 protected:
149 
151  const G4ParticleDefinition*) = 0;
152 
153  //------------------------------------------------------------------------
154  // Methods with standard implementation; may be overwritten if needed
155  //------------------------------------------------------------------------
156 
158  const G4Material*, G4double cut);
159 
160  //------------------------------------------------------------------------
161  // Virtual methods implementation common to all EM ContinuousDiscrete
162  // processes. Further inheritance is not assumed
163  //------------------------------------------------------------------------
164 
165 public:
166 
167  // prepare all tables
168  virtual void PreparePhysicsTable(const G4ParticleDefinition&);
169 
170  // build all tables
171  virtual void BuildPhysicsTable(const G4ParticleDefinition&);
172 
173  // build a table
175 
176  // build a table
178 
179  // summary printout after initialisation
180  void PrintInfoDefinition(const G4ParticleDefinition& part);
181 
182  // Called before tracking of each new G4Track
183  virtual void StartTracking(G4Track*);
184 
185  // Step limit from AlongStep
187  const G4Track&,
188  G4double previousStepSize,
189  G4double currentMinimumStep,
190  G4double& currentSafety,
191  G4GPILSelection* selection);
192 
193  // Step limit from cross section
195  const G4Track& track,
196  G4double previousStepSize,
198 
199  // AlongStep computations
200  virtual G4VParticleChange* AlongStepDoIt(const G4Track&, const G4Step&);
201 
202  // Sampling of secondaries in vicinity of geometrical boundary
203  // Return sum of secodaries energy
204  G4double SampleSubCutSecondaries(std::vector<G4Track*>&, const G4Step&,
205  G4VEmModel* model, G4int matIdx);
206 
207  // PostStep sampling of secondaries
208  virtual G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step&);
209 
210  // Store all PhysicsTable in files.
211  // Return false in case of any fatal failure at I/O
213  const G4String& directory,
214  G4bool ascii = false);
215 
216  // Retrieve all Physics from a files.
217  // Return true if all the Physics Table are built.
218  // Return false if any fatal failure.
220  const G4String& directory,
221  G4bool ascii);
222 
223 private:
224  // store a table
226  G4PhysicsTable*, G4bool ascii,
227  const G4String& directory,
228  const G4String& tname);
229 
230  // retrieve a table
232  G4PhysicsTable*, G4bool ascii,
233  const G4String& directory,
234  const G4String& tname,
235  G4bool mandatory);
236 
237  //------------------------------------------------------------------------
238  // Public interface to cross section, mfp and sampling of fluctuations
239  // These methods are not used in run time
240  //------------------------------------------------------------------------
241 
242 public:
243 
244  // access to dispersion of restricted energy loss
246  const G4DynamicParticle* dp,
247  G4double length);
248 
249  // Access to cross section table
251  const G4MaterialCutsCouple* couple);
252 
253  // access to cross section
254  G4double MeanFreePath(const G4Track& track);
255 
256  // access to step limit
257  G4double ContinuousStepLimit(const G4Track& track,
258  G4double previousStepSize,
259  G4double currentMinimumStep,
260  G4double& currentSafety);
261 
262 protected:
263 
264  // implementation of the pure virtual method
265  virtual G4double GetMeanFreePath(const G4Track& track,
266  G4double previousStepSize,
267  G4ForceCondition* condition);
268 
269  // implementation of the pure virtual method
270  virtual G4double GetContinuousStepLimit(const G4Track& track,
271  G4double previousStepSize,
272  G4double currentMinimumStep,
273  G4double& currentSafety);
274 
275  //------------------------------------------------------------------------
276  // Run time method which may be also used by derived processes
277  //------------------------------------------------------------------------
278 
279  // creeation of an empty vector for cross section
281  G4double cut);
282 
283  inline size_t CurrentMaterialCutsCoupleIndex() const;
284 
285  //------------------------------------------------------------------------
286  // Specific methods to set, access, modify models
287  //------------------------------------------------------------------------
288 
289  // Select model in run time
290  inline void SelectModel(G4double kinEnergy);
291 
292 public:
293  // Select model by energy and region index
294  inline G4VEmModel* SelectModelForMaterial(G4double kinEnergy,
295  size_t& idx) const;
296 
297  // Add EM model coupled with fluctuation model for region, smaller value
298  // of order defines which pair of models will be selected for a given
299  // energy interval
300  void AddEmModel(G4int, G4VEmModel*,
301  G4VEmFluctuationModel* fluc = 0,
302  const G4Region* region = 0);
303 
304  // Define new energy range for the model identified by the name
305  void UpdateEmModel(const G4String&, G4double, G4double);
306 
307  // Assign a model to a process
308  void SetEmModel(G4VEmModel*, G4int index=1);
309 
310  // return the assigned model
311  G4VEmModel* EmModel(G4int index=1) const;
312 
313  // Access to models
314  G4VEmModel* GetModelByIndex(G4int idx = 0, G4bool ver = false) const;
315 
316  G4int NumberOfModels() const;
317 
318  // Assign a fluctuation model to a process
320 
321  // return the assigned fluctuation model
323 
324  //------------------------------------------------------------------------
325  // Define and access particle type
326  //------------------------------------------------------------------------
327 
328 protected:
329  inline void SetParticle(const G4ParticleDefinition* p);
330  inline void SetSecondaryParticle(const G4ParticleDefinition* p);
331 
332 public:
333  inline void SetBaseParticle(const G4ParticleDefinition* p);
334  inline const G4ParticleDefinition* Particle() const;
335  inline const G4ParticleDefinition* BaseParticle() const;
336  inline const G4ParticleDefinition* SecondaryParticle() const;
337 
338  //------------------------------------------------------------------------
339  // Get/set parameters to configure the process at initialisation time
340  //------------------------------------------------------------------------
341 
342  // Add subcutoff option for the region
343  void ActivateSubCutoff(G4bool val, const G4Region* region = 0);
344 
345  // Activate biasing
346  void SetCrossSectionBiasingFactor(G4double f, G4bool flag = true);
347 
348  void ActivateForcedInteraction(G4double length = 0.0,
349  const G4String& region = "",
350  G4bool flag = true);
351 
352  void ActivateSecondaryBiasing(const G4String& region, G4double factor,
353  G4double energyLimit);
354 
355  // Add subcutoff process (bremsstrahlung) to sample secondary
356  // particle production in vicinity of the geometry boundary
358 
359  inline void SetLossFluctuations(G4bool val);
360 
361  inline void SetIntegral(G4bool val);
362  inline G4bool IsIntegral() const;
363 
364  // Set/Get flag "isIonisation"
365  void SetIonisation(G4bool val);
366  inline G4bool IsIonisationProcess() const;
367 
368  // Redefine parameteters for stepping control
369  void SetLinearLossLimit(G4double val);
370  void SetStepFunction(G4double v1, G4double v2);
372 
373  inline G4int NumberOfSubCutoffRegions() const;
374 
375  //------------------------------------------------------------------------
376  // Specific methods to path Physics Tables to the process
377  //------------------------------------------------------------------------
378 
379  void SetDEDXTable(G4PhysicsTable* p, G4EmTableType tType);
380  void SetCSDARangeTable(G4PhysicsTable* pRange);
384 
387 
388  // Binning for dEdx, range, inverse range and labda tables
389  void SetDEDXBinning(G4int nbins);
390 
391  // Min kinetic energy for tables
392  void SetMinKinEnergy(G4double e);
393  inline G4double MinKinEnergy() const;
394 
395  // Max kinetic energy for tables
396  void SetMaxKinEnergy(G4double e);
397  inline G4double MaxKinEnergy() const;
398 
399  // Biasing parameters
400  inline G4double CrossSectionBiasingFactor() const;
401 
402  // Return values for given G4MaterialCutsCouple
403  inline G4double GetDEDX(G4double& kineticEnergy,
404  const G4MaterialCutsCouple*);
405  inline G4double GetDEDXForSubsec(G4double& kineticEnergy,
406  const G4MaterialCutsCouple*);
407  inline G4double GetRange(G4double& kineticEnergy,
408  const G4MaterialCutsCouple*);
409  inline G4double GetCSDARange(G4double& kineticEnergy,
410  const G4MaterialCutsCouple*);
411  inline G4double GetRangeForLoss(G4double& kineticEnergy,
412  const G4MaterialCutsCouple*);
413  inline G4double GetKineticEnergy(G4double& range,
414  const G4MaterialCutsCouple*);
415  inline G4double GetLambda(G4double& kineticEnergy,
416  const G4MaterialCutsCouple*);
417 
418  inline G4bool TablesAreBuilt() const;
419 
420  // Access to specific tables
421  inline G4PhysicsTable* DEDXTable() const;
422  inline G4PhysicsTable* DEDXTableForSubsec() const;
423  inline G4PhysicsTable* DEDXunRestrictedTable() const;
424  inline G4PhysicsTable* IonisationTable() const;
426  inline G4PhysicsTable* CSDARangeTable() const;
427  inline G4PhysicsTable* SecondaryRangeTable() const;
428  inline G4PhysicsTable* RangeTableForLoss() const;
429  inline G4PhysicsTable* InverseRangeTable() const;
430  inline G4PhysicsTable* LambdaTable() const;
431  inline G4PhysicsTable* SubLambdaTable() const;
432 
433  //------------------------------------------------------------------------
434  // Run time method for simulation of ionisation
435  //------------------------------------------------------------------------
436 
437  // access atom on which interaction happens
438  const G4Element* GetCurrentElement() const;
439 
440  // Set scaling parameters for ions is needed to G4EmCalculator
441  inline void SetDynamicMassCharge(G4double massratio, G4double charge2ratio);
442 
443 private:
444 
445  void FillSecondariesAlongStep(G4double& eloss, G4double& weight);
446 
447  void PrintWarning(G4String, G4double val);
448 
449  // define material and indexes
450  inline void DefineMaterial(const G4MaterialCutsCouple* couple);
451 
452  //------------------------------------------------------------------------
453  // Compute values using scaling relation, mass and charge of based particle
454  //------------------------------------------------------------------------
455 
456  inline G4double GetDEDXForScaledEnergy(G4double scaledKinEnergy);
457  inline G4double GetSubDEDXForScaledEnergy(G4double scaledKinEnergy);
458  inline G4double GetIonisationForScaledEnergy(G4double scaledKinEnergy);
459  inline G4double GetSubIonisationForScaledEnergy(G4double scaledKinEnergy);
460  inline G4double GetScaledRangeForScaledEnergy(G4double scaledKinEnergy);
461  inline G4double GetLimitScaledRangeForScaledEnergy(G4double scaledKinEnergy);
463  inline G4double GetLambdaForScaledEnergy(G4double scaledKinEnergy);
464  inline void ComputeLambdaForScaledEnergy(G4double scaledKinEnergy);
465 
466  // hide assignment operator
469 
470  // ======== Parameters of the class fixed at construction =========
471 
477 
483 
484  // ======== Parameters of the class fixed at initialisation =======
485 
486  std::vector<G4VEmModel*> emModels;
490  std::vector<const G4Region*> scoffRegions;
493 
494  std::vector<G4VEnergyLossProcess*> scProcesses;
496 
497  // tables and vectors
509 
510  size_t idxDEDX;
511  size_t idxDEDXSub;
515  size_t idxRange;
516  size_t idxCSDA;
517  size_t idxSecRange;
519  size_t idxLambda;
520  size_t idxSubLambda;
521 
522  std::vector<G4double> theDEDXAtMaxEnergy;
523  std::vector<G4double> theRangeAtMaxEnergy;
524  std::vector<G4double> theEnergyOfCrossSectionMax;
525  std::vector<G4double> theCrossSectionMax;
526 
527  const std::vector<G4double>* theDensityFactor;
528  const std::vector<G4int>* theDensityIdx;
529 
532 
534 
537 
542 
548 
565 
566 protected:
567 
569 
570  // ======== Cached values - may be state dependent ================
571 
572 private:
573 
574  std::vector<G4DynamicParticle*> secParticles;
575  std::vector<G4Track*> scTracks;
576 
578 
584  size_t lastIdx;
585 
587 
592 
600 
602 
606 };
607 
608 // ======== Run time inline methods ================
609 
611 {
612  return currentCoupleIndex;
613 }
614 
615 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
616 
618 {
621 }
622 
623 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
624 
626  G4double kinEnergy, size_t& idx) const
627 {
628  return modelManager->SelectModel(kinEnergy, idx);
629 }
630 
631 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
632 
633 inline void
635 {
636  if(couple != currentCouple) {
637  currentCouple = couple;
638  currentMaterial = couple->GetMaterial();
639  currentCoupleIndex = couple->GetIndex();
640  basedCoupleIndex = (*theDensityIdx)[currentCoupleIndex];
641  fFactor = chargeSqRatio*biasFactor*(*theDensityFactor)[currentCoupleIndex];
644  idxLambda = idxSubLambda = 0;
645  }
646 }
647 
648 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
649 
651  G4double charge2ratio)
652 {
653  massRatio = massratio;
654  fFactor = charge2ratio*biasFactor*(*theDensityFactor)[currentCoupleIndex];
655  chargeSqRatio = charge2ratio;
657 }
658 
659 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
660 
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]->Value(e, idxDEDX);
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;
713  preStepRangeEnergy = e;
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 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
726 
727 inline G4double
729 {
730  G4double x;
731  if (e < maxKinEnergyCSDA) {
732  x = ((*theCSDARangeTable)[basedCoupleIndex])->Value(e, idxCSDA);
733  if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
734  } else {
737  }
738  return x;
739 }
740 
741 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
742 
744 {
745  //G4cout << "G4VEnergyLossProcess::GetEnergy: Idx= "
746  // << basedCoupleIndex << " R(mm)= " << r << " "
747  // << theInverseRangeTable << G4endl;
748  G4PhysicsVector* v = (*theInverseRangeTable)[basedCoupleIndex];
749  G4double rmin = v->Energy(0);
750  G4double e = 0.0;
751  if(r >= rmin) { e = v->Value(r, idxInverseRange); }
752  else if(r > 0.0) {
753  G4double x = r/rmin;
754  e = minKinEnergy*x*x;
755  }
756  return e;
757 }
758 
759 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
760 
762 {
763  return fFactor*((*theLambdaTable)[basedCoupleIndex])->Value(e, idxLambda);
764 }
765 
766 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
767 
768 inline G4double
770  const G4MaterialCutsCouple* couple)
771 {
772  DefineMaterial(couple);
773  return GetDEDXForScaledEnergy(kineticEnergy*massRatio);
774 }
775 
776 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
777 
778 inline G4double
780  const G4MaterialCutsCouple* couple)
781 {
782  DefineMaterial(couple);
783  return GetSubDEDXForScaledEnergy(kineticEnergy*massRatio);
784 }
785 
786 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
787 
788 inline G4double
790  const G4MaterialCutsCouple* couple)
791 {
792  G4double x = fRange;
793  DefineMaterial(couple);
794  if(theCSDARangeTable) {
796  * reduceFactor;
797  } else if(theRangeTableForLoss) {
799  }
800  return x;
801 }
802 
803 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
804 
805 inline G4double
807  const G4MaterialCutsCouple* couple)
808 {
809  DefineMaterial(couple);
810  G4double x = DBL_MAX;
811  if(theCSDARangeTable) {
813  }
814  return x;
815 }
816 
817 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
818 
819 inline G4double
821  const G4MaterialCutsCouple* couple)
822 {
823  // G4cout << "GetRangeForLoss: Range from " << GetProcessName() << G4endl;
824  DefineMaterial(couple);
825  G4double x =
827  //G4cout << "GetRangeForLoss: Range from " << GetProcessName()
828  // << " e= " << kineticEnergy << " r= " << x << G4endl;
829  return x;
830 }
831 
832 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
833 
834 inline G4double
836  const G4MaterialCutsCouple* couple)
837 {
838  DefineMaterial(couple);
840 }
841 
842 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
843 
844 inline G4double
846  const G4MaterialCutsCouple* couple)
847 {
848  DefineMaterial(couple);
849  G4double x = 0.0;
850  if(theLambdaTable) {
851  x = GetLambdaForScaledEnergy(kineticEnergy*massRatio);
852  }
853  return x;
854 }
855 
856 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
857 
859 {
861  if (e <= mfpKinEnergy) {
863 
864  } else {
866  if(e1 > mfpKinEnergy) {
868  G4double preStepLambda1 = GetLambdaForScaledEnergy(e1);
869  if(preStepLambda1 > preStepLambda) {
870  mfpKinEnergy = e1;
871  preStepLambda = preStepLambda1;
872  }
873  } else {
875  }
876  }
877 }
878 
879 // ======== Get/Set inline methods used at initialisation ================
880 
882 {
883  fluctModel = p;
884 }
885 
886 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
887 
889 {
890  return fluctModel;
891 }
892 
893 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
894 
896 {
897  particle = p;
898 }
899 
900 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
901 
902 inline void
904 {
905  secondaryParticle = p;
906 }
907 
908 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
909 
910 inline void
912 {
913  baseParticle = p;
914 }
915 
916 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
917 
919 {
920  return particle;
921 }
922 
923 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
924 
926 {
927  return baseParticle;
928 }
929 
930 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
931 
932 inline const G4ParticleDefinition*
934 {
935  return secondaryParticle;
936 }
937 
938 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
939 
941 {
942  lossFluctuationFlag = val;
943 }
944 
945 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
946 
948 {
949  integral = val;
950 }
951 
952 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
953 
955 {
956  return integral;
957 }
958 
959 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
960 
962 {
963  return isIonisation;
964 }
965 
966 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
967 
969 {
970  return nSCoffRegions;
971 }
972 
973 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
974 
976 {
977  return minKinEnergy;
978 }
979 
980 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
981 
983 {
984  return maxKinEnergy;
985 }
986 
987 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
988 
990 {
991  return biasFactor;
992 }
993 
994 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
995 
997 {
998  return tablesAreBuilt;
999 }
1000 
1001 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1002 
1004 {
1005  return theDEDXTable;
1006 }
1007 
1008 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1009 
1011 {
1012  return theDEDXSubTable;
1013 }
1014 
1015 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1016 
1018 {
1019  return theDEDXunRestrictedTable;
1020 }
1021 
1022 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1023 
1025 {
1026  return theIonisationTable;
1027 }
1028 
1029 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1030 
1032 {
1033  return theIonisationSubTable;
1034 }
1035 
1036 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1037 
1039 {
1040  return theCSDARangeTable;
1041 }
1042 
1043 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1044 
1046 {
1047  return theSecondaryRangeTable;
1048 }
1049 
1050 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1051 
1053 {
1054  return theRangeTableForLoss;
1055 }
1056 
1057 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1058 
1060 {
1061  return theInverseRangeTable;
1062 }
1063 
1064 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1065 
1067 {
1068  return theLambdaTable;
1069 }
1070 
1071 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1072 
1074 {
1075  return theSubLambdaTable;
1076 }
1077 
1078 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1079 
1080 #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)
G4GPILSelection aGPILSelection
void SetDynamicMassCharge(G4double massratio, G4double charge2ratio)
void PrintInfoDefinition(const G4ParticleDefinition &part)
void SetIonisation(G4bool val)
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
G4double GetDEDX(G4double &kineticEnergy, const G4MaterialCutsCouple *)
G4PhysicsTable * SubLambdaTable() const
G4VEnergyLossProcess & operator=(const G4VEnergyLossProcess &right)
G4PhysicsTable * RangeTableForLoss() const
virtual 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
virtual 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
virtual G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
void SetStepFunction(G4double v1, G4double v2)
G4PhysicsTable * CSDARangeTable() const
G4PhysicsTable * theSecondaryRangeTable
virtual G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false)
G4VEmFluctuationModel * fluctModel
G4ParticleChangeForLoss fParticleChange
G4double GetCSDARange(G4double &kineticEnergy, const G4MaterialCutsCouple *)
G4VSubCutProducer * subcutProducer
void SelectModel(G4double kinEnergy)
G4PhysicsTable * theIonisationTable
G4VEmModel * GetModelByIndex(G4int idx=0, G4bool ver=false) const
virtual 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
void PrintWarning(G4String, G4double val)
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)
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
virtual 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
virtual void ProcessDescription(std::ostream &outFile) const
G4PhysicsTable * DEDXTableForSubsec() const
G4double Energy(size_t index) const
G4double Value(G4double theEnergy, size_t &lastidx) const
virtual 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
virtual G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii)
virtual void PreparePhysicsTable(const G4ParticleDefinition &)
const G4ParticleDefinition * theElectron
const G4Material * currentMaterial
void SetLambdaTable(G4PhysicsTable *p)
G4PhysicsTable * theInverseRangeTable
G4SafetyHelper * safetyHelper
std::vector< G4double > theCrossSectionMax
static const G4double factor
virtual void PrintInfo()=0
G4PhysicsTable * InverseRangeTable() const
G4double CrossSectionBiasingFactor() const
virtual void StartTracking(G4Track *)
const G4ParticleDefinition * Particle() const
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:445
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)
G4bool TablesAreBuilt() const
const G4double x[NPOINTSGL]
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 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
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
G4EmParameters * theParameters
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
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