source: trunk/source/processes/electromagnetic/adjoint/src/G4ContinuousGainOfEnergy.cc @ 1317

Last change on this file since 1317 was 1228, checked in by garnier, 15 years ago

update geant4.9.3 tag

File size: 10.3 KB
Line 
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,v 1.4 2009/11/20 10:31:20 ldesorgh Exp $
27// GEANT4 tag $Name: geant4-09-03 $
28//
29#include "G4ContinuousGainOfEnergy.hh"
30#include "G4Step.hh"
31#include "G4ParticleDefinition.hh"
32#include "G4VEmModel.hh"
33#include "G4VEmFluctuationModel.hh"
34#include "G4VParticleChange.hh"
35#include "G4UnitsTable.hh"
36#include "G4AdjointCSManager.hh"
37#include "G4LossTableManager.hh"
38
39
40///////////////////////////////////////////////////////
41//
42G4ContinuousGainOfEnergy::G4ContinuousGainOfEnergy(const G4String& name, 
43  G4ProcessType type): G4VContinuousProcess(name, type)
44{
45
46
47  linLossLimit=0.05;
48  lossFluctuationArePossible =true;
49  lossFluctuationFlag=true;
50  is_integral = false;
51 
52  //Will be properly set in SetDirectParticle()
53  IsIon=false;
54  massRatio =1.;
55  chargeSqRatio=1.;
56  preStepChargeSqRatio=1.;
57
58 
59
60 
61 
62}
63
64///////////////////////////////////////////////////////
65//
66G4ContinuousGainOfEnergy::~G4ContinuousGainOfEnergy()
67{
68 
69}
70///////////////////////////////////////////////////////
71//
72
73void G4ContinuousGainOfEnergy::PreparePhysicsTable(
74     const G4ParticleDefinition& )
75{//theDirectEnergyLossProcess->PreparePhysicsTable(part);
76
77; 
78}
79
80///////////////////////////////////////////////////////
81//
82
83void G4ContinuousGainOfEnergy::BuildPhysicsTable(const G4ParticleDefinition&)
84{//theDirectEnergyLossProcess->BuildPhysicsTable(part);
85;
86}
87
88///////////////////////////////////////////////////////
89//
90void  G4ContinuousGainOfEnergy::SetDirectParticle(G4ParticleDefinition* p)
91{theDirectPartDef=p;
92 if (theDirectPartDef->GetParticleType()== "nucleus") {
93         IsIon=true;
94         massRatio = proton_mass_c2/theDirectPartDef->GetPDGMass();
95         G4double q=theDirectPartDef->GetPDGCharge();
96         chargeSqRatio=q*q;
97       
98         
99 }
100 
101}
102
103///////////////////////////////////////////////////////
104//
105//
106G4VParticleChange* G4ContinuousGainOfEnergy::AlongStepDoIt(const G4Track& track,
107                                                       const G4Step& step)
108{
109   
110  //Caution in this method the  step length should be the true step length
111  // 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
112  //Forward sim. Nothing we can really do against that at this time. This is inherent to the MS method
113  //
114 
115 
116 
117  aParticleChange.Initialize(track);
118 
119  // Get the actual (true) Step length
120  //----------------------------------
121  G4double length = step.GetStepLength();
122  G4double degain  = 0.0;
123 
124 
125 
126  // Compute this for weight change after continuous energy loss
127  //-------------------------------------------------------------
128  G4double DEDX_before = theDirectEnergyLossProcess->GetDEDX(preStepKinEnergy, currentCouple);
129   
130 
131 
132  // For the fluctuation we generate a new dynamic particle with energy =preEnergy+egain
133  // and then compute the fluctuation given in  the direct case.
134  //-----------------------------------------------------------------------
135  G4DynamicParticle* dynParticle = new G4DynamicParticle();
136  *dynParticle = *(track.GetDynamicParticle());
137  dynParticle->SetDefinition(theDirectPartDef);
138  G4double Tkin = dynParticle->GetKineticEnergy(); 
139  G4double Tkin1=Tkin*0.001;
140
141  size_t n=1;
142  if (is_integral ) n=10;
143  n=1;
144  G4double dlength= length/n; 
145  for (size_t i=0;i<n;i++) {
146        G4double factor_dE=1.;
147        if (Tkin != preStepKinEnergy && IsIon) {
148                chargeSqRatio =  currentModel->GetChargeSquareRatio(theDirectPartDef,currentMaterial,Tkin);
149                theDirectEnergyLossProcess->SetDynamicMassCharge(massRatio,chargeSqRatio); 
150       
151        }
152 
153        G4double r = theDirectEnergyLossProcess->GetRange(Tkin, currentCouple);
154        if( dlength <= linLossLimit * r ) {
155                degain = DEDX_before*dlength;
156                G4double degain1 = dlength*theDirectEnergyLossProcess->GetDEDX(Tkin1, currentCouple);
157                factor_dE=1.+(degain1-degain)/(Tkin1-Tkin);
158        } 
159        else {
160                G4double x = r + dlength;
161                //degain = theDirectEnergyLossProcess->GetKineticEnergy(x,currentCouple) - theDirectEnergyLossProcess->GetKineticEnergy(r,currentCouple);
162                G4double E = theDirectEnergyLossProcess->GetKineticEnergy(x,currentCouple);
163                if (IsIon){
164                        chargeSqRatio =  currentModel->GetChargeSquareRatio(theDirectPartDef,currentMaterial,E);
165                        theDirectEnergyLossProcess->SetDynamicMassCharge(massRatio,chargeSqRatio);
166                        G4double x1= theDirectEnergyLossProcess->GetRange(E, currentCouple);
167                        while (std::abs(x-x1)>0.01*x) {
168                                E = theDirectEnergyLossProcess->GetKineticEnergy(x,currentCouple);
169                                chargeSqRatio =  currentModel->GetChargeSquareRatio(theDirectPartDef,currentMaterial,E);
170                                theDirectEnergyLossProcess->SetDynamicMassCharge(massRatio,chargeSqRatio);
171                                x1= theDirectEnergyLossProcess->GetRange(E, currentCouple);
172                       
173                        } 
174                }
175                G4double r1 = theDirectEnergyLossProcess->GetRange(Tkin1, currentCouple);
176                G4double x1 = r1 + dlength;
177                G4double E1 = theDirectEnergyLossProcess->GetKineticEnergy(x1,currentCouple);
178                factor_dE=(E1-E)/(Tkin1-Tkin);
179                degain=E-Tkin; 
180               
181               
182               
183        }
184        //G4cout<<degain<<G4endl;
185        G4double tmax = currentModel->MaxSecondaryKinEnergy(dynParticle);
186        tmax = std::min(tmax,currentTcut);
187       
188       
189        dynParticle->SetKineticEnergy(Tkin+degain);
190
191        // Corrections, which cannot be tabulated for ions
192        //----------------------------------------
193        G4double esecdep=0;//not used in most models
194        currentModel->CorrectionsAlongStep(currentCouple, dynParticle, degain,esecdep, dlength); 
195
196        // Sample fluctuations
197        //-------------------
198       
199       
200        G4double deltaE =0.;
201        if (lossFluctuationFlag ) {
202                deltaE = currentModel->GetModelOfFluctuations()->
203                                                SampleFluctuations(currentMaterial,dynParticle,tmax,dlength,degain)-degain;
204        }
205       
206        G4double egain=degain+deltaE;
207        if (egain <=0) egain=degain;
208        Tkin+=egain;
209        dynParticle->SetKineticEnergy(Tkin);
210 }
211 
212 
213 
214
215 
216  delete dynParticle;
217 
218  if (IsIon){
219        chargeSqRatio =  currentModel->GetChargeSquareRatio(theDirectPartDef,currentMaterial,Tkin);
220        theDirectEnergyLossProcess->SetDynamicMassCharge(massRatio,chargeSqRatio);
221               
222  }
223 
224  G4double DEDX_after = theDirectEnergyLossProcess->GetDEDX(Tkin, currentCouple);
225 
226 
227  G4double weight_correction=DEDX_after/DEDX_before;
228 
229 
230  aParticleChange.ProposeEnergy(Tkin);
231 
232  //we still need to register in the particleChange the modification of the weight of the particle
233  G4double new_weight=weight_correction*track.GetWeight();
234  aParticleChange.SetParentWeightByProcess(false);
235  aParticleChange.ProposeParentWeight(new_weight);
236 
237
238  return &aParticleChange;
239
240}
241///////////////////////////////////////////////////////
242//
243void G4ContinuousGainOfEnergy::SetLossFluctuations(G4bool val)
244{
245  if(val && !lossFluctuationArePossible) return;
246  lossFluctuationFlag = val;
247}
248///////////////////////////////////////////////////////
249//
250
251
252
253G4double G4ContinuousGainOfEnergy::GetContinuousStepLimit(const G4Track& track,
254                G4double , G4double , G4double& )
255{ 
256  G4double x = DBL_MAX;
257  x=.1*mm;
258 
259 
260  DefineMaterial(track.GetMaterialCutsCouple());
261 
262  preStepKinEnergy = track.GetKineticEnergy(); 
263  preStepScaledKinEnergy = track.GetKineticEnergy()*massRatio;
264  currentModel = theDirectEnergyLossProcess->SelectModelForMaterial(preStepScaledKinEnergy,currentCoupleIndex);
265  G4double emax_model=currentModel->HighEnergyLimit();
266  if (IsIon) {
267        chargeSqRatio =  currentModel->GetChargeSquareRatio(theDirectPartDef,currentMaterial,preStepKinEnergy);
268        preStepChargeSqRatio = chargeSqRatio;
269        theDirectEnergyLossProcess->SetDynamicMassCharge(massRatio,preStepChargeSqRatio);
270  } 
271 
272 
273  G4double maxE =1.1*preStepKinEnergy;
274  /*if (preStepKinEnergy< 0.05*MeV) maxE =2.*preStepKinEnergy;
275  else if (preStepKinEnergy< 0.1*MeV) maxE =1.5*preStepKinEnergy;
276  else if (preStepKinEnergy< 0.5*MeV) maxE =1.25*preStepKinEnergy;*/
277   
278  if (preStepKinEnergy < currentTcut) maxE = std::min(currentTcut,maxE);
279 
280  maxE=std::min(emax_model*1.001,maxE);
281       
282  G4double r = theDirectEnergyLossProcess->GetRange(preStepKinEnergy, currentCouple);
283 
284  if (IsIon) {
285        G4double chargeSqRatioAtEmax = currentModel->GetChargeSquareRatio(theDirectPartDef,currentMaterial,maxE);
286        theDirectEnergyLossProcess->SetDynamicMassCharge(massRatio,chargeSqRatioAtEmax);
287  }     
288 
289  G4double r1 = theDirectEnergyLossProcess->GetRange(maxE, currentCouple);
290 
291  if (IsIon) theDirectEnergyLossProcess->SetDynamicMassCharge(massRatio,preStepChargeSqRatio);
292 
293 
294
295  x=r1-r;
296  x=std::max(r1-r,0.001*mm);
297 
298  return x;
299 
300 
301}
302#include "G4EmCorrections.hh"
303///////////////////////////////////////////////////////
304//
305
306void G4ContinuousGainOfEnergy::SetDynamicMassCharge(const G4Track& ,G4double energy)
307{ 
308
309  G4double ChargeSqRatio= G4LossTableManager::Instance()->EmCorrections()->EffectiveChargeSquareRatio(theDirectPartDef,currentMaterial,energy); 
310  if (theDirectEnergyLossProcess) theDirectEnergyLossProcess->SetDynamicMassCharge(massRatio,ChargeSqRatio);
311}
Note: See TracBrowser for help on using the repository browser.