Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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// $Id$
27//
28
30
32#include "G4SystemOfUnits.hh"
33#include "G4Step.hh"
35#include "G4VEmModel.hh"
37#include "G4VParticleChange.hh"
38#include "G4UnitsTable.hh"
39#include "G4AdjointCSManager.hh"
40#include "G4LossTableManager.hh"
41
42///////////////////////////////////////////////////////
43//
45 G4ProcessType type): G4VContinuousProcess(name, type)
46{
47
48
49 linLossLimit=0.05;
50 lossFluctuationArePossible =true;
51 lossFluctuationFlag=true;
52 is_integral = false;
53
54 //Will be properly set in SetDirectParticle()
55 IsIon=false;
56 massRatio =1.;
57 chargeSqRatio=1.;
58 preStepChargeSqRatio=1.;
59
60 //Some initialization
61 currentCoupleIndex=9999999;
62 currentCutInRange=0.;
63 currentMaterialIndex=9999999;
64 currentTcut=0.;
65 preStepKinEnergy=0.;
66 preStepRange=0.;
67 preStepScaledKinEnergy=0.;
68
69 currentCouple=0;
70}
71
72///////////////////////////////////////////////////////
73//
75{
76
77}
78///////////////////////////////////////////////////////
79//
80
83{//theDirectEnergyLossProcess->PreparePhysicsTable(part);
84
85;
86}
87
88///////////////////////////////////////////////////////
89//
90
92{//theDirectEnergyLossProcess->BuildPhysicsTable(part);
93;
94}
95
96///////////////////////////////////////////////////////
97//
99{theDirectPartDef=p;
100 if (theDirectPartDef->GetParticleType()== "nucleus") {
101 IsIon=true;
102 massRatio = proton_mass_c2/theDirectPartDef->GetPDGMass();
103 G4double q=theDirectPartDef->GetPDGCharge();
104 chargeSqRatio=q*q;
105
106
107 }
108
109}
110
111///////////////////////////////////////////////////////
112//
113//
115 const G4Step& step)
116{
117
118 //Caution in this method the step length should be the true step length
119 // A problem is that this is compute by the multiple scattering that does not know the energy at the end of the adjoint step. This energy is used during the
120 //Forward sim. Nothing we can really do against that at this time. This is inherent to the MS method
121 //
122
123
124
126
127 // Get the actual (true) Step length
128 //----------------------------------
129 G4double length = step.GetStepLength();
130 G4double degain = 0.0;
131
132
133
134 // Compute this for weight change after continuous energy loss
135 //-------------------------------------------------------------
136 G4double DEDX_before = theDirectEnergyLossProcess->GetDEDX(preStepKinEnergy, currentCouple);
137
138
139
140 // For the fluctuation we generate a new dynamic particle with energy =preEnergy+egain
141 // and then compute the fluctuation given in the direct case.
142 //-----------------------------------------------------------------------
143 G4DynamicParticle* dynParticle = new G4DynamicParticle();
144 *dynParticle = *(track.GetDynamicParticle());
145 dynParticle->SetDefinition(theDirectPartDef);
146 G4double Tkin = dynParticle->GetKineticEnergy();
147
148
149 size_t n=1;
150 if (is_integral ) n=10;
151 n=1;
152 G4double dlength= length/n;
153 for (size_t i=0;i<n;i++) {
154 if (Tkin != preStepKinEnergy && IsIon) {
155 chargeSqRatio = currentModel->GetChargeSquareRatio(theDirectPartDef,currentMaterial,Tkin);
156 theDirectEnergyLossProcess->SetDynamicMassCharge(massRatio,chargeSqRatio);
157
158 }
159
160 G4double r = theDirectEnergyLossProcess->GetRange(Tkin, currentCouple);
161 if( dlength <= linLossLimit * r ) {
162 degain = DEDX_before*dlength;
163 }
164 else {
165 G4double x = r + dlength;
166 //degain = theDirectEnergyLossProcess->GetKineticEnergy(x,currentCouple) - theDirectEnergyLossProcess->GetKineticEnergy(r,currentCouple);
167 G4double E = theDirectEnergyLossProcess->GetKineticEnergy(x,currentCouple);
168 if (IsIon){
169 chargeSqRatio = currentModel->GetChargeSquareRatio(theDirectPartDef,currentMaterial,E);
170 theDirectEnergyLossProcess->SetDynamicMassCharge(massRatio,chargeSqRatio);
171 G4double x1= theDirectEnergyLossProcess->GetRange(E, currentCouple);
172 while (std::abs(x-x1)>0.01*x) {
173 E = theDirectEnergyLossProcess->GetKineticEnergy(x,currentCouple);
174 chargeSqRatio = currentModel->GetChargeSquareRatio(theDirectPartDef,currentMaterial,E);
175 theDirectEnergyLossProcess->SetDynamicMassCharge(massRatio,chargeSqRatio);
176 x1= theDirectEnergyLossProcess->GetRange(E, currentCouple);
177
178 }
179 }
180
181 degain=E-Tkin;
182
183
184
185 }
186 //G4cout<<degain<<G4endl;
187 G4double tmax = currentModel->MaxSecondaryKinEnergy(dynParticle);
188 tmax = std::min(tmax,currentTcut);
189
190
191 dynParticle->SetKineticEnergy(Tkin+degain);
192
193 // Corrections, which cannot be tabulated for ions
194 //----------------------------------------
195 G4double esecdep=0;//not used in most models
196 currentModel->CorrectionsAlongStep(currentCouple, dynParticle, degain,esecdep, dlength);
197
198 // Sample fluctuations
199 //-------------------
200
201
202 G4double deltaE =0.;
203 if (lossFluctuationFlag ) {
204 deltaE = currentModel->GetModelOfFluctuations()->
205 SampleFluctuations(currentMaterial,dynParticle,tmax,dlength,degain)-degain;
206 }
207
208 G4double egain=degain+deltaE;
209 if (egain <=0) egain=degain;
210 Tkin+=egain;
211 dynParticle->SetKineticEnergy(Tkin);
212 }
213
214
215
216
217
218 delete dynParticle;
219
220 if (IsIon){
221 chargeSqRatio = currentModel->GetChargeSquareRatio(theDirectPartDef,currentMaterial,Tkin);
222 theDirectEnergyLossProcess->SetDynamicMassCharge(massRatio,chargeSqRatio);
223
224 }
225
226 G4double DEDX_after = theDirectEnergyLossProcess->GetDEDX(Tkin, currentCouple);
227
228
229 G4double weight_correction=DEDX_after/DEDX_before;
230
231
233
234
235 //Caution!!!
236 // It is important to select the weight of the post_step_point
237 // as the current weight and not the weight of the track, as t
238 // the weight of the track is changed after having applied all
239 // the along_step_do_it.
240
241 // G4double new_weight=weight_correction*track.GetWeight(); //old
242 G4double new_weight=weight_correction*step.GetPostStepPoint()->GetWeight();
245
246
247 return &aParticleChange;
248
249}
250///////////////////////////////////////////////////////
251//
253{
254 if(val && !lossFluctuationArePossible) return;
255 lossFluctuationFlag = val;
256}
257///////////////////////////////////////////////////////
258//
259
260
261
264{
265 G4double x = DBL_MAX;
266 x=.1*mm;
267
268
269 DefineMaterial(track.GetMaterialCutsCouple());
270
271 preStepKinEnergy = track.GetKineticEnergy();
272 preStepScaledKinEnergy = track.GetKineticEnergy()*massRatio;
273 currentModel = theDirectEnergyLossProcess->SelectModelForMaterial(preStepScaledKinEnergy,currentCoupleIndex);
274 G4double emax_model=currentModel->HighEnergyLimit();
275 if (IsIon) {
276 chargeSqRatio = currentModel->GetChargeSquareRatio(theDirectPartDef,currentMaterial,preStepKinEnergy);
277 preStepChargeSqRatio = chargeSqRatio;
278 theDirectEnergyLossProcess->SetDynamicMassCharge(massRatio,preStepChargeSqRatio);
279 }
280
281
282 G4double maxE =1.1*preStepKinEnergy;
283 /*if (preStepKinEnergy< 0.05*MeV) maxE =2.*preStepKinEnergy;
284 else if (preStepKinEnergy< 0.1*MeV) maxE =1.5*preStepKinEnergy;
285 else if (preStepKinEnergy< 0.5*MeV) maxE =1.25*preStepKinEnergy;*/
286
287 if (preStepKinEnergy < currentTcut) maxE = std::min(currentTcut,maxE);
288
289 maxE=std::min(emax_model*1.001,maxE);
290
291 preStepRange = theDirectEnergyLossProcess->GetRange(preStepKinEnergy, currentCouple);
292
293 if (IsIon) {
294 G4double chargeSqRatioAtEmax = currentModel->GetChargeSquareRatio(theDirectPartDef,currentMaterial,maxE);
295 theDirectEnergyLossProcess->SetDynamicMassCharge(massRatio,chargeSqRatioAtEmax);
296 }
297
298 G4double r1 = theDirectEnergyLossProcess->GetRange(maxE, currentCouple);
299
300 if (IsIon) theDirectEnergyLossProcess->SetDynamicMassCharge(massRatio,preStepChargeSqRatio);
301
302
303
304 x=r1-preStepRange;
305 x=std::max(r1-preStepRange,0.001*mm);
306
307 return x;
308
309
310}
311#include "G4EmCorrections.hh"
312///////////////////////////////////////////////////////
313//
314
315void G4ContinuousGainOfEnergy::SetDynamicMassCharge(const G4Track& ,G4double energy)
316{
317
318 G4double ChargeSqRatio= G4LossTableManager::Instance()->EmCorrections()->EffectiveChargeSquareRatio(theDirectPartDef,currentMaterial,energy);
319 if (theDirectEnergyLossProcess) theDirectEnergyLossProcess->SetDynamicMassCharge(massRatio,ChargeSqRatio);
320}
G4ProcessType
double G4double
Definition: G4Types.hh:64
bool G4bool
Definition: G4Types.hh:67
virtual G4double GetContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)
void BuildPhysicsTable(const G4ParticleDefinition &)
G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &)
void SetDirectParticle(G4ParticleDefinition *p)
void PreparePhysicsTable(const G4ParticleDefinition &)
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 ProposeEnergy(G4double finalEnergy)
virtual void Initialize(const G4Track &)
const G4String & GetParticleType() const
G4double GetPDGCharge() const
G4double GetWeight() const
Definition: G4Step.hh:78
G4double GetStepLength() const
G4StepPoint * GetPostStepPoint() const
const G4DynamicParticle * GetDynamicParticle() const
G4double GetKineticEnergy() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
virtual void CorrectionsAlongStep(const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
Definition: G4VEmModel.cc:279
G4VEmFluctuationModel * GetModelOfFluctuations()
Definition: G4VEmModel.hh:501
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:522
virtual G4double GetChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:262
G4double MaxSecondaryKinEnergy(const G4DynamicParticle *dynParticle)
Definition: G4VEmModel.hh:399
G4double GetKineticEnergy(G4double &range, const G4MaterialCutsCouple *)
G4VEmModel * SelectModelForMaterial(G4double kinEnergy, size_t &idx) const
void SetDynamicMassCharge(G4double massratio, G4double charge2ratio)
G4double GetDEDX(G4double &kineticEnergy, const G4MaterialCutsCouple *)
G4double GetRange(G4double &kineticEnergy, const G4MaterialCutsCouple *)
void SetParentWeightByProcess(G4bool)
void ProposeParentWeight(G4double finalWeight)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
#define DBL_MAX
Definition: templates.hh:83