Geant4_10
G4BetheBlochModel.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: G4BetheBlochModel.cc 74790 2013-10-22 07:31:37Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class header file
31 //
32 //
33 // File name: G4BetheBlochModel
34 //
35 // Author: Vladimir Ivanchenko on base of Laszlo Urban code
36 //
37 // Creation date: 03.01.2002
38 //
39 // Modifications:
40 //
41 // 04-12-02 Fix problem of G4DynamicParticle constructor (V.Ivanchenko)
42 // 23-12-02 Change interface in order to move to cut per region (V.Ivanchenko)
43 // 27-01-03 Make models region aware (V.Ivanchenko)
44 // 13-02-03 Add name (V.Ivanchenko)
45 // 24-03-05 Add G4EmCorrections (V.Ivanchenko)
46 // 11-04-05 Major optimisation of internal interfaces (V.Ivanchenko)
47 // 11-02-06 ComputeCrossSectionPerElectron, ComputeCrossSectionPerAtom (mma)
48 // 12-02-06 move G4LossTableManager::Instance()->EmCorrections()
49 // in constructor (mma)
50 // 12-08-08 Added methods GetParticleCharge, GetChargeSquareRatio,
51 // CorrectionsAlongStep needed for ions(V.Ivanchenko)
52 //
53 // -------------------------------------------------------------------
54 //
55 
56 
57 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
58 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
59 
60 #include "G4BetheBlochModel.hh"
61 #include "Randomize.hh"
62 #include "G4PhysicalConstants.hh"
63 #include "G4SystemOfUnits.hh"
64 #include "G4Electron.hh"
65 #include "G4LossTableManager.hh"
66 #include "G4EmCorrections.hh"
68 #include "G4Log.hh"
69 #include "G4DeltaAngle.hh"
70 
71 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
72 
73 using namespace std;
74 
76  const G4String& nam)
77  : G4VEmModel(nam),
78  particle(0),
79  tlimit(DBL_MAX),
80  twoln10(2.0*G4Log(10.0)),
81  bg2lim(0.0169),
82  taulim(8.4146e-3),
83  isIon(false),
84  isInitialised(false)
85 {
86  fParticleChange = 0;
87  theElectron = G4Electron::Electron();
88  if(p) {
89  SetGenericIon(p);
90  SetParticle(p);
91  } else {
92  SetParticle(theElectron);
93  }
95  nist = G4NistManager::Instance();
97 }
98 
99 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
100 
102 {}
103 
104 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
105 
107  const G4DataVector&)
108 {
109  SetGenericIon(p);
110  SetParticle(p);
111 
112  //G4cout << "G4BetheBlochModel::Initialise for " << p->GetParticleName()
113  // << " isIon= " << isIon
114  // << G4endl;
115 
116  // always false before the run
117  SetDeexcitationFlag(false);
118 
119  if(!isInitialised) {
120  isInitialised = true;
121  fParticleChange = GetParticleChangeForLoss();
124  }
125  }
126 }
127 
128 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
129 
131  const G4Material* mat,
132  G4double kineticEnergy)
133 {
134  // this method is called only for ions
135  G4double q2 = corr->EffectiveChargeSquareRatio(p,mat,kineticEnergy);
136  corrFactor = q2*corr->EffectiveChargeCorrection(p,mat,kineticEnergy);
137  return corrFactor;
138 }
139 
140 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
141 
143  const G4Material* mat,
144  G4double kineticEnergy)
145 {
146  // this method is called only for ions, so no check if it is an ion
147  return corr->GetParticleCharge(p,mat,kineticEnergy);
148 }
149 
150 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
151 
152 void G4BetheBlochModel::SetupParameters()
153 {
154  mass = particle->GetPDGMass();
155  spin = particle->GetPDGSpin();
156  G4double q = particle->GetPDGCharge()/eplus;
157  chargeSquare = q*q;
158  corrFactor = chargeSquare;
159  ratio = electron_mass_c2/mass;
160  G4double magmom =
161  particle->GetPDGMagneticMoment()*mass/(0.5*eplus*hbar_Planck*c_squared);
162  magMoment2 = magmom*magmom - 1.0;
163  formfact = 0.0;
164  if(particle->GetLeptonNumber() == 0) {
165  G4double x = 0.8426*GeV;
166  if(spin == 0.0 && mass < GeV) {x = 0.736*GeV;}
167  else if(mass > GeV) {
168  x /= nist->GetZ13(mass/proton_mass_c2);
169  // tlimit = 51.2*GeV*A13[iz]*A13[iz];
170  }
171  formfact = 2.0*electron_mass_c2/(x*x);
172  tlimit = 2.0/formfact;
173  }
174 }
175 
176 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
177 
179  const G4MaterialCutsCouple* couple)
180 {
181  return couple->GetMaterial()->GetIonisation()->GetMeanExcitationEnergy();
182 }
183 
184 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
185 
186 G4double
188  G4double kineticEnergy,
189  G4double cutEnergy,
190  G4double maxKinEnergy)
191 {
192  G4double cross = 0.0;
193  G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
194  G4double maxEnergy = min(tmax,maxKinEnergy);
195  if(cutEnergy < maxEnergy) {
196 
197  G4double totEnergy = kineticEnergy + mass;
198  G4double energy2 = totEnergy*totEnergy;
199  G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
200 
201  cross = 1.0/cutEnergy - 1.0/maxEnergy
202  - beta2*G4Log(maxEnergy/cutEnergy)/tmax;
203 
204  // +term for spin=1/2 particle
205  if( 0.5 == spin ) { cross += 0.5*(maxEnergy - cutEnergy)/energy2; }
206 
207  // High order correction different for hadrons and ions
208  // nevetheless they are applied to reduce high energy transfers
209  // if(!isIon)
210  //cross += corr->FiniteSizeCorrectionXS(p,currentMaterial,
211  // kineticEnergy,cutEnergy);
212 
213  cross *= twopi_mc2_rcl2*chargeSquare/beta2;
214  }
215 
216  // G4cout << "BB: e= " << kineticEnergy << " tmin= " << cutEnergy
217  // << " tmax= " << tmax << " cross= " << cross << G4endl;
218 
219  return cross;
220 }
221 
222 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
223 
225  const G4ParticleDefinition* p,
226  G4double kineticEnergy,
228  G4double cutEnergy,
229  G4double maxEnergy)
230 {
232  (p,kineticEnergy,cutEnergy,maxEnergy);
233  return cross;
234 }
235 
236 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
237 
239  const G4Material* material,
240  const G4ParticleDefinition* p,
241  G4double kineticEnergy,
242  G4double cutEnergy,
243  G4double maxEnergy)
244 {
245  G4double eDensity = material->GetElectronDensity();
247  (p,kineticEnergy,cutEnergy,maxEnergy);
248  return cross;
249 }
250 
251 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
252 
254  const G4ParticleDefinition* p,
255  G4double kineticEnergy,
256  G4double cut)
257 {
258  G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
259  G4double cutEnergy = std::min(cut,tmax);
260 
261  G4double tau = kineticEnergy/mass;
262  G4double gam = tau + 1.0;
263  G4double bg2 = tau * (tau+2.0);
264  G4double beta2 = bg2/(gam*gam);
265 
266  G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy();
267  G4double eexc2 = eexc*eexc;
268 
269  G4double eDensity = material->GetElectronDensity();
270 
271  G4double dedx = G4Log(2.0*electron_mass_c2*bg2*cutEnergy/eexc2)
272  - (1.0 + cutEnergy/tmax)*beta2;
273 
274  if(0.5 == spin) {
275  G4double del = 0.5*cutEnergy/(kineticEnergy + mass);
276  dedx += del*del;
277  }
278 
279  // density correction
280  G4double x = G4Log(bg2)/twoln10;
281  dedx -= material->GetIonisation()->DensityCorrection(x);
282 
283  // shell correction
284  dedx -= 2.0*corr->ShellCorrection(p,material,kineticEnergy);
285 
286  // now compute the total ionization loss
287  dedx *= twopi_mc2_rcl2*chargeSquare*eDensity/beta2;
288 
289  //High order correction different for hadrons and ions
290  if(isIon) {
291  dedx += corr->IonBarkasCorrection(p,material,kineticEnergy);
292  } else {
293  dedx += corr->HighOrderCorrections(p,material,kineticEnergy,cutEnergy);
294  }
295 
296  if (dedx < 0.0) { dedx = 0.0; }
297 
298  //G4cout << "E(MeV)= " << kineticEnergy/MeV << " dedx= " << dedx
299  // << " " << material->GetName() << G4endl;
300 
301  return dedx;
302 }
303 
304 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
305 
307  const G4DynamicParticle* dp,
308  G4double& eloss,
309  G4double&,
310  G4double length)
311 {
312  if(isIon) {
313  const G4ParticleDefinition* p = dp->GetDefinition();
314  const G4Material* mat = couple->GetMaterial();
315  G4double preKinEnergy = dp->GetKineticEnergy();
316  G4double e = preKinEnergy - eloss*0.5;
317  if(e < preKinEnergy*0.75) { e = preKinEnergy*0.75; }
318 
319  G4double q2 = corr->EffectiveChargeSquareRatio(p,mat,e);
321  G4double qfactor = q2*corr->EffectiveChargeCorrection(p,mat,e)/corrFactor;
322  G4double highOrder = length*corr->IonHighOrderCorrections(p,couple,e);
323  G4double elossnew = eloss*qfactor + highOrder;
324  if(elossnew > preKinEnergy) { elossnew = preKinEnergy; }
325  else if(elossnew < eloss*0.5) { elossnew = eloss*0.5; }
326  eloss = elossnew;
327  //G4cout << "G4BetheBlochModel::CorrectionsAlongStep: e= " << preKinEnergy
328  // << " qfactor= " << qfactor
329  // << " highOrder= " << highOrder << " ("
330  // << highOrder/eloss << ")" << G4endl;
331  }
332 }
333 
334 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
335 
336 void G4BetheBlochModel::SampleSecondaries(vector<G4DynamicParticle*>* vdp,
337  const G4MaterialCutsCouple* couple,
338  const G4DynamicParticle* dp,
339  G4double minKinEnergy,
340  G4double maxEnergy)
341 {
342  G4double kineticEnergy = dp->GetKineticEnergy();
343  G4double tmax = MaxSecondaryEnergy(dp->GetDefinition(),kineticEnergy);
344 
345  G4double maxKinEnergy = std::min(maxEnergy,tmax);
346  if(minKinEnergy >= maxKinEnergy) { return; }
347 
348  G4double totEnergy = kineticEnergy + mass;
349  G4double etot2 = totEnergy*totEnergy;
350  G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/etot2;
351 
352  G4double deltaKinEnergy, f;
353  G4double f1 = 0.0;
354  G4double fmax = 1.0;
355  if( 0.5 == spin ) { fmax += 0.5*maxKinEnergy*maxKinEnergy/etot2; }
356 
357  // sampling without nuclear size effect
358  do {
359  G4double q = G4UniformRand();
360  deltaKinEnergy = minKinEnergy*maxKinEnergy
361  /(minKinEnergy*(1.0 - q) + maxKinEnergy*q);
362 
363  f = 1.0 - beta2*deltaKinEnergy/tmax;
364  if( 0.5 == spin ) {
365  f1 = 0.5*deltaKinEnergy*deltaKinEnergy/etot2;
366  f += f1;
367  }
368 
369  } while( fmax*G4UniformRand() > f);
370 
371  // projectile formfactor - suppresion of high energy
372  // delta-electron production at high energy
373 
374  G4double x = formfact*deltaKinEnergy;
375  if(x > 1.e-6) {
376 
377  G4double x1 = 1.0 + x;
378  G4double grej = 1.0/(x1*x1);
379  if( 0.5 == spin ) {
380  G4double x2 = 0.5*electron_mass_c2*deltaKinEnergy/(mass*mass);
381  grej *= (1.0 + magMoment2*(x2 - f1/f)/(1.0 + x2));
382  }
383  if(grej > 1.1) {
384  G4cout << "### G4BetheBlochModel WARNING: grej= " << grej
385  << " " << dp->GetDefinition()->GetParticleName()
386  << " Ekin(MeV)= " << kineticEnergy
387  << " delEkin(MeV)= " << deltaKinEnergy
388  << G4endl;
389  }
390  if(G4UniformRand() > grej) return;
391  }
392 
393  G4ThreeVector deltaDirection;
394 
396 
397  const G4Material* mat = couple->GetMaterial();
399 
400  deltaDirection =
401  GetAngularDistribution()->SampleDirection(dp, deltaKinEnergy, Z, mat);
402 
403  } else {
404 
405  G4double deltaMomentum =
406  sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
407  G4double cost = deltaKinEnergy * (totEnergy + electron_mass_c2) /
408  (deltaMomentum * dp->GetTotalMomentum());
409  if(cost > 1.0) { cost = 1.0; }
410  G4double sint = sqrt((1.0 - cost)*(1.0 + cost));
411 
412  G4double phi = twopi * G4UniformRand() ;
413 
414  deltaDirection.set(sint*cos(phi),sint*sin(phi), cost) ;
415  deltaDirection.rotateUz(dp->GetMomentumDirection());
416  }
417 
418  /*
419  G4cout << "### G4BetheBlochModel "
420  << dp->GetDefinition()->GetParticleName()
421  << " Ekin(MeV)= " << kineticEnergy
422  << " delEkin(MeV)= " << deltaKinEnergy
423  << " tmin(MeV)= " << minKinEnergy
424  << " tmax(MeV)= " << maxKinEnergy
425  << " dir= " << dp->GetMomentumDirection()
426  << " dirDelta= " << deltaDirection
427  << G4endl;
428  }
429  */
430 
431  // create G4DynamicParticle object for delta ray
432  G4DynamicParticle* delta =
433  new G4DynamicParticle(theElectron,deltaDirection,deltaKinEnergy);
434 
435  vdp->push_back(delta);
436 
437  // Change kinematics of primary particle
438  kineticEnergy -= deltaKinEnergy;
439  G4ThreeVector finalP = dp->GetMomentum() - delta->GetMomentum();
440  finalP = finalP.unit();
441 
442  fParticleChange->SetProposedKineticEnergy(kineticEnergy);
443  fParticleChange->SetProposedMomentumDirection(finalP);
444 }
445 
446 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
447 
449  G4double kinEnergy)
450 {
451  // here particle type is checked for any method
452  SetParticle(pd);
453  G4double tau = kinEnergy/mass;
454  G4double tmax = 2.0*electron_mass_c2*tau*(tau + 2.) /
455  (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
456  return std::min(tmax,tlimit);
457 }
458 
459 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
void set(double x, double y, double z)
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:224
G4double HighOrderCorrections(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy, G4double cutEnergy)
static G4LossTableManager * Instance()
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:107
Double_t x2[nxs]
Definition: Style.C:19
G4double GetKineticEnergy() const
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
G4double EffectiveChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
const char * p
Definition: xmltok.h:285
G4double GetZ13(G4double Z)
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:578
G4double GetParticleCharge(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double EffectiveChargeCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4ParticleDefinition * GetDefinition() const
tuple x
Definition: test.py:50
G4VEmFluctuationModel * GetModelOfFluctuations()
Definition: G4VEmModel.hh:571
int G4int
Definition: G4Types.hh:78
static G4NistManager * Instance()
const G4String & GetParticleName() const
TFile f
Definition: plotHisto.C:6
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy)
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *couple)
Float_t mat
Definition: plot.C:40
G4double GetTotalMomentum() const
G4BetheBlochModel(const G4ParticleDefinition *p=0, const G4String &nam="BetheBloch")
string material
Definition: eplot.py:19
Float_t f1
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy)
virtual G4double GetParticleCharge(const G4ParticleDefinition *p, const G4Material *mat, G4double kineticEnergy)
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
Float_t Z
Definition: plot.C:39
G4double GetElectronDensity() const
Definition: G4Material.hh:215
G4bool UseAngularGeneratorFlag() const
Definition: G4VEmModel.hh:655
G4double IonBarkasCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4EmCorrections * EmCorrections()
const G4ThreeVector & GetMomentumDirection() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
Double_t x1[nxs]
Definition: Style.C:18
void SetProposedKineticEnergy(G4double proposedKinEnergy)
G4double IonHighOrderCorrections(const G4ParticleDefinition *, const G4MaterialCutsCouple *, G4double kineticEnergy)
float proton_mass_c2
Definition: hepunit.py:275
float electron_mass_c2
Definition: hepunit.py:274
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
void SetProposedMomentumDirection(const G4ThreeVector &dir)
G4double G4Log(G4double x)
Definition: G4Log.hh:227
virtual void CorrectionsAlongStep(const G4MaterialCutsCouple *couple, const G4DynamicParticle *dp, G4double &eloss, G4double &, G4double length)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
G4double ShellCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double DensityCorrection(G4double x)
G4double GetPDGMass() const
Hep3Vector unit() const
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:585
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4double GetPDGSpin() const
G4double GetMeanExcitationEnergy() const
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
G4double GetPDGMagneticMoment() const
double G4double
Definition: G4Types.hh:76
G4double GetChargeSquareRatio() const
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:690
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:739
G4double GetPDGCharge() const
virtual void SetParticleAndCharge(const G4ParticleDefinition *, G4double q2)
#define DBL_MAX
Definition: templates.hh:83
G4ThreeVector GetMomentum() const
const G4Material * GetMaterial() const
G4int SelectRandomAtomNumber(const G4Material *)
Definition: G4VEmModel.hh:529