Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ErrorEnergyLoss.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26
27#include "G4ErrorEnergyLoss.hh"
30
31//-------------------------------------------------------------------
33 G4ProcessType type)
34 : G4VContinuousProcess(processName, type)
35{
36 if (verboseLevel>2) {
37 G4cout << GetProcessName() << " is created " << G4endl;
38 }
39
40 theELossForExtrapolator = new G4EnergyLossForExtrapolator;
41
42 theStepLimit = 1.;
43}
44
45//-------------------------------------------------------------------
46void G4ErrorEnergyLoss::InstantiateEforExtrapolator()
47{}
48
49
50//-------------------------------------------------------------------
52{
53 delete theELossForExtrapolator;
54}
55
56//-------------------------------------------------------------------
59{
61
63
64 G4double kinEnergyStart = aTrack.GetKineticEnergy();
65 G4double step_length = aStep.GetStepLength();
66
67 const G4Material* aMaterial = aTrack.GetMaterial();
68 const G4ParticleDefinition* aParticleDef = aTrack.GetDynamicParticle()->GetDefinition();
69 G4double kinEnergyEnd = kinEnergyStart;
70
71 if( g4edata->GetMode() == G4ErrorMode(G4ErrorMode_PropBackwards) ) {
72 kinEnergyEnd = theELossForExtrapolator->EnergyBeforeStep( kinEnergyStart,
73 step_length,
74 aMaterial,
75 aParticleDef );
76 G4double kinEnergyHalfStep = kinEnergyStart - (kinEnergyStart-kinEnergyEnd)/2.;
77
78#ifdef G4VERBOSE
80 G4cout << " G4ErrorEnergyLoss FWD end " << kinEnergyEnd
81 << " halfstep " << kinEnergyHalfStep << G4endl;
82#endif
83
84 //--- rescale to energy lost at 1/2 step
85 kinEnergyEnd = theELossForExtrapolator->EnergyBeforeStep( kinEnergyHalfStep,
86 step_length,
87 aMaterial,
88 aParticleDef );
89 kinEnergyEnd = kinEnergyStart - (kinEnergyHalfStep - kinEnergyEnd );
90 }else if( g4edata->GetMode() == G4ErrorMode(G4ErrorMode_PropForwards) ) {
91 kinEnergyEnd = theELossForExtrapolator->EnergyAfterStep( kinEnergyStart,
92 step_length,
93 aMaterial,
94 aParticleDef );
95 G4double kinEnergyHalfStep = kinEnergyStart - (kinEnergyStart-kinEnergyEnd)/2.;
96#ifdef G4VERBOSE
98 G4cout << " G4ErrorEnergyLoss BCKD end " << kinEnergyEnd
99 << " halfstep " << kinEnergyHalfStep << G4endl;
100#endif
101
102 //--- rescale to energy lost at 1/2 step
103 kinEnergyEnd = theELossForExtrapolator->EnergyAfterStep( kinEnergyHalfStep,
104 step_length,
105 aMaterial,
106 aParticleDef );
107 kinEnergyEnd = kinEnergyStart - (kinEnergyHalfStep - kinEnergyEnd );
108 }
109
110 G4double edepo = kinEnergyEnd - kinEnergyStart;
111
112#ifdef G4VERBOSE
114 G4cout << "AlongStepDoIt Estart= " << kinEnergyStart << " Eend " << kinEnergyEnd
115 << " Ediff " << kinEnergyStart-kinEnergyEnd << " step= " << step_length
116 << " mate= " << aMaterial->GetName()
117 << " particle= " << aParticleDef->GetParticleName() << G4endl;
118#endif
119
123
124 aParticleChange.ProposeEnergy( kinEnergyEnd );
125
126 return &aParticleChange;
127}
128
129
130//-------------------------------------------------------------------
132 G4double ,
133 G4double currentMinimumStep,
134 G4double& )
135{
136 G4double Step = DBL_MAX;
137 if( theStepLimit != 1. ) {
138 G4double kinEnergyStart = aTrack.GetKineticEnergy();
139 G4double kinEnergyLoss = kinEnergyStart;
140 const G4Material* aMaterial = aTrack.GetMaterial();
141 const G4ParticleDefinition* aParticleDef = aTrack.GetDynamicParticle()->GetDefinition();
143 if( g4edata->GetMode() == G4ErrorMode(G4ErrorMode_PropBackwards) ) {
144 kinEnergyLoss = - kinEnergyStart +
145 theELossForExtrapolator->EnergyBeforeStep( kinEnergyStart, currentMinimumStep,
146 aMaterial, aParticleDef );
147 }else if( g4edata->GetMode() == G4ErrorMode(G4ErrorMode_PropForwards) ) {
148 kinEnergyLoss = kinEnergyStart -
149 theELossForExtrapolator->EnergyAfterStep( kinEnergyStart, currentMinimumStep,
150 aMaterial, aParticleDef );
151 }
152#ifdef G4VERBOSE
154 G4cout << " G4ErrorEnergyLoss: currentMinimumStep " <<currentMinimumStep
155 << " kinEnergyLoss " << kinEnergyLoss
156 << " kinEnergyStart " << kinEnergyStart << G4endl;
157#endif
158 if( kinEnergyLoss / kinEnergyStart > theStepLimit ) {
159 Step = theStepLimit / (kinEnergyLoss / kinEnergyStart) * currentMinimumStep;
160#ifdef G4VERBOSE
162 G4cout << " G4ErrorEnergyLoss: limiting Step " << Step
163 << " energy loss fraction " << kinEnergyLoss / kinEnergyStart
164 << " > " << theStepLimit << G4endl;
165#endif
166 }
167 }
168
169 return Step;
170
171}
@ G4ErrorMode_PropForwards
@ G4ErrorMode_PropBackwards
G4ProcessType
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4ParticleDefinition * GetDefinition() const
G4double EnergyBeforeStep(G4double kinEnergy, G4double step, const G4Material *, const G4ParticleDefinition *)
G4double EnergyAfterStep(G4double kinEnergy, G4double step, const G4Material *, const G4ParticleDefinition *)
G4double GetContinuousStepLimit(const G4Track &aTrack, G4double, G4double currentMinimumStep, G4double &) override
G4VParticleChange * AlongStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
G4ErrorEnergyLoss(const G4String &processName="G4ErrorEnergyLoss", G4ProcessType type=fElectromagnetic)
virtual ~G4ErrorEnergyLoss()
static G4ErrorPropagatorData * GetErrorPropagatorData()
G4ErrorMode GetMode() const
const G4String & GetName() const
Definition: G4Material.hh:175
void ProposeEnergy(G4double finalEnergy)
virtual void Initialize(const G4Track &)
const G4String & GetParticleName() const
Definition: G4Step.hh:62
G4double GetStepLength() const
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
G4double GetKineticEnergy() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetNumberOfSecondaries(G4int totSecondaries)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:327
G4int verboseLevel
Definition: G4VProcess.hh:356
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382
#define DBL_MAX
Definition: templates.hh:62