Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VEnergyLossProcess Class Referenceabstract

#include <G4VEnergyLossProcess.hh>

+ Inheritance diagram for G4VEnergyLossProcess:

Public Member Functions

 G4VEnergyLossProcess (const G4String &name="EnergyLoss", G4ProcessType type=fElectromagnetic)
 
virtual ~G4VEnergyLossProcess ()
 
virtual G4bool IsApplicable (const G4ParticleDefinition &p) override=0
 
virtual void PrintInfo ()
 
virtual void ProcessDescription (std::ostream &outFile) const override
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &) override
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &) override
 
G4PhysicsTableBuildDEDXTable (G4EmTableType tType=fRestricted)
 
G4PhysicsTableBuildLambdaTable (G4EmTableType tType=fRestricted)
 
virtual void StartTracking (G4Track *) override
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection) override
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &) override
 
G4double SampleSubCutSecondaries (std::vector< G4Track * > &, const G4Step &, G4VEmModel *model, G4int matIdx)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &) override
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false) override
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii) override
 
G4double GetDEDXDispersion (const G4MaterialCutsCouple *couple, const G4DynamicParticle *dp, G4double length)
 
G4double CrossSectionPerVolume (G4double kineticEnergy, const G4MaterialCutsCouple *couple)
 
G4double CrossSectionPerVolume (G4double kineticEnergy, const G4MaterialCutsCouple *couple, G4double logKineticEnergy)
 
G4double MeanFreePath (const G4Track &track)
 
G4double ContinuousStepLimit (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)
 
G4VEmModelSelectModelForMaterial (G4double kinEnergy, size_t &idx) const
 
void AddEmModel (G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=nullptr)
 
void UpdateEmModel (const G4String &, G4double, G4double)
 
void SetEmModel (G4VEmModel *, G4int index=0)
 
G4VEmModelEmModel (size_t index=0) const
 
G4VEmModelGetModelByIndex (G4int idx=0, G4bool ver=false) const
 
G4int NumberOfModels () const
 
void SetFluctModel (G4VEmFluctuationModel *)
 
G4VEmFluctuationModelFluctModel ()
 
void SetBaseParticle (const G4ParticleDefinition *p)
 
const G4ParticleDefinitionParticle () const
 
const G4ParticleDefinitionBaseParticle () const
 
const G4ParticleDefinitionSecondaryParticle () const
 
void ActivateSubCutoff (G4bool val, const G4Region *region=nullptr)
 
void SetCrossSectionBiasingFactor (G4double f, G4bool flag=true)
 
void ActivateForcedInteraction (G4double length, const G4String &region, G4bool flag=true)
 
void ActivateSecondaryBiasing (const G4String &region, G4double factor, G4double energyLimit)
 
void AddCollaborativeProcess (G4VEnergyLossProcess *)
 
void SetLossFluctuations (G4bool val)
 
void SetIntegral (G4bool val)
 
G4bool IsIntegral () const
 
void SetIonisation (G4bool val)
 
G4bool IsIonisationProcess () const
 
void SetLinearLossLimit (G4double val)
 
void SetStepFunction (G4double v1, G4double v2)
 
void SetLowestEnergyLimit (G4double)
 
G4int NumberOfSubCutoffRegions () const
 
void SetDEDXTable (G4PhysicsTable *p, G4EmTableType tType)
 
void SetCSDARangeTable (G4PhysicsTable *pRange)
 
void SetRangeTableForLoss (G4PhysicsTable *p)
 
void SetSecondaryRangeTable (G4PhysicsTable *p)
 
void SetInverseRangeTable (G4PhysicsTable *p)
 
void SetLambdaTable (G4PhysicsTable *p)
 
void SetSubLambdaTable (G4PhysicsTable *p)
 
void SetDEDXBinning (G4int nbins)
 
void SetMinKinEnergy (G4double e)
 
G4double MinKinEnergy () const
 
void SetMaxKinEnergy (G4double e)
 
G4double MaxKinEnergy () const
 
G4double CrossSectionBiasingFactor () const
 
G4double GetDEDX (G4double kineticEnergy, const G4MaterialCutsCouple *)
 
G4double GetDEDX (G4double kineticEnergy, const G4MaterialCutsCouple *, G4double logKineticEnergy)
 
G4double GetDEDXForSubsec (G4double kineticEnergy, const G4MaterialCutsCouple *)
 
G4double GetRange (G4double kineticEnergy, const G4MaterialCutsCouple *)
 
G4double GetRange (G4double kineticEnergy, const G4MaterialCutsCouple *, G4double logKineticEnergy)
 
G4double GetCSDARange (G4double kineticEnergy, const G4MaterialCutsCouple *)
 
G4double GetRangeForLoss (G4double kineticEnergy, const G4MaterialCutsCouple *)
 
G4double GetRangeForLoss (G4double kineticEnergy, const G4MaterialCutsCouple *, G4double logKineticEnergy)
 
G4double GetKineticEnergy (G4double range, const G4MaterialCutsCouple *)
 
G4double GetLambda (G4double kineticEnergy, const G4MaterialCutsCouple *)
 
G4double GetLambda (G4double kineticEnergy, const G4MaterialCutsCouple *, G4double logKineticEnergy)
 
G4bool TablesAreBuilt () const
 
G4PhysicsTableDEDXTable () const
 
G4PhysicsTableDEDXTableForSubsec () const
 
G4PhysicsTableDEDXunRestrictedTable () const
 
G4PhysicsTableIonisationTable () const
 
G4PhysicsTableIonisationTableForSubsec () const
 
G4PhysicsTableCSDARangeTable () const
 
G4PhysicsTableSecondaryRangeTable () const
 
G4PhysicsTableRangeTableForLoss () const
 
G4PhysicsTableInverseRangeTable () const
 
G4PhysicsTableLambdaTable () const
 
G4PhysicsTableSubLambdaTable () const
 
const G4ElementGetCurrentElement () const
 
void SetDynamicMassCharge (G4double massratio, G4double charge2ratio)
 
- Public Member Functions inherited from G4VContinuousDiscreteProcess
 G4VContinuousDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VContinuousDiscreteProcess (G4VContinuousDiscreteProcess &)
 
virtual ~G4VContinuousDiscreteProcess ()
 
G4VContinuousDiscreteProcessoperator= (const G4VContinuousDiscreteProcess &)=delete
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4VProcessoperator= (const G4VProcess &)=delete
 
G4bool operator== (const G4VProcess &right) const
 
G4bool operator!= (const G4VProcess &right) const
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual void StartTracking (G4Track *)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
virtual void ProcessDescription (std::ostream &outfile) const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Protected Member Functions

virtual void StreamProcessInfo (std::ostream &) const
 
virtual void InitialiseEnergyLossProcess (const G4ParticleDefinition *, const G4ParticleDefinition *)=0
 
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *, G4double cut)
 
virtual G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
virtual G4double GetContinuousStepLimit (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety) override
 
G4PhysicsVectorLambdaPhysicsVector (const G4MaterialCutsCouple *, G4double cut)
 
size_t CurrentMaterialCutsCoupleIndex () const
 
void SelectModel (G4double kinEnergy)
 
void SetParticle (const G4ParticleDefinition *p)
 
void SetSecondaryParticle (const G4ParticleDefinition *p)
 
- Protected Member Functions inherited from G4VContinuousDiscreteProcess
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
virtual G4double GetContinuousStepLimit (const G4Track &aTrack, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)=0
 
void SetGPILSelection (G4GPILSelection selection)
 
G4GPILSelection GetGPILSelection () const
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Protected Attributes

G4ParticleChangeForLoss fParticleChange
 
const G4MaterialcurrentMaterial
 
const G4MaterialCutsCouplecurrentCouple
 
size_t currentCoupleIndex
 
G4double preStepLambda
 
G4double fRange
 
G4double computedRange
 
G4double preStepKinEnergy
 
G4double preStepLogKinEnergy
 
G4double preStepScaledEnergy
 
G4double preStepLogScaledEnergy
 
G4double preStepRangeEnergy
 
G4double mfpKinEnergy
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft = -1.0
 
G4double currentInteractionLength = -1.0
 
G4double theInitialNumberOfInteractionLength = -1.0
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType = fNotDefined
 
G4int theProcessSubType = -1
 
G4double thePILfactor = 1.0
 
G4int verboseLevel = 0
 
G4bool enableAtRestDoIt = true
 
G4bool enableAlongStepDoIt = true
 
G4bool enablePostStepDoIt = true
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 

Detailed Description

Definition at line 81 of file G4VEnergyLossProcess.hh.

Constructor & Destructor Documentation

◆ G4VEnergyLossProcess()

G4VEnergyLossProcess::G4VEnergyLossProcess ( const G4String name = "EnergyLoss",
G4ProcessType  type = fElectromagnetic 
)

Definition at line 86 of file G4VEnergyLossProcess.cc.

87 :
89 secondaryParticle(nullptr),
90 nSCoffRegions(0),
91 idxSCoffRegions(nullptr),
92 nProcesses(0),
93 theDEDXTable(nullptr),
94 theDEDXSubTable(nullptr),
95 theDEDXunRestrictedTable(nullptr),
96 theIonisationTable(nullptr),
97 theIonisationSubTable(nullptr),
98 theRangeTableForLoss(nullptr),
99 theCSDARangeTable(nullptr),
100 theSecondaryRangeTable(nullptr),
101 theInverseRangeTable(nullptr),
102 theLambdaTable(nullptr),
103 theSubLambdaTable(nullptr),
104 baseParticle(nullptr),
105 lossFluctuationFlag(true),
106 rndmStepFlag(false),
107 tablesAreBuilt(false),
108 integral(true),
109 isIon(false),
110 isIonisation(true),
111 useSubCutoff(false),
112 useDeexcitation(false),
113 currentCouple(nullptr),
114 mfpKinEnergy(0.0),
115 particle(nullptr)
116{
117 theParameters = G4EmParameters::Instance();
119
120 // low energy limit
121 lowestKinEnergy = theParameters->LowestElectronEnergy();
122 preStepKinEnergy = 0.0;
124 preStepRangeEnergy = 0.0;
126
127 // Size of tables assuming spline
128 minKinEnergy = 0.1*keV;
129 maxKinEnergy = 100.0*TeV;
130 nBins = 84;
131 maxKinEnergyCSDA = 1.0*GeV;
132 nBinsCSDA = 35;
133 actMinKinEnergy = actMaxKinEnergy = actBinning = actLinLossLimit
134 = actLossFluc = actIntegral = false;
135
136 // default linear loss limit for spline
137 linLossLimit = 0.01;
138 dRoverRange = 0.2;
139 finalRange = CLHEP::mm;
140
141 // default lambda factor
142 lambdaFactor = 0.8;
143 logLambdafactor = G4Log(lambdaFactor);
144
145 // cross section biasing
146 biasFactor = 1.0;
147
148 // particle types
149 theElectron = G4Electron::Electron();
150 thePositron = G4Positron::Positron();
151 theGamma = G4Gamma::Gamma();
152 theGenericIon = nullptr;
153
154 // run time objects
157 modelManager = new G4EmModelManager();
159 ->GetSafetyHelper();
160 aGPILSelection = CandidateForSelection;
161
162 // initialise model
163 lManager = G4LossTableManager::Instance();
164 lManager->Register(this);
165 G4LossTableBuilder* bld = lManager->GetTableBuilder();
166 theDensityFactor = bld->GetDensityFactors();
167 theDensityIdx = bld->GetCoupleIndexes();
168
169 fluctModel = nullptr;
170 currentModel = nullptr;
171 atomDeexcitation = nullptr;
172 subcutProducer = nullptr;
173
174 biasManager = nullptr;
175 biasFlag = false;
176 weightFlag = false;
177 isMaster = true;
178 lastIdx = 0;
179
180 scTracks.reserve(5);
181 secParticles.reserve(5);
182
183 theCuts = theSubCuts = nullptr;
184 currentMaterial = nullptr;
185 currentCoupleIndex = basedCoupleIndex = 0;
186 massRatio = fFactor = reduceFactor = chargeSqRatio = 1.0;
187 preStepLambda = preStepScaledEnergy = fRange = logMassRatio = 0.0;
189
190 secID = biasID = subsecID = -1;
191}
@ CandidateForSelection
G4double G4Log(G4double x)
Definition: G4Log.hh:226
static G4Electron * Electron()
Definition: G4Electron.cc:93
static G4EmParameters * Instance()
G4double LowestElectronEnergy() const
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
const std::vector< G4double > * GetDensityFactors() const
const std::vector< G4int > * GetCoupleIndexes() const
static G4LossTableManager * Instance()
G4LossTableBuilder * GetTableBuilder()
void Register(G4VEnergyLossProcess *p)
static G4Positron * Positron()
Definition: G4Positron.cc:93
static G4TransportationManager * GetTransportationManager()
G4SafetyHelper * GetSafetyHelper() const
G4ParticleChangeForLoss fParticleChange
const G4MaterialCutsCouple * currentCouple
const G4Material * currentMaterial
void SetSecondaryWeightByProcess(G4bool)
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:412
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:321
#define LOG_EKIN_MIN
Definition: templates.hh:98
#define DBL_MAX
Definition: templates.hh:62

◆ ~G4VEnergyLossProcess()

G4VEnergyLossProcess::~G4VEnergyLossProcess ( )
virtual

Definition at line 195 of file G4VEnergyLossProcess.cc.

196{
197 /*
198 G4cout << "** G4VEnergyLossProcess::~G4VEnergyLossProcess() for "
199 << GetProcessName() << " isMaster: " << isMaster
200 << " basePart: " << baseParticle
201 << G4endl;
202 */
203 Clean();
204
205 // G4cout << " isIonisation " << isIonisation << " "
206 // << theDEDXTable << " " << theIonisationTable << G4endl;
207
208 if (isMaster && !baseParticle) {
209 if(theDEDXTable) {
210
211 //G4cout << " theIonisationTable " << theIonisationTable << G4endl;
212 if(theIonisationTable == theDEDXTable) { theIonisationTable = nullptr; }
213 //G4cout << " delete theDEDXTable " << theDEDXTable << G4endl;
214 theDEDXTable->clearAndDestroy();
215 delete theDEDXTable;
216 theDEDXTable = nullptr;
217 if(theDEDXSubTable) {
218 if(theIonisationSubTable == theDEDXSubTable)
219 { theIonisationSubTable = nullptr; }
220 theDEDXSubTable->clearAndDestroy();
221 delete theDEDXSubTable;
222 theDEDXSubTable = nullptr;
223 }
224 }
225 //G4cout << " theIonisationTable " << theIonisationTable << G4endl;
226 if(theIonisationTable) {
227 //G4cout << " delete theIonisationTable " << theIonisationTable << G4endl;
228 theIonisationTable->clearAndDestroy();
229 delete theIonisationTable;
230 theIonisationTable = nullptr;
231 }
232 if(theIonisationSubTable) {
233 theIonisationSubTable->clearAndDestroy();
234 delete theIonisationSubTable;
235 theIonisationSubTable = nullptr;
236 }
237 if(theDEDXunRestrictedTable && isIonisation) {
238 theDEDXunRestrictedTable->clearAndDestroy();
239 delete theDEDXunRestrictedTable;
240 theDEDXunRestrictedTable = nullptr;
241 }
242 if(theCSDARangeTable && isIonisation) {
243 theCSDARangeTable->clearAndDestroy();
244 delete theCSDARangeTable;
245 theCSDARangeTable = nullptr;
246 }
247 //G4cout << "delete RangeTable: " << theRangeTableForLoss << G4endl;
248 if(theRangeTableForLoss && isIonisation) {
249 theRangeTableForLoss->clearAndDestroy();
250 delete theRangeTableForLoss;
251 theRangeTableForLoss = nullptr;
252 }
253 //G4cout << "delete InvRangeTable: " << theInverseRangeTable << G4endl;
254 if(theInverseRangeTable && isIonisation /*&& !isIon*/) {
255 theInverseRangeTable->clearAndDestroy();
256 delete theInverseRangeTable;
257 theInverseRangeTable = nullptr;
258 }
259 //G4cout << "delete LambdaTable: " << theLambdaTable << G4endl;
260 if(theLambdaTable) {
261 theLambdaTable->clearAndDestroy();
262 delete theLambdaTable;
263 theLambdaTable = nullptr;
264 }
265 if(theSubLambdaTable) {
266 theSubLambdaTable->clearAndDestroy();
267 delete theSubLambdaTable;
268 theSubLambdaTable = nullptr;
269 }
270 }
271
272 delete modelManager;
273 delete biasManager;
274 lManager->DeRegister(this);
275 //G4cout << "** all removed" << G4endl;
276}
void DeRegister(G4VEnergyLossProcess *p)
void clearAndDestroy()

Member Function Documentation

◆ ActivateForcedInteraction()

void G4VEnergyLossProcess::ActivateForcedInteraction ( G4double  length,
const G4String region,
G4bool  flag = true 
)

Definition at line 2236 of file G4VEnergyLossProcess.cc.

2239{
2240 if(!biasManager) { biasManager = new G4EmBiasingManager(); }
2241 if(1 < verboseLevel) {
2242 G4cout << "### ActivateForcedInteraction: for "
2243 << " process " << GetProcessName()
2244 << " length(mm)= " << length/mm
2245 << " in G4Region <" << region
2246 << "> weightFlag= " << flag
2247 << G4endl;
2248 }
2249 weightFlag = flag;
2250 biasManager->ActivateForcedInteraction(length, region);
2251}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void ActivateForcedInteraction(G4double length=0.0, const G4String &r="")
G4int verboseLevel
Definition: G4VProcess.hh:356
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382

Referenced by G4EmExtraParameters::DefineRegParamForLoss().

◆ ActivateSecondaryBiasing()

void G4VEnergyLossProcess::ActivateSecondaryBiasing ( const G4String region,
G4double  factor,
G4double  energyLimit 
)

Definition at line 2256 of file G4VEnergyLossProcess.cc.

2259{
2260 if (0.0 <= factor) {
2261
2262 // Range cut can be applied only for e-
2263 if(0.0 == factor && secondaryParticle != G4Electron::Electron())
2264 { return; }
2265
2266 if(!biasManager) { biasManager = new G4EmBiasingManager(); }
2267 biasManager->ActivateSecondaryBiasing(region, factor, energyLimit);
2268 if(1 < verboseLevel) {
2269 G4cout << "### ActivateSecondaryBiasing: for "
2270 << " process " << GetProcessName()
2271 << " factor= " << factor
2272 << " in G4Region <" << region
2273 << "> energyLimit(MeV)= " << energyLimit/MeV
2274 << G4endl;
2275 }
2276 }
2277}
void ActivateSecondaryBiasing(const G4String &region, G4double factor, G4double energyLimit)

Referenced by G4EmExtraParameters::DefineRegParamForLoss().

◆ ActivateSubCutoff()

void G4VEnergyLossProcess::ActivateSubCutoff ( G4bool  val,
const G4Region region = nullptr 
)

Definition at line 956 of file G4VEnergyLossProcess.cc.

957{
959 const G4Region* reg = r;
960 if (!reg) {
961 reg = regionStore->GetRegion("DefaultRegionForTheWorld", false);
962 }
963
964 // the region is in the list
965 if (nSCoffRegions > 0) {
966 for (G4int i=0; i<nSCoffRegions; ++i) {
967 if (reg == scoffRegions[i]) {
968 return;
969 }
970 }
971 }
972 // new region
973 if(val) {
974 scoffRegions.push_back(reg);
975 ++nSCoffRegions;
976 }
977}
int G4int
Definition: G4Types.hh:85
static G4RegionStore * GetInstance()
G4Region * GetRegion(const G4String &name, G4bool verbose=true) const

Referenced by G4EmExtraParameters::DefineRegParamForLoss().

◆ AddCollaborativeProcess()

void G4VEnergyLossProcess::AddCollaborativeProcess ( G4VEnergyLossProcess p)

Definition at line 1985 of file G4VEnergyLossProcess.cc.

1987{
1988 G4bool add = true;
1989 if(p->GetProcessName() != "eBrem") { add = false; }
1990 if(add && nProcesses > 0) {
1991 for(G4int i=0; i<nProcesses; ++i) {
1992 if(p == scProcesses[i]) {
1993 add = false;
1994 break;
1995 }
1996 }
1997 }
1998 if(add) {
1999 scProcesses.push_back(p);
2000 ++nProcesses;
2001 if (1 < verboseLevel) {
2002 G4cout << "### The process " << p->GetProcessName()
2003 << " is added to the list of collaborative processes of "
2004 << GetProcessName() << G4endl;
2005 }
2006 }
2007}
bool G4bool
Definition: G4Types.hh:86

◆ AddEmModel()

void G4VEnergyLossProcess::AddEmModel ( G4int  order,
G4VEmModel p,
G4VEmFluctuationModel fluc = 0,
const G4Region region = nullptr 
)

◆ AlongStepDoIt()

G4VParticleChange * G4VEnergyLossProcess::AlongStepDoIt ( const G4Track track,
const G4Step step 
)
overridevirtual

Reimplemented from G4VContinuousDiscreteProcess.

Definition at line 1186 of file G4VEnergyLossProcess.cc.

1188{
1190 // The process has range table - calculate energy loss
1191 if(!isIonisation || !currentModel->IsActive(preStepScaledEnergy)) {
1192 return &fParticleChange;
1193 }
1194
1195 // Get the actual (true) Step length
1196 G4double length = step.GetStepLength();
1197 if(length <= 0.0) { return &fParticleChange; }
1198 G4double eloss = 0.0;
1199
1200 /*
1201 if(-1 < verboseLevel) {
1202 const G4ParticleDefinition* d = track.GetParticleDefinition();
1203 G4cout << "AlongStepDoIt for "
1204 << GetProcessName() << " and particle "
1205 << d->GetParticleName()
1206 << " eScaled(MeV)= " << preStepScaledEnergy/MeV
1207 << " range(mm)= " << fRange/mm
1208 << " s(mm)= " << length/mm
1209 << " rf= " << reduceFactor
1210 << " q^2= " << chargeSqRatio
1211 << " md= " << d->GetPDGMass()
1212 << " status= " << track.GetTrackStatus()
1213 << " " << track.GetMaterial()->GetName()
1214 << G4endl;
1215 }
1216 */
1217
1218 const G4DynamicParticle* dynParticle = track.GetDynamicParticle();
1219
1220 // define new weight for primary and secondaries
1222 if(weightFlag) {
1223 weight /= biasFactor;
1225 }
1226
1227 // stopping
1228 if (length >= fRange || preStepKinEnergy <= lowestKinEnergy) {
1229 eloss = preStepKinEnergy;
1230 if (useDeexcitation) {
1231 atomDeexcitation->AlongStepDeexcitation(scTracks, step,
1232 eloss, currentCoupleIndex);
1233 if(scTracks.size() > 0) { FillSecondariesAlongStep(eloss, weight); }
1234 eloss = std::max(eloss, 0.0);
1235 }
1238 return &fParticleChange;
1239 }
1240 //G4cout << theDEDXTable << " idx= " << basedCoupleIndex
1241 // << " " << GetProcessName() << " "<< currentMaterial->GetName()<<G4endl;
1242 //if(particle->GetParticleName() == "e-")G4cout << (*theDEDXTable) <<G4endl;
1243 // Short step
1244 eloss = GetDEDXForScaledEnergy(preStepScaledEnergy, preStepLogScaledEnergy);
1245 eloss *= length;
1246
1247 //G4cout << "eloss= " << eloss << G4endl;
1248
1249 // Long step
1250 if(eloss > preStepKinEnergy*linLossLimit) {
1251
1252 G4double x = (fRange - length)/reduceFactor;
1253 //G4cout << "x= " << x << " " << theInverseRangeTable << G4endl;
1254 eloss = preStepKinEnergy - ScaledKinEnergyForLoss(x)/massRatio;
1255
1256 /*
1257 if(-1 < verboseLevel)
1258 G4cout << "Long STEP: rPre(mm)= "
1259 << GetScaledRangeForScaledEnergy(preStepScaledEnergy)/mm
1260 << " rPost(mm)= " << x/mm
1261 << " ePre(MeV)= " << preStepScaledEnergy/MeV
1262 << " eloss(MeV)= " << eloss/MeV
1263 << " eloss0(MeV)= "
1264 << GetDEDXForScaledEnergy(preStepScaledEnergy)*length/MeV
1265 << " lim(MeV)= " << preStepKinEnergy*linLossLimit/MeV
1266 << G4endl;
1267 */
1268 }
1269
1270 /*
1271 G4double eloss0 = eloss;
1272 if(-1 < verboseLevel ) {
1273 G4cout << "Before fluct: eloss(MeV)= " << eloss/MeV
1274 << " e-eloss= " << preStepKinEnergy-eloss
1275 << " step(mm)= " << length/mm
1276 << " range(mm)= " << fRange/mm
1277 << " fluct= " << lossFluctuationFlag
1278 << G4endl;
1279 }
1280 */
1281
1282 G4double cut = (*theCuts)[currentCoupleIndex];
1283 G4double esec = 0.0;
1284
1285 //G4cout << "cut= " << cut << " useSubCut= " << useSubCutoff << G4endl;
1286
1287 // SubCutOff
1288 if(useSubCutoff && !subcutProducer) {
1289 if(idxSCoffRegions[currentCoupleIndex]) {
1290
1291 G4bool yes = false;
1292 const G4StepPoint* prePoint = step.GetPreStepPoint();
1293
1294 // Check boundary
1295 if(prePoint->GetStepStatus() == fGeomBoundary) { yes = true; }
1296
1297 // Check PrePoint
1298 else {
1299 G4double preSafety = prePoint->GetSafety();
1300 G4double rcut =
1302
1303 // recompute presafety
1304 if(preSafety < rcut) {
1305 preSafety = safetyHelper->ComputeSafety(prePoint->GetPosition(),
1306 rcut);
1307 }
1308
1309 if(preSafety < rcut) { yes = true; }
1310
1311 // Check PostPoint
1312 else {
1313 G4double postSafety = preSafety - length;
1314 if(postSafety < rcut) {
1315 postSafety = safetyHelper->ComputeSafety(
1316 step.GetPostStepPoint()->GetPosition(), rcut);
1317 if(postSafety < rcut) { yes = true; }
1318 }
1319 }
1320 }
1321
1322 // Decided to start subcut sampling
1323 if(yes) {
1324
1325 cut = (*theSubCuts)[currentCoupleIndex];
1326 eloss -= GetSubDEDXForScaledEnergy(preStepScaledEnergy)*length;
1327 esec = SampleSubCutSecondaries(scTracks, step,
1328 currentModel,currentCoupleIndex);
1329 // add bremsstrahlung sampling
1330 /*
1331 if(nProcesses > 0) {
1332 for(G4int i=0; i<nProcesses; ++i) {
1333 (scProcesses[i])->SampleSubCutSecondaries(
1334 scTracks, step, (scProcesses[i])->
1335 SelectModelForMaterial(preStepKinEnergy, currentCoupleIndex),
1336 currentCoupleIndex);
1337 }
1338 }
1339 */
1340 }
1341 }
1342 }
1343
1344 // Corrections, which cannot be tabulated
1345 if(isIon) {
1346 G4double eadd = 0.0;
1347 G4double eloss_before = eloss;
1348 currentModel->CorrectionsAlongStep(currentCouple, dynParticle,
1349 eloss, eadd, length);
1350 if(eloss < 0.0) { eloss = 0.5*eloss_before; }
1351 }
1352
1353 // Sample fluctuations
1354 if (lossFluctuationFlag) {
1355 G4VEmFluctuationModel* fluc = currentModel->GetModelOfFluctuations();
1356 if(eloss + esec < preStepKinEnergy) {
1357
1358 G4double tmax =
1359 std::min(currentModel->MaxSecondaryKinEnergy(dynParticle),cut);
1360 eloss = fluc->SampleFluctuations(currentCouple,dynParticle,
1361 tmax,length,eloss);
1362 /*
1363 if(-1 < verboseLevel)
1364 G4cout << "After fluct: eloss(MeV)= " << eloss/MeV
1365 << " fluc= " << (eloss-eloss0)/MeV
1366 << " ChargeSqRatio= " << chargeSqRatio
1367 << " massRatio= " << massRatio
1368 << " tmax= " << tmax
1369 << G4endl;
1370 */
1371 }
1372 }
1373
1374 // deexcitation
1375 if (useDeexcitation) {
1376 G4double esecfluo = preStepKinEnergy - esec;
1377 G4double de = esecfluo;
1378 //G4double eloss0 = eloss;
1379 /*
1380 G4cout << "### 1: E(keV)= " << preStepKinEnergy/keV
1381 << " Efluomax(keV)= " << de/keV
1382 << " Eloss(keV)= " << eloss/keV << G4endl;
1383 */
1384 atomDeexcitation->AlongStepDeexcitation(scTracks, step,
1385 de, currentCoupleIndex);
1386
1387 // sum of de-excitation energies
1388 esecfluo -= de;
1389
1390 // subtracted from energy loss
1391 if(eloss >= esecfluo) {
1392 esec += esecfluo;
1393 eloss -= esecfluo;
1394 } else {
1395 esec += esecfluo;
1396 eloss = 0.0;
1397 }
1398 /*
1399 if(esecfluo > 0.0) {
1400 G4cout << "### 2: E(keV)= " << preStepKinEnergy/keV
1401 << " Esec(keV)= " << esec/keV
1402 << " Esecf(kV)= " << esecfluo/keV
1403 << " Eloss0(kV)= " << eloss0/keV
1404 << " Eloss(keV)= " << eloss/keV
1405 << G4endl;
1406 }
1407 */
1408 }
1409 if(subcutProducer && idxSCoffRegions[currentCoupleIndex]) {
1410 subcutProducer->SampleSecondaries(step, scTracks, eloss, cut);
1411 }
1412 if(scTracks.size() > 0) { FillSecondariesAlongStep(eloss, weight); }
1413
1414 // Energy balance
1415 G4double finalT = preStepKinEnergy - eloss - esec;
1416 if (finalT <= lowestKinEnergy) {
1417 eloss += finalT;
1418 finalT = 0.0;
1419 } else if(isIon) {
1421 currentModel->GetParticleCharge(track.GetParticleDefinition(),
1422 currentMaterial,finalT));
1423 }
1424
1425 eloss = std::max(eloss, 0.0);
1426
1429 /*
1430 if(-1 < verboseLevel) {
1431 G4double del = finalT + eloss + esec - preStepKinEnergy;
1432 G4cout << "Final value eloss(MeV)= " << eloss/MeV
1433 << " preStepKinEnergy= " << preStepKinEnergy
1434 << " postStepKinEnergy= " << finalT
1435 << " de(keV)= " << del/keV
1436 << " lossFlag= " << lossFluctuationFlag
1437 << " status= " << track.GetTrackStatus()
1438 << G4endl;
1439 }
1440 */
1441 return &fParticleChange;
1442}
@ fGeomBoundary
Definition: G4StepStatus.hh:43
double G4double
Definition: G4Types.hh:83
G4ProductionCuts * GetProductionCuts() const
void InitializeForAlongStep(const G4Track &)
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void SetProposedCharge(G4double theCharge)
G4double GetProductionCut(G4int index) const
G4double ComputeSafety(const G4ThreeVector &pGlobalPoint, G4double maxRadius=DBL_MAX)
G4StepStatus GetStepStatus() const
G4double GetSafety() const
const G4ThreeVector & GetPosition() const
G4StepPoint * GetPreStepPoint() const
G4double GetStepLength() const
G4StepPoint * GetPostStepPoint() const
const G4ParticleDefinition * GetParticleDefinition() const
const G4DynamicParticle * GetDynamicParticle() const
void AlongStepDeexcitation(std::vector< G4Track * > &tracks, const G4Step &step, G4double &eLoss, G4int coupleIndex)
virtual G4double SampleFluctuations(const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmax, G4double length, G4double meanLoss)=0
virtual void CorrectionsAlongStep(const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
Definition: G4VEmModel.cc:408
G4VEmFluctuationModel * GetModelOfFluctuations()
Definition: G4VEmModel.hh:604
virtual G4double GetParticleCharge(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:400
G4bool IsActive(G4double kinEnergy) const
Definition: G4VEmModel.hh:785
G4double MaxSecondaryKinEnergy(const G4DynamicParticle *dynParticle)
Definition: G4VEmModel.hh:510
G4double SampleSubCutSecondaries(std::vector< G4Track * > &, const G4Step &, G4VEmModel *model, G4int matIdx)
G4double GetParentWeight() const
void ProposeWeight(G4double finalWeight)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
virtual void SampleSecondaries(const G4Step &step, std::vector< G4Track * > &tracks, G4double &eloss, G4double cut) const =0

◆ AlongStepGetPhysicalInteractionLength()

G4double G4VEnergyLossProcess::AlongStepGetPhysicalInteractionLength ( const G4Track ,
G4double  previousStepSize,
G4double  currentMinimumStep,
G4double currentSafety,
G4GPILSelection selection 
)
overridevirtual

Reimplemented from G4VContinuousDiscreteProcess.

Definition at line 1023 of file G4VEnergyLossProcess.cc.

1026{
1027 G4double x = DBL_MAX;
1028 *selection = aGPILSelection;
1029 if(isIonisation && currentModel->IsActive(preStepScaledEnergy)) {
1030 fRange = reduceFactor*GetScaledRangeForScaledEnergy(preStepScaledEnergy,
1032 G4double finR = (rndmStepFlag) ? std::min(finalRange,
1034 x = (fRange > finR) ?
1035 fRange*dRoverRange + finR*(1.0-dRoverRange)*(2.0-finR/fRange) : fRange;
1036 // if(particle->GetPDGMass() > 0.9*GeV)
1037 /*
1038 G4cout<<GetProcessName()<<": e= "<<preStepKinEnergy
1039 <<" range= "<<fRange << " idx= " << basedCoupleIndex
1040 << " finR= " << finR
1041 << " limit= " << x <<G4endl;
1042 G4cout << "massRatio= " << massRatio << " Q^2= " << chargeSqRatio
1043 << " finR= " << finR << " dRoverRange= " << dRoverRange
1044 << " finalRange= " << finalRange << G4endl;
1045 */
1046 }
1047 //G4cout<<GetProcessName()<<": e= "<<preStepKinEnergy
1048 //<<" stepLimit= "<<x<<G4endl;
1049 return x;
1050}

Referenced by ContinuousStepLimit().

◆ BaseParticle()

const G4ParticleDefinition * G4VEnergyLossProcess::BaseParticle ( ) const
inline

◆ BuildDEDXTable()

G4PhysicsTable * G4VEnergyLossProcess::BuildDEDXTable ( G4EmTableType  tType = fRestricted)

Definition at line 713 of file G4VEnergyLossProcess.cc.

714{
715 if(1 < verboseLevel ) {
716 G4cout << "G4VEnergyLossProcess::BuildDEDXTable() of type " << tType
717 << " for " << GetProcessName()
718 << " and particle " << particle->GetParticleName()
719 << G4endl;
720 }
721 G4PhysicsTable* table = nullptr;
722 G4double emax = maxKinEnergy;
723 G4int bin = nBins;
724
725 if(fTotal == tType) {
726 emax = maxKinEnergyCSDA;
727 bin = nBinsCSDA;
728 table = theDEDXunRestrictedTable;
729 } else if(fRestricted == tType) {
730 table = theDEDXTable;
731 } else if(fSubRestricted == tType) {
732 table = theDEDXSubTable;
733 } else {
734 G4cout << "G4VEnergyLossProcess::BuildDEDXTable WARNING: wrong type "
735 << tType << G4endl;
736 }
737
738 // Access to materials
739 const G4ProductionCutsTable* theCoupleTable=
741 size_t numOfCouples = theCoupleTable->GetTableSize();
742
743 if(1 < verboseLevel) {
744 G4cout << numOfCouples << " materials"
745 << " minKinEnergy= " << minKinEnergy
746 << " maxKinEnergy= " << emax
747 << " nbin= " << bin
748 << " EmTableType= " << tType
749 << " table= " << table << " " << this
750 << G4endl;
751 }
752 if(!table) { return table; }
753
754 G4LossTableBuilder* bld = lManager->GetTableBuilder();
755 G4bool splineFlag = theParameters->Spline();
756 G4PhysicsLogVector* aVector = nullptr;
757 G4PhysicsLogVector* bVector = nullptr;
758
759 for(size_t i=0; i<numOfCouples; ++i) {
760
761 if(1 < verboseLevel) {
762 G4cout << "G4VEnergyLossProcess::BuildDEDXVector Idx= " << i
763 << " flagTable= " << table->GetFlag(i)
764 << " Flag= " << bld->GetFlag(i) << G4endl;
765 }
766 if(bld->GetFlag(i)) {
767
768 // create physics vector and fill it
769 const G4MaterialCutsCouple* couple =
770 theCoupleTable->GetMaterialCutsCouple(i);
771 if((*table)[i]) { delete (*table)[i]; }
772 if(bVector) {
773 aVector = new G4PhysicsLogVector(*bVector);
774 } else {
775 bVector = new G4PhysicsLogVector(minKinEnergy, emax, bin);
776 aVector = bVector;
777 }
778 aVector->SetSpline(splineFlag);
779
780 modelManager->FillDEDXVector(aVector, couple, tType);
781 if(splineFlag) { aVector->FillSecondDerivatives(); }
782
783 // Insert vector for this material into the table
785 }
786 }
787
788 if(1 < verboseLevel) {
789 G4cout << "G4VEnergyLossProcess::BuildDEDXTable(): table is built for "
790 << particle->GetParticleName()
791 << " and process " << GetProcessName()
792 << G4endl;
793 if(2 < verboseLevel) G4cout << (*table) << G4endl;
794 }
795
796 return table;
797}
@ fSubRestricted
@ fTotal
@ fRestricted
void FillDEDXVector(G4PhysicsVector *, const G4MaterialCutsCouple *, G4EmTableType t=fRestricted)
G4bool Spline() const
G4bool GetFlag(size_t idx)
const G4String & GetParticleName() const
static void SetPhysicsVector(G4PhysicsTable *physTable, std::size_t idx, G4PhysicsVector *vec)
G4bool GetFlag(std::size_t i) const
void FillSecondDerivatives()
void SetSpline(G4bool)
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()

◆ BuildLambdaTable()

G4PhysicsTable * G4VEnergyLossProcess::BuildLambdaTable ( G4EmTableType  tType = fRestricted)

Definition at line 801 of file G4VEnergyLossProcess.cc.

802{
803 G4PhysicsTable* table = nullptr;
804
805 if(fRestricted == tType) {
806 table = theLambdaTable;
807 } else if(fSubRestricted == tType) {
808 table = theSubLambdaTable;
809 } else {
810 G4cout << "G4VEnergyLossProcess::BuildLambdaTable WARNING: wrong type "
811 << tType << G4endl;
812 }
813
814 if(1 < verboseLevel) {
815 G4cout << "G4VEnergyLossProcess::BuildLambdaTable() of type "
816 << tType << " for process "
817 << GetProcessName() << " and particle "
818 << particle->GetParticleName()
819 << " EmTableType= " << tType
820 << " table= " << table
821 << G4endl;
822 }
823 if(!table) {return table;}
824
825 // Access to materials
826 const G4ProductionCutsTable* theCoupleTable=
828 size_t numOfCouples = theCoupleTable->GetTableSize();
829
830 G4LossTableBuilder* bld = lManager->GetTableBuilder();
831 theDensityFactor = bld->GetDensityFactors();
832 theDensityIdx = bld->GetCoupleIndexes();
833
834 G4bool splineFlag = theParameters->Spline();
835 G4PhysicsLogVector* aVector = nullptr;
836 G4double scale = G4Log(maxKinEnergy/minKinEnergy);
837
838 for(size_t i=0; i<numOfCouples; ++i) {
839
840 if (bld->GetFlag(i)) {
841
842 // create physics vector and fill it
843 const G4MaterialCutsCouple* couple =
844 theCoupleTable->GetMaterialCutsCouple(i);
845 delete (*table)[i];
846
847 G4bool startNull = true;
848 G4double emin =
849 MinPrimaryEnergy(particle,couple->GetMaterial(),(*theCuts)[i]);
850 if(minKinEnergy > emin) {
851 emin = minKinEnergy;
852 startNull = false;
853 }
854
855 G4double emax = maxKinEnergy;
856 if(emax <= emin) { emax = 2*emin; }
857 G4int bin = G4lrint(nBins*G4Log(emax/emin)/scale);
858 bin = std::max(bin, 3);
859 aVector = new G4PhysicsLogVector(emin, emax, bin);
860 aVector->SetSpline(splineFlag);
861
862 modelManager->FillLambdaVector(aVector, couple, startNull, tType);
863 if(splineFlag) { aVector->FillSecondDerivatives(); }
864
865 // Insert vector for this material into the table
867 }
868 }
869
870 if(1 < verboseLevel) {
871 G4cout << "Lambda table is built for "
872 << particle->GetParticleName()
873 << G4endl;
874 }
875
876 return table;
877}
void FillLambdaVector(G4PhysicsVector *, const G4MaterialCutsCouple *, G4bool startFromNull=true, G4EmTableType t=fRestricted)
const G4Material * GetMaterial() const
virtual G4double MinPrimaryEnergy(const G4ParticleDefinition *, const G4Material *, G4double cut)
int G4lrint(double ad)
Definition: templates.hh:134

◆ BuildPhysicsTable()

void G4VEnergyLossProcess::BuildPhysicsTable ( const G4ParticleDefinition part)
overridevirtual

Reimplemented from G4VProcess.

Reimplemented in G4ePolarizedIonisation.

Definition at line 613 of file G4VEnergyLossProcess.cc.

614{
615 if(1 < verboseLevel) {
616 G4cout << "### G4VEnergyLossProcess::BuildPhysicsTable() for "
617 << GetProcessName()
618 << " and particle " << part.GetParticleName()
619 << "; local: " << particle->GetParticleName();
620 if(baseParticle) {
621 G4cout << "; base: " << baseParticle->GetParticleName();
622 }
623 G4cout << " TablesAreBuilt= " << tablesAreBuilt
624 << " isIon= " << isIon << " " << this << G4endl;
625 }
626
627 if(&part == particle) {
628
629 if(isMaster) {
630 lManager->BuildPhysicsTable(particle, this);
631
632 } else {
633
634 const G4VEnergyLossProcess* masterProcess =
635 static_cast<const G4VEnergyLossProcess*>(GetMasterProcess());
636
637 // copy table pointers from master thread
638 SetDEDXTable(masterProcess->DEDXTable(),fRestricted);
640 SetDEDXTable(masterProcess->DEDXunRestrictedTable(),fTotal);
643 SetRangeTableForLoss(masterProcess->RangeTableForLoss());
644 SetCSDARangeTable(masterProcess->CSDARangeTable());
646 SetInverseRangeTable(masterProcess->InverseRangeTable());
647 SetLambdaTable(masterProcess->LambdaTable());
648 SetSubLambdaTable(masterProcess->SubLambdaTable());
649 isIonisation = masterProcess->IsIonisationProcess();
650
651 tablesAreBuilt = true;
652 // local initialisation of models
653 G4bool printing = true;
654 G4int numberOfModels = modelManager->NumberOfModels();
655 for(G4int i=0; i<numberOfModels; ++i) {
656 G4VEmModel* mod = GetModelByIndex(i, printing);
657 G4VEmModel* mod0= masterProcess->GetModelByIndex(i,printing);
658 mod->InitialiseLocal(particle, mod0);
659 }
660
661 lManager->LocalPhysicsTables(particle, this);
662 }
663
664 // needs to be done only once
665 safetyHelper->InitialiseHelper();
666 }
667 // explicitly defined printout by particle name
668 G4String num = part.GetParticleName();
669 if(1 < verboseLevel ||
670 (0 < verboseLevel && (num == "e-" ||
671 num == "e+" || num == "mu+" ||
672 num == "mu-" || num == "proton"||
673 num == "pi+" || num == "pi-" ||
674 num == "kaon+" || num == "kaon-" ||
675 num == "alpha" || num == "anti_proton" ||
676 num == "GenericIon"|| num == "alpha++" ||
677 num == "alpha+" )))
678 {
679 StreamInfo(G4cout, part);
680 }
681
682 // Added tracking cut to avoid tracking artifacts
683 // identify deexcitation flag
684 if(isIonisation) {
685 atomDeexcitation = lManager->AtomDeexcitation();
686 if(nSCoffRegions > 0) { subcutProducer = lManager->SubCutProducer(); }
687 if(atomDeexcitation) {
688 if(atomDeexcitation->IsPIXEActive()) { useDeexcitation = true; }
689 }
690 }
691 /*
692 G4cout << "** G4VEnergyLossProcess::BuildPhysicsTable() for "
693 << GetProcessName() << " and " << particle->GetParticleName()
694 << " isMaster: " << isMaster << " isIonisation: "
695 << isIonisation << G4endl;
696 G4cout << " theDEDX: " << theDEDXTable
697 << " theRange: " << theRangeTableForLoss
698 << " theInverse: " << theInverseRangeTable
699 << " theLambda: " << theLambdaTable << G4endl;
700 */
701 //if(1 < verboseLevel || verb) {
702 if(1 < verboseLevel) {
703 G4cout << "### G4VEnergyLossProcess::BuildPhysicsTable() done for "
704 << GetProcessName()
705 << " and particle " << part.GetParticleName();
706 if(isIonisation) { G4cout << " isIonisation flag = 1"; }
707 G4cout << G4endl;
708 }
709}
@ fIsSubIonisation
@ fIsIonisation
G4int NumberOfModels() const
void LocalPhysicsTables(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p)
void BuildPhysicsTable(const G4ParticleDefinition *aParticle)
G4VSubCutProducer * SubCutProducer()
G4VAtomDeexcitation * AtomDeexcitation()
void InitialiseHelper()
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
Definition: G4VEmModel.cc:220
G4PhysicsTable * RangeTableForLoss() const
G4PhysicsTable * InverseRangeTable() const
G4PhysicsTable * CSDARangeTable() const
G4VEmModel * GetModelByIndex(G4int idx=0, G4bool ver=false) const
void SetRangeTableForLoss(G4PhysicsTable *p)
G4PhysicsTable * DEDXTableForSubsec() const
G4PhysicsTable * IonisationTableForSubsec() const
G4PhysicsTable * SecondaryRangeTable() const
void SetInverseRangeTable(G4PhysicsTable *p)
G4PhysicsTable * SubLambdaTable() const
G4bool IsIonisationProcess() const
void SetDEDXTable(G4PhysicsTable *p, G4EmTableType tType)
void SetSecondaryRangeTable(G4PhysicsTable *p)
void SetSubLambdaTable(G4PhysicsTable *p)
void SetLambdaTable(G4PhysicsTable *p)
G4PhysicsTable * IonisationTable() const
G4PhysicsTable * LambdaTable() const
void SetCSDARangeTable(G4PhysicsTable *pRange)
G4PhysicsTable * DEDXunRestrictedTable() const
G4PhysicsTable * DEDXTable() const
const G4VProcess * GetMasterProcess() const
Definition: G4VProcess.hh:518

Referenced by G4ePolarizedIonisation::BuildPhysicsTable().

◆ ContinuousStepLimit()

G4double G4VEnergyLossProcess::ContinuousStepLimit ( const G4Track track,
G4double  previousStepSize,
G4double  currentMinimumStep,
G4double currentSafety 
)

Definition at line 1942 of file G4VEnergyLossProcess.cc.

1945{
1946 G4GPILSelection sel;
1947 return AlongStepGetPhysicalInteractionLength(track, x, y, z, &sel);
1948}
G4GPILSelection
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection) override

◆ CrossSectionBiasingFactor()

G4double G4VEnergyLossProcess::CrossSectionBiasingFactor ( ) const
inline

Definition at line 1046 of file G4VEnergyLossProcess.hh.

1047{
1048 return biasFactor;
1049}

◆ CrossSectionPerVolume() [1/2]

G4double G4VEnergyLossProcess::CrossSectionPerVolume ( G4double  kineticEnergy,
const G4MaterialCutsCouple couple 
)

◆ CrossSectionPerVolume() [2/2]

G4double G4VEnergyLossProcess::CrossSectionPerVolume ( G4double  kineticEnergy,
const G4MaterialCutsCouple couple,
G4double  logKineticEnergy 
)

Definition at line 1908 of file G4VEnergyLossProcess.cc.

1911{
1912 // Cross section per volume is calculated
1913 DefineMaterial(couple);
1914 G4double cross = 0.0;
1915 if (theLambdaTable) {
1916 cross = GetLambdaForScaledEnergy(kineticEnergy * massRatio,
1917 logKineticEnergy + logMassRatio);
1918 } else {
1919 SelectModel(kineticEnergy*massRatio);
1920 cross = biasFactor*(*theDensityFactor)[currentCoupleIndex]
1921 *(currentModel->CrossSectionPerVolume(currentMaterial,
1922 particle, kineticEnergy,
1923 (*theCuts)[currentCoupleIndex]));
1924 }
1925 return std::max(cross, 0.0);
1926}
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:254
void SelectModel(G4double kinEnergy)

◆ CSDARangeTable()

G4PhysicsTable * G4VEnergyLossProcess::CSDARangeTable ( ) const
inline

Definition at line 1095 of file G4VEnergyLossProcess.hh.

1096{
1097 return theCSDARangeTable;
1098}

Referenced by BuildPhysicsTable().

◆ CurrentMaterialCutsCoupleIndex()

size_t G4VEnergyLossProcess::CurrentMaterialCutsCoupleIndex ( ) const
inlineprotected

Definition at line 595 of file G4VEnergyLossProcess.hh.

596{
597 return currentCoupleIndex;
598}

◆ DEDXTable()

G4PhysicsTable * G4VEnergyLossProcess::DEDXTable ( ) const
inline

Definition at line 1060 of file G4VEnergyLossProcess.hh.

1061{
1062 return theDEDXTable;
1063}

Referenced by BuildPhysicsTable(), G4LossTableManager::LocalPhysicsTables(), and G4EmCalculator::PrintDEDXTable().

◆ DEDXTableForSubsec()

G4PhysicsTable * G4VEnergyLossProcess::DEDXTableForSubsec ( ) const
inline

Definition at line 1067 of file G4VEnergyLossProcess.hh.

1068{
1069 return theDEDXSubTable;
1070}

Referenced by BuildPhysicsTable().

◆ DEDXunRestrictedTable()

G4PhysicsTable * G4VEnergyLossProcess::DEDXunRestrictedTable ( ) const
inline

Definition at line 1074 of file G4VEnergyLossProcess.hh.

1075{
1076 return theDEDXunRestrictedTable;
1077}

Referenced by BuildPhysicsTable().

◆ EmModel()

◆ FluctModel()

◆ GetContinuousStepLimit()

G4double G4VEnergyLossProcess::GetContinuousStepLimit ( const G4Track track,
G4double  previousStepSize,
G4double  currentMinimumStep,
G4double currentSafety 
)
overrideprotectedvirtual

Implements G4VContinuousDiscreteProcess.

Definition at line 1964 of file G4VEnergyLossProcess.cc.

1967{
1968 return DBL_MAX;
1969}

◆ GetCSDARange()

G4double G4VEnergyLossProcess::GetCSDARange ( G4double  kineticEnergy,
const G4MaterialCutsCouple couple 
)
inline

Definition at line 871 of file G4VEnergyLossProcess.hh.

873{
874 DefineMaterial(couple);
875 return (theCSDARangeTable) ?
876 GetLimitScaledRangeForScaledEnergy(kineticEnergy*massRatio)*reduceFactor
877 : DBL_MAX;
878}

Referenced by G4LossTableManager::GetCSDARange().

◆ GetCurrentElement()

const G4Element * G4VEnergyLossProcess::GetCurrentElement ( ) const

Definition at line 2209 of file G4VEnergyLossProcess.cc.

2210{
2211 const G4Element* elm = nullptr;
2212 if(currentModel) { elm = currentModel->GetCurrentElement(); }
2213 return elm;
2214}
const G4Element * GetCurrentElement() const
Definition: G4VEmModel.hh:495

◆ GetDEDX() [1/2]

G4double G4VEnergyLossProcess::GetDEDX ( G4double  kineticEnergy,
const G4MaterialCutsCouple couple 
)
inline

Definition at line 809 of file G4VEnergyLossProcess.hh.

811{
812 DefineMaterial(couple);
813 return GetDEDXForScaledEnergy(kinEnergy*massRatio);
814}

Referenced by G4ContinuousGainOfEnergy::AlongStepDoIt(), G4LossTableManager::GetDEDX(), and G4VMscModel::GetDEDX().

◆ GetDEDX() [2/2]

G4double G4VEnergyLossProcess::GetDEDX ( G4double  kineticEnergy,
const G4MaterialCutsCouple couple,
G4double  logKineticEnergy 
)
inline

Definition at line 817 of file G4VEnergyLossProcess.hh.

820{
821 DefineMaterial(couple);
822 return GetDEDXForScaledEnergy(kinEnergy*massRatio, logKinEnergy+logMassRatio);
823}

◆ GetDEDXDispersion()

G4double G4VEnergyLossProcess::GetDEDXDispersion ( const G4MaterialCutsCouple couple,
const G4DynamicParticle dp,
G4double  length 
)

Definition at line 1889 of file G4VEnergyLossProcess.cc.

1893{
1894 DefineMaterial(couple);
1895 G4double ekin = dp->GetKineticEnergy();
1896 SelectModel(ekin*massRatio);
1897 G4double tmax = currentModel->MaxSecondaryKinEnergy(dp);
1898 tmax = std::min(tmax,(*theCuts)[currentCoupleIndex]);
1899 G4double d = 0.0;
1900 G4VEmFluctuationModel* fm = currentModel->GetModelOfFluctuations();
1901 if(fm) { d = fm->Dispersion(currentMaterial,dp,tmax,length); }
1902 return d;
1903}
G4double GetKineticEnergy() const
virtual G4double Dispersion(const G4Material *, const G4DynamicParticle *, G4double tmax, G4double length)=0

Referenced by G4LossTableManager::GetDEDXDispersion().

◆ GetDEDXForSubsec()

G4double G4VEnergyLossProcess::GetDEDXForSubsec ( G4double  kineticEnergy,
const G4MaterialCutsCouple couple 
)
inline

Definition at line 828 of file G4VEnergyLossProcess.hh.

830{
831 DefineMaterial(couple);
832 return GetSubDEDXForScaledEnergy(kineticEnergy*massRatio);
833}

Referenced by G4LossTableManager::GetSubDEDX().

◆ GetKineticEnergy()

G4double G4VEnergyLossProcess::GetKineticEnergy ( G4double  range,
const G4MaterialCutsCouple couple 
)
inline

Definition at line 906 of file G4VEnergyLossProcess.hh.

908{
909 DefineMaterial(couple);
910 return ScaledKinEnergyForLoss(range/reduceFactor)/massRatio;
911}

Referenced by G4ContinuousGainOfEnergy::AlongStepDoIt(), G4LossTableManager::GetEnergy(), and G4VMscModel::GetEnergy().

◆ GetLambda() [1/2]

G4double G4VEnergyLossProcess::GetLambda ( G4double  kineticEnergy,
const G4MaterialCutsCouple couple 
)
inline

Definition at line 916 of file G4VEnergyLossProcess.hh.

918{
919 DefineMaterial(couple);
920 return theLambdaTable ? GetLambdaForScaledEnergy(kinEnergy*massRatio) : 0.0;
921}

◆ GetLambda() [2/2]

G4double G4VEnergyLossProcess::GetLambda ( G4double  kineticEnergy,
const G4MaterialCutsCouple couple,
G4double  logKineticEnergy 
)
inline

Definition at line 924 of file G4VEnergyLossProcess.hh.

927{
928 DefineMaterial(couple);
929 return theLambdaTable
930 ? GetLambdaForScaledEnergy(kinEnergy*massRatio, logKinEnergy+logMassRatio)
931 : 0.0;
932}

◆ GetMeanFreePath()

G4double G4VEnergyLossProcess::GetMeanFreePath ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
overrideprotectedvirtual

Implements G4VContinuousDiscreteProcess.

Reimplemented in G4ePolarizedIonisation.

Definition at line 1952 of file G4VEnergyLossProcess.cc.

1957{
1959 return MeanFreePath(track);
1960}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced
G4double MeanFreePath(const G4Track &track)

Referenced by G4ePolarizedIonisation::GetMeanFreePath().

◆ GetModelByIndex()

G4VEmModel * G4VEnergyLossProcess::GetModelByIndex ( G4int  idx = 0,
G4bool  ver = false 
) const

Definition at line 345 of file G4VEnergyLossProcess.cc.

346{
347 return modelManager->GetModel(idx, ver);
348}
G4VEmModel * GetModel(G4int idx, G4bool ver=false)

Referenced by BuildPhysicsTable().

◆ GetRange() [1/2]

G4double G4VEnergyLossProcess::GetRange ( G4double  kineticEnergy,
const G4MaterialCutsCouple couple 
)
inline

Definition at line 838 of file G4VEnergyLossProcess.hh.

840{
841 G4double x = fRange;
842 DefineMaterial(couple);
843 if(theCSDARangeTable) {
844 x = reduceFactor * GetLimitScaledRangeForScaledEnergy(kinEnergy*massRatio);
845 } else if(theRangeTableForLoss) {
846 x = reduceFactor * GetScaledRangeForScaledEnergy(kinEnergy*massRatio);
847 }
848 return x;
849}

Referenced by G4ContinuousGainOfEnergy::AlongStepDoIt(), G4ContinuousGainOfEnergy::GetContinuousStepLimit(), and G4LossTableManager::GetRange().

◆ GetRange() [2/2]

G4double G4VEnergyLossProcess::GetRange ( G4double  kineticEnergy,
const G4MaterialCutsCouple couple,
G4double  logKineticEnergy 
)
inline

Definition at line 852 of file G4VEnergyLossProcess.hh.

855{
856 G4double x = fRange;
857 DefineMaterial(couple);
858 if(theCSDARangeTable) {
859 x = reduceFactor * GetLimitScaledRangeForScaledEnergy(kinEnergy*massRatio,
860 logKinEnergy+logMassRatio);
861 } else if(theRangeTableForLoss) {
862 x = reduceFactor * GetScaledRangeForScaledEnergy(kinEnergy*massRatio,
863 logKinEnergy+logMassRatio);
864 }
865 return x;
866}

◆ GetRangeForLoss() [1/2]

G4double G4VEnergyLossProcess::GetRangeForLoss ( G4double  kineticEnergy,
const G4MaterialCutsCouple couple 
)
inline

Definition at line 883 of file G4VEnergyLossProcess.hh.

885{
886 // G4cout << "GetRangeForLoss: Range from " << GetProcessName() << G4endl;
887 DefineMaterial(couple);
888 return reduceFactor * GetScaledRangeForScaledEnergy(kinEnergy*massRatio);
889}

Referenced by G4VMscModel::GetRange(), and G4LossTableManager::GetRangeFromRestricteDEDX().

◆ GetRangeForLoss() [2/2]

G4double G4VEnergyLossProcess::GetRangeForLoss ( G4double  kineticEnergy,
const G4MaterialCutsCouple couple,
G4double  logKineticEnergy 
)
inline

Definition at line 893 of file G4VEnergyLossProcess.hh.

896{
897 // G4cout << "GetRangeForLoss: Range from " << GetProcessName() << G4endl;
898 DefineMaterial(couple);
899 return reduceFactor * GetScaledRangeForScaledEnergy(kinEnergy*massRatio,
900 logKinEnergy+logMassRatio);
901}

◆ InitialiseEnergyLossProcess()

◆ InverseRangeTable()

G4PhysicsTable * G4VEnergyLossProcess::InverseRangeTable ( ) const
inline

Definition at line 1116 of file G4VEnergyLossProcess.hh.

1117{
1118 return theInverseRangeTable;
1119}

Referenced by BuildPhysicsTable(), G4LossTableManager::LocalPhysicsTables(), and G4EmCalculator::PrintInverseRangeTable().

◆ IonisationTable()

G4PhysicsTable * G4VEnergyLossProcess::IonisationTable ( ) const
inline

Definition at line 1081 of file G4VEnergyLossProcess.hh.

1082{
1083 return theIonisationTable;
1084}

Referenced by BuildPhysicsTable().

◆ IonisationTableForSubsec()

G4PhysicsTable * G4VEnergyLossProcess::IonisationTableForSubsec ( ) const
inline

Definition at line 1088 of file G4VEnergyLossProcess.hh.

1089{
1090 return theIonisationSubTable;
1091}

Referenced by BuildPhysicsTable().

◆ IsApplicable()

◆ IsIntegral()

G4bool G4VEnergyLossProcess::IsIntegral ( ) const
inline

Definition at line 1011 of file G4VEnergyLossProcess.hh.

1012{
1013 return integral;
1014}

◆ IsIonisationProcess()

G4bool G4VEnergyLossProcess::IsIonisationProcess ( ) const
inline

Definition at line 1018 of file G4VEnergyLossProcess.hh.

1019{
1020 return isIonisation;
1021}

Referenced by BuildPhysicsTable(), G4LossTableManager::BuildPhysicsTable(), and G4LossTableManager::LocalPhysicsTables().

◆ LambdaPhysicsVector()

G4PhysicsVector * G4VEnergyLossProcess::LambdaPhysicsVector ( const G4MaterialCutsCouple ,
G4double  cut 
)
protected

Definition at line 1974 of file G4VEnergyLossProcess.cc.

1976{
1977 G4PhysicsVector* v =
1978 new G4PhysicsLogVector(minKinEnergy, maxKinEnergy, nBins);
1979 v->SetSpline(theParameters->Spline());
1980 return v;
1981}

◆ LambdaTable()

G4PhysicsTable * G4VEnergyLossProcess::LambdaTable ( ) const
inline

Definition at line 1123 of file G4VEnergyLossProcess.hh.

1124{
1125 return theLambdaTable;
1126}

Referenced by BuildPhysicsTable().

◆ MaxKinEnergy()

G4double G4VEnergyLossProcess::MaxKinEnergy ( ) const
inline

Definition at line 1039 of file G4VEnergyLossProcess.hh.

1040{
1041 return maxKinEnergy;
1042}

◆ MeanFreePath()

G4double G4VEnergyLossProcess::MeanFreePath ( const G4Track track)

Definition at line 1930 of file G4VEnergyLossProcess.cc.

1931{
1932 DefineMaterial(track.GetMaterialCutsCouple());
1933 const G4double kinEnergy = track.GetKineticEnergy();
1934 const G4double logKinEnergy = track.GetDynamicParticle()->GetLogKineticEnergy();
1935 const G4double cs = GetLambdaForScaledEnergy(kinEnergy * massRatio,
1936 logKinEnergy + logMassRatio);
1937 return (0.0 < cs) ? 1.0/cs : DBL_MAX;
1938}
G4double GetLogKineticEnergy() const
G4double GetKineticEnergy() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const

Referenced by GetMeanFreePath().

◆ MinKinEnergy()

G4double G4VEnergyLossProcess::MinKinEnergy ( ) const
inline

Definition at line 1032 of file G4VEnergyLossProcess.hh.

1033{
1034 return minKinEnergy;
1035}

◆ MinPrimaryEnergy()

G4double G4VEnergyLossProcess::MinPrimaryEnergy ( const G4ParticleDefinition ,
const G4Material ,
G4double  cut 
)
protectedvirtual

◆ NumberOfModels()

G4int G4VEnergyLossProcess::NumberOfModels ( ) const

Definition at line 352 of file G4VEnergyLossProcess.cc.

353{
354 return modelManager->NumberOfModels();
355}

◆ NumberOfSubCutoffRegions()

G4int G4VEnergyLossProcess::NumberOfSubCutoffRegions ( ) const
inline

Definition at line 1025 of file G4VEnergyLossProcess.hh.

1026{
1027 return nSCoffRegions;
1028}

◆ Particle()

const G4ParticleDefinition * G4VEnergyLossProcess::Particle ( ) const
inline

Definition at line 973 of file G4VEnergyLossProcess.hh.

974{
975 return particle;
976}

Referenced by G4LossTableManager::BuildPhysicsTable(), and G4LossTableManager::LocalPhysicsTables().

◆ PostStepDoIt()

G4VParticleChange * G4VEnergyLossProcess::PostStepDoIt ( const G4Track track,
const G4Step step 
)
overridevirtual

Reimplemented from G4VContinuousDiscreteProcess.

Definition at line 1557 of file G4VEnergyLossProcess.cc.

1559{
1560 // In all cases clear number of interaction lengths
1563
1565 G4double finalT = track.GetKineticEnergy();
1566 if(finalT <= lowestKinEnergy) { return &fParticleChange; }
1567
1568 G4double postStepScaledEnergy = finalT*massRatio;
1569 SelectModel(postStepScaledEnergy);
1570
1571 if(!currentModel->IsActive(postStepScaledEnergy)) {
1572 return &fParticleChange;
1573 }
1574 /*
1575 if(-1 < verboseLevel) {
1576 G4cout << GetProcessName()
1577 << "::PostStepDoIt: E(MeV)= " << finalT/MeV
1578 << G4endl;
1579 }
1580 */
1581
1582 // forced process - should happen only once per track
1583 if(biasFlag) {
1585 biasFlag = false;
1586 }
1587 }
1588
1589 const G4DynamicParticle* dp = track.GetDynamicParticle();
1590 const G4double logFinalT = dp->GetLogKineticEnergy();
1591 // postStepLogScaledEnergy = logFinalT + logMassRatio;
1592
1593 // Integral approach
1594 if (integral) {
1595 const G4double lx = GetLambdaForScaledEnergy(postStepScaledEnergy,
1596 logFinalT + logMassRatio);
1597 /*
1598 if(preStepLambda<lx && 1 < verboseLevel) {
1599 G4cout << "WARNING: for " << particle->GetParticleName()
1600 << " and " << GetProcessName()
1601 << " E(MeV)= " << finalT/MeV
1602 << " preLambda= " << preStepLambda
1603 << " < " << lx << " (postLambda) "
1604 << G4endl;
1605 }
1606 */
1607 if(lx <= 0.0 || preStepLambda*G4UniformRand() > lx) {
1608 return &fParticleChange;
1609 }
1610 }
1611
1612 SelectModel(postStepScaledEnergy);
1613
1614 // define new weight for primary and secondaries
1616 if(weightFlag) {
1617 weight /= biasFactor;
1619 }
1620
1621 G4double tcut = (*theCuts)[currentCoupleIndex];
1622
1623 // sample secondaries
1624 secParticles.clear();
1625 //G4cout<< "@@@ Eprimary= "<<dynParticle->GetKineticEnergy()/MeV
1626 // << " cut= " << tcut/MeV << G4endl;
1627 currentModel->SampleSecondaries(&secParticles, currentCouple, dp, tcut);
1628
1629 G4int num0 = secParticles.size();
1630
1631 // bremsstrahlung splitting or Russian roulette
1632 if(biasManager) {
1633 if(biasManager->SecondaryBiasingRegion(currentCoupleIndex)) {
1634 G4double eloss = 0.0;
1635 weight *= biasManager->ApplySecondaryBiasing(
1636 secParticles,
1637 track, currentModel,
1638 &fParticleChange, eloss,
1639 currentCoupleIndex, tcut,
1640 step.GetPostStepPoint()->GetSafety());
1641 if(eloss > 0.0) {
1644 }
1645 }
1646 }
1647
1648 // save secondaries
1649 G4int num = secParticles.size();
1650 if(num > 0) {
1651
1653 G4double time = track.GetGlobalTime();
1654
1655 for (G4int i=0; i<num; ++i) {
1656 if(secParticles[i]) {
1657 G4Track* t = new G4Track(secParticles[i], time, track.GetPosition());
1659 if (biasManager) {
1660 t->SetWeight(weight * biasManager->GetWeight(i));
1661 } else {
1662 t->SetWeight(weight);
1663 }
1664 if(i < num0) { t->SetCreatorModelIndex(secID); }
1665 else { t->SetCreatorModelIndex(biasID); }
1666
1667 //G4cout << "Secondary(post step) has weight " << t->GetWeight()
1668 // << ", kenergy " << t->GetKineticEnergy()/MeV << " MeV"
1669 // << " time= " << time/ns << " ns " << G4endl;
1671 }
1672 }
1673 }
1674
1677 if(particle->GetProcessManager()->GetAtRestProcessVector()->size() > 0)
1680 }
1681
1682 /*
1683 if(-1 < verboseLevel) {
1684 G4cout << "::PostStepDoIt: Sample secondary; Efin= "
1685 << fParticleChange.GetProposedKineticEnergy()/MeV
1686 << " MeV; model= (" << currentModel->LowEnergyLimit()
1687 << ", " << currentModel->HighEnergyLimit() << ")"
1688 << " preStepLambda= " << preStepLambda
1689 << " dir= " << track.GetMomentumDirection()
1690 << " status= " << track.GetTrackStatus()
1691 << G4endl;
1692 }
1693 */
1694 return &fParticleChange;
1695}
@ fAlive
@ fStopAndKill
@ fStopButAlive
#define G4UniformRand()
Definition: Randomize.hh:52
G4bool ForcedInteractionRegion(G4int coupleIdx)
G4double ApplySecondaryBiasing(std::vector< G4DynamicParticle * > &, const G4Track &track, G4VEmModel *currentModel, G4ParticleChangeForGamma *pParticleChange, G4double &eloss, G4int coupleIdx, G4double tcut, G4double safety=0.0)
G4double GetWeight(G4int i)
G4bool SecondaryBiasingRegion(G4int coupleIdx)
void InitializeForPostStep(const G4Track &)
G4double GetProposedKineticEnergy() const
G4ProcessManager * GetProcessManager() const
G4ProcessVector * GetAtRestProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
std::size_t size() const
void SetWeight(G4double aValue)
const G4ThreeVector & GetPosition() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4double GetGlobalTime() const
const G4TouchableHandle & GetTouchableHandle() const
void SetCreatorModelIndex(G4int idx)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
void ProposeTrackStatus(G4TrackStatus status)
G4double GetLocalEnergyDeposit() const
void AddSecondary(G4Track *aSecondary)
void SetNumberOfSecondaries(G4int totSecondaries)
G4TrackStatus GetTrackStatus() const
G4double currentInteractionLength
Definition: G4VProcess.hh:335
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:331

◆ PostStepGetPhysicalInteractionLength()

G4double G4VEnergyLossProcess::PostStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
overridevirtual

Reimplemented from G4VContinuousDiscreteProcess.

Reimplemented in G4ePolarizedIonisation.

Definition at line 1054 of file G4VEnergyLossProcess.cc.

1058{
1059 // condition is set to "Not Forced"
1061 G4double x = DBL_MAX;
1062
1063 // initialisation of material, mass, charge, model
1064 // at the beginning of the step
1065 DefineMaterial(track.GetMaterialCutsCouple());
1071
1072 if(!currentModel->IsActive(preStepScaledEnergy)) {
1075 return x;
1076 }
1077
1078 // change effective charge of an ion on fly
1079 if(isIon) {
1080 G4double q2 = currentModel->ChargeSquareRatio(track);
1081 if(q2 != chargeSqRatio && q2 > 0.0) {
1082 chargeSqRatio = q2;
1083 fFactor = q2*biasFactor*(*theDensityFactor)[currentCoupleIndex];
1084 reduceFactor = 1.0/(fFactor*massRatio);
1085 }
1086 }
1087 // if(particle->GetPDGMass() > 0.9*GeV)
1088 //G4cout << "q2= "<<chargeSqRatio << " massRatio= " << massRatio << G4endl;
1089
1090 // forced biasing only for primary particles
1091 if(biasManager) {
1092 if(0 == track.GetParentID() && biasFlag &&
1094 return biasManager->GetStepLimit(currentCoupleIndex, previousStepSize);
1095 }
1096 }
1097
1098 // compute mean free path
1100 if (integral) {
1101 ComputeLambdaForScaledEnergy(preStepScaledEnergy, preStepLogScaledEnergy);
1102 } else {
1103 preStepLambda =
1104 GetLambdaForScaledEnergy(preStepScaledEnergy, preStepLogScaledEnergy);
1105 }
1106
1107 // zero cross section
1108 if(preStepLambda <= 0.0) {
1111 }
1112 }
1113
1114 // non-zero cross section
1115 if(preStepLambda > 0.0) {
1117
1118 // beggining of tracking (or just after DoIt of this process)
1121
1122 } else if(currentInteractionLength < DBL_MAX) {
1123
1124 // subtract NumberOfInteractionLengthLeft using previous step
1126 previousStepSize/currentInteractionLength;
1127
1129 std::max(theNumberOfInteractionLengthLeft, 0.0);
1130 }
1131
1132 // new mean free path and step limit
1135 }
1136#ifdef G4VERBOSE
1137 if (verboseLevel>2){
1138 // if(particle->GetPDGMass() > 0.9*GeV){
1139 G4cout << "G4VEnergyLossProcess::PostStepGetPhysicalInteractionLength ";
1140 G4cout << "[ " << GetProcessName() << "]" << G4endl;
1141 G4cout << " for " << track.GetDefinition()->GetParticleName()
1142 << " in Material " << currentMaterial->GetName()
1143 << " Ekin(MeV)= " << preStepKinEnergy/MeV
1144 << " " << track.GetMaterial()->GetName()
1145 <<G4endl;
1146 G4cout << "MeanFreePath = " << currentInteractionLength/cm << "[cm]"
1147 << "InteractionLength= " << x/cm <<"[cm] " <<G4endl;
1148 }
1149#endif
1150 return x;
1151}
G4double GetStepLimit(G4int coupleIdx, G4double previousStep)
const G4String & GetName() const
Definition: G4Material.hh:175
G4Material * GetMaterial() const
G4ParticleDefinition * GetDefinition() const
G4int GetParentID() const
virtual G4double ChargeSquareRatio(const G4Track &)
Definition: G4VEmModel.cc:383
G4double theInitialNumberOfInteractionLength
Definition: G4VProcess.hh:338

Referenced by G4ePolarizedIonisation::PostStepGetPhysicalInteractionLength().

◆ PreparePhysicsTable()

void G4VEnergyLossProcess::PreparePhysicsTable ( const G4ParticleDefinition part)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 360 of file G4VEnergyLossProcess.cc.

361{
362 if(1 < verboseLevel) {
363 G4cout << "G4VEnergyLossProcess::PreparePhysicsTable for "
364 << GetProcessName() << " for " << part.GetParticleName()
365 << " " << this << G4endl;
366 }
367 isMaster = lManager->IsMaster();
368
369 currentCouple = nullptr;
370 preStepLambda = 0.0;
372 fRange = DBL_MAX;
373 preStepKinEnergy = 0.0;
375 preStepRangeEnergy = 0.0;
376 chargeSqRatio = 1.0;
377 massRatio = 1.0;
378 logMassRatio = 0.;
379 reduceFactor = 1.0;
380 fFactor = 1.0;
381 lastIdx = 0;
382
383 // Are particle defined?
384 if( !particle ) { particle = &part; }
385
386 if(part.GetParticleType() == "nucleus") {
387
388 G4String pname = part.GetParticleName();
389 if(pname != "deuteron" && pname != "triton" &&
390 pname != "alpha+" && pname != "helium" &&
391 pname != "hydrogen") {
392
393 if(!theGenericIon) {
394 theGenericIon =
396 }
397 isIon = true;
398 if(theGenericIon && particle != theGenericIon) {
399 G4ProcessManager* pm = theGenericIon->GetProcessManager();
401 size_t n = v->size();
402 for(size_t j=0; j<n; ++j) {
403 if((*v)[j] == this) {
404 particle = theGenericIon;
405 break;
406 }
407 }
408 }
409 }
410 }
411
412 if( particle != &part ) {
413 if(!isIon) {
414 lManager->RegisterExtraParticle(&part, this);
415 }
416 if(1 < verboseLevel) {
417 G4cout << "### G4VEnergyLossProcess::PreparePhysicsTable()"
418 << " interrupted for "
419 << part.GetParticleName() << " isIon= " << isIon
420 << " particle " << particle << " GenericIon " << theGenericIon
421 << G4endl;
422 }
423 return;
424 }
425
426 Clean();
427 lManager->PreparePhysicsTable(&part, this, isMaster);
428 G4LossTableBuilder* bld = lManager->GetTableBuilder();
429
430 // Base particle and set of models can be defined here
431 InitialiseEnergyLossProcess(particle, baseParticle);
432
433 const G4ProductionCutsTable* theCoupleTable=
435 size_t n = theCoupleTable->GetTableSize();
436
437 theDEDXAtMaxEnergy.resize(n, 0.0);
438 theRangeAtMaxEnergy.resize(n, 0.0);
439 theEnergyOfCrossSectionMax.resize(n, 0.0);
440 theCrossSectionMax.resize(n, DBL_MAX);
441
442 // parameters of the process
443 if(!actIntegral) { integral = theParameters->Integral(); }
444 if(!actLossFluc) { lossFluctuationFlag = theParameters->LossFluctuation(); }
445 rndmStepFlag = theParameters->UseCutAsFinalRange();
446 if(!actMinKinEnergy) { minKinEnergy = theParameters->MinKinEnergy(); }
447 if(!actMaxKinEnergy) { maxKinEnergy = theParameters->MaxKinEnergy(); }
448 if(!actBinning) {
449 nBins = theParameters->NumberOfBinsPerDecade()
450 *G4lrint(std::log10(maxKinEnergy/minKinEnergy));
451 }
452 maxKinEnergyCSDA = theParameters->MaxEnergyForCSDARange();
453 nBinsCSDA = theParameters->NumberOfBinsPerDecade()
454 *G4lrint(std::log10(maxKinEnergyCSDA/minKinEnergy));
455 if(!actLinLossLimit) { linLossLimit = theParameters->LinearLossLimit(); }
456 lambdaFactor = theParameters->LambdaFactor();
457 logLambdafactor = G4Log(lambdaFactor);
458 if(isMaster) { SetVerboseLevel(theParameters->Verbose()); }
459 else { SetVerboseLevel(theParameters->WorkerVerbose()); }
460
461 theParameters->DefineRegParamForLoss(this);
462
463 G4double initialCharge = particle->GetPDGCharge();
464 G4double initialMass = particle->GetPDGMass();
465
466 theParameters->FillStepFunction(particle, this);
467
468 if (baseParticle) {
469 massRatio = (baseParticle->GetPDGMass())/initialMass;
470 logMassRatio = G4Log(massRatio);
471 G4double q = initialCharge/baseParticle->GetPDGCharge();
472 chargeSqRatio = q*q;
473 if(chargeSqRatio > 0.0) { reduceFactor = 1.0/(chargeSqRatio*massRatio); }
474 }
475 lowestKinEnergy = (initialMass < MeV) ? theParameters->LowestElectronEnergy()
476 : theParameters->LowestMuHadEnergy();
477
478 // Tables preparation
479 if (isMaster && !baseParticle) {
480
481 if(theDEDXTable && isIonisation) {
482 if(theIonisationTable && theDEDXTable != theIonisationTable) {
483 theDEDXTable->clearAndDestroy();
484 delete theDEDXTable;
485 theDEDXTable = theIonisationTable;
486 }
487 if(theDEDXSubTable && theIonisationSubTable &&
488 theDEDXSubTable != theIonisationSubTable) {
489 theDEDXSubTable->clearAndDestroy();
490 delete theDEDXSubTable;
491 theDEDXSubTable = theIonisationSubTable;
492 }
493 }
494
495 theDEDXTable = G4PhysicsTableHelper::PreparePhysicsTable(theDEDXTable);
496 bld->InitialiseBaseMaterials(theDEDXTable);
497
498 if(theDEDXSubTable) {
499 theDEDXSubTable =
501 }
502
503 if (theParameters->BuildCSDARange()) {
504 theDEDXunRestrictedTable =
505 G4PhysicsTableHelper::PreparePhysicsTable(theDEDXunRestrictedTable);
506 theCSDARangeTable =
508 }
509
510 theLambdaTable = G4PhysicsTableHelper::PreparePhysicsTable(theLambdaTable);
511
512 if(isIonisation) {
513 theRangeTableForLoss =
514 G4PhysicsTableHelper::PreparePhysicsTable(theRangeTableForLoss);
515 theInverseRangeTable =
516 G4PhysicsTableHelper::PreparePhysicsTable(theInverseRangeTable);
517 }
518
519 if (nSCoffRegions && !lManager->SubCutProducer()) {
520 theDEDXSubTable =
522 theSubLambdaTable =
524 }
525 }
526 /*
527 G4cout << "** G4VEnergyLossProcess::PreparePhysicsTable() for "
528 << GetProcessName() << " and " << particle->GetParticleName()
529 << " isMaster: " << isMaster << " isIonisation: "
530 << isIonisation << G4endl;
531 G4cout << " theDEDX: " << theDEDXTable
532 << " theRange: " << theRangeTableForLoss
533 << " theInverse: " << theInverseRangeTable
534 << " theLambda: " << theLambdaTable << G4endl;
535 */
536 // forced biasing
537 if(biasManager) {
538 biasManager->Initialise(part,GetProcessName(),verboseLevel);
539 biasFlag = false;
540 }
541
542 // defined ID of secondary particles
543 if(isMaster) {
544 G4String nam1 = GetProcessName();
545 G4String nam4 = nam1 + "_split";
546 G4String nam5 = nam1 + "_subcut";
547 secID = G4PhysicsModelCatalog::Register(nam1);
548 biasID = G4PhysicsModelCatalog::Register(nam4);
549 subsecID= G4PhysicsModelCatalog::Register(nam5);
550 }
551
552 // initialisation of models
553 G4int nmod = modelManager->NumberOfModels();
554 for(G4int i=0; i<nmod; ++i) {
555 G4VEmModel* mod = modelManager->GetModel(i);
556 mod->SetMasterThread(isMaster);
558 theParameters->UseAngularGeneratorForIonisation());
559 if(mod->HighEnergyLimit() > maxKinEnergy) {
560 mod->SetHighEnergyLimit(maxKinEnergy);
561 }
562 }
563 theCuts = modelManager->Initialise(particle, secondaryParticle,
564 theParameters->MinSubRange(),
566
567 // Sub Cutoff
568 if(nSCoffRegions > 0) {
569 if(theParameters->MinSubRange() < 1.0) { useSubCutoff = true; }
570
571 theSubCuts = modelManager->SubCutoff();
572
573 idxSCoffRegions = new G4bool[n];
574 for (size_t j=0; j<n; ++j) {
575
576 const G4MaterialCutsCouple* couple =
577 theCoupleTable->GetMaterialCutsCouple(j);
578 const G4ProductionCuts* pcuts = couple->GetProductionCuts();
579
580 G4bool reg = false;
581 for(G4int i=0; i<nSCoffRegions; ++i) {
582 if( pcuts == scoffRegions[i]->GetProductionCuts()) {
583 reg = true;
584 break;
585 }
586 }
587 idxSCoffRegions[j] = reg;
588 }
589 }
590
591 if(1 < verboseLevel) {
592 G4cout << "G4VEnergyLossProcess::PrepearPhysicsTable() is done "
593 << " for local " << particle->GetParticleName()
594 << " isIon= " << isIon;
595 if(baseParticle) {
596 G4cout << "; base: " << baseParticle->GetParticleName();
597 }
598 G4cout << " chargeSqRatio= " << chargeSqRatio
599 << " massRatio= " << massRatio
600 << " reduceFactor= " << reduceFactor << G4endl;
601 if (nSCoffRegions) {
602 G4cout << " SubCutoff Regime is ON for regions: " << G4endl;
603 for (G4int i=0; i<nSCoffRegions; ++i) {
604 const G4Region* r = scoffRegions[i];
605 G4cout << " " << r->GetName() << G4endl;
606 }
607 }
608 }
609}
void Initialise(const G4ParticleDefinition &part, const G4String &procName, G4int verbose)
const G4DataVector * Initialise(const G4ParticleDefinition *part, const G4ParticleDefinition *secPart, G4double minSubRange, G4int verb)
const G4DataVector * SubCutoff() const
void DefineRegParamForLoss(G4VEnergyLossProcess *) const
void FillStepFunction(const G4ParticleDefinition *, G4VEnergyLossProcess *) const
G4double MinKinEnergy() const
G4int NumberOfBinsPerDecade() const
G4bool BuildCSDARange() const
G4bool LossFluctuation() const
G4bool UseCutAsFinalRange() const
G4int Verbose() const
G4double MinSubRange() const
G4int WorkerVerbose() const
G4double MaxKinEnergy() const
G4bool Integral() const
G4double MaxEnergyForCSDARange() const
G4bool UseAngularGeneratorForIonisation() const
G4double LinearLossLimit() const
G4double LowestMuHadEnergy() const
G4double LambdaFactor() const
void InitialiseBaseMaterials(const G4PhysicsTable *table=nullptr)
void PreparePhysicsTable(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p, G4bool theMaster)
void RegisterExtraParticle(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p)
const G4String & GetParticleType() const
G4double GetPDGCharge() const
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
static G4int Register(const G4String &)
static G4PhysicsTable * PreparePhysicsTable(G4PhysicsTable *physTable)
G4ProcessVector * GetAlongStepProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
const G4String & GetName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:757
void SetMasterThread(G4bool val)
Definition: G4VEmModel.hh:729
void SetAngularGeneratorFlag(G4bool)
Definition: G4VEmModel.hh:715
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:645
virtual void InitialiseEnergyLossProcess(const G4ParticleDefinition *, const G4ParticleDefinition *)=0

◆ PrintInfo()

virtual void G4VEnergyLossProcess::PrintInfo ( )
inlinevirtual

◆ ProcessDescription()

◆ RangeTableForLoss()

G4PhysicsTable * G4VEnergyLossProcess::RangeTableForLoss ( ) const
inline

Definition at line 1109 of file G4VEnergyLossProcess.hh.

1110{
1111 return theRangeTableForLoss;
1112}

Referenced by BuildPhysicsTable(), G4LossTableManager::LocalPhysicsTables(), and G4EmCalculator::PrintRangeTable().

◆ RetrievePhysicsTable()

G4bool G4VEnergyLossProcess::RetrievePhysicsTable ( const G4ParticleDefinition part,
const G4String directory,
G4bool  ascii 
)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 1762 of file G4VEnergyLossProcess.cc.

1765{
1766 G4bool res = true;
1767 if (!isMaster) return res;
1768 const G4String particleName = part->GetParticleName();
1769
1770 if(1 < verboseLevel) {
1771 G4cout << "G4VEnergyLossProcess::RetrievePhysicsTable() for "
1772 << particleName << " and process " << GetProcessName()
1773 << "; tables_are_built= " << tablesAreBuilt
1774 << G4endl;
1775 }
1776 if(particle == part) {
1777
1778 if ( !baseParticle ) {
1779
1780 G4bool fpi = true;
1781 if(!RetrieveTable(part,theDEDXTable,ascii,directory,"DEDX",fpi))
1782 {fpi = false;}
1783
1784 // ionisation table keeps individual dEdx and not sum of sub-processes
1785 if(!RetrieveTable(part,theDEDXTable,ascii,directory,"Ionisation",false))
1786 {fpi = false;}
1787
1788 if(!RetrieveTable(part,theRangeTableForLoss,ascii,directory,"Range",fpi))
1789 {res = false;}
1790
1791 if(!RetrieveTable(part,theDEDXunRestrictedTable,ascii,directory,
1792 "DEDXnr",false))
1793 {res = false;}
1794
1795 if(!RetrieveTable(part,theCSDARangeTable,ascii,directory,
1796 "CSDARange",false))
1797 {res = false;}
1798
1799 if(!RetrieveTable(part,theInverseRangeTable,ascii,directory,
1800 "InverseRange",fpi))
1801 {res = false;}
1802
1803 if(!RetrieveTable(part,theLambdaTable,ascii,directory,"Lambda",true))
1804 {res = false;}
1805
1806 G4bool yes = false;
1807 if(nSCoffRegions > 0) {yes = true;}
1808
1809 if(!RetrieveTable(part,theDEDXSubTable,ascii,directory,"SubDEDX",yes))
1810 {res = false;}
1811
1812 if(!RetrieveTable(part,theSubLambdaTable,ascii,directory,
1813 "SubLambda",yes))
1814 {res = false;}
1815
1816 if(!fpi) yes = false;
1817 if(!RetrieveTable(part,theIonisationSubTable,ascii,directory,
1818 "SubIonisation",yes))
1819 {res = false;}
1820 }
1821 }
1822
1823 return res;
1824}

◆ SampleSubCutSecondaries()

G4double G4VEnergyLossProcess::SampleSubCutSecondaries ( std::vector< G4Track * > &  tracks,
const G4Step step,
G4VEmModel model,
G4int  matIdx 
)

Definition at line 1479 of file G4VEnergyLossProcess.cc.

1483{
1484 // Fast check weather subcutoff can work
1485 G4double esec = 0.0;
1486 G4double subcut = (*theSubCuts)[idx];
1487 G4double cut = (*theCuts)[idx];
1488 if(cut <= subcut) { return esec; }
1489
1490 const G4Track* track = step.GetTrack();
1491 const G4DynamicParticle* dp = track->GetDynamicParticle();
1492 G4double e = dp->GetKineticEnergy()*massRatio;
1493 G4double cross = (*theDensityFactor)[idx]*chargeSqRatio
1494 *(((*theSubLambdaTable)[(*theDensityIdx)[idx]])->Value(e, idxSubLambda));
1495 G4double length = step.GetStepLength();
1496
1497 // negligible probability to get any interaction
1498 if(length*cross < perMillion) { return esec; }
1499 /*
1500 if(-1 < verboseLevel)
1501 G4cout << "<<< Subcutoff for " << GetProcessName()
1502 << " cross(1/mm)= " << cross*mm << ">>>"
1503 << " e(MeV)= " << preStepScaledEnergy
1504 << " matIdx= " << currentCoupleIndex
1505 << G4endl;
1506 */
1507
1508 // Sample subcutoff secondaries
1509 G4StepPoint* preStepPoint = step.GetPreStepPoint();
1510 G4StepPoint* postStepPoint = step.GetPostStepPoint();
1511 G4ThreeVector prepoint = preStepPoint->GetPosition();
1512 G4ThreeVector dr = postStepPoint->GetPosition() - prepoint;
1513 G4double pretime = preStepPoint->GetGlobalTime();
1514 G4double dt = postStepPoint->GetGlobalTime() - pretime;
1515 G4double fragment = 0.0;
1516
1517 do {
1518 G4double del = -G4Log(G4UniformRand())/cross;
1519 fragment += del/length;
1520 if (fragment > 1.0) { break; }
1521
1522 // sample secondaries
1523 secParticles.clear();
1524 model->SampleSecondaries(&secParticles,track->GetMaterialCutsCouple(),
1525 dp,subcut,cut);
1526
1527 // position of subcutoff particles
1528 G4ThreeVector r = prepoint + fragment*dr;
1529 std::vector<G4DynamicParticle*>::iterator it;
1530 for(it=secParticles.begin(); it!=secParticles.end(); ++it) {
1531
1532 G4Track* t = new G4Track((*it), pretime + fragment*dt, r);
1534 t->SetCreatorModelIndex(subsecID);
1535 tracks.push_back(t);
1536 esec += t->GetKineticEnergy();
1537 if (t->GetParticleDefinition() == thePositron) {
1538 esec += 2.0*electron_mass_c2;
1539 }
1540
1541 /*
1542 if(-1 < verboseLevel)
1543 G4cout << "New track "
1544 << t->GetParticleDefinition()->GetParticleName()
1545 << " e(keV)= " << t->GetKineticEnergy()/keV
1546 << " fragment= " << fragment
1547 << G4endl;
1548 */
1549 }
1550 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
1551 } while (fragment <= 1.0);
1552 return esec;
1553}
G4double GetGlobalTime() const
G4Track * GetTrack() const

Referenced by AlongStepDoIt().

◆ SecondaryParticle()

const G4ParticleDefinition * G4VEnergyLossProcess::SecondaryParticle ( ) const
inline

Definition at line 988 of file G4VEnergyLossProcess.hh.

989{
990 return secondaryParticle;
991}

◆ SecondaryRangeTable()

G4PhysicsTable * G4VEnergyLossProcess::SecondaryRangeTable ( ) const
inline

Definition at line 1102 of file G4VEnergyLossProcess.hh.

1103{
1104 return theSecondaryRangeTable;
1105}

Referenced by BuildPhysicsTable().

◆ SelectModel()

void G4VEnergyLossProcess::SelectModel ( G4double  kinEnergy)
inlineprotected

Definition at line 602 of file G4VEnergyLossProcess.hh.

603{
604 currentModel = modelManager->SelectModel(kinEnergy, currentCoupleIndex);
605 currentModel->SetCurrentCouple(currentCouple);
606}
G4VEmModel * SelectModel(G4double &energy, size_t &index)
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:465

Referenced by CrossSectionPerVolume(), GetDEDXDispersion(), PostStepDoIt(), and PostStepGetPhysicalInteractionLength().

◆ SelectModelForMaterial()

G4VEmModel * G4VEnergyLossProcess::SelectModelForMaterial ( G4double  kinEnergy,
size_t &  idx 
) const
inline

Definition at line 610 of file G4VEnergyLossProcess.hh.

612{
613 return modelManager->SelectModel(kinEnergy, idx);
614}

Referenced by G4ContinuousGainOfEnergy::GetContinuousStepLimit().

◆ SetBaseParticle()

◆ SetCrossSectionBiasingFactor()

void G4VEnergyLossProcess::SetCrossSectionBiasingFactor ( G4double  f,
G4bool  flag = true 
)

Definition at line 2218 of file G4VEnergyLossProcess.cc.

2220{
2221 if(f > 0.0) {
2222 biasFactor = f;
2223 weightFlag = flag;
2224 if(1 < verboseLevel) {
2225 G4cout << "### SetCrossSectionBiasingFactor: for "
2226 << " process " << GetProcessName()
2227 << " biasFactor= " << f << " weightFlag= " << flag
2228 << G4endl;
2229 }
2230 }
2231}

Referenced by G4EmExtraParameters::DefineRegParamForLoss().

◆ SetCSDARangeTable()

void G4VEnergyLossProcess::SetCSDARangeTable ( G4PhysicsTable pRange)

Definition at line 2073 of file G4VEnergyLossProcess.cc.

2074{
2075 theCSDARangeTable = p;
2076
2077 if(p) {
2078 size_t n = p->length();
2079 G4PhysicsVector* pv;
2080 G4double emax = maxKinEnergyCSDA;
2081
2082 for (size_t i=0; i<n; ++i) {
2083 pv = (*p)[i];
2084 G4double rmax = 0.0;
2085 if(pv) { rmax = pv->Value(emax, idxCSDA); }
2086 else {
2087 pv = (*p)[(*theDensityIdx)[i]];
2088 if(pv) { rmax = pv->Value(emax, idxCSDA)/(*theDensityFactor)[i]; }
2089 }
2090 theRangeAtMaxEnergy[i] = rmax;
2091 //G4cout << "i= " << i << " Emax(MeV)= " << emax/MeV << " Rmax= "
2092 //<< rmax<< G4endl;
2093 }
2094 }
2095}
G4double Value(G4double theEnergy, std::size_t &lastidx) const

Referenced by BuildPhysicsTable().

◆ SetDEDXBinning()

void G4VEnergyLossProcess::SetDEDXBinning ( G4int  nbins)

Definition at line 2320 of file G4VEnergyLossProcess.cc.

2321{
2322 if(2 < n && n < 1000000000) {
2323 nBins = n;
2324 actBinning = true;
2325 } else {
2326 G4double e = (G4double)n;
2327 PrintWarning("SetDEDXBinning", e);
2328 }
2329}

Referenced by G4hhIonisation::InitialiseEnergyLossProcess(), G4mplIonisation::InitialiseEnergyLossProcess(), and G4hIonisation::InitialiseEnergyLossProcess().

◆ SetDEDXTable()

void G4VEnergyLossProcess::SetDEDXTable ( G4PhysicsTable p,
G4EmTableType  tType 
)

Definition at line 2012 of file G4VEnergyLossProcess.cc.

2013{
2014 if(fTotal == tType) {
2015 theDEDXunRestrictedTable = p;
2016 if(p) {
2017 size_t n = p->length();
2018 G4PhysicsVector* pv = (*p)[0];
2019 G4double emax = maxKinEnergyCSDA;
2020
2021 G4LossTableBuilder* bld = lManager->GetTableBuilder();
2022 theDensityFactor = bld->GetDensityFactors();
2023 theDensityIdx = bld->GetCoupleIndexes();
2024
2025 for (size_t i=0; i<n; ++i) {
2026 G4double dedx = 0.0;
2027 pv = (*p)[i];
2028 if(pv) {
2029 dedx = pv->Value(emax, idxDEDXunRestricted);
2030 } else {
2031 pv = (*p)[(*theDensityIdx)[i]];
2032 if(pv) {
2033 dedx =
2034 pv->Value(emax, idxDEDXunRestricted)*(*theDensityFactor)[i];
2035 }
2036 }
2037 theDEDXAtMaxEnergy[i] = dedx;
2038 //G4cout << "i= " << i << " emax(MeV)= " << emax/MeV<< " dedx= "
2039 // << dedx << G4endl;
2040 }
2041 }
2042
2043 } else if(fRestricted == tType) {
2044 /*
2045 G4cout<< "G4VEnergyLossProcess::SetDEDXTable "
2046 << particle->GetParticleName()
2047 << " oldTable " << theDEDXTable << " newTable " << p
2048 << " ion " << theIonisationTable
2049 << " IsMaster " << isMaster
2050 << " " << GetProcessName() << G4endl;
2051 G4cout << (*p) << G4endl;
2052 */
2053 theDEDXTable = p;
2054 } else if(fSubRestricted == tType) {
2055 theDEDXSubTable = p;
2056 } else if(fIsIonisation == tType) {
2057 /*
2058 G4cout<< "G4VEnergyLossProcess::SetIonisationTable "
2059 << particle->GetParticleName()
2060 << " oldTable " << theDEDXTable << " newTable " << p
2061 << " ion " << theIonisationTable
2062 << " IsMaster " << isMaster
2063 << " " << GetProcessName() << G4endl;
2064 */
2065 theIonisationTable = p;
2066 } else if(fIsSubIonisation == tType) {
2067 theIonisationSubTable = p;
2068 }
2069}
std::size_t length() const

Referenced by BuildPhysicsTable().

◆ SetDynamicMassCharge()

void G4VEnergyLossProcess::SetDynamicMassCharge ( G4double  massratio,
G4double  charge2ratio 
)
inline

Definition at line 635 of file G4VEnergyLossProcess.hh.

637{
638 massRatio = massratio;
639 logMassRatio = G4Log(massRatio);
640 fFactor = charge2ratio*biasFactor*(*theDensityFactor)[currentCoupleIndex];
641 chargeSqRatio = charge2ratio;
642 reduceFactor = 1.0/(fFactor*massRatio);
643}

Referenced by G4ContinuousGainOfEnergy::AlongStepDoIt(), and G4ContinuousGainOfEnergy::GetContinuousStepLimit().

◆ SetEmModel()

◆ SetFluctModel()

◆ SetIntegral()

void G4VEnergyLossProcess::SetIntegral ( G4bool  val)
inline

Definition at line 1003 of file G4VEnergyLossProcess.hh.

1004{
1005 integral = val;
1006 actIntegral = true;
1007}

◆ SetInverseRangeTable()

void G4VEnergyLossProcess::SetInverseRangeTable ( G4PhysicsTable p)

Definition at line 2123 of file G4VEnergyLossProcess.cc.

2124{
2125 theInverseRangeTable = p;
2126 if(1 < verboseLevel) {
2127 G4cout << "### Set InverseRange table " << p
2128 << " for " << particle->GetParticleName()
2129 << " and process " << GetProcessName() << G4endl;
2130 }
2131}

Referenced by BuildPhysicsTable().

◆ SetIonisation()

◆ SetLambdaTable()

void G4VEnergyLossProcess::SetLambdaTable ( G4PhysicsTable p)

Definition at line 2135 of file G4VEnergyLossProcess.cc.

2136{
2137 if(1 < verboseLevel) {
2138 G4cout << "### Set Lambda table " << p
2139 << " for " << particle->GetParticleName()
2140 << " and process " << GetProcessName() << G4endl;
2141 //G4cout << *p << G4endl;
2142 }
2143 theLambdaTable = p;
2144 tablesAreBuilt = true;
2145
2146 G4LossTableBuilder* bld = lManager->GetTableBuilder();
2147 theDensityFactor = bld->GetDensityFactors();
2148 theDensityIdx = bld->GetCoupleIndexes();
2149
2150 if(theLambdaTable) {
2151 size_t n = theLambdaTable->length();
2152 G4PhysicsVector* pv = (*theLambdaTable)[0];
2153 G4double e, ss, smax, emax;
2154
2155 size_t i;
2156
2157 // first loop on existing vectors
2158 for (i=0; i<n; ++i) {
2159 pv = (*theLambdaTable)[i];
2160 if(pv) {
2161 size_t nb = pv->GetVectorLength();
2162 emax = DBL_MAX;
2163 smax = 0.0;
2164 if(nb > 0) {
2165 for (size_t j=0; j<nb; ++j) {
2166 e = pv->Energy(j);
2167 ss = (*pv)(j);
2168 if(ss > smax) {
2169 smax = ss;
2170 emax = e;
2171 }
2172 }
2173 }
2174 theEnergyOfCrossSectionMax[i] = emax;
2175 theCrossSectionMax[i] = smax;
2176 if(1 < verboseLevel) {
2177 G4cout << "For " << particle->GetParticleName()
2178 << " Max CS at i= " << i << " emax(MeV)= " << emax/MeV
2179 << " lambda= " << smax << G4endl;
2180 }
2181 }
2182 }
2183 // second loop using base materials
2184 for (i=0; i<n; ++i) {
2185 pv = (*theLambdaTable)[i];
2186 if(!pv){
2187 G4int j = (*theDensityIdx)[i];
2188 theEnergyOfCrossSectionMax[i] = theEnergyOfCrossSectionMax[j];
2189 theCrossSectionMax[i] = (*theDensityFactor)[i]*theCrossSectionMax[j];
2190 }
2191 }
2192 }
2193}
G4double Energy(std::size_t index) const
std::size_t GetVectorLength() const

Referenced by BuildPhysicsTable().

◆ SetLinearLossLimit()

void G4VEnergyLossProcess::SetLinearLossLimit ( G4double  val)

Definition at line 2289 of file G4VEnergyLossProcess.cc.

2290{
2291 if(0.0 < val && val < 1.0) {
2292 linLossLimit = val;
2293 actLinLossLimit = true;
2294 } else { PrintWarning("SetLinearLossLimit", val); }
2295}

Referenced by G4alphaIonisation::G4alphaIonisation(), and G4ionIonisation::G4ionIonisation().

◆ SetLossFluctuations()

void G4VEnergyLossProcess::SetLossFluctuations ( G4bool  val)
inline

Definition at line 995 of file G4VEnergyLossProcess.hh.

996{
997 lossFluctuationFlag = val;
998 actLossFluc = true;
999}

◆ SetLowestEnergyLimit()

void G4VEnergyLossProcess::SetLowestEnergyLimit ( G4double  val)

Definition at line 2312 of file G4VEnergyLossProcess.cc.

2313{
2314 if(1.e-18 < val && val < 1.e+50) { lowestKinEnergy = val; }
2315 else { PrintWarning("SetLowestEnergyLimit", val); }
2316}

◆ SetMaxKinEnergy()

void G4VEnergyLossProcess::SetMaxKinEnergy ( G4double  e)

Definition at line 2343 of file G4VEnergyLossProcess.cc.

2344{
2345 if(minKinEnergy < e && e < 1.e+50) {
2346 maxKinEnergy = e;
2347 actMaxKinEnergy = true;
2348 if(e < maxKinEnergyCSDA) { maxKinEnergyCSDA = e; }
2349 } else { PrintWarning("SetMaxKinEnergy", e); }
2350}

Referenced by G4hhIonisation::InitialiseEnergyLossProcess(), G4mplIonisation::InitialiseEnergyLossProcess(), and G4hIonisation::InitialiseEnergyLossProcess().

◆ SetMinKinEnergy()

void G4VEnergyLossProcess::SetMinKinEnergy ( G4double  e)

Definition at line 2333 of file G4VEnergyLossProcess.cc.

2334{
2335 if(1.e-18 < e && e < maxKinEnergy) {
2336 minKinEnergy = e;
2337 actMinKinEnergy = true;
2338 } else { PrintWarning("SetMinKinEnergy", e); }
2339}

Referenced by G4hhIonisation::InitialiseEnergyLossProcess(), G4mplIonisation::InitialiseEnergyLossProcess(), and G4hIonisation::InitialiseEnergyLossProcess().

◆ SetParticle()

void G4VEnergyLossProcess::SetParticle ( const G4ParticleDefinition p)
inlineprotected

Definition at line 950 of file G4VEnergyLossProcess.hh.

951{
952 particle = p;
953}

◆ SetRangeTableForLoss()

void G4VEnergyLossProcess::SetRangeTableForLoss ( G4PhysicsTable p)

Definition at line 2099 of file G4VEnergyLossProcess.cc.

2100{
2101 theRangeTableForLoss = p;
2102 if(1 < verboseLevel) {
2103 G4cout << "### Set Range table " << p
2104 << " for " << particle->GetParticleName()
2105 << " and process " << GetProcessName() << G4endl;
2106 }
2107}

Referenced by BuildPhysicsTable().

◆ SetSecondaryParticle()

◆ SetSecondaryRangeTable()

void G4VEnergyLossProcess::SetSecondaryRangeTable ( G4PhysicsTable p)

Definition at line 2111 of file G4VEnergyLossProcess.cc.

2112{
2113 theSecondaryRangeTable = p;
2114 if(1 < verboseLevel) {
2115 G4cout << "### Set SecondaryRange table " << p
2116 << " for " << particle->GetParticleName()
2117 << " and process " << GetProcessName() << G4endl;
2118 }
2119}

Referenced by BuildPhysicsTable().

◆ SetStepFunction()

◆ SetSubLambdaTable()

void G4VEnergyLossProcess::SetSubLambdaTable ( G4PhysicsTable p)

Definition at line 2197 of file G4VEnergyLossProcess.cc.

2198{
2199 theSubLambdaTable = p;
2200 if(1 < verboseLevel) {
2201 G4cout << "### Set SebLambda table " << p
2202 << " for " << particle->GetParticleName()
2203 << " and process " << GetProcessName() << G4endl;
2204 }
2205}

Referenced by BuildPhysicsTable().

◆ StartTracking()

void G4VEnergyLossProcess::StartTracking ( G4Track track)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 981 of file G4VEnergyLossProcess.cc.

982{
983 /*
984 G4cout << track->GetDefinition()->GetParticleName()
985 << " e(MeV)= " << track->GetKineticEnergy()
986 << " baseParticle " << baseParticle << " proc " << this;
987 if(particle) G4cout << " " << particle->GetParticleName();
988 G4cout << " isIon= " << isIon << " dedx " << theDEDXTable <<G4endl;
989 */
990 // reset parameters for the new track
993 preStepRangeEnergy = 0.0;
994
995 // reset ion
996 if(isIon) {
997 chargeSqRatio = 0.5;
998
999 G4double newmass = track->GetDefinition()->GetPDGMass();
1000 if(baseParticle) {
1001 massRatio = baseParticle->GetPDGMass()/newmass;
1002 logMassRatio = G4Log(massRatio);
1003 } else if(theGenericIon) {
1004 massRatio = proton_mass_c2/newmass;
1005 logMassRatio = G4Log(massRatio);
1006 } else {
1007 massRatio = 1.0;
1008 logMassRatio = 0.0;
1009 }
1010 }
1011 // forced biasing only for primary particles
1012 if(biasManager) {
1013 if(0 == track->GetParentID()) {
1014 // primary particle
1015 biasFlag = true;
1016 biasManager->ResetForcedInteraction();
1017 }
1018 }
1019}

◆ StorePhysicsTable()

G4bool G4VEnergyLossProcess::StorePhysicsTable ( const G4ParticleDefinition part,
const G4String directory,
G4bool  ascii = false 
)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 1699 of file G4VEnergyLossProcess.cc.

1702{
1703 G4bool res = true;
1704 //G4cout << "G4VEnergyLossProcess::StorePhysicsTable: " << part->GetParticleName()
1705 // << " " << directory << " " << ascii << G4endl;
1706 if (!isMaster || baseParticle || part != particle ) return res;
1707
1708 if(!StoreTable(part,theDEDXTable,ascii,directory,"DEDX"))
1709 {res = false;}
1710
1711 if(!StoreTable(part,theDEDXunRestrictedTable,ascii,directory,"DEDXnr"))
1712 {res = false;}
1713
1714 if(!StoreTable(part,theDEDXSubTable,ascii,directory,"SubDEDX"))
1715 {res = false;}
1716
1717 if(!StoreTable(part,theIonisationTable,ascii,directory,"Ionisation"))
1718 {res = false;}
1719
1720 if(!StoreTable(part,theIonisationSubTable,ascii,directory,"SubIonisation"))
1721 {res = false;}
1722
1723 if(isIonisation &&
1724 !StoreTable(part,theCSDARangeTable,ascii,directory,"CSDARange"))
1725 {res = false;}
1726
1727 if(isIonisation &&
1728 !StoreTable(part,theRangeTableForLoss,ascii,directory,"Range"))
1729 {res = false;}
1730
1731 if(isIonisation &&
1732 !StoreTable(part,theInverseRangeTable,ascii,directory,"InverseRange"))
1733 {res = false;}
1734
1735 if(!StoreTable(part,theLambdaTable,ascii,directory,"Lambda"))
1736 {res = false;}
1737
1738 if(!StoreTable(part,theSubLambdaTable,ascii,directory,"SubLambda"))
1739 {res = false;}
1740
1741 if ( !res ) {
1742 if(1 < verboseLevel) {
1743 G4cout << "Physics tables are stored for "
1744 << particle->GetParticleName()
1745 << " and process " << GetProcessName()
1746 << " in the directory <" << directory
1747 << "> " << G4endl;
1748 }
1749 } else {
1750 G4cout << "Fail to store Physics Tables for "
1751 << particle->GetParticleName()
1752 << " and process " << GetProcessName()
1753 << " in the directory <" << directory
1754 << "> " << G4endl;
1755 }
1756 return res;
1757}

◆ StreamProcessInfo()

virtual void G4VEnergyLossProcess::StreamProcessInfo ( std::ostream &  ) const
inlineprotectedvirtual

Reimplemented in G4ePairProduction, G4MuPairProduction, G4eBremsstrahlung, and G4ionIonisation.

Definition at line 108 of file G4VEnergyLossProcess.hh.

108{};

◆ SubLambdaTable()

G4PhysicsTable * G4VEnergyLossProcess::SubLambdaTable ( ) const
inline

Definition at line 1130 of file G4VEnergyLossProcess.hh.

1131{
1132 return theSubLambdaTable;
1133}

Referenced by BuildPhysicsTable().

◆ TablesAreBuilt()

G4bool G4VEnergyLossProcess::TablesAreBuilt ( ) const
inline

Definition at line 1053 of file G4VEnergyLossProcess.hh.

1054{
1055 return tablesAreBuilt;
1056}

◆ UpdateEmModel()

void G4VEnergyLossProcess::UpdateEmModel ( const G4String nam,
G4double  emin,
G4double  emax 
)

Definition at line 322 of file G4VEnergyLossProcess.cc.

324{
325 modelManager->UpdateEmModel(nam, emin, emax);
326}
void UpdateEmModel(const G4String &model_name, G4double emin, G4double emax)

Member Data Documentation

◆ computedRange

G4double G4VEnergyLossProcess::computedRange
protected

Definition at line 559 of file G4VEnergyLossProcess.hh.

Referenced by G4VEnergyLossProcess().

◆ currentCouple

const G4MaterialCutsCouple* G4VEnergyLossProcess::currentCouple
protected

◆ currentCoupleIndex

◆ currentMaterial

const G4Material* G4VEnergyLossProcess::currentMaterial
protected

◆ fParticleChange

G4ParticleChangeForLoss G4VEnergyLossProcess::fParticleChange
protected

Definition at line 552 of file G4VEnergyLossProcess.hh.

Referenced by AlongStepDoIt(), G4VEnergyLossProcess(), and PostStepDoIt().

◆ fRange

G4double G4VEnergyLossProcess::fRange
protected

◆ mfpKinEnergy

G4double G4VEnergyLossProcess::mfpKinEnergy
protected

◆ preStepKinEnergy

G4double G4VEnergyLossProcess::preStepKinEnergy
protected

◆ preStepLambda

G4double G4VEnergyLossProcess::preStepLambda
protected

◆ preStepLogKinEnergy

G4double G4VEnergyLossProcess::preStepLogKinEnergy
protected

◆ preStepLogScaledEnergy

G4double G4VEnergyLossProcess::preStepLogScaledEnergy
protected

◆ preStepRangeEnergy

G4double G4VEnergyLossProcess::preStepRangeEnergy
protected

◆ preStepScaledEnergy

G4double G4VEnergyLossProcess::preStepScaledEnergy
protected

The documentation for this class was generated from the following files: