Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4VEnergyLossProcess.hh
이 파일의 문서화 페이지로 가기
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // $Id: G4VEnergyLossProcess.hh 106714 2017-10-20 09:38:06Z 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  // obsolete to be removed
145  virtual void PrintInfo() {};
146 
147  virtual void ProcessDescription(std::ostream& outFile) const override;
148 
149 protected:
150 
151  virtual void StreamProcessInfo(std::ostream&, G4String) const {};
152 
154  const G4ParticleDefinition*) = 0;
155 
156  //------------------------------------------------------------------------
157  // Methods with standard implementation; may be overwritten if needed
158  //------------------------------------------------------------------------
159 
161  const G4Material*, G4double cut);
162 
163  //------------------------------------------------------------------------
164  // Virtual methods implementation common to all EM ContinuousDiscrete
165  // processes. Further inheritance is not assumed
166  //------------------------------------------------------------------------
167 
168 public:
169 
170  // prepare all tables
171  virtual void PreparePhysicsTable(const G4ParticleDefinition&) override;
172 
173  // build all tables
174  virtual void BuildPhysicsTable(const G4ParticleDefinition&) override;
175 
176  // build a table
178 
179  // build a table
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 
227  // summary printout after initialisation
228  void StreamInfo(std::ostream& out, const G4ParticleDefinition& part,
229  G4String endOfLine=G4String("\n")) const;
230 
231  // store a table
233  G4PhysicsTable*, G4bool ascii,
234  const G4String& directory,
235  const G4String& tname);
236 
237  // retrieve a table
239  G4PhysicsTable*, G4bool ascii,
240  const G4String& directory,
241  const G4String& tname,
242  G4bool mandatory);
243 
244  //------------------------------------------------------------------------
245  // Public interface to cross section, mfp and sampling of fluctuations
246  // These methods are not used in run time
247  //------------------------------------------------------------------------
248 
249 public:
250 
251  // access to dispersion of restricted energy loss
253  const G4DynamicParticle* dp,
254  G4double length);
255 
256  // Access to cross section table
258  const G4MaterialCutsCouple* couple);
259 
260  // access to cross section
262 
263  // access to step limit
265  G4double previousStepSize,
266  G4double currentMinimumStep,
267  G4double& currentSafety);
268 
269 protected:
270 
271  // implementation of the pure virtual method
272  virtual G4double GetMeanFreePath(const G4Track& track,
273  G4double previousStepSize,
274  G4ForceCondition* condition) override;
275 
276  // implementation of the pure virtual method
278  G4double previousStepSize,
279  G4double currentMinimumStep,
280  G4double& currentSafety) override;
281 
282  //------------------------------------------------------------------------
283  // Run time method which may be also used by derived processes
284  //------------------------------------------------------------------------
285 
286  // creeation of an empty vector for cross section
288  G4double cut);
289 
290  inline size_t CurrentMaterialCutsCoupleIndex() const;
291 
292  //------------------------------------------------------------------------
293  // Specific methods to set, access, modify models
294  //------------------------------------------------------------------------
295 
296  // Select model in run time
297  inline void SelectModel(G4double kinEnergy);
298 
299 public:
300  // Select model by energy and region index
301  inline G4VEmModel* SelectModelForMaterial(G4double kinEnergy,
302  size_t& idx) const;
303 
304  // Add EM model coupled with fluctuation model for region, smaller value
305  // of order defines which pair of models will be selected for a given
306  // energy interval
307  void AddEmModel(G4int, G4VEmModel*,
308  G4VEmFluctuationModel* fluc = 0,
309  const G4Region* region = nullptr);
310 
311  // Define new energy range for the model identified by the name
312  void UpdateEmModel(const G4String&, G4double, G4double);
313 
314  // Assign a model to a process local list, to enable the list in run time
315  // the derived process should execute AddEmModel(..) for all such models
316  void SetEmModel(G4VEmModel*, G4int index=0);
317 
318  // return a model from the local list
319  G4VEmModel* EmModel(size_t index=0) const;
320 
321  // Access to models
322  G4VEmModel* GetModelByIndex(G4int idx = 0, G4bool ver = false) const;
323 
324  G4int NumberOfModels() const;
325 
326  // Assign a fluctuation model to a process
327  inline void SetFluctModel(G4VEmFluctuationModel*);
328 
329  // return the assigned fluctuation model
331 
332  //------------------------------------------------------------------------
333  // Define and access particle type
334  //------------------------------------------------------------------------
335 
336 protected:
337  inline void SetParticle(const G4ParticleDefinition* p);
338  inline void SetSecondaryParticle(const G4ParticleDefinition* p);
339 
340 public:
341  inline void SetBaseParticle(const G4ParticleDefinition* p);
342  inline const G4ParticleDefinition* Particle() const;
343  inline const G4ParticleDefinition* BaseParticle() const;
344  inline const G4ParticleDefinition* SecondaryParticle() const;
345 
346  //------------------------------------------------------------------------
347  // Get/set parameters to configure the process at initialisation time
348  //------------------------------------------------------------------------
349 
350  // Add subcutoff option for the region
351  void ActivateSubCutoff(G4bool val, const G4Region* region = nullptr);
352 
353  // Activate biasing
354  void SetCrossSectionBiasingFactor(G4double f, G4bool flag = true);
355 
356  void ActivateForcedInteraction(G4double length,
357  const G4String& region,
358  G4bool flag = true);
359 
360  void ActivateSecondaryBiasing(const G4String& region, G4double factor,
361  G4double energyLimit);
362 
363  // Add subcutoff process (bremsstrahlung) to sample secondary
364  // particle production in vicinity of the geometry boundary
366 
367  inline void SetLossFluctuations(G4bool val);
368 
369  inline void SetIntegral(G4bool val);
370  inline G4bool IsIntegral() const;
371 
372  // Set/Get flag "isIonisation"
373  void SetIonisation(G4bool val);
374  inline G4bool IsIonisationProcess() const;
375 
376  // Redefine parameteters for stepping control
377  void SetLinearLossLimit(G4double val);
378  void SetStepFunction(G4double v1, G4double v2, G4bool lock=true);
380 
381  inline G4int NumberOfSubCutoffRegions() const;
382 
383  //------------------------------------------------------------------------
384  // Specific methods to path Physics Tables to the process
385  //------------------------------------------------------------------------
386 
387  void SetDEDXTable(G4PhysicsTable* p, G4EmTableType tType);
388  void SetCSDARangeTable(G4PhysicsTable* pRange);
392 
395 
396  // Binning for dEdx, range, inverse range and labda tables
397  void SetDEDXBinning(G4int nbins);
398 
399  // Min kinetic energy for tables
400  void SetMinKinEnergy(G4double e);
401  inline G4double MinKinEnergy() const;
402 
403  // Max kinetic energy for tables
404  void SetMaxKinEnergy(G4double e);
405  inline G4double MaxKinEnergy() const;
406 
407  // Biasing parameters
408  inline G4double CrossSectionBiasingFactor() const;
409 
410  // Return values for given G4MaterialCutsCouple
412  const G4MaterialCutsCouple*);
414  const G4MaterialCutsCouple*);
416  const G4MaterialCutsCouple*);
418  const G4MaterialCutsCouple*);
420  const G4MaterialCutsCouple*);
422  const G4MaterialCutsCouple*);
424  const G4MaterialCutsCouple*);
425 
426  inline G4bool TablesAreBuilt() const;
427 
428  // Access to specific tables
429  inline G4PhysicsTable* DEDXTable() const;
430  inline G4PhysicsTable* DEDXTableForSubsec() const;
431  inline G4PhysicsTable* DEDXunRestrictedTable() const;
432  inline G4PhysicsTable* IonisationTable() const;
434  inline G4PhysicsTable* CSDARangeTable() const;
435  inline G4PhysicsTable* SecondaryRangeTable() const;
436  inline G4PhysicsTable* RangeTableForLoss() const;
437  inline G4PhysicsTable* InverseRangeTable() const;
438  inline G4PhysicsTable* LambdaTable() const;
439  inline G4PhysicsTable* SubLambdaTable() const;
440 
441  //------------------------------------------------------------------------
442  // Run time method for simulation of ionisation
443  //------------------------------------------------------------------------
444 
445  // access atom on which interaction happens
446  const G4Element* GetCurrentElement() const;
447 
448  // Set scaling parameters for ions is needed to G4EmCalculator
449  inline void SetDynamicMassCharge(G4double massratio, G4double charge2ratio);
450 
451 private:
452 
454 
455  void PrintWarning(G4String, G4double val);
456 
457  // define material and indexes
458  inline void DefineMaterial(const G4MaterialCutsCouple* couple);
459 
460  //------------------------------------------------------------------------
461  // Compute values using scaling relation, mass and charge of based particle
462  //------------------------------------------------------------------------
463 
464  inline G4double GetDEDXForScaledEnergy(G4double scaledKinEnergy);
465  inline G4double GetSubDEDXForScaledEnergy(G4double scaledKinEnergy);
466  inline G4double GetIonisationForScaledEnergy(G4double scaledKinEnergy);
467  inline G4double GetSubIonisationForScaledEnergy(G4double scaledKinEnergy);
468  inline G4double GetScaledRangeForScaledEnergy(G4double scaledKinEnergy);
469  inline G4double GetLimitScaledRangeForScaledEnergy(G4double scaledKinEnergy);
471  inline G4double GetLambdaForScaledEnergy(G4double scaledKinEnergy);
472  inline void ComputeLambdaForScaledEnergy(G4double scaledKinEnergy);
473 
474  // hide assignment operator
477 
478  // ======== Parameters of the class fixed at construction =========
479 
485 
491 
492  // ======== Parameters of the class fixed at initialisation =======
493 
494  std::vector<G4VEmModel*> emModels;
498  std::vector<const G4Region*> scoffRegions;
501 
502  std::vector<G4VEnergyLossProcess*> scProcesses;
504 
505  // tables and vectors
517 
518  size_t idxDEDX;
519  size_t idxDEDXSub;
523  size_t idxRange;
524  size_t idxCSDA;
525  size_t idxSecRange;
527  size_t idxLambda;
528  size_t idxSubLambda;
529 
530  std::vector<G4double> theDEDXAtMaxEnergy;
531  std::vector<G4double> theRangeAtMaxEnergy;
532  std::vector<G4double> theEnergyOfCrossSectionMax;
533  std::vector<G4double> theCrossSectionMax;
534 
535  const std::vector<G4double>* theDensityFactor;
536  const std::vector<G4int>* theDensityIdx;
537 
540 
542 
545 
550 
556 
575 
576 protected:
577 
579 
580  // ======== Cached values - may be state dependent ================
581 
582 private:
583 
584  std::vector<G4DynamicParticle*> secParticles;
585  std::vector<G4Track*> scTracks;
586 
588 
594  size_t lastIdx;
595 
600 
608 
610 
614 };
615 
616 // ======== Run time inline methods ================
617 
619 {
620  return currentCoupleIndex;
621 }
622 
623 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
624 
626 {
629 }
630 
631 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
632 
634  G4double kinEnergy, size_t& idx) const
635 {
636  return modelManager->SelectModel(kinEnergy, idx);
637 }
638 
639 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
640 
641 inline void
643 {
644  if(couple != currentCouple) {
645  currentCouple = couple;
646  currentMaterial = couple->GetMaterial();
647  currentCoupleIndex = couple->GetIndex();
648  basedCoupleIndex = (*theDensityIdx)[currentCoupleIndex];
649  fFactor = chargeSqRatio*biasFactor*(*theDensityFactor)[currentCoupleIndex];
652  idxLambda = idxSubLambda = 0;
653  }
654 }
655 
656 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
657 
659  G4double charge2ratio)
660 {
661  massRatio = massratio;
662  fFactor = charge2ratio*biasFactor*(*theDensityFactor)[currentCoupleIndex];
663  chargeSqRatio = charge2ratio;
665 }
666 
667 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
668 
670 {
671  /*
672  G4cout << "G4VEnergyLossProcess::GetDEDX: Idx= "
673  << basedCoupleIndex << " E(MeV)= " << e
674  << " Emin= " << minKinEnergy << " Factor= " << fFactor
675  << " " << theDEDXTable << G4endl; */
676  G4double x = fFactor*(*theDEDXTable)[basedCoupleIndex]->Value(e, idxDEDX);
677  if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
678  return x;
679 }
680 
681 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
682 
684 {
685  G4double x =
686  fFactor*(*theDEDXSubTable)[basedCoupleIndex]->Value(e, idxDEDXSub);
687  if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
688  return x;
689 }
690 
691 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
692 
694 {
695  G4double x =
696  fFactor*(*theIonisationTable)[basedCoupleIndex]->Value(e, idxIonisation);
697  if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
698  return x;
699 }
700 
701 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
702 
703 inline
705 {
706  G4double x = fFactor*
707  (*theIonisationSubTable)[basedCoupleIndex]->Value(e, idxIonisationSub);
708  if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
709  return x;
710 }
711 
712 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
713 
715 {
716  //G4cout << "G4VEnergyLossProcess::GetScaledRange: Idx= "
717  // << basedCoupleIndex << " E(MeV)= " << e
718  // << " lastIdx= " << lastIdx << " " << theRangeTableForLoss << G4endl;
722  computedRange =
723  ((*theRangeTableForLoss)[basedCoupleIndex])->Value(e, idxRange);
724  if(e < minKinEnergy) { computedRange *= std::sqrt(e/minKinEnergy); }
725  }
726  //G4cout << "G4VEnergyLossProcess::GetScaledRange: Idx= "
727  // << basedCoupleIndex << " E(MeV)= " << e
728  // << " R= " << fRange << " " << theRangeTableForLoss << G4endl;
729 
730  return computedRange;
731 }
732 
733 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
734 
735 inline G4double
737 {
738  G4double x;
739  if (e < maxKinEnergyCSDA) {
740  x = ((*theCSDARangeTable)[basedCoupleIndex])->Value(e, idxCSDA);
741  if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
742  } else {
745  }
746  return x;
747 }
748 
749 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
750 
752 {
753  //G4cout << "G4VEnergyLossProcess::GetEnergy: Idx= "
754  // << basedCoupleIndex << " R(mm)= " << r << " "
755  // << theInverseRangeTable << G4endl;
756  G4PhysicsVector* v = (*theInverseRangeTable)[basedCoupleIndex];
757  G4double rmin = v->Energy(0);
758  G4double e = 0.0;
759  if(r >= rmin) { e = v->Value(r, idxInverseRange); }
760  else if(r > 0.0) {
761  G4double x = r/rmin;
762  e = minKinEnergy*x*x;
763  }
764  return e;
765 }
766 
767 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
768 
770 {
771  return fFactor*((*theLambdaTable)[basedCoupleIndex])->Value(e, idxLambda);
772 }
773 
774 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
775 
776 inline G4double
778  const G4MaterialCutsCouple* couple)
779 {
780  DefineMaterial(couple);
781  return GetDEDXForScaledEnergy(kineticEnergy*massRatio);
782 }
783 
784 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
785 
786 inline G4double
788  const G4MaterialCutsCouple* couple)
789 {
790  DefineMaterial(couple);
791  return GetSubDEDXForScaledEnergy(kineticEnergy*massRatio);
792 }
793 
794 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
795 
796 inline G4double
798  const G4MaterialCutsCouple* couple)
799 {
800  G4double x = fRange;
801  DefineMaterial(couple);
802  if(theCSDARangeTable) {
804  * reduceFactor;
805  } else if(theRangeTableForLoss) {
807  }
808  return x;
809 }
810 
811 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
812 
813 inline G4double
815  const G4MaterialCutsCouple* couple)
816 {
817  DefineMaterial(couple);
818  return (theCSDARangeTable) ?
820  : DBL_MAX;
821 }
822 
823 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
824 
825 inline G4double
827  const G4MaterialCutsCouple* couple)
828 {
829  // G4cout << "GetRangeForLoss: Range from " << GetProcessName() << G4endl;
830  DefineMaterial(couple);
832 }
833 
834 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
835 
836 inline G4double
838  const G4MaterialCutsCouple* couple)
839 {
840  DefineMaterial(couple);
842 }
843 
844 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
845 
846 inline G4double
848  const G4MaterialCutsCouple* couple)
849 {
850  DefineMaterial(couple);
851  return theLambdaTable ? GetLambdaForScaledEnergy(kineticEnergy*massRatio) : 0.0;
852 }
853 
854 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
855 
857 {
859  if (e <= mfpKinEnergy) {
861 
862  } else {
863  G4double e1 = e*lambdaFactor;
864  if(e1 > mfpKinEnergy) {
866  G4double preStepLambda1 = GetLambdaForScaledEnergy(e1);
867  if(preStepLambda1 > preStepLambda) {
868  mfpKinEnergy = e1;
869  preStepLambda = preStepLambda1;
870  }
871  } else {
873  }
874  }
875 }
876 
877 // ======== Get/Set inline methods used at initialisation ================
878 
880 {
881  fluctModel = p;
882 }
883 
884 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
885 
887 {
888  return fluctModel;
889 }
890 
891 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
892 
894 {
895  particle = p;
896 }
897 
898 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
899 
900 inline void
902 {
904 }
905 
906 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
907 
908 inline void
910 {
911  baseParticle = p;
912 }
913 
914 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
915 
917 {
918  return particle;
919 }
920 
921 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
922 
924 {
925  return baseParticle;
926 }
927 
928 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
929 
930 inline const G4ParticleDefinition*
932 {
933  return secondaryParticle;
934 }
935 
936 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
937 
939 {
940  lossFluctuationFlag = val;
941  actLossFluc = true;
942 }
943 
944 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
945 
947 {
948  integral = val;
949  actIntegral = true;
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
virtual G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii) override
Float_t x
Definition: compare.C:6
void SetDEDXTable(G4PhysicsTable *p, G4EmTableType tType)
const std::vector< G4double > * theDensityFactor
G4double Energy(size_t index) const
G4double CrossSectionBiasingFactor() const
const G4Element * GetCurrentElement() const
const XML_Char * name
Definition: expat.h:151
G4PhysicsTable * theCSDARangeTable
G4VEmFluctuationModel * fluctModel
std::vector< G4double > theCrossSectionMax
G4double GetCSDARange(G4double &kineticEnergy, const G4MaterialCutsCouple *)
G4int NumberOfSubCutoffRegions() const
std::vector< G4DynamicParticle * > secParticles
void SetMinKinEnergy(G4double e)
G4PhysicsTable * DEDXTable() const
const G4ParticleDefinition * theElectron
void AddCollaborativeProcess(G4VEnergyLossProcess *)
void SetStepFunction(G4double v1, G4double v2, G4bool lock=true)
G4PhysicsTable * BuildDEDXTable(G4EmTableType tType=fRestricted)
G4PhysicsTable * theLambdaTable
void SetEmModel(G4VEmModel *, G4int index=0)
G4double GetRange(G4double &kineticEnergy, const G4MaterialCutsCouple *)
std::vector< G4double > theRangeAtMaxEnergy
G4double GetDEDX(G4double &kineticEnergy, const G4MaterialCutsCouple *)
G4PhysicsTable * IonisationTableForSubsec() const
G4PhysicsTable * theIonisationTable
G4double GetLimitScaledRangeForScaledEnergy(G4double scaledKinEnergy)
G4PhysicsTable * theDEDXTable
G4PhysicsTable * theRangeTableForLoss
const char * p
Definition: xmltok.h:285
G4double SampleSubCutSecondaries(std::vector< G4Track * > &, const G4Step &, G4VEmModel *model, G4int matIdx)
std::vector< G4double > theEnergyOfCrossSectionMax
G4PhysicsTable * RangeTableForLoss() const
std::vector< G4double > theDEDXAtMaxEnergy
G4double GetScaledRangeForScaledEnergy(G4double scaledKinEnergy)
const G4ParticleDefinition * thePositron
virtual void StartTracking(G4Track *) override
const G4DataVector * theSubCuts
const G4ParticleDefinition * BaseParticle() const
virtual void PreparePhysicsTable(const G4ParticleDefinition &) override
void ActivateForcedInteraction(G4double length, const G4String &region, G4bool flag=true)
const G4ParticleDefinition * theGenericIon
void SetSubLambdaTable(G4PhysicsTable *p)
G4PhysicsTable * SubLambdaTable() const
G4double Value(G4double theEnergy, size_t &lastidx) const
G4double condition(const G4ErrorSymMatrix &m)
const std::vector< G4int > * theDensityIdx
G4PhysicsTable * theInverseRangeTable
void SetDynamicMassCharge(G4double massratio, G4double charge2ratio)
void SetDEDXBinning(G4int nbins)
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:455
std::vector< G4VEmModel * > emModels
G4EmTableType
G4PhysicsTable * theDEDXSubTable
G4EmParameters * theParameters
std::vector< G4Track * > scTracks
const G4ParticleDefinition * theGamma
const G4MaterialCutsCouple * currentCouple
void SetLowestEnergyLimit(G4double)
void SetRangeTableForLoss(G4PhysicsTable *p)
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
TString part[npart]
Definition: Style.C:32
G4EmBiasingManager * biasManager
size_t CurrentMaterialCutsCoupleIndex() const
G4double GetRangeForLoss(G4double &kineticEnergy, const G4MaterialCutsCouple *)
G4VEmModel * GetModelByIndex(G4int idx=0, G4bool ver=false) const
void DefineMaterial(const G4MaterialCutsCouple *couple)
void SetIonisation(G4bool val)
virtual G4double GetContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety) override
G4PhysicsTable * InverseRangeTable() const
const G4ParticleDefinition const G4Material *G4double range
G4ProcessType
void SetCSDARangeTable(G4PhysicsTable *pRange)
G4double GetSubDEDXForScaledEnergy(G4double scaledKinEnergy)
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=nullptr)
G4double ScaledKinEnergyForLoss(G4double range)
const G4DataVector * theCuts
virtual void InitialiseEnergyLossProcess(const G4ParticleDefinition *, const G4ParticleDefinition *)=0
void SetParticle(const G4ParticleDefinition *p)
virtual G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &) override
G4GPILSelection
G4PhysicsTable * theIonisationSubTable
G4VEmModel * SelectModel(G4double &energy, size_t &index)
G4VEnergyLossProcess & operator=(const G4VEnergyLossProcess &right)=delete
G4PhysicsTable * theDEDXunRestrictedTable
void SetLambdaTable(G4PhysicsTable *p)
G4PhysicsTable * BuildLambdaTable(G4EmTableType tType=fRestricted)
G4PhysicsVector * LambdaPhysicsVector(const G4MaterialCutsCouple *, G4double cut)
G4GPILSelection aGPILSelection
void SetMaxKinEnergy(G4double e)
G4EmModelManager * modelManager
G4double GetKineticEnergy(G4double &range, const G4MaterialCutsCouple *)
G4PhysicsTable * DEDXTableForSubsec() const
G4double GetSubIonisationForScaledEnergy(G4double scaledKinEnergy)
Definition: G4Step.hh:76
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &) override
G4ParticleChangeForLoss fParticleChange
virtual G4bool IsApplicable(const G4ParticleDefinition &p) override=0
const G4Material * currentMaterial
G4LossTableManager * lManager
const G4ParticleDefinition * secondaryParticle
void ComputeLambdaForScaledEnergy(G4double scaledKinEnergy)
const G4ParticleDefinition * SecondaryParticle() const
G4SafetyHelper * safetyHelper
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection) override
std::vector< G4VEnergyLossProcess * > scProcesses
void SetFluctModel(G4VEmFluctuationModel *)
G4VAtomDeexcitation * atomDeexcitation
void SetSecondaryParticle(const G4ParticleDefinition *p)
G4bool StoreTable(const G4ParticleDefinition *p, G4PhysicsTable *, G4bool ascii, const G4String &directory, const G4String &tname)
void SetBaseParticle(const G4ParticleDefinition *p)
G4bool TablesAreBuilt() const
void SetSecondaryRangeTable(G4PhysicsTable *p)
G4PhysicsTable * theSubLambdaTable
void SetIntegral(G4bool val)
void PrintWarning(G4String, G4double val)
G4double GetDEDXForScaledEnergy(G4double scaledKinEnergy)
const G4ParticleDefinition * particle
G4double GetDEDXForSubsec(G4double &kineticEnergy, const G4MaterialCutsCouple *)
G4bool IsIonisationProcess() const
virtual G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false) override
virtual G4double MinPrimaryEnergy(const G4ParticleDefinition *, const G4Material *, G4double cut)
G4double GetLambda(G4double &kineticEnergy, const G4MaterialCutsCouple *)
virtual G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
G4double ContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)
G4PhysicsTable * SecondaryRangeTable() const
double weight
Definition: plottest35.C:25
int G4int
Definition: G4Types.hh:78
G4VSubCutProducer * subcutProducer
void ActivateSecondaryBiasing(const G4String &region, G4double factor, G4double energyLimit)
G4VEmModel * EmModel(size_t index=0) const
void SelectModel(G4double kinEnergy)
G4PhysicsTable * LambdaTable() const
G4bool RetrieveTable(const G4ParticleDefinition *p, G4PhysicsTable *, G4bool ascii, const G4String &directory, const G4String &tname, G4bool mandatory)
void FillSecondariesAlongStep(G4double &eloss, G4double &weight)
G4VEnergyLossProcess(const G4String &name="EnergyLoss", G4ProcessType type=fElectromagnetic)
G4double CrossSectionPerVolume(G4double kineticEnergy, const G4MaterialCutsCouple *couple)
void SetLossFluctuations(G4bool val)
G4ForceCondition
virtual void BuildPhysicsTable(const G4ParticleDefinition &) override
std::vector< const G4Region * > scoffRegions
G4PhysicsTable * IonisationTable() const
G4double MinKinEnergy() const
G4PhysicsTable * CSDARangeTable() const
const G4Material * GetMaterial() const
G4VEmModel * SelectModelForMaterial(G4double kinEnergy, size_t &idx) const
void SetInverseRangeTable(G4PhysicsTable *p)
G4PhysicsTable * DEDXunRestrictedTable() const
void SetCrossSectionBiasingFactor(G4double f, G4bool flag=true)
G4double GetIonisationForScaledEnergy(G4double scaledKinEnergy)
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
G4double MeanFreePath(const G4Track &track)
void StreamInfo(std::ostream &out, const G4ParticleDefinition &part, G4String endOfLine=G4String("\n")) const
virtual void StreamProcessInfo(std::ostream &, G4String) const
#define DBL_MAX
Definition: templates.hh:83
const XML_Char XML_Content * model
Definition: expat.h:151
G4double GetDEDXDispersion(const G4MaterialCutsCouple *couple, const G4DynamicParticle *dp, G4double length)
G4double MaxKinEnergy() const
G4double GetLambdaForScaledEnergy(G4double scaledKinEnergy)
void UpdateEmModel(const G4String &, G4double, G4double)
G4PhysicsTable * theSecondaryRangeTable
void ActivateSubCutoff(G4bool val, const G4Region *region=nullptr)
G4VEmFluctuationModel * FluctModel()
const G4ParticleDefinition * Particle() const
void SetLinearLossLimit(G4double val)
virtual void ProcessDescription(std::ostream &outFile) const override
const G4ParticleDefinition * baseParticle