Geant4  10.01.p03
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 85011 2014-10-23 09:41:59Z 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 protected:
147 
149  const G4ParticleDefinition*) = 0;
150 
151  //------------------------------------------------------------------------
152  // Methods with standard implementation; may be overwritten if needed
153  //------------------------------------------------------------------------
154 
156  const G4Material*, G4double cut);
157 
158  //------------------------------------------------------------------------
159  // Virtual methods implementation common to all EM ContinuousDiscrete
160  // processes. Further inheritance is not assumed
161  //------------------------------------------------------------------------
162 
163 public:
164 
165  // prepare all tables
167 
168  // build all tables
170 
171  // build a table
173 
174  // build a table
176 
177  // summary printout after initialisation
178  void PrintInfoDefinition(const G4ParticleDefinition& part);
179 
180  // Called before tracking of each new G4Track
181  void StartTracking(G4Track*);
182 
183  // Step limit from AlongStep
185  G4double previousStepSize,
186  G4double currentMinimumStep,
187  G4double& currentSafety,
188  G4GPILSelection* selection);
189 
190  // Step limit from cross section
192  G4double previousStepSize,
194 
195  // AlongStep computations
197 
198  // Sampling of secondaries in vicinity of geometrical boundary
199  // Return sum of secodaries energy
200  G4double SampleSubCutSecondaries(std::vector<G4Track*>&, const G4Step&,
201  G4VEmModel* model, G4int matIdx);
202 
203  // PostStep sampling of secondaries
204  G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step&);
205 
206  // Store all PhysicsTable in files.
207  // Return false in case of any fatal failure at I/O
209  const G4String& directory,
210  G4bool ascii = false);
211 
212  // Retrieve all Physics from a files.
213  // Return true if all the Physics Table are built.
214  // Return false if any fatal failure.
216  const G4String& directory,
217  G4bool ascii);
218 
219 private:
220  // store a table
222  G4PhysicsTable*, G4bool ascii,
223  const G4String& directory,
224  const G4String& tname);
225 
226  // retrieve a table
228  G4PhysicsTable*, G4bool ascii,
229  const G4String& directory,
230  const G4String& tname,
231  G4bool mandatory);
232 
233  //------------------------------------------------------------------------
234  // Public interface to cross section, mfp and sampling of fluctuations
235  // These methods are not used in run time
236  //------------------------------------------------------------------------
237 
238 public:
239 
240  // access to dispersion of restricted energy loss
242  const G4DynamicParticle* dp,
243  G4double length);
244 
245  // Access to cross section table
247  const G4MaterialCutsCouple* couple);
248 
249  // access to cross section
250  G4double MeanFreePath(const G4Track& track);
251 
252  // access to step limit
253  G4double ContinuousStepLimit(const G4Track& track,
254  G4double previousStepSize,
255  G4double currentMinimumStep,
256  G4double& currentSafety);
257 
258 protected:
259 
260  // implementation of the pure virtual method
261  G4double GetMeanFreePath(const G4Track& track,
262  G4double previousStepSize,
263  G4ForceCondition* condition);
264 
265  // implementation of the pure virtual method
267  G4double previousStepSize,
268  G4double currentMinimumStep,
269  G4double& currentSafety);
270 
271  //------------------------------------------------------------------------
272  // Run time method which may be also used by derived processes
273  //------------------------------------------------------------------------
274 
275  // creeation of an empty vector for cross section
277  G4double cut);
278 
279  inline size_t CurrentMaterialCutsCoupleIndex() const;
280 
281  //------------------------------------------------------------------------
282  // Specific methods to set, access, modify models
283  //------------------------------------------------------------------------
284 
285  // Select model in run time
286  inline void SelectModel(G4double kinEnergy);
287 
288 public:
289  // Select model by energy and region index
290  inline G4VEmModel* SelectModelForMaterial(G4double kinEnergy,
291  size_t& idx) const;
292 
293  // Add EM model coupled with fluctuation model for region, smaller value
294  // of order defines which pair of models will be selected for a given
295  // energy interval
296  void AddEmModel(G4int, G4VEmModel*,
297  G4VEmFluctuationModel* fluc = 0,
298  const G4Region* region = 0);
299 
300  // Define new energy range for the model identified by the name
301  void UpdateEmModel(const G4String&, G4double, G4double);
302 
303  // Assign a model to a process
304  void SetEmModel(G4VEmModel*, G4int index=1);
305 
306  // return the assigned model
307  G4VEmModel* EmModel(G4int index=1) const;
308 
309  // Access to models
310  G4VEmModel* GetModelByIndex(G4int idx = 0, G4bool ver = false) const;
311 
312  G4int NumberOfModels() const;
313 
314  // Assign a fluctuation model to a process
316 
317  // return the assigned fluctuation model
319 
320  //------------------------------------------------------------------------
321  // Define and access particle type
322  //------------------------------------------------------------------------
323 
324 protected:
325  inline void SetParticle(const G4ParticleDefinition* p);
326  inline void SetSecondaryParticle(const G4ParticleDefinition* p);
327 
328 public:
329  inline void SetBaseParticle(const G4ParticleDefinition* p);
330  inline const G4ParticleDefinition* Particle() const;
331  inline const G4ParticleDefinition* BaseParticle() const;
332  inline const G4ParticleDefinition* SecondaryParticle() const;
333 
334  //------------------------------------------------------------------------
335  // Get/set parameters to configure the process at initialisation time
336  //------------------------------------------------------------------------
337 
338  // Add subcutoff option for the region
339  void ActivateSubCutoff(G4bool val, const G4Region* region = 0);
340 
341  // Activate biasing
342  void SetCrossSectionBiasingFactor(G4double f, G4bool flag = true);
343 
344  void ActivateForcedInteraction(G4double length = 0.0,
345  const G4String& region = "",
346  G4bool flag = true);
347 
348  void ActivateSecondaryBiasing(const G4String& region, G4double factor,
349  G4double energyLimit);
350 
351  // Add subcutoff process (bremsstrahlung) to sample secondary
352  // particle production in vicinity of the geometry boundary
354 
355  inline void SetLossFluctuations(G4bool val);
356  //inline void SetRandomStep(G4bool val);
357 
358  inline void SetIntegral(G4bool val);
359  inline G4bool IsIntegral() const;
360 
361  // Set/Get flag "isIonisation"
362  void SetIonisation(G4bool val);
363  inline G4bool IsIonisationProcess() const;
364 
365  // Redefine parameteters for stepping control
366  void SetLinearLossLimit(G4double val);
367  // void SetMinSubRange(G4double val);
368  // void SetLambdaFactor(G4double val);
369  void SetStepFunction(G4double v1, G4double v2);
371 
372  inline G4int NumberOfSubCutoffRegions() const;
373 
374  //------------------------------------------------------------------------
375  // Specific methods to path Physics Tables to the process
376  //------------------------------------------------------------------------
377 
378  void SetDEDXTable(G4PhysicsTable* p, G4EmTableType tType);
379  void SetCSDARangeTable(G4PhysicsTable* pRange);
383 
386 
387  // Binning for dEdx, range, inverse range and labda tables
388  void SetDEDXBinning(G4int nbins);
389 
390  // Min kinetic energy for tables
391  void SetMinKinEnergy(G4double e);
392  inline G4double MinKinEnergy() const;
393 
394  // Max kinetic energy for tables
395  void SetMaxKinEnergy(G4double e);
396  inline G4double MaxKinEnergy() const;
397 
398  // Biasing parameters
399  inline G4double CrossSectionBiasingFactor() const;
400 
401  // Return values for given G4MaterialCutsCouple
402  inline G4double GetDEDX(G4double& kineticEnergy,
403  const G4MaterialCutsCouple*);
404  inline G4double GetDEDXForSubsec(G4double& kineticEnergy,
405  const G4MaterialCutsCouple*);
406  inline G4double GetRange(G4double& kineticEnergy,
407  const G4MaterialCutsCouple*);
408  inline G4double GetCSDARange(G4double& kineticEnergy,
409  const G4MaterialCutsCouple*);
410  inline G4double GetRangeForLoss(G4double& kineticEnergy,
411  const G4MaterialCutsCouple*);
412  inline G4double GetKineticEnergy(G4double& range,
413  const G4MaterialCutsCouple*);
414  inline G4double GetLambda(G4double& kineticEnergy,
415  const G4MaterialCutsCouple*);
416 
417  inline G4bool TablesAreBuilt() const;
418 
419  // Access to specific tables
420  inline G4PhysicsTable* DEDXTable() const;
421  inline G4PhysicsTable* DEDXTableForSubsec() const;
422  inline G4PhysicsTable* DEDXunRestrictedTable() const;
423  inline G4PhysicsTable* IonisationTable() const;
425  inline G4PhysicsTable* CSDARangeTable() const;
426  inline G4PhysicsTable* SecondaryRangeTable() const;
427  inline G4PhysicsTable* RangeTableForLoss() const;
428  inline G4PhysicsTable* InverseRangeTable() const;
429  inline G4PhysicsTable* LambdaTable() const;
430  inline G4PhysicsTable* SubLambdaTable() const;
431 
432  //------------------------------------------------------------------------
433  // Run time method for simulation of ionisation
434  //------------------------------------------------------------------------
435 
436  // access atom on which interaction happens
437  const G4Element* GetCurrentElement() const;
438 
439  // sample range at the end of a step
440  // inline G4double SampleRange();
441 
442  // Set scaling parameters for ions is needed to G4EmCalculator
443  inline void SetDynamicMassCharge(G4double massratio, G4double charge2ratio);
444 
445 private:
446 
447  void FillSecondariesAlongStep(G4double& eloss, G4double& weight);
448 
449  void PrintWarning(G4String, G4double val);
450 
451  // define material and indexes
452  inline void DefineMaterial(const G4MaterialCutsCouple* couple);
453 
454  //------------------------------------------------------------------------
455  // Compute values using scaling relation, mass and charge of based particle
456  //------------------------------------------------------------------------
457 
458  inline G4double GetDEDXForScaledEnergy(G4double scaledKinEnergy);
459  inline G4double GetSubDEDXForScaledEnergy(G4double scaledKinEnergy);
460  inline G4double GetIonisationForScaledEnergy(G4double scaledKinEnergy);
461  inline G4double GetSubIonisationForScaledEnergy(G4double scaledKinEnergy);
462  inline G4double GetScaledRangeForScaledEnergy(G4double scaledKinEnergy);
463  inline G4double GetLimitScaledRangeForScaledEnergy(G4double scaledKinEnergy);
465  inline G4double GetLambdaForScaledEnergy(G4double scaledKinEnergy);
466  inline void ComputeLambdaForScaledEnergy(G4double scaledKinEnergy);
467 
468  // hide assignment operator
471 
472  // ======== Parameters of the class fixed at construction =========
473 
479 
485 
486  // ======== Parameters of the class fixed at initialisation =======
487 
488  std::vector<G4VEmModel*> emModels;
492  std::vector<const G4Region*> scoffRegions;
495 
496  std::vector<G4VEnergyLossProcess*> scProcesses;
498 
499  // tables and vectors
511 
512  size_t idxDEDX;
513  size_t idxDEDXSub;
517  size_t idxRange;
518  size_t idxCSDA;
519  size_t idxSecRange;
521  size_t idxLambda;
522  size_t idxSubLambda;
523 
524  std::vector<G4double> theDEDXAtMaxEnergy;
525  std::vector<G4double> theRangeAtMaxEnergy;
526  std::vector<G4double> theEnergyOfCrossSectionMax;
527  std::vector<G4double> theCrossSectionMax;
528 
529  const std::vector<G4double>* theDensityFactor;
530  const std::vector<G4int>* theDensityIdx;
531 
534 
536 
539 
544 
546  // G4double minSubRange;
551 
568 
569 protected:
570 
572 
573  // ======== Cached values - may be state dependent ================
574 
575 private:
576 
577  std::vector<G4DynamicParticle*> secParticles;
578  std::vector<G4Track*> scTracks;
579 
581 
587  size_t lastIdx;
588 
590 
595 
603 
605 
609 };
610 
611 // ======== Run time inline methods ================
612 
614 {
615  return currentCoupleIndex;
616 }
617 
618 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
619 
621 {
624 }
625 
626 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
627 
629  G4double kinEnergy, size_t& idx) const
630 {
631  return modelManager->SelectModel(kinEnergy, idx);
632 }
633 
634 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
635 
636 inline void
638 {
639  if(couple != currentCouple) {
640  currentCouple = couple;
641  currentMaterial = couple->GetMaterial();
642  currentCoupleIndex = couple->GetIndex();
643  basedCoupleIndex = (*theDensityIdx)[currentCoupleIndex];
644  fFactor = chargeSqRatio*biasFactor*(*theDensityFactor)[currentCoupleIndex];
647  idxLambda = idxSubLambda = 0;
648  }
649 }
650 
651 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
652 
654  G4double charge2ratio)
655 {
656  massRatio = massratio;
657  fFactor = charge2ratio*biasFactor*(*theDensityFactor)[currentCoupleIndex];
658  chargeSqRatio = charge2ratio;
660 }
661 
662 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
663 
665 {
666  /*
667  G4cout << "G4VEnergyLossProcess::GetDEDX: Idx= "
668  << basedCoupleIndex << " E(MeV)= " << e
669  << " Emin= " << minKinEnergy << " Factor= " << fFactor
670  << " " << theDEDXTable << G4endl; */
671  G4double x = fFactor*(*theDEDXTable)[basedCoupleIndex]->Value(e, idxDEDX);
672  if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
673  return x;
674 }
675 
676 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
677 
679 {
680  G4double x =
681  fFactor*(*theDEDXSubTable)[basedCoupleIndex]->Value(e, idxDEDXSub);
682  if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
683  return x;
684 }
685 
686 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
687 
689 {
690  G4double x =
691  fFactor*(*theIonisationTable)[basedCoupleIndex]->Value(e, idxIonisation);
692  if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
693  return x;
694 }
695 
696 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
697 
698 inline
700 {
701  G4double x = fFactor*
702  (*theIonisationSubTable)[basedCoupleIndex]->Value(e, idxIonisationSub);
703  if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
704  return x;
705 }
706 
707 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
708 
710 {
711  //G4cout << "G4VEnergyLossProcess::GetScaledRange: Idx= "
712  // << basedCoupleIndex << " E(MeV)= " << e
713  // << " lastIdx= " << lastIdx << " " << theRangeTableForLoss << G4endl;
716  preStepRangeEnergy = e;
717  computedRange =
718  ((*theRangeTableForLoss)[basedCoupleIndex])->Value(e, idxRange);
719  if(e < minKinEnergy) { computedRange *= std::sqrt(e/minKinEnergy); }
720  }
721  //G4cout << "G4VEnergyLossProcess::GetScaledRange: Idx= "
722  // << basedCoupleIndex << " E(MeV)= " << e
723  // << " R= " << fRange << " " << theRangeTableForLoss << G4endl;
724 
725  return computedRange;
726 }
727 
728 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
729 
730 inline G4double
732 {
733  G4double x;
734  if (e < maxKinEnergyCSDA) {
735  x = ((*theCSDARangeTable)[basedCoupleIndex])->Value(e, idxCSDA);
736  if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
737  } else {
740  }
741  return x;
742 }
743 
744 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
745 
747 {
748  //G4cout << "G4VEnergyLossProcess::GetEnergy: Idx= "
749  // << basedCoupleIndex << " R(mm)= " << r << " "
750  // << theInverseRangeTable << G4endl;
751  G4PhysicsVector* v = (*theInverseRangeTable)[basedCoupleIndex];
752  G4double rmin = v->Energy(0);
753  G4double e = 0.0;
754  if(r >= rmin) { e = v->Value(r, idxInverseRange); }
755  else if(r > 0.0) {
756  G4double x = r/rmin;
757  e = minKinEnergy*x*x;
758  }
759  return e;
760 }
761 
762 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
763 
765 {
766  return fFactor*((*theLambdaTable)[basedCoupleIndex])->Value(e, idxLambda);
767 }
768 
769 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
770 
771 inline G4double
773  const G4MaterialCutsCouple* couple)
774 {
775  DefineMaterial(couple);
776  return GetDEDXForScaledEnergy(kineticEnergy*massRatio);
777 }
778 
779 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
780 
781 inline G4double
783  const G4MaterialCutsCouple* couple)
784 {
785  DefineMaterial(couple);
786  return GetSubDEDXForScaledEnergy(kineticEnergy*massRatio);
787 }
788 
789 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
790 
791 inline G4double
793  const G4MaterialCutsCouple* couple)
794 {
795  G4double x = fRange;
796  DefineMaterial(couple);
797  if(theCSDARangeTable) {
799  * reduceFactor;
800  } else if(theRangeTableForLoss) {
802  }
803  return x;
804 }
805 
806 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
807 
808 inline G4double
810  const G4MaterialCutsCouple* couple)
811 {
812  DefineMaterial(couple);
813  G4double x = DBL_MAX;
814  if(theCSDARangeTable) {
816  }
817  return x;
818 }
819 
820 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
821 
822 inline G4double
824  const G4MaterialCutsCouple* couple)
825 {
826  // G4cout << "GetRangeForLoss: Range from " << GetProcessName() << G4endl;
827  DefineMaterial(couple);
828  G4double x =
830  //G4cout << "GetRangeForLoss: Range from " << GetProcessName()
831  // << " e= " << kineticEnergy << " r= " << x << G4endl;
832  return x;
833 }
834 
835 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
836 
837 inline G4double
839  const G4MaterialCutsCouple* couple)
840 {
841  DefineMaterial(couple);
843 }
844 
845 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
846 
847 inline G4double
849  const G4MaterialCutsCouple* couple)
850 {
851  DefineMaterial(couple);
852  G4double x = 0.0;
853  if(theLambdaTable) {
854  x = GetLambdaForScaledEnergy(kineticEnergy*massRatio);
855  }
856  return x;
857 }
858 
859 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
860 
862 {
864  if (e <= mfpKinEnergy) {
866 
867  } else {
869  if(e1 > mfpKinEnergy) {
871  G4double preStepLambda1 = GetLambdaForScaledEnergy(e1);
872  if(preStepLambda1 > preStepLambda) {
873  mfpKinEnergy = e1;
874  preStepLambda = preStepLambda1;
875  }
876  } else {
878  }
879  }
880 }
881 
882 // ======== Get/Set inline methods used at initialisation ================
883 
885 {
886  fluctModel = p;
887 }
888 
889 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
890 
892 {
893  return fluctModel;
894 }
895 
896 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
897 
899 {
900  particle = p;
901 }
902 
903 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
904 
905 inline void
907 {
908  secondaryParticle = p;
909 }
910 
911 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
912 
913 inline void
915 {
916  baseParticle = p;
917 }
918 
919 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
920 
922 {
923  return particle;
924 }
925 
926 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
927 
929 {
930  return baseParticle;
931 }
932 
933 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
934 
935 inline const G4ParticleDefinition*
937 {
938  return secondaryParticle;
939 }
940 
941 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
942 
944 {
945  lossFluctuationFlag = val;
946 }
947 
948 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
949 
951 {
952  integral = val;
953 }
954 
955 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
956 
958 {
959  return integral;
960 }
961 
962 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
963 
965 {
966  return isIonisation;
967 }
968 
969 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
970 
972 {
973  return nSCoffRegions;
974 }
975 
976 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
977 
979 {
980  return minKinEnergy;
981 }
982 
983 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
984 
986 {
987  return maxKinEnergy;
988 }
989 
990 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
991 
993 {
994  return biasFactor;
995 }
996 
997 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
998 
1000 {
1001  return tablesAreBuilt;
1002 }
1003 
1004 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1005 
1007 {
1008  return theDEDXTable;
1009 }
1010 
1011 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1012 
1014 {
1015  return theDEDXSubTable;
1016 }
1017 
1018 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1019 
1021 {
1022  return theDEDXunRestrictedTable;
1023 }
1024 
1025 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1026 
1028 {
1029  return theIonisationTable;
1030 }
1031 
1032 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1033 
1035 {
1036  return theIonisationSubTable;
1037 }
1038 
1039 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1040 
1042 {
1043  return theCSDARangeTable;
1044 }
1045 
1046 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1047 
1049 {
1050  return theSecondaryRangeTable;
1051 }
1052 
1053 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1054 
1056 {
1057  return theRangeTableForLoss;
1058 }
1059 
1060 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1061 
1063 {
1064  return theInverseRangeTable;
1065 }
1066 
1067 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1068 
1070 {
1071  return theLambdaTable;
1072 }
1073 
1074 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1075 
1077 {
1078  return theSubLambdaTable;
1079 }
1080 
1081 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1082 
1083 #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)
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 *)
G4VSubCutProducer * subcutProducer
void SelectModel(G4double kinEnergy)
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
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
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
static const G4double factor
virtual void PrintInfo()=0
G4PhysicsTable * InverseRangeTable() const
G4double CrossSectionBiasingFactor() const
const G4ParticleDefinition * Particle() const
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:426
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
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