Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4ContinuousGainOfEnergy.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
29
30#include "G4EmCorrections.hh"
31#include "G4LossTableManager.hh"
32#include "G4Material.hh"
34#include "G4ParticleChange.hh"
37#include "G4Step.hh"
38#include "G4SystemOfUnits.hh"
40#include "G4VEmModel.hh"
42#include "G4VParticleChange.hh"
43
44///////////////////////////////////////////////////////
46 G4ProcessType type)
47 : G4VContinuousProcess(name, type)
48{}
49
50///////////////////////////////////////////////////////
52
53///////////////////////////////////////////////////////
55{
56 out << "Continuous process acting on adjoint particles to compute the "
57 "continuous gain of energy of charged particles when they are "
58 "tracked back.\n";
59}
60
61///////////////////////////////////////////////////////
63{
64 fDirectPartDef = p;
65 if(fDirectPartDef->GetParticleType() == "nucleus")
66 {
67 fIsIon = true;
68 fMassRatio = proton_mass_c2 / fDirectPartDef->GetPDGMass();
69 }
70}
71
72///////////////////////////////////////////////////////
74 const G4Step& step)
75{
76 // Caution in this method the step length should be the true step length
77 // A problem is that this is computed by the multiple scattering that does
78 // not know the energy at the end of the adjoint step. This energy is used
79 // during the forward sim. Nothing we can really do against that at this
80 // time. This is inherent to the MS method
81
83
84 // Get the actual (true) Step length
85 G4double length = step.GetStepLength();
86 G4double degain = 0.0;
87
88 // Compute this for weight change after continuous energy loss
89 G4double DEDX_before =
90 fDirectEnergyLossProcess->GetDEDX(fPreStepKinEnergy, fCurrentCouple);
91
92 // For the fluctuation we generate a new dynamic particle with energy
93 // = preEnergy+egain and then compute the fluctuation given in the direct
94 // case.
95 G4DynamicParticle* dynParticle = new G4DynamicParticle();
96 *dynParticle = *(track.GetDynamicParticle());
97 dynParticle->SetDefinition(fDirectPartDef);
98 G4double Tkin = dynParticle->GetKineticEnergy();
99
100 G4double dlength = length;
101 if(Tkin != fPreStepKinEnergy && fIsIon)
102 {
103 G4double chargeSqRatio = fCurrentModel->GetChargeSquareRatio(
104 fDirectPartDef, fCurrentMaterial, Tkin);
105 fDirectEnergyLossProcess->SetDynamicMassCharge(fMassRatio, chargeSqRatio);
106 }
107
108 G4double r = fDirectEnergyLossProcess->GetRange(Tkin, fCurrentCouple);
109 if(dlength <= fLinLossLimit * r)
110 {
111 degain = DEDX_before * dlength;
112 }
113 else
114 {
115 G4double x = r + dlength;
116 G4double E = fDirectEnergyLossProcess->GetKineticEnergy(x, fCurrentCouple);
117 if(fIsIon)
118 {
119 G4double chargeSqRatio = fCurrentModel->GetChargeSquareRatio(
120 fDirectPartDef, fCurrentMaterial, E);
121 fDirectEnergyLossProcess->SetDynamicMassCharge(fMassRatio, chargeSqRatio);
122 G4double x1 = fDirectEnergyLossProcess->GetRange(E, fCurrentCouple);
123
124 G4int ii = 0;
125 constexpr G4int iimax = 100;
126 while(std::abs(x - x1) > 0.01 * x)
127 {
128 E = fDirectEnergyLossProcess->GetKineticEnergy(x, fCurrentCouple);
129 chargeSqRatio = fCurrentModel->GetChargeSquareRatio(
130 fDirectPartDef, fCurrentMaterial, E);
131 fDirectEnergyLossProcess->SetDynamicMassCharge(fMassRatio,
132 chargeSqRatio);
133 x1 = fDirectEnergyLossProcess->GetRange(E, fCurrentCouple);
134 ++ii;
135 if(ii >= iimax)
136 {
137 break;
138 }
139 }
140 }
141
142 degain = E - Tkin;
143 }
144 G4double tmax = fCurrentModel->MaxSecondaryKinEnergy(dynParticle);
145 fCurrentTcut = std::min(fCurrentTcut, tmax);
146
147 dynParticle->SetKineticEnergy(Tkin + degain);
148
149 // Corrections, which cannot be tabulated for ions
150 fCurrentModel->CorrectionsAlongStep(fCurrentCouple, dynParticle, dlength, degain);
151
152 // Sample fluctuations
153 G4double deltaE = 0.;
154 if(fLossFluctuationFlag)
155 {
156 deltaE = fCurrentModel->GetModelOfFluctuations()->SampleFluctuations(
157 fCurrentCouple, dynParticle, fCurrentTcut, tmax, dlength, degain)
158 - degain;
159 }
160
161 G4double egain = degain + deltaE;
162 if(egain <= 0.)
163 egain = degain;
164 Tkin += egain;
165 dynParticle->SetKineticEnergy(Tkin);
166
167 delete dynParticle;
168
169 if(fIsIon)
170 {
171 G4double chargeSqRatio = fCurrentModel->GetChargeSquareRatio(
172 fDirectPartDef, fCurrentMaterial, Tkin);
173 fDirectEnergyLossProcess->SetDynamicMassCharge(fMassRatio, chargeSqRatio);
174 }
175
176 G4double DEDX_after = fDirectEnergyLossProcess->GetDEDX(Tkin, fCurrentCouple);
177 G4double weight_correction = DEDX_after / DEDX_before;
178
180
181 // Caution!!! It is important to select the weight of the post_step_point
182 // as the current weight and not the weight of the track, as the weight of
183 // the track is changed after having applied all the along_step_do_it.
184
185 G4double new_weight =
186 weight_correction * step.GetPostStepPoint()->GetWeight();
189
190 return &aParticleChange;
191}
192
193///////////////////////////////////////////////////////
195{
196 if(val && !fLossFluctuationArePossible)
197 return;
198 fLossFluctuationFlag = val;
199}
200
201///////////////////////////////////////////////////////
204 G4double&)
205{
206 DefineMaterial(track.GetMaterialCutsCouple());
207
208 fPreStepKinEnergy = track.GetKineticEnergy();
209 fCurrentModel = fDirectEnergyLossProcess->SelectModelForMaterial(
210 track.GetKineticEnergy() * fMassRatio, fCurrentCoupleIndex);
211 G4double emax_model = fCurrentModel->HighEnergyLimit();
212 G4double preStepChargeSqRatio = 0.;
213 if(fIsIon)
214 {
215 G4double chargeSqRatio = fCurrentModel->GetChargeSquareRatio(
216 fDirectPartDef, fCurrentMaterial, fPreStepKinEnergy);
217 preStepChargeSqRatio = chargeSqRatio;
218 fDirectEnergyLossProcess->SetDynamicMassCharge(fMassRatio,
219 preStepChargeSqRatio);
220 }
221
222 G4double maxE = 1.1 * fPreStepKinEnergy;
223
224 if(fPreStepKinEnergy < fCurrentTcut)
225 maxE = std::min(fCurrentTcut, maxE);
226
227 maxE = std::min(emax_model * 1.001, maxE);
228
229 G4double preStepRange =
230 fDirectEnergyLossProcess->GetRange(fPreStepKinEnergy, fCurrentCouple);
231
232 if(fIsIon)
233 {
234 G4double chargeSqRatioAtEmax = fCurrentModel->GetChargeSquareRatio(
235 fDirectPartDef, fCurrentMaterial, maxE);
236 fDirectEnergyLossProcess->SetDynamicMassCharge(fMassRatio,
237 chargeSqRatioAtEmax);
238 }
239
240 G4double r1 = fDirectEnergyLossProcess->GetRange(maxE, fCurrentCouple);
241
242 if(fIsIon)
243 fDirectEnergyLossProcess->SetDynamicMassCharge(fMassRatio,
244 preStepChargeSqRatio);
245
246 return std::max(r1 - preStepRange, 0.001 * mm);
247}
248
249///////////////////////////////////////////////////////
250void G4ContinuousGainOfEnergy::SetDynamicMassCharge(const G4Track&,
251 G4double energy)
252{
253 G4double ChargeSqRatio =
255 fDirectPartDef, fCurrentMaterial, energy);
256 if(fDirectEnergyLossProcess)
257 fDirectEnergyLossProcess->SetDynamicMassCharge(fMassRatio, ChargeSqRatio);
258}
G4ProcessType
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
void SetDirectParticle(G4ParticleDefinition *p)
G4double GetContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety) override
G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &) override
void ProcessDescription(std::ostream &) const override
G4ContinuousGainOfEnergy(const G4String &name="EnergyGain", G4ProcessType type=fElectromagnetic)
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
G4double GetKineticEnergy() const
void SetKineticEnergy(G4double aEnergy)
G4double EffectiveChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
static G4LossTableManager * Instance()
G4EmCorrections * EmCorrections()
void Initialize(const G4Track &) override
void ProposeEnergy(G4double finalEnergy)
const G4String & GetParticleType() const
G4double GetWeight() const
Definition: G4Step.hh:62
G4double GetStepLength() const
G4StepPoint * GetPostStepPoint() const
const G4DynamicParticle * GetDynamicParticle() const
G4double GetKineticEnergy() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
virtual G4double SampleFluctuations(const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double tcut, const G4double tmax, const G4double length, const G4double meanLoss)=0
G4VEmFluctuationModel * GetModelOfFluctuations()
Definition: G4VEmModel.hh:593
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
virtual void CorrectionsAlongStep(const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double &length, G4double &eloss)
Definition: G4VEmModel.cc:342
virtual G4double GetChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:325
G4double MaxSecondaryKinEnergy(const G4DynamicParticle *dynParticle)
Definition: G4VEmModel.hh:501
G4double GetKineticEnergy(G4double range, const G4MaterialCutsCouple *)
void SetDynamicMassCharge(G4double massratio, G4double charge2ratio)
G4double GetDEDX(G4double kineticEnergy, const G4MaterialCutsCouple *)
G4double GetRange(G4double kineticEnergy, const G4MaterialCutsCouple *)
G4VEmModel * SelectModelForMaterial(G4double kinEnergy, std::size_t &idxCouple) const
void SetParentWeightByProcess(G4bool)
void ProposeParentWeight(G4double finalWeight)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:331