Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 96139 2016-03-17 14:10:31Z 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) override = 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&) override;
169 
170  // build all tables
171  virtual void BuildPhysicsTable(const G4ParticleDefinition&) override;
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*) override;
184 
185  // Step limit from AlongStep
187  const G4Track&,
188  G4double previousStepSize,
189  G4double currentMinimumStep,
190  G4double& currentSafety,
191  G4GPILSelection* selection) override;
192 
193  // Step limit from cross section
195  const G4Track& track,
196  G4double previousStepSize,
197  G4ForceCondition* condition) override;
198 
199  // AlongStep computations
200  virtual G4VParticleChange* AlongStepDoIt(const G4Track&,
201  const G4Step&) override;
202 
203  // Sampling of secondaries in vicinity of geometrical boundary
204  // Return sum of secodaries energy
205  G4double SampleSubCutSecondaries(std::vector<G4Track*>&, const G4Step&,
206  G4VEmModel* model, G4int matIdx);
207 
208  // PostStep sampling of secondaries
209  virtual G4VParticleChange* PostStepDoIt(const G4Track&,
210  const G4Step&) override;
211 
212  // Store all PhysicsTable in files.
213  // Return false in case of any fatal failure at I/O
215  const G4String& directory,
216  G4bool ascii = false) override;
217 
218  // Retrieve all Physics from a files.
219  // Return true if all the Physics Table are built.
220  // Return false if any fatal failure.
222  const G4String& directory,
223  G4bool ascii) override;
224 
225 private:
226  // store a table
227  G4bool StoreTable(const G4ParticleDefinition* p,
228  G4PhysicsTable*, G4bool ascii,
229  const G4String& directory,
230  const G4String& tname);
231 
232  // retrieve a table
233  G4bool RetrieveTable(const G4ParticleDefinition* p,
234  G4PhysicsTable*, G4bool ascii,
235  const G4String& directory,
236  const G4String& tname,
237  G4bool mandatory);
238 
239  //------------------------------------------------------------------------
240  // Public interface to cross section, mfp and sampling of fluctuations
241  // These methods are not used in run time
242  //------------------------------------------------------------------------
243 
244 public:
245 
246  // access to dispersion of restricted energy loss
248  const G4DynamicParticle* dp,
249  G4double length);
250 
251  // Access to cross section table
253  const G4MaterialCutsCouple* couple);
254 
255  // access to cross section
256  G4double MeanFreePath(const G4Track& track);
257 
258  // access to step limit
259  G4double ContinuousStepLimit(const G4Track& track,
260  G4double previousStepSize,
261  G4double currentMinimumStep,
262  G4double& currentSafety);
263 
264 protected:
265 
266  // implementation of the pure virtual method
267  virtual G4double GetMeanFreePath(const G4Track& track,
268  G4double previousStepSize,
269  G4ForceCondition* condition) override;
270 
271  // implementation of the pure virtual method
272  virtual G4double GetContinuousStepLimit(const G4Track& track,
273  G4double previousStepSize,
274  G4double currentMinimumStep,
275  G4double& currentSafety) override;
276 
277  //------------------------------------------------------------------------
278  // Run time method which may be also used by derived processes
279  //------------------------------------------------------------------------
280 
281  // creeation of an empty vector for cross section
283  G4double cut);
284 
285  inline size_t CurrentMaterialCutsCoupleIndex() const;
286 
287  //------------------------------------------------------------------------
288  // Specific methods to set, access, modify models
289  //------------------------------------------------------------------------
290 
291  // Select model in run time
292  inline void SelectModel(G4double kinEnergy);
293 
294 public:
295  // Select model by energy and region index
296  inline G4VEmModel* SelectModelForMaterial(G4double kinEnergy,
297  size_t& idx) const;
298 
299  // Add EM model coupled with fluctuation model for region, smaller value
300  // of order defines which pair of models will be selected for a given
301  // energy interval
302  void AddEmModel(G4int, G4VEmModel*,
303  G4VEmFluctuationModel* fluc = 0,
304  const G4Region* region = nullptr);
305 
306  // Define new energy range for the model identified by the name
307  void UpdateEmModel(const G4String&, G4double, G4double);
308 
309  // Assign a model to a process
310  void SetEmModel(G4VEmModel*, G4int index=1);
311 
312  // return the assigned model
313  G4VEmModel* EmModel(G4int index=1) const;
314 
315  // Access to models
316  G4VEmModel* GetModelByIndex(G4int idx = 0, G4bool ver = false) const;
317 
318  G4int NumberOfModels() const;
319 
320  // Assign a fluctuation model to a process
322 
323  // return the assigned fluctuation model
325 
326  //------------------------------------------------------------------------
327  // Define and access particle type
328  //------------------------------------------------------------------------
329 
330 protected:
331  inline void SetParticle(const G4ParticleDefinition* p);
332  inline void SetSecondaryParticle(const G4ParticleDefinition* p);
333 
334 public:
335  inline void SetBaseParticle(const G4ParticleDefinition* p);
336  inline const G4ParticleDefinition* Particle() const;
337  inline const G4ParticleDefinition* BaseParticle() const;
338  inline const G4ParticleDefinition* SecondaryParticle() const;
339 
340  //------------------------------------------------------------------------
341  // Get/set parameters to configure the process at initialisation time
342  //------------------------------------------------------------------------
343 
344  // Add subcutoff option for the region
345  void ActivateSubCutoff(G4bool val, const G4Region* region = nullptr);
346 
347  // Activate biasing
348  void SetCrossSectionBiasingFactor(G4double f, G4bool flag = true);
349 
350  void ActivateForcedInteraction(G4double length,
351  const G4String& region,
352  G4bool flag = true);
353 
354  void ActivateSecondaryBiasing(const G4String& region, G4double factor,
355  G4double energyLimit);
356 
357  // Add subcutoff process (bremsstrahlung) to sample secondary
358  // particle production in vicinity of the geometry boundary
360 
361  inline void SetLossFluctuations(G4bool val);
362 
363  inline void SetIntegral(G4bool val);
364  inline G4bool IsIntegral() const;
365 
366  // Set/Get flag "isIonisation"
367  void SetIonisation(G4bool val);
368  inline G4bool IsIonisationProcess() const;
369 
370  // Redefine parameteters for stepping control
371  void SetLinearLossLimit(G4double val);
372  void SetStepFunction(G4double v1, G4double v2, G4bool lock=true);
374 
375  inline G4int NumberOfSubCutoffRegions() const;
376 
377  //------------------------------------------------------------------------
378  // Specific methods to path Physics Tables to the process
379  //------------------------------------------------------------------------
380 
381  void SetDEDXTable(G4PhysicsTable* p, G4EmTableType tType);
382  void SetCSDARangeTable(G4PhysicsTable* pRange);
386 
389 
390  // Binning for dEdx, range, inverse range and labda tables
391  void SetDEDXBinning(G4int nbins);
392 
393  // Min kinetic energy for tables
394  void SetMinKinEnergy(G4double e);
395  inline G4double MinKinEnergy() const;
396 
397  // Max kinetic energy for tables
398  void SetMaxKinEnergy(G4double e);
399  inline G4double MaxKinEnergy() const;
400 
401  // Biasing parameters
402  inline G4double CrossSectionBiasingFactor() const;
403 
404  // Return values for given G4MaterialCutsCouple
405  inline G4double GetDEDX(G4double& kineticEnergy,
406  const G4MaterialCutsCouple*);
407  inline G4double GetDEDXForSubsec(G4double& kineticEnergy,
408  const G4MaterialCutsCouple*);
409  inline G4double GetRange(G4double& kineticEnergy,
410  const G4MaterialCutsCouple*);
411  inline G4double GetCSDARange(G4double& kineticEnergy,
412  const G4MaterialCutsCouple*);
413  inline G4double GetRangeForLoss(G4double& kineticEnergy,
414  const G4MaterialCutsCouple*);
415  inline G4double GetKineticEnergy(G4double& range,
416  const G4MaterialCutsCouple*);
417  inline G4double GetLambda(G4double& kineticEnergy,
418  const G4MaterialCutsCouple*);
419 
420  inline G4bool TablesAreBuilt() const;
421 
422  // Access to specific tables
423  inline G4PhysicsTable* DEDXTable() const;
424  inline G4PhysicsTable* DEDXTableForSubsec() const;
425  inline G4PhysicsTable* DEDXunRestrictedTable() const;
426  inline G4PhysicsTable* IonisationTable() const;
428  inline G4PhysicsTable* CSDARangeTable() const;
429  inline G4PhysicsTable* SecondaryRangeTable() const;
430  inline G4PhysicsTable* RangeTableForLoss() const;
431  inline G4PhysicsTable* InverseRangeTable() const;
432  inline G4PhysicsTable* LambdaTable() const;
433  inline G4PhysicsTable* SubLambdaTable() const;
434 
435  //------------------------------------------------------------------------
436  // Run time method for simulation of ionisation
437  //------------------------------------------------------------------------
438 
439  // access atom on which interaction happens
440  const G4Element* GetCurrentElement() const;
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);
464  inline G4double ScaledKinEnergyForLoss(G4double range);
465  inline G4double GetLambdaForScaledEnergy(G4double scaledKinEnergy);
466  inline void ComputeLambdaForScaledEnergy(G4double scaledKinEnergy);
467 
468  // hide assignment operator
470  G4VEnergyLossProcess & operator=(const G4VEnergyLossProcess &right) = delete;
471 
472  // ======== Parameters of the class fixed at construction =========
473 
474  G4LossTableManager* lManager;
475  G4EmModelManager* modelManager;
476  G4EmBiasingManager* biasManager;
477  G4SafetyHelper* safetyHelper;
478  G4EmParameters* theParameters;
479 
480  const G4ParticleDefinition* secondaryParticle;
481  const G4ParticleDefinition* theElectron;
482  const G4ParticleDefinition* thePositron;
483  const G4ParticleDefinition* theGamma;
484  const G4ParticleDefinition* theGenericIon;
485 
486  // ======== Parameters of the class fixed at initialisation =======
487 
488  std::vector<G4VEmModel*> emModels;
489  G4VEmFluctuationModel* fluctModel;
490  G4VAtomDeexcitation* atomDeexcitation;
491  G4VSubCutProducer* subcutProducer;
492  std::vector<const G4Region*> scoffRegions;
493  G4int nSCoffRegions;
494  G4bool* idxSCoffRegions;
495 
496  std::vector<G4VEnergyLossProcess*> scProcesses;
497  G4int nProcesses;
498 
499  // tables and vectors
500  G4PhysicsTable* theDEDXTable;
501  G4PhysicsTable* theDEDXSubTable;
502  G4PhysicsTable* theDEDXunRestrictedTable;
503  G4PhysicsTable* theIonisationTable;
504  G4PhysicsTable* theIonisationSubTable;
505  G4PhysicsTable* theRangeTableForLoss;
506  G4PhysicsTable* theCSDARangeTable;
507  G4PhysicsTable* theSecondaryRangeTable;
508  G4PhysicsTable* theInverseRangeTable;
509  G4PhysicsTable* theLambdaTable;
510  G4PhysicsTable* theSubLambdaTable;
511 
512  size_t idxDEDX;
513  size_t idxDEDXSub;
514  size_t idxDEDXunRestricted;
515  size_t idxIonisation;
516  size_t idxIonisationSub;
517  size_t idxRange;
518  size_t idxCSDA;
519  size_t idxSecRange;
520  size_t idxInverseRange;
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 
532  const G4DataVector* theCuts;
533  const G4DataVector* theSubCuts;
534 
535  const G4ParticleDefinition* baseParticle;
536 
537  G4int nBins;
538  G4int nBinsCSDA;
539 
540  G4double lowestKinEnergy;
541  G4double minKinEnergy;
542  G4double maxKinEnergy;
543  G4double maxKinEnergyCSDA;
544 
545  G4double linLossLimit;
546  G4double dRoverRange;
547  G4double finalRange;
548  G4double lambdaFactor;
549  G4double biasFactor;
550 
551  G4bool lossFluctuationFlag;
552  G4bool rndmStepFlag;
553  G4bool tablesAreBuilt;
554  G4bool integral;
555  G4bool isIon;
556  G4bool isIonisation;
557  G4bool useSubCutoff;
558  G4bool useDeexcitation;
559  G4bool biasFlag;
560  G4bool weightFlag;
561  G4bool isMaster;
562  G4bool actIntegral;
563  G4bool actStepFunc;
564  G4bool actLinLossLimit;
565  G4bool actLossFluc;
566  G4bool actBinning;
567  G4bool actMinKinEnergy;
568  G4bool actMaxKinEnergy;
569 
570 protected:
571 
573 
574  // ======== Cached values - may be state dependent ================
575 
576 private:
577 
578  std::vector<G4DynamicParticle*> secParticles;
579  std::vector<G4Track*> scTracks;
580 
581  const G4ParticleDefinition* particle;
582 
583  G4VEmModel* currentModel;
584  const G4Material* currentMaterial;
585  const G4MaterialCutsCouple* currentCouple;
586  size_t currentCoupleIndex;
587  size_t basedCoupleIndex;
588  size_t lastIdx;
589 
590  G4double massRatio;
591  G4double fFactor;
592  G4double reduceFactor;
593  G4double chargeSqRatio;
594 
595  G4double preStepLambda;
596  G4double fRange;
597  G4double computedRange;
598  G4double preStepKinEnergy;
599  G4double preStepScaledEnergy;
600  G4double preStepRangeEnergy;
601  G4double mfpKinEnergy;
602 
603  G4GPILSelection aGPILSelection;
604 
605  G4int secID;
606  G4int subsecID;
607  G4int biasID;
608 };
609 
610 // ======== Run time inline methods ================
611 
613 {
614  return currentCoupleIndex;
615 }
616 
617 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
618 
620 {
621  currentModel = modelManager->SelectModel(kinEnergy, currentCoupleIndex);
622  currentModel->SetCurrentCouple(currentCouple);
623 }
624 
625 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
626 
628  G4double kinEnergy, size_t& idx) const
629 {
630  return modelManager->SelectModel(kinEnergy, idx);
631 }
632 
633 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
634 
635 inline void
636 G4VEnergyLossProcess::DefineMaterial(const G4MaterialCutsCouple* couple)
637 {
638  if(couple != currentCouple) {
639  currentCouple = couple;
640  currentMaterial = couple->GetMaterial();
641  currentCoupleIndex = couple->GetIndex();
642  basedCoupleIndex = (*theDensityIdx)[currentCoupleIndex];
643  fFactor = chargeSqRatio*biasFactor*(*theDensityFactor)[currentCoupleIndex];
644  reduceFactor = 1.0/(fFactor*massRatio);
645  mfpKinEnergy = DBL_MAX;
646  idxLambda = idxSubLambda = 0;
647  }
648 }
649 
650 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
651 
653  G4double charge2ratio)
654 {
655  massRatio = massratio;
656  fFactor = charge2ratio*biasFactor*(*theDensityFactor)[currentCoupleIndex];
657  chargeSqRatio = charge2ratio;
658  reduceFactor = 1.0/(fFactor*massRatio);
659 }
660 
661 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
662 
663 inline G4double G4VEnergyLossProcess::GetDEDXForScaledEnergy(G4double e)
664 {
665  /*
666  G4cout << "G4VEnergyLossProcess::GetDEDX: Idx= "
667  << basedCoupleIndex << " E(MeV)= " << e
668  << " Emin= " << minKinEnergy << " Factor= " << fFactor
669  << " " << theDEDXTable << G4endl; */
670  G4double x = fFactor*(*theDEDXTable)[basedCoupleIndex]->Value(e, idxDEDX);
671  if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
672  return x;
673 }
674 
675 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
676 
677 inline G4double G4VEnergyLossProcess::GetSubDEDXForScaledEnergy(G4double e)
678 {
679  G4double x =
680  fFactor*(*theDEDXSubTable)[basedCoupleIndex]->Value(e, idxDEDXSub);
681  if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
682  return x;
683 }
684 
685 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
686 
687 inline G4double G4VEnergyLossProcess::GetIonisationForScaledEnergy(G4double e)
688 {
689  G4double x =
690  fFactor*(*theIonisationTable)[basedCoupleIndex]->Value(e, idxIonisation);
691  if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
692  return x;
693 }
694 
695 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
696 
697 inline
698 G4double G4VEnergyLossProcess::GetSubIonisationForScaledEnergy(G4double e)
699 {
700  G4double x = fFactor*
701  (*theIonisationSubTable)[basedCoupleIndex]->Value(e, idxIonisationSub);
702  if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
703  return x;
704 }
705 
706 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
707 
708 inline G4double G4VEnergyLossProcess::GetScaledRangeForScaledEnergy(G4double e)
709 {
710  //G4cout << "G4VEnergyLossProcess::GetScaledRange: Idx= "
711  // << basedCoupleIndex << " E(MeV)= " << e
712  // << " lastIdx= " << lastIdx << " " << theRangeTableForLoss << G4endl;
713  if(basedCoupleIndex != lastIdx || preStepRangeEnergy != e) {
714  lastIdx = basedCoupleIndex;
715  preStepRangeEnergy = e;
716  computedRange =
717  ((*theRangeTableForLoss)[basedCoupleIndex])->Value(e, idxRange);
718  if(e < minKinEnergy) { computedRange *= std::sqrt(e/minKinEnergy); }
719  }
720  //G4cout << "G4VEnergyLossProcess::GetScaledRange: Idx= "
721  // << basedCoupleIndex << " E(MeV)= " << e
722  // << " R= " << fRange << " " << theRangeTableForLoss << G4endl;
723 
724  return computedRange;
725 }
726 
727 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
728 
729 inline G4double
730 G4VEnergyLossProcess::GetLimitScaledRangeForScaledEnergy(G4double e)
731 {
732  G4double x;
733  if (e < maxKinEnergyCSDA) {
734  x = ((*theCSDARangeTable)[basedCoupleIndex])->Value(e, idxCSDA);
735  if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
736  } else {
737  x = theRangeAtMaxEnergy[basedCoupleIndex] +
738  (e - maxKinEnergyCSDA)/theDEDXAtMaxEnergy[basedCoupleIndex];
739  }
740  return x;
741 }
742 
743 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
744 
745 inline G4double G4VEnergyLossProcess::ScaledKinEnergyForLoss(G4double r)
746 {
747  //G4cout << "G4VEnergyLossProcess::GetEnergy: Idx= "
748  // << basedCoupleIndex << " R(mm)= " << r << " "
749  // << theInverseRangeTable << G4endl;
750  G4PhysicsVector* v = (*theInverseRangeTable)[basedCoupleIndex];
751  G4double rmin = v->Energy(0);
752  G4double e = 0.0;
753  if(r >= rmin) { e = v->Value(r, idxInverseRange); }
754  else if(r > 0.0) {
755  G4double x = r/rmin;
756  e = minKinEnergy*x*x;
757  }
758  return e;
759 }
760 
761 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
762 
763 inline G4double G4VEnergyLossProcess::GetLambdaForScaledEnergy(G4double e)
764 {
765  return fFactor*((*theLambdaTable)[basedCoupleIndex])->Value(e, idxLambda);
766 }
767 
768 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
769 
770 inline G4double
772  const G4MaterialCutsCouple* couple)
773 {
774  DefineMaterial(couple);
775  return GetDEDXForScaledEnergy(kineticEnergy*massRatio);
776 }
777 
778 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
779 
780 inline G4double
782  const G4MaterialCutsCouple* couple)
783 {
784  DefineMaterial(couple);
785  return GetSubDEDXForScaledEnergy(kineticEnergy*massRatio);
786 }
787 
788 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
789 
790 inline G4double
792  const G4MaterialCutsCouple* couple)
793 {
794  G4double x = fRange;
795  DefineMaterial(couple);
796  if(theCSDARangeTable) {
797  x = GetLimitScaledRangeForScaledEnergy(kineticEnergy*massRatio)
798  * reduceFactor;
799  } else if(theRangeTableForLoss) {
800  x = GetScaledRangeForScaledEnergy(kineticEnergy*massRatio)*reduceFactor;
801  }
802  return x;
803 }
804 
805 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
806 
807 inline G4double
809  const G4MaterialCutsCouple* couple)
810 {
811  DefineMaterial(couple);
812  G4double x = DBL_MAX;
813  if(theCSDARangeTable) {
814  x=GetLimitScaledRangeForScaledEnergy(kineticEnergy*massRatio)*reduceFactor;
815  }
816  return x;
817 }
818 
819 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
820 
821 inline G4double
823  const G4MaterialCutsCouple* couple)
824 {
825  // G4cout << "GetRangeForLoss: Range from " << GetProcessName() << G4endl;
826  DefineMaterial(couple);
827  G4double x =
828  GetScaledRangeForScaledEnergy(kineticEnergy*massRatio)*reduceFactor;
829  //G4cout << "GetRangeForLoss: Range from " << GetProcessName()
830  // << " e= " << kineticEnergy << " r= " << x << G4endl;
831  return x;
832 }
833 
834 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
835 
836 inline G4double
838  const G4MaterialCutsCouple* couple)
839 {
840  DefineMaterial(couple);
841  return ScaledKinEnergyForLoss(range/reduceFactor)/massRatio;
842 }
843 
844 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
845 
846 inline G4double
848  const G4MaterialCutsCouple* couple)
849 {
850  DefineMaterial(couple);
851  G4double x = 0.0;
852  if(theLambdaTable) {
853  x = GetLambdaForScaledEnergy(kineticEnergy*massRatio);
854  }
855  return x;
856 }
857 
858 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
859 
860 inline void G4VEnergyLossProcess::ComputeLambdaForScaledEnergy(G4double e)
861 {
862  mfpKinEnergy = theEnergyOfCrossSectionMax[currentCoupleIndex];
863  if (e <= mfpKinEnergy) {
864  preStepLambda = GetLambdaForScaledEnergy(e);
865 
866  } else {
867  G4double e1 = e*lambdaFactor;
868  if(e1 > mfpKinEnergy) {
869  preStepLambda = GetLambdaForScaledEnergy(e);
870  G4double preStepLambda1 = GetLambdaForScaledEnergy(e1);
871  if(preStepLambda1 > preStepLambda) {
872  mfpKinEnergy = e1;
873  preStepLambda = preStepLambda1;
874  }
875  } else {
876  preStepLambda = fFactor*theCrossSectionMax[currentCoupleIndex];
877  }
878  }
879 }
880 
881 // ======== Get/Set inline methods used at initialisation ================
882 
884 {
885  fluctModel = p;
886 }
887 
888 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
889 
891 {
892  return fluctModel;
893 }
894 
895 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
896 
898 {
899  particle = p;
900 }
901 
902 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
903 
904 inline void
906 {
907  secondaryParticle = p;
908 }
909 
910 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
911 
912 inline void
914 {
915  baseParticle = p;
916 }
917 
918 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
919 
921 {
922  return particle;
923 }
924 
925 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
926 
928 {
929  return baseParticle;
930 }
931 
932 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
933 
934 inline const G4ParticleDefinition*
936 {
937  return secondaryParticle;
938 }
939 
940 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
941 
943 {
944  lossFluctuationFlag = val;
945  actLossFluc = true;
946 }
947 
948 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
949 
951 {
952  integral = val;
953  actIntegral = true;
954 }
955 
956 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
957 
959 {
960  return integral;
961 }
962 
963 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
964 
966 {
967  return isIonisation;
968 }
969 
970 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
971 
973 {
974  return nSCoffRegions;
975 }
976 
977 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
978 
980 {
981  return minKinEnergy;
982 }
983 
984 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
985 
987 {
988  return maxKinEnergy;
989 }
990 
991 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
992 
994 {
995  return biasFactor;
996 }
997 
998 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
999 
1001 {
1002  return tablesAreBuilt;
1003 }
1004 
1005 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1006 
1008 {
1009  return theDEDXTable;
1010 }
1011 
1012 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1013 
1015 {
1016  return theDEDXSubTable;
1017 }
1018 
1019 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1020 
1022 {
1023  return theDEDXunRestrictedTable;
1024 }
1025 
1026 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1027 
1029 {
1030  return theIonisationTable;
1031 }
1032 
1033 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1034 
1036 {
1037  return theIonisationSubTable;
1038 }
1039 
1040 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1041 
1043 {
1044  return theCSDARangeTable;
1045 }
1046 
1047 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1048 
1050 {
1051  return theSecondaryRangeTable;
1052 }
1053 
1054 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1055 
1057 {
1058  return theRangeTableForLoss;
1059 }
1060 
1061 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1062 
1064 {
1065  return theInverseRangeTable;
1066 }
1067 
1068 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1069 
1071 {
1072  return theLambdaTable;
1073 }
1074 
1075 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1076 
1078 {
1079  return theSubLambdaTable;
1080 }
1081 
1082 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1083 
1084 #endif
G4double GetKineticEnergy(G4double &range, const G4MaterialCutsCouple *)
G4VEnergyLossProcess(const G4String &name="EnergyLoss", G4ProcessType type=fElectromagnetic)
G4double condition(const G4ErrorSymMatrix &m)
const XML_Char * name
Definition: expat.h:151
void SetDEDXTable(G4PhysicsTable *p, G4EmTableType tType)
void SetIntegral(G4bool val)
void SetDynamicMassCharge(G4double massratio, G4double charge2ratio)
void PrintInfoDefinition(const G4ParticleDefinition &part)
void SetIonisation(G4bool val)
virtual void BuildPhysicsTable(const G4ParticleDefinition &) override
G4double GetDEDX(G4double &kineticEnergy, const G4MaterialCutsCouple *)
G4PhysicsTable * SubLambdaTable() const
virtual G4double GetContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety) override
G4PhysicsTable * RangeTableForLoss() const
const char * p
Definition: xmltok.h:285
G4PhysicsTable * SecondaryRangeTable() const
G4double GetRangeForLoss(G4double &kineticEnergy, const G4MaterialCutsCouple *)
void AddCollaborativeProcess(G4VEnergyLossProcess *)
void SetLinearLossLimit(G4double val)
G4PhysicsTable * CSDARangeTable() const
G4ParticleChangeForLoss fParticleChange
G4double GetCSDARange(G4double &kineticEnergy, const G4MaterialCutsCouple *)
virtual G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
void SelectModel(G4double kinEnergy)
G4VEmModel * GetModelByIndex(G4int idx=0, G4bool ver=false) const
virtual G4bool IsApplicable(const G4ParticleDefinition &p) override=0
G4PhysicsTable * IonisationTableForSubsec() const
G4double GetDEDXForSubsec(G4double &kineticEnergy, const G4MaterialCutsCouple *)
virtual G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &) override
void ActivateSecondaryBiasing(const G4String &region, G4double factor, G4double energyLimit)
void SetFluctModel(G4VEmFluctuationModel *)
G4PhysicsVector * LambdaPhysicsVector(const G4MaterialCutsCouple *, G4double cut)
G4PhysicsTable * IonisationTable() const
int G4int
Definition: G4Types.hh:78
G4PhysicsTable * BuildLambdaTable(G4EmTableType tType=fRestricted)
void SetStepFunction(G4double v1, G4double v2, G4bool lock=true)
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=nullptr)
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)
G4PhysicsTable * DEDXTable() const
const G4Element * GetCurrentElement() const
bool G4bool
Definition: G4Types.hh:79
void SetLossFluctuations(G4bool val)
void ActivateForcedInteraction(G4double length, const G4String &region, G4bool flag=true)
G4int NumberOfSubCutoffRegions() const
void SetMaxKinEnergy(G4double e)
G4double GetLambda(G4double &kineticEnergy, const G4MaterialCutsCouple *)
const G4ParticleDefinition * BaseParticle() const
Definition: G4Step.hh:76
virtual void ProcessDescription(std::ostream &outFile) const
G4PhysicsTable * DEDXTableForSubsec() const
G4double Energy(size_t index) const
G4double Value(G4double theEnergy, size_t &lastidx) const
void SetSecondaryParticle(const G4ParticleDefinition *p)
void SetLambdaTable(G4PhysicsTable *p)
virtual void PrintInfo()=0
G4PhysicsTable * InverseRangeTable() const
G4double CrossSectionBiasingFactor() const
virtual void StartTracking(G4Track *) override
const G4ParticleDefinition * Particle() const
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:444
virtual G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii) override
G4double MeanFreePath(const G4Track &track)
void UpdateEmModel(const G4String &, G4double, G4double)
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection) override
virtual void PreparePhysicsTable(const G4ParticleDefinition &) override
void SetEmModel(G4VEmModel *, G4int index=1)
G4double ContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)
G4double MinKinEnergy() const
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
void SetSecondaryRangeTable(G4PhysicsTable *p)
G4PhysicsTable * DEDXunRestrictedTable() const
void ActivateSubCutoff(G4bool val, const G4Region *region=nullptr)
double G4double
Definition: G4Types.hh:76
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &) override
void SetLowestEnergyLimit(G4double)
void SetSubLambdaTable(G4PhysicsTable *p)
virtual G4double MinPrimaryEnergy(const G4ParticleDefinition *, const G4Material *, G4double cut)
G4ForceCondition
void SetRangeTableForLoss(G4PhysicsTable *p)
#define DBL_MAX
Definition: templates.hh:83
virtual G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false) override
const XML_Char XML_Content * model
Definition: expat.h:151
void SetParticle(const G4ParticleDefinition *p)
void SetBaseParticle(const G4ParticleDefinition *p)
G4PhysicsTable * BuildDEDXTable(G4EmTableType tType=fRestricted)
G4VEmModel * EmModel(G4int index=1) const
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
size_t CurrentMaterialCutsCoupleIndex() const
G4double CrossSectionPerVolume(G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4ProcessType
G4bool IsIonisationProcess() const
void SetMinKinEnergy(G4double e)