Geant4  10.01
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: G4ContinuousGainOfEnergy.cc 75591 2013-11-04 12:33:11Z gcosmo $
27 //
28 
30 
31 #include "G4PhysicalConstants.hh"
32 #include "G4SystemOfUnits.hh"
33 #include "G4Step.hh"
34 #include "G4ParticleDefinition.hh"
35 #include "G4VEmModel.hh"
36 #include "G4VEmFluctuationModel.hh"
37 #include "G4VParticleChange.hh"
38 #include "G4AdjointCSManager.hh"
39 #include "G4LossTableManager.hh"
40 #include "G4SystemOfUnits.hh"
41 #include "G4PhysicalConstants.hh"
42 
44 //
46  G4ProcessType type): G4VContinuousProcess(name, type)
47 {
48 
49 
50  linLossLimit=0.05;
53  is_integral = false;
54 
55  //Will be properly set in SetDirectParticle()
56  IsIon=false;
57  massRatio =1.;
58  chargeSqRatio=1.;
60 
61  //Some initialization
62  currentCoupleIndex=9999999;
64  currentMaterialIndex=9999999;
65  currentTcut=0.;
67  preStepRange=0.;
69 
70  currentCouple=0;
71 }
72 
74 //
76 {
77 
78 }
80 //
81 
83  const G4ParticleDefinition& )
84 {//theDirectEnergyLossProcess->PreparePhysicsTable(part);
85 
86 ;
87 }
88 
90 //
91 
93 {//theDirectEnergyLossProcess->BuildPhysicsTable(part);
94 ;
95 }
96 
98 //
101  if (theDirectPartDef->GetParticleType()== "nucleus") {
102  IsIon=true;
103  massRatio = proton_mass_c2/theDirectPartDef->GetPDGMass();
105  chargeSqRatio=q*q;
106 
107 
108  }
109 
110 }
111 
113 //
114 //
116  const G4Step& step)
117 {
118 
119  //Caution in this method the step length should be the true step length
120  // 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
121  //Forward sim. Nothing we can really do against that at this time. This is inherent to the MS method
122  //
123 
124 
125 
127 
128  // Get the actual (true) Step length
129  //----------------------------------
130  G4double length = step.GetStepLength();
131  G4double degain = 0.0;
132 
133 
134 
135  // Compute this for weight change after continuous energy loss
136  //-------------------------------------------------------------
138 
139 
140 
141  // For the fluctuation we generate a new dynamic particle with energy =preEnergy+egain
142  // and then compute the fluctuation given in the direct case.
143  //-----------------------------------------------------------------------
144  G4DynamicParticle* dynParticle = new G4DynamicParticle();
145  *dynParticle = *(track.GetDynamicParticle());
146  dynParticle->SetDefinition(theDirectPartDef);
147  G4double Tkin = dynParticle->GetKineticEnergy();
148 
149 
150  size_t n=1;
151  if (is_integral ) n=10;
152  n=1;
153  G4double dlength= length/n;
154  for (size_t i=0;i<n;i++) {
155  if (Tkin != preStepKinEnergy && IsIon) {
158 
159  }
160 
162  if( dlength <= linLossLimit * r ) {
163  degain = DEDX_before*dlength;
164  }
165  else {
166  G4double x = r + dlength;
167  //degain = theDirectEnergyLossProcess->GetKineticEnergy(x,currentCouple) - theDirectEnergyLossProcess->GetKineticEnergy(r,currentCouple);
169  if (IsIon){
173  while (std::abs(x-x1)>0.01*x) {
178 
179  }
180  }
181 
182  degain=E-Tkin;
183 
184 
185 
186  }
187  //G4cout<<degain<<G4endl;
188  G4double tmax = currentModel->MaxSecondaryKinEnergy(dynParticle);
189  tmax = std::min(tmax,currentTcut);
190 
191 
192  dynParticle->SetKineticEnergy(Tkin+degain);
193 
194  // Corrections, which cannot be tabulated for ions
195  //----------------------------------------
196  G4double esecdep=0;//not used in most models
197  currentModel->CorrectionsAlongStep(currentCouple, dynParticle, degain,esecdep, dlength);
198 
199  // Sample fluctuations
200  //-------------------
201 
202 
203  G4double deltaE =0.;
204  if (lossFluctuationFlag ) {
206  SampleFluctuations(currentCouple,dynParticle,tmax,dlength,degain)-degain;
207  }
208 
209  G4double egain=degain+deltaE;
210  if (egain <=0) egain=degain;
211  Tkin+=egain;
212  dynParticle->SetKineticEnergy(Tkin);
213  }
214 
215 
216 
217 
218 
219  delete dynParticle;
220 
221  if (IsIon){
224 
225  }
226 
228 
229 
230  G4double weight_correction=DEDX_after/DEDX_before;
231 
232 
234 
235 
236  //Caution!!!
237  // It is important to select the weight of the post_step_point
238  // as the current weight and not the weight of the track, as t
239  // the weight of the track is changed after having applied all
240  // the along_step_do_it.
241 
242  // G4double new_weight=weight_correction*track.GetWeight(); //old
243  G4double new_weight=weight_correction*step.GetPostStepPoint()->GetWeight();
246 
247 
248  return &aParticleChange;
249 
250 }
252 //
254 {
255  if(val && !lossFluctuationArePossible) return;
256  lossFluctuationFlag = val;
257 }
259 //
260 
261 
262 
264  G4double , G4double , G4double& )
265 {
266  G4double x = DBL_MAX;
267  x=.1*mm;
268 
269 
271 
275  G4double emax_model=currentModel->HighEnergyLimit();
276  if (IsIon) {
280  }
281 
282 
283  G4double maxE =1.1*preStepKinEnergy;
284  /*if (preStepKinEnergy< 0.05*MeV) maxE =2.*preStepKinEnergy;
285  else if (preStepKinEnergy< 0.1*MeV) maxE =1.5*preStepKinEnergy;
286  else if (preStepKinEnergy< 0.5*MeV) maxE =1.25*preStepKinEnergy;*/
287 
288  if (preStepKinEnergy < currentTcut) maxE = std::min(currentTcut,maxE);
289 
290  maxE=std::min(emax_model*1.001,maxE);
291 
293 
294  if (IsIon) {
297  }
298 
300 
302 
303 
304 
305  x=r1-preStepRange;
306  x=std::max(r1-preStepRange,0.001*mm);
307 
308  return x;
309 
310 
311 }
312 #include "G4EmCorrections.hh"
314 //
315 
317 {
318 
321 }
G4double GetKineticEnergy(G4double &range, const G4MaterialCutsCouple *)
G4double MaxSecondaryKinEnergy(const G4DynamicParticle *dynParticle)
Definition: G4VEmModel.hh:463
static G4LossTableManager * Instance()
G4double GetWeight() const
G4ContinuousGainOfEnergy(const G4String &name="EnergyGain", G4ProcessType type=fElectromagnetic)
virtual void CorrectionsAlongStep(const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
Definition: G4VEmModel.cc:365
void SetDynamicMassCharge(G4double massratio, G4double charge2ratio)
G4double GetKineticEnergy() const
G4double GetStepLength() const
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:615
const G4DynamicParticle * GetDynamicParticle() const
G4double GetDEDX(G4double &kineticEnergy, const G4MaterialCutsCouple *)
G4double EffectiveChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4String name
Definition: TRTMaterials.hh:40
void ProposeParentWeight(G4double finalWeight)
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
void BuildPhysicsTable(const G4ParticleDefinition &)
G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &)
G4VEmFluctuationModel * GetModelOfFluctuations()
Definition: G4VEmModel.hh:591
void SetDynamicMassCharge(const G4Track &track, G4double energy)
G4ParticleDefinition * theDirectPartDef
G4VEnergyLossProcess * theDirectEnergyLossProcess
void SetParentWeightByProcess(G4bool)
G4double GetKineticEnergy() const
const G4MaterialCutsCouple * currentCouple
bool G4bool
Definition: G4Types.hh:79
G4EmCorrections * EmCorrections()
void SetDirectParticle(G4ParticleDefinition *p)
const G4String & GetParticleType() const
Definition: G4Step.hh:76
const G4int n
void SetKineticEnergy(G4double aEnergy)
virtual G4double GetChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:348
virtual G4double GetContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)
virtual void Initialize(const G4Track &)
G4double GetPDGMass() const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4double energy(const ThreeVector &p, const G4double m)
G4StepPoint * GetPostStepPoint() const
void PreparePhysicsTable(const G4ParticleDefinition &)
void ProposeEnergy(G4double finalEnergy)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
void DefineMaterial(const G4MaterialCutsCouple *couple)
double G4double
Definition: G4Types.hh:76
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
G4double GetPDGCharge() const
#define DBL_MAX
Definition: templates.hh:83
static const double mm
Definition: G4SIunits.hh:102
G4VEmModel * SelectModelForMaterial(G4double kinEnergy, size_t &idx) const
G4double GetRange(G4double &kineticEnergy, const G4MaterialCutsCouple *)
G4ProcessType