Geant4-11
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//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class header file
30//
31//
32// File name: G4VEnergyLossProcess
33//
34// Author: Vladimir Ivanchenko on base of Laszlo Urban code
35//
36// Creation date: 03.01.2002
37//
38// Modifications: Vladimir Ivanchenko
39//
40// Class Description:
41//
42// It is the unified energy loss process it calculates the continuous
43// energy loss for charged particles using a set of Energy Loss
44// models valid for different energy regions. There are a possibility
45// to create and access to dE/dx and range tables, or to calculate
46// that information on fly.
47
48// -------------------------------------------------------------------
49//
50
51#ifndef G4VEnergyLossProcess_h
52#define G4VEnergyLossProcess_h 1
53
55#include "globals.hh"
56#include "G4Material.hh"
58#include "G4Track.hh"
59#include "G4EmModelManager.hh"
61#include "G4EmTableType.hh"
63#include "G4PhysicsTable.hh"
64#include "G4PhysicsVector.hh"
65
66class G4Step;
68class G4EmParameters;
69class G4VEmModel;
71class G4DataVector;
72class G4Region;
73class G4SafetyHelper;
78
79//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
80
82{
83public:
84
85 G4VEnergyLossProcess(const G4String& name = "EnergyLoss",
87
88 ~G4VEnergyLossProcess() override;
89
90 //------------------------------------------------------------------------
91 // Virtual methods to be implemented in concrete processes
92 //------------------------------------------------------------------------
93
94protected:
95
96 // description of specific process parameters
97 virtual void StreamProcessInfo(std::ostream&) const {};
98
100 const G4ParticleDefinition*) = 0;
101
102 //------------------------------------------------------------------------
103 // Methods with standard implementation; may be overwritten if needed
104 //------------------------------------------------------------------------
105
107 const G4Material*, G4double cut);
108
109 //------------------------------------------------------------------------
110 // Virtual methods implementation common to all EM ContinuousDiscrete
111 // processes. Further inheritance is not assumed
112 //------------------------------------------------------------------------
113
114public:
115
116 // print documentation in html format
117 void ProcessDescription(std::ostream& outFile) const override;
118
119 // prepare all tables
120 void PreparePhysicsTable(const G4ParticleDefinition&) override;
121
122 // build all tables
123 void BuildPhysicsTable(const G4ParticleDefinition&) override;
124
125 // build a table
127
128 // build a table
130
131 // Called before tracking of each new G4Track
132 void StartTracking(G4Track*) override;
133
134 // Step limit from AlongStep
136 const G4Track&,
137 G4double previousStepSize,
138 G4double currentMinimumStep,
139 G4double& currentSafety,
140 G4GPILSelection* selection) override;
141
142 // Step limit from cross section
144 const G4Track& track,
145 G4double previousStepSize,
146 G4ForceCondition* condition) override;
147
148 // AlongStep computations
149 G4VParticleChange* AlongStepDoIt(const G4Track&, const G4Step&) override;
150
151 // PostStep sampling of secondaries
152 G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step&) override;
153
154 // Store all PhysicsTable in files.
155 // Return false in case of any fatal failure at I/O
157 const G4String& directory,
158 G4bool ascii = false) override;
159
160 // Retrieve all Physics from a files.
161 // Return true if all the Physics Table are built.
162 // Return false if any fatal failure.
164 const G4String& directory,
165 G4bool ascii) override;
166
167private:
168
169 // summary printout after initialisation
170 void StreamInfo(std::ostream& out, const G4ParticleDefinition& part,
171 G4bool rst=false) const;
172
173 // store a table
175 G4PhysicsTable*, G4bool ascii,
176 const G4String& directory,
177 const G4String& tname);
178
179 // retrieve a table
181 G4PhysicsTable*, G4bool ascii,
182 const G4String& directory,
183 const G4String& tname,
184 G4bool mandatory);
185
186 //------------------------------------------------------------------------
187 // Public interface to cross section, mfp and sampling of fluctuations
188 // These methods are not used in run time
189 //------------------------------------------------------------------------
190
191public:
192
193 // access to dispersion of restricted energy loss
195 const G4DynamicParticle* dp,
196 G4double length);
197
198 // Access to cross section table
200 const G4MaterialCutsCouple* couple);
202 const G4MaterialCutsCouple* couple,
203 G4double logKineticEnergy);
204
205 // access to cross section
206 G4double MeanFreePath(const G4Track& track);
207
208 // access to step limit
210 G4double previousStepSize,
211 G4double currentMinimumStep,
212 G4double& currentSafety);
213
214protected:
215
216 // implementation of the pure virtual method
217 G4double GetMeanFreePath(const G4Track& track,
218 G4double previousStepSize,
219 G4ForceCondition* condition) override;
220
221 // implementation of the pure virtual method
223 G4double previousStepSize,
224 G4double currentMinimumStep,
225 G4double& currentSafety) override;
226
227 // creation of an empty vector for cross sections for derived processes
229 G4double cut);
230
231 inline size_t CurrentMaterialCutsCoupleIndex() const;
232
233 //------------------------------------------------------------------------
234 // Specific methods to set, access, modify models
235 //------------------------------------------------------------------------
236
237 // Select model in run time
238 inline void SelectModel(G4double kinEnergy);
239
240public:
241 // Select model by energy and couple index
242 // Not for run time processing
243 inline G4VEmModel* SelectModelForMaterial(G4double kinEnergy,
244 size_t& idxCouple) const;
245
246 // Add EM model coupled with fluctuation model for region, smaller value
247 // of order defines which pair of models will be selected for a given
248 // energy interval
250 G4VEmFluctuationModel* fluc = nullptr,
251 const G4Region* region = nullptr);
252
253 // Assign a model to a process local list, to enable the list in run time
254 // the derived process should execute AddEmModel(..) for all such models
255 void SetEmModel(G4VEmModel*, G4int index=0);
256
257 // Access to models
258 inline size_t NumberOfModels() const;
259
260 // Return a model from the local list
261 inline G4VEmModel* EmModel(size_t index=0) const;
262
263 // Access to models from G4EmModelManager list
264 inline G4VEmModel* GetModelByIndex(size_t idx = 0, G4bool ver = false) const;
265
266 // Assign a fluctuation model to a process
268
269 // Return the assigned fluctuation model
270 inline G4VEmFluctuationModel* FluctModel() const;
271
272 //------------------------------------------------------------------------
273 // Define and access particle type
274 //------------------------------------------------------------------------
275
276protected:
277 inline void SetParticle(const G4ParticleDefinition* p);
278 inline void SetSecondaryParticle(const G4ParticleDefinition* p);
279
280public:
281 inline void SetBaseParticle(const G4ParticleDefinition* p);
282 inline const G4ParticleDefinition* Particle() const;
283 inline const G4ParticleDefinition* BaseParticle() const;
284 inline const G4ParticleDefinition* SecondaryParticle() const;
285
286 // hide assignment operator
289
290 //------------------------------------------------------------------------
291 // Get/set parameters to configure the process at initialisation time
292 //------------------------------------------------------------------------
293
294 // Add subcut processor for the region
295 void ActivateSubCutoff(const G4Region* region);
296
297 // Activate biasing
298 void SetCrossSectionBiasingFactor(G4double f, G4bool flag = true);
299
301 const G4String& region,
302 G4bool flag = true);
303
304 void ActivateSecondaryBiasing(const G4String& region, G4double factor,
305 G4double energyLimit);
306
307 inline void SetLossFluctuations(G4bool val);
308
309 inline void SetSpline(G4bool val);
312
313 // Set/Get flag "isIonisation"
314 void SetIonisation(G4bool val);
315 inline G4bool IsIonisationProcess() const;
316
317 // Redefine parameteters for stepping control
319 void SetStepFunction(G4double v1, G4double v2);
321
322 inline G4int NumberOfSubCutoffRegions() const;
323
324 //------------------------------------------------------------------------
325 // Specific methods to path Physics Tables to the process
326 //------------------------------------------------------------------------
327
329 void SetCSDARangeTable(G4PhysicsTable* pRange);
334 void SetTwoPeaksXS(std::vector<G4TwoPeaksXS*>* p);
335
336 //------------------------------------------------------------------------
337 // Specific methods to define custom Physics Tables to the process
338 //------------------------------------------------------------------------
339
340 // Binning for dEdx, range, inverse range and labda tables
341 void SetDEDXBinning(G4int nbins);
342
343 // Min kinetic energy for tables
345 inline G4double MinKinEnergy() const;
346
347 // Max kinetic energy for tables
349 inline G4double MaxKinEnergy() const;
350
351 // Biasing parameters
352 inline G4double CrossSectionBiasingFactor() const;
353
354 // Return values for given G4MaterialCutsCouple
355 inline G4double GetDEDX(G4double kineticEnergy, const G4MaterialCutsCouple*);
356 inline G4double GetCSDADEDX(G4double kineticEnergy,
357 const G4MaterialCutsCouple*);
358 inline G4double GetDEDX(G4double kineticEnergy, const G4MaterialCutsCouple*,
359 G4double logKineticEnergy);
360 inline G4double GetRange(G4double kineticEnergy, const G4MaterialCutsCouple*);
361 inline G4double GetRange(G4double kineticEnergy, const G4MaterialCutsCouple*,
362 G4double logKineticEnergy);
363 inline G4double GetCSDARange(G4double kineticEnergy,
364 const G4MaterialCutsCouple*);
365 inline G4double GetKineticEnergy(G4double range,
366 const G4MaterialCutsCouple*);
367 inline G4double GetLambda(G4double kineticEnergy,const G4MaterialCutsCouple*);
368 inline G4double GetLambda(G4double kineticEnergy,const G4MaterialCutsCouple*,
369 G4double logKineticEnergy);
370
371 inline G4bool TablesAreBuilt() const;
372
373 // Access to specific tables
374 inline G4PhysicsTable* DEDXTable() const;
376 inline G4PhysicsTable* IonisationTable() const;
377 inline G4PhysicsTable* CSDARangeTable() const;
378 inline G4PhysicsTable* SecondaryRangeTable() const;
379 inline G4PhysicsTable* RangeTableForLoss() const;
380 inline G4PhysicsTable* InverseRangeTable() const;
381 inline G4PhysicsTable* LambdaTable() const;
382 inline std::vector<G4TwoPeaksXS*>* TwoPeaksXS() const;
383
384 inline G4bool UseBaseMaterial() const;
385
386 //------------------------------------------------------------------------
387 // Run time method for simulation of ionisation
388 //------------------------------------------------------------------------
389
390 // access atom on which interaction happens
391 const G4Element* GetCurrentElement() const;
392
393 // Set scaling parameters for ions is needed to G4EmCalculator
394 void SetDynamicMassCharge(G4double massratio, G4double charge2ratio);
395
396private:
397
399
400 void PrintWarning(const G4String&, G4double val) const;
401
402 // define material and indexes
403 inline void DefineMaterial(const G4MaterialCutsCouple* couple);
404
405 //------------------------------------------------------------------------
406 // Compute values using scaling relation, mass and charge of based particle
407 //------------------------------------------------------------------------
408 inline G4double GetDEDXForScaledEnergy(G4double scaledKinE);
409 inline G4double GetDEDXForScaledEnergy(G4double scaledKinE,
410 G4double logScaledKinE);
414 G4double logScaledKinE);
415
418 G4double logScaledKinE);
419
421 inline G4double GetLambdaForScaledEnergy(G4double scaledKinE);
422 inline G4double GetLambdaForScaledEnergy(G4double scaledKinE,
423 G4double logScaledKinE);
424
426 G4double logScaledKinE);
427
429
430protected:
431
433 const G4Material* currentMaterial = nullptr;
435
436private:
437
447
455
465
466 std::vector<const G4Region*>* scoffRegions = nullptr;
467 std::vector<G4VEmModel*>* emModels = nullptr;
468 const std::vector<G4int>* theDensityIdx = nullptr;
469 const std::vector<G4double>* theDensityFactor = nullptr;
470 const G4DataVector* theCuts = nullptr;
471
472 std::vector<G4TwoPeaksXS*>* fXSpeaks = nullptr;
473
478
485
495
496protected:
497
504
506
507private:
508
517
519 size_t coupleIdxRange = 0;
520 size_t coupleIdxLambda = 0;
521 size_t idxDEDX = 0;
523 size_t idxIonisation = 0;
524 size_t idxRange = 0;
525 size_t idxCSDA = 0;
526 size_t idxSecRange = 0;
527 size_t idxInverseRange = 0;
528 size_t idxLambda = 0;
529
532
537 G4bool isIon = false;
549
550 std::vector<G4DynamicParticle*> secParticles;
551 std::vector<G4Track*> scTracks;
552
553};
554
555// ======== Run time inline methods ================
556
558{
559 return currentCoupleIndex;
560}
561
562//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
563
565{
568}
569
570//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
571
573 G4double kinEnergy, size_t& idx) const
574{
575 return modelManager->SelectModel(kinEnergy, idx);
576}
577
578//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
579
580inline void
582{
583 if(couple != currentCouple) {
584 currentCouple = couple;
585 currentMaterial = couple->GetMaterial();
589 idxLambda = 0;
590 if(baseMat) {
591 basedCoupleIndex = (*theDensityIdx)[currentCoupleIndex];
592 fFactor *= (*theDensityFactor)[currentCoupleIndex];
593 }
595 }
596}
597
598//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
599
601{
602 /*
603 G4cout << "G4VEnergyLossProcess::GetDEDX: Idx= "
604 << basedCoupleIndex << " E(MeV)= " << e
605 << " Emin= " << minKinEnergy << " Factor= " << fFactor
606 << " " << theDEDXTable << G4endl; */
607 G4double x = fFactor*(*theDEDXTable)[basedCoupleIndex]->Value(e, idxDEDX);
608 if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
609 return x;
610}
611
612inline
614{
615 /*
616 G4cout << "G4VEnergyLossProcess::GetDEDX: Idx= "
617 << basedCoupleIndex << " E(MeV)= " << e
618 << " Emin= " << minKinEnergy << " Factor= " << fFactor
619 << " " << theDEDXTable << G4endl; */
620 G4double x = fFactor*(*theDEDXTable)[basedCoupleIndex]->LogVectorValue(e,loge);
621 if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
622 return x;
623}
624
625//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
626
628{
629 G4double x =
630 fFactor*(*theIonisationTable)[basedCoupleIndex]->Value(e, idxIonisation);
631 if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
632 return x;
633}
634
635//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
636
638{
639 //G4cout << "G4VEnergyLossProcess::GetScaledRange: Idx= "
640 // << basedCoupleIndex << " E(MeV)= " << e
641 // << " lastIdx= " << lastIdx << " " << theRangeTableForLoss << G4endl;
644 fRangeEnergy = e;
645 fRange = reduceFactor*((*theRangeTableForLoss)[basedCoupleIndex])->Value(e, idxRange);
646 if(e < minKinEnergy) { fRange *= std::sqrt(e/minKinEnergy); }
647 }
648 //G4cout << "G4VEnergyLossProcess::GetScaledRange: Idx= "
649 // << basedCoupleIndex << " E(MeV)= " << e
650 // << " R= " << computedRange << " " << theRangeTableForLoss << G4endl;
651 return fRange;
652}
653
654inline G4double
656{
657 //G4cout << "G4VEnergyLossProcess::GetScaledRange: Idx= "
658 // << basedCoupleIndex << " E(MeV)= " << e
659 // << " lastIdx= " << lastIdx << " " << theRangeTableForLoss << G4endl;
662 fRangeEnergy = e;
663 fRange = reduceFactor*((*theRangeTableForLoss)[basedCoupleIndex])->LogVectorValue(e, loge);
664 if(e < minKinEnergy) { fRange *= std::sqrt(e/minKinEnergy); }
665 }
666 //G4cout << "G4VEnergyLossProcess::GetScaledRange: Idx= "
667 // << basedCoupleIndex << " E(MeV)= " << e
668 // << " R= " << fRange << " " << theRangeTableForLoss << G4endl;
669 return fRange;
670}
671
672//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
673
674inline G4double
676{
677 G4double x = ((*theCSDARangeTable)[basedCoupleIndex])->Value(e, idxCSDA);
678 if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
679 return x;
680}
681
682//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
683
684inline G4double
686 G4double loge)
687{
688 G4double x = ((*theCSDARangeTable)[basedCoupleIndex])->LogVectorValue(e, loge);
689 if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
690 return x;
691}
692
693//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
694
696{
697 //G4cout << "G4VEnergyLossProcess::GetEnergy: Idx= "
698 // << basedCoupleIndex << " R(mm)= " << r << " "
699 // << theInverseRangeTable << G4endl;
700 G4PhysicsVector* v = (*theInverseRangeTable)[basedCoupleIndex];
701 G4double rmin = v->Energy(0);
702 G4double e = 0.0;
703 if(r >= rmin) { e = v->Value(r, idxInverseRange); }
704 else if(r > 0.0) {
705 G4double x = r/rmin;
706 e = minKinEnergy*x*x;
707 }
708 return e;
709}
710
711//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
712
714{
717 fLambdaEnergy = e;
718 fLambda = fFactor*((*theLambdaTable)[basedCoupleIndex])->Value(e, idxLambda);
719 }
720 return fLambda;
721}
722
723//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
724
725inline G4double
727{
730 fLambdaEnergy = e;
732 ((*theLambdaTable)[basedCoupleIndex])->LogVectorValue(e, loge);
733 }
734 return fLambda;
735}
736
737//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
738
739inline G4double
741 const G4MaterialCutsCouple* couple)
742{
743 DefineMaterial(couple);
744 return GetDEDXForScaledEnergy(kinEnergy*massRatio);
745}
746
747//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
748
749inline G4double
751 const G4MaterialCutsCouple* couple,
752 G4double logKinEnergy)
753{
754 DefineMaterial(couple);
755 return GetDEDXForScaledEnergy(kinEnergy*massRatio, logKinEnergy+logMassRatio);
756}
757
758//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
759
760inline G4double
762 const G4MaterialCutsCouple* couple)
763{
764 DefineMaterial(couple);
766}
767
768//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
769
770inline G4double
772 const G4MaterialCutsCouple* couple,
773 G4double logKinEnergy)
774{
775 DefineMaterial(couple);
776 return GetScaledRangeForScaledEnergy(kinEnergy*massRatio, logKinEnergy+logMassRatio);
777}
778
779//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
780
781inline G4double
783 const G4MaterialCutsCouple* couple)
784{
785 DefineMaterial(couple);
786 return (nullptr == theCSDARangeTable) ? DBL_MAX :
788}
789
790//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
791
792inline G4double
794 const G4MaterialCutsCouple* couple)
795{
796 DefineMaterial(couple);
798}
799
800//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
801
802inline G4double
804 const G4MaterialCutsCouple* couple)
805{
806 DefineMaterial(couple);
807 return (nullptr != theLambdaTable) ?
808 GetLambdaForScaledEnergy(kinEnergy*massRatio) : 0.0;
809}
810
811//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
812
813inline G4double
815 const G4MaterialCutsCouple* couple,
816 G4double logKinEnergy)
817{
818 DefineMaterial(couple);
819 return (nullptr != theLambdaTable) ?
820 GetLambdaForScaledEnergy(kinEnergy*massRatio, logKinEnergy+logMassRatio)
821 : 0.0;
822}
823
824// ======== Get/Set inline methods used at initialisation ================
825
827{
828 fluctModel = p;
829}
830
831//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
832
834{
835 return fluctModel;
836}
837
838//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
839
841{
842 particle = p;
843}
844
845//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
846
847inline void
849{
851}
852
853//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
854
855inline void
857{
858 baseParticle = p;
859}
860
861//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
862
864{
865 return particle;
866}
867
868//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
869
871{
872 return baseParticle;
873}
874
875//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
876
877inline const G4ParticleDefinition*
879{
880 return secondaryParticle;
881}
882
883//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
884
886{
888 actLossFluc = true;
889}
890
891//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
892
894{
895 spline = val;
896}
897
898//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
899
901{
902 fXSType = val;
903}
904
905//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
906
908{
909 return fXSType;
910}
911
912//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
913
915{
916 return isIonisation;
917}
918
919//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
920
922{
923 return nSCoffRegions;
924}
925
926//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
927
929{
930 return minKinEnergy;
931}
932
933//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
934
936{
937 return maxKinEnergy;
938}
939
940//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
941
943{
944 return biasFactor;
945}
946
947//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
948
950{
951 return tablesAreBuilt;
952}
953
954//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
955
957{
958 return theDEDXTable;
959}
960
961//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
962
964{
966}
967
968//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
969
971{
972 return theIonisationTable;
973}
974
975//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
976
978{
979 return theCSDARangeTable;
980}
981
982//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
983
985{
987}
988
989//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
990
992{
994}
995
996//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
997
999{
1000 return theInverseRangeTable;
1001}
1002
1003//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1004
1006{
1007 return theLambdaTable;
1008}
1009
1010//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1011
1013{
1014 return baseMat;
1015}
1016
1017//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1018
1019inline std::vector<G4TwoPeaksXS*>* G4VEnergyLossProcess::TwoPeaksXS() const
1020{
1021 return fXSpeaks;
1022}
1023
1024//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1025
1027{
1028 return numberOfModels;
1029}
1030
1031//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1032
1034{
1035 return (nullptr != emModels && index < emModels->size())
1036 ? (*emModels)[index] : nullptr;
1037}
1038
1039//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1040
1041inline G4VEmModel*
1043{
1044 return modelManager->GetModel(idx, ver);
1045}
1046
1047//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1048
1049#endif
G4EmTableType
@ fRestricted
G4CrossSectionType
@ fEmIncreasing
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
G4GPILSelection
G4ProcessType
@ fElectromagnetic
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4VEmModel * GetModel(G4int idx, G4bool ver=false)
G4VEmModel * SelectModel(G4double energy, size_t index)
const G4Material * GetMaterial() const
G4double Energy(const std::size_t index) const
G4double Value(const G4double energy, std::size_t &lastidx) const
Definition: G4Step.hh:62
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:472
const G4ParticleDefinition * BaseParticle() const
G4double GetDEDXForScaledEnergy(G4double scaledKinE)
G4double ScaledKinEnergyForLoss(G4double range)
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=nullptr, const G4Region *region=nullptr)
G4CrossSectionType fXSType
G4PhysicsTable * RangeTableForLoss() const
G4double GetLimitScaledRangeForScaledEnergy(G4double scaledKinE)
const G4ParticleDefinition * theGamma
G4double GetContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety) override
void SetMaxKinEnergy(G4double e)
G4ParticleChangeForLoss fParticleChange
void PreparePhysicsTable(const G4ParticleDefinition &) override
G4PhysicsTable * theDEDXTable
G4double GetCSDADEDX(G4double kineticEnergy, const G4MaterialCutsCouple *)
G4PhysicsTable * InverseRangeTable() const
G4double MeanFreePath(const G4Track &track)
void SetFluctModel(G4VEmFluctuationModel *)
const G4DataVector * theCuts
G4VEnergyLossProcess(G4VEnergyLossProcess &)=delete
G4double GetKineticEnergy(G4double range, const G4MaterialCutsCouple *)
G4PhysicsTable * CSDARangeTable() const
G4PhysicsTable * theIonisationTable
virtual void InitialiseEnergyLossProcess(const G4ParticleDefinition *, const G4ParticleDefinition *)=0
void SelectModel(G4double kinEnergy)
void SetRangeTableForLoss(G4PhysicsTable *p)
const G4ParticleDefinition * baseParticle
G4int NumberOfSubCutoffRegions() const
G4PhysicsTable * theInverseRangeTable
G4PhysicsVector * LambdaPhysicsVector(const G4MaterialCutsCouple *, G4double cut)
void ProcessDescription(std::ostream &outFile) const override
const G4ParticleDefinition * thePositron
G4double GetScaledRangeForScaledEnergy(G4double scaledKinE)
G4double GetIonisationForScaledEnergy(G4double scaledKinE)
G4PhysicsTable * theIonisationSubTable
G4PhysicsTable * BuildDEDXTable(G4EmTableType tType=fRestricted)
G4VEmModel * GetModelByIndex(size_t idx=0, G4bool ver=false) const
const G4MaterialCutsCouple * currentCouple
G4double MaxKinEnergy() const
const G4ParticleDefinition * theGenericIon
virtual G4double MinPrimaryEnergy(const G4ParticleDefinition *, const G4Material *, G4double cut)
void ActivateSecondaryBiasing(const G4String &region, G4double factor, G4double energyLimit)
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
const G4Element * GetCurrentElement() const
void SetDEDXBinning(G4int nbins)
void SetStepFunction(G4double v1, G4double v2)
G4double GetLambda(G4double kineticEnergy, const G4MaterialCutsCouple *)
G4VSubCutProducer * subcutProducer
void SetEmModel(G4VEmModel *, G4int index=0)
void SetCrossSectionBiasingFactor(G4double f, G4bool flag=true)
const G4Material * currentMaterial
G4VEnergyLossProcess & operator=(const G4VEnergyLossProcess &right)=delete
G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
void FillSecondariesAlongStep(G4double weight)
G4EmModelManager * modelManager
G4double ContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)
std::vector< G4TwoPeaksXS * > * TwoPeaksXS() const
G4VEmFluctuationModel * fluctModel
void SetParticle(const G4ParticleDefinition *p)
G4EmParameters * theParameters
G4VEmModel * EmModel(size_t index=0) const
G4PhysicsTable * SecondaryRangeTable() const
const std::vector< G4double > * theDensityFactor
void StreamInfo(std::ostream &out, const G4ParticleDefinition &part, G4bool rst=false) const
G4PhysicsTable * theRangeTableForLoss
void SetLossFluctuations(G4bool val)
void SetCrossSectionType(G4CrossSectionType val)
const G4ParticleDefinition * Particle() const
void SetInverseRangeTable(G4PhysicsTable *p)
G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii) override
G4VEmModel * SelectModelForMaterial(G4double kinEnergy, size_t &idxCouple) const
G4bool StoreTable(const G4ParticleDefinition *p, G4PhysicsTable *, G4bool ascii, const G4String &directory, const G4String &tname)
void ActivateForcedInteraction(G4double length, const G4String &region, G4bool flag=true)
G4CrossSectionType CrossSectionType() const
void SetDynamicMassCharge(G4double massratio, G4double charge2ratio)
G4double MinKinEnergy() const
G4GPILSelection aGPILSelection
void SetBaseParticle(const G4ParticleDefinition *p)
G4double CrossSectionBiasingFactor() const
std::vector< const G4Region * > * scoffRegions
std::vector< G4TwoPeaksXS * > * fXSpeaks
G4bool IsIonisationProcess() const
G4PhysicsTable * theLambdaTable
G4double CrossSectionPerVolume(G4double kineticEnergy, const G4MaterialCutsCouple *couple)
const std::vector< G4int > * theDensityIdx
void SetDEDXTable(G4PhysicsTable *p, G4EmTableType tType)
void SetTwoPeaksXS(std::vector< G4TwoPeaksXS * > *p)
void SetSecondaryRangeTable(G4PhysicsTable *p)
G4PhysicsTable * BuildLambdaTable(G4EmTableType tType=fRestricted)
G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false) override
G4double GetCSDARange(G4double kineticEnergy, const G4MaterialCutsCouple *)
void SetIonisation(G4bool val)
G4bool IsRegionForCubcutProcessor(const G4Track &aTrack)
void DefineMaterial(const G4MaterialCutsCouple *couple)
std::vector< G4Track * > scTracks
G4double GetDEDX(G4double kineticEnergy, const G4MaterialCutsCouple *)
G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection) override
G4VEmFluctuationModel * FluctModel() const
void SetLinearLossLimit(G4double val)
void SetLowestEnergyLimit(G4double)
G4bool RetrieveTable(const G4ParticleDefinition *p, G4PhysicsTable *, G4bool ascii, const G4String &directory, const G4String &tname, G4bool mandatory)
G4EmBiasingManager * biasManager
void SetLambdaTable(G4PhysicsTable *p)
const G4ParticleDefinition * secondaryParticle
void BuildPhysicsTable(const G4ParticleDefinition &) override
G4VEnergyLossProcess(const G4String &name="EnergyLoss", G4ProcessType type=fElectromagnetic)
G4PhysicsTable * theSecondaryRangeTable
G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &) override
G4VAtomDeexcitation * atomDeexcitation
G4PhysicsTable * theDEDXunRestrictedTable
G4double GetRange(G4double kineticEnergy, const G4MaterialCutsCouple *)
std::vector< G4VEmModel * > * emModels
G4double GetDEDXDispersion(const G4MaterialCutsCouple *couple, const G4DynamicParticle *dp, G4double length)
size_t CurrentMaterialCutsCoupleIndex() const
G4PhysicsTable * IonisationTable() const
void ActivateSubCutoff(const G4Region *region)
void PrintWarning(const G4String &, G4double val) const
void SetSecondaryParticle(const G4ParticleDefinition *p)
G4PhysicsTable * LambdaTable() const
void SetCSDARangeTable(G4PhysicsTable *pRange)
const G4ParticleDefinition * particle
void ComputeLambdaForScaledEnergy(G4double scaledKinE, G4double logScaledKinE)
G4PhysicsTable * DEDXunRestrictedTable() const
const G4ParticleDefinition * theElectron
G4double GetLambdaForScaledEnergy(G4double scaledKinE)
G4PhysicsTable * theCSDARangeTable
G4LossTableManager * lManager
virtual void StreamProcessInfo(std::ostream &) const
G4SafetyHelper * safetyHelper
G4PhysicsTable * DEDXTable() const
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &) override
void StartTracking(G4Track *) override
void SetMinKinEnergy(G4double e)
std::vector< G4DynamicParticle * > secParticles
const G4ParticleDefinition * SecondaryParticle() const
const char * name(G4int ptype)
#define LOG_EKIN_MIN
Definition: templates.hh:98
#define DBL_MAX
Definition: templates.hh:62