Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4MuBetheBlochModel.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: G4MuBetheBlochModel.cc 74019 2013-09-19 13:35:34Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class header file
31 //
32 //
33 // File name: G4MuBetheBlochModel
34 //
35 // Author: Vladimir Ivanchenko on base of Laszlo Urban code
36 //
37 // Creation date: 09.08.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 // 10-02-04 Calculation of radiative corrections using R.Kokoulin model (V.I)
46 // 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
47 // 12-04-05 Add usage of G4EmCorrections (V.Ivanchenko)
48 // 13-02-06 ComputeCrossSectionPerElectron, ComputeCrossSectionPerAtom (mma)
49 //
50 
51 //
52 // -------------------------------------------------------------------
53 //
54 
55 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
56 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
57 
58 #include "G4MuBetheBlochModel.hh"
59 #include "G4PhysicalConstants.hh"
60 #include "G4SystemOfUnits.hh"
61 #include "Randomize.hh"
62 #include "G4Electron.hh"
63 #include "G4LossTableManager.hh"
64 #include "G4EmCorrections.hh"
66 #include "G4Log.hh"
67 #include "G4Exp.hh"
68 
69 G4double G4MuBetheBlochModel::xgi[]={ 0.0199, 0.1017, 0.2372, 0.4083, 0.5917,
70  0.7628, 0.8983, 0.9801 };
71 
72 G4double G4MuBetheBlochModel::wgi[]={ 0.0506, 0.1112, 0.1569, 0.1813, 0.1813,
73  0.1569, 0.1112, 0.0506 };
74 
75 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
76 
77 using namespace std;
78 
80  const G4String& nam)
81  : G4VEmModel(nam),
82  particle(0),
83  limitKinEnergy(100.*keV),
84  logLimitKinEnergy(G4Log(limitKinEnergy)),
85  twoln10(2.0*G4Log(10.0)),
86  bg2lim(0.0169),
87  taulim(8.4146e-3),
88  alphaprime(fine_structure_const/twopi)
89 {
90  theElectron = G4Electron::Electron();
92  fParticleChange = 0;
93 
94  // initial initialisation of memeber should be overwritten
95  // by SetParticle
96  mass = massSquare = ratio = 1.0;
97 
98  if(p) { SetParticle(p); }
99 }
100 
101 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
102 
104 {}
105 
106 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
107 
109  const G4MaterialCutsCouple* couple)
110 {
111  return couple->GetMaterial()->GetIonisation()->GetMeanExcitationEnergy();
112 }
113 
114 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
115 
117  G4double kinEnergy)
118 {
119  G4double tau = kinEnergy/mass;
120  G4double tmax = 2.0*electron_mass_c2*tau*(tau + 2.) /
121  (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
122  return tmax;
123 }
124 
125 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
126 
128  const G4DataVector&)
129 {
130  if(p) { SetParticle(p); }
131  if(!fParticleChange) { fParticleChange = GetParticleChangeForLoss(); }
132 }
133 
134 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
135 
137  const G4ParticleDefinition* p,
138  G4double kineticEnergy,
139  G4double cutEnergy,
140  G4double maxKinEnergy)
141 {
142  G4double cross = 0.0;
143  G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
144  G4double maxEnergy = min(tmax,maxKinEnergy);
145  if(cutEnergy < maxEnergy) {
146 
147  G4double totEnergy = kineticEnergy + mass;
148  G4double energy2 = totEnergy*totEnergy;
149  G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
150 
151  cross = 1.0/cutEnergy - 1.0/maxEnergy - beta2*G4Log(maxEnergy/cutEnergy)/tmax
152  + 0.5*(maxEnergy - cutEnergy)/energy2;
153 
154  // radiative corrections of R. Kokoulin
155  if (maxEnergy > limitKinEnergy) {
156 
157  G4double logtmax = G4Log(maxEnergy);
158  G4double logtmin = G4Log(max(cutEnergy,limitKinEnergy));
159  G4double logstep = logtmax - logtmin;
160  G4double dcross = 0.0;
161 
162  for (G4int ll=0; ll<8; ll++)
163  {
164  G4double ep = G4Exp(logtmin + xgi[ll]*logstep);
165  G4double a1 = G4Log(1.0 + 2.0*ep/electron_mass_c2);
166  G4double a3 = G4Log(4.0*totEnergy*(totEnergy - ep)/massSquare);
167  dcross += wgi[ll]*(1.0/ep - beta2/tmax + 0.5*ep/energy2)*a1*(a3 - a1);
168  }
169 
170  cross += dcross*logstep*alphaprime;
171  }
172 
173  cross *= twopi_mc2_rcl2/beta2;
174 
175  }
176 
177  // G4cout << "tmin= " << cutEnergy << " tmax= " << tmax
178  // << " cross= " << cross << G4endl;
179 
180  return cross;
181 }
182 
183 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
184 
186  const G4ParticleDefinition* p,
187  G4double kineticEnergy,
188  G4double Z, G4double,
189  G4double cutEnergy,
190  G4double maxEnergy)
191 {
193  (p,kineticEnergy,cutEnergy,maxEnergy);
194  return cross;
195 }
196 
197 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
198 
200  const G4Material* material,
201  const G4ParticleDefinition* p,
202  G4double kineticEnergy,
203  G4double cutEnergy,
204  G4double maxEnergy)
205 {
206  G4double eDensity = material->GetElectronDensity();
208  (p,kineticEnergy,cutEnergy,maxEnergy);
209  return cross;
210 }
211 
212 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
213 
215  const G4ParticleDefinition* p,
216  G4double kineticEnergy,
217  G4double cut)
218 {
219  G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
220  G4double tau = kineticEnergy/mass;
221  G4double cutEnergy = min(cut,tmax);
222  G4double gam = tau + 1.0;
223  G4double bg2 = tau * (tau+2.0);
224  G4double beta2 = bg2/(gam*gam);
225 
226  G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy();
227  G4double eexc2 = eexc*eexc;
228  //G4double cden = material->GetIonisation()->GetCdensity();
229  //G4double mden = material->GetIonisation()->GetMdensity();
230  //G4double aden = material->GetIonisation()->GetAdensity();
231  //G4double x0den = material->GetIonisation()->GetX0density();
232  //G4double x1den = material->GetIonisation()->GetX1density();
233 
234  G4double eDensity = material->GetElectronDensity();
235 
236  G4double dedx = G4Log(2.0*electron_mass_c2*bg2*cutEnergy/eexc2)
237  -(1.0 + cutEnergy/tmax)*beta2;
238 
239  G4double totEnergy = kineticEnergy + mass;
240  G4double del = 0.5*cutEnergy/totEnergy;
241  dedx += del*del;
242 
243  // density correction
244  G4double x = G4Log(bg2)/twoln10;
245  //if ( x >= x0den ) {
246  // dedx -= twoln10*x - cden ;
247  // if ( x < x1den ) dedx -= aden*pow((x1den-x),mden) ;
248  //}
249  dedx -= material->GetIonisation()->DensityCorrection(x);
250 
251  // shell correction
252  dedx -= 2.0*corr->ShellCorrection(p,material,kineticEnergy);
253 
254  // now compute the total ionization loss
255 
256  if (dedx < 0.0) dedx = 0.0 ;
257 
258  // radiative corrections of R. Kokoulin
259  if (cutEnergy > limitKinEnergy) {
260 
261  G4double logtmax = G4Log(cutEnergy);
262  G4double logstep = logtmax - logLimitKinEnergy;
263  G4double dloss = 0.0;
264  G4double ftot2= 0.5/(totEnergy*totEnergy);
265 
266  for (G4int ll=0; ll<8; ll++)
267  {
268  G4double ep = G4Exp(logLimitKinEnergy + xgi[ll]*logstep);
269  G4double a1 = G4Log(1.0 + 2.0*ep/electron_mass_c2);
270  G4double a3 = G4Log(4.0*totEnergy*(totEnergy - ep)/massSquare);
271  dloss += wgi[ll]*(1.0 - beta2*ep/tmax + ep*ep*ftot2)*a1*(a3 - a1);
272  }
273  dedx += dloss*logstep*alphaprime;
274  }
275 
276  dedx *= twopi_mc2_rcl2*eDensity/beta2;
277 
278  //High order corrections
279  dedx += corr->HighOrderCorrections(p,material,kineticEnergy,cutEnergy);
280 
281  return dedx;
282 }
283 
284 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
285 
286 void G4MuBetheBlochModel::SampleSecondaries(vector<G4DynamicParticle*>* vdp,
287  const G4MaterialCutsCouple*,
288  const G4DynamicParticle* dp,
289  G4double minKinEnergy,
290  G4double maxEnergy)
291 {
292  G4double tmax = MaxSecondaryKinEnergy(dp);
293  G4double maxKinEnergy = min(maxEnergy,tmax);
294  if(minKinEnergy >= maxKinEnergy) { return; }
295 
296  G4double kineticEnergy = dp->GetKineticEnergy();
297  G4double totEnergy = kineticEnergy + mass;
298  G4double etot2 = totEnergy*totEnergy;
299  G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/etot2;
300 
301  G4double grej = 1.;
302  if(tmax > limitKinEnergy) {
303  G4double a0 = G4Log(2.*totEnergy/mass);
304  grej += alphaprime*a0*a0;
305  }
306 
307  G4double deltaKinEnergy, f;
308 
309  // sampling follows ...
310  do {
311  G4double q = G4UniformRand();
312  deltaKinEnergy = minKinEnergy*maxKinEnergy
313  /(minKinEnergy*(1.0 - q) + maxKinEnergy*q);
314 
315 
316  f = 1.0 - beta2*deltaKinEnergy/tmax
317  + 0.5*deltaKinEnergy*deltaKinEnergy/etot2;
318 
319  if(deltaKinEnergy > limitKinEnergy) {
320  G4double a1 = G4Log(1.0 + 2.0*deltaKinEnergy/electron_mass_c2);
321  G4double a3 = G4Log(4.0*totEnergy*(totEnergy - deltaKinEnergy)/massSquare);
322  f *= (1. + alphaprime*a1*(a3 - a1));
323  }
324 
325  if(f > grej) {
326  G4cout << "G4MuBetheBlochModel::SampleSecondary Warning! "
327  << "Majorant " << grej << " < "
328  << f << " for edelta= " << deltaKinEnergy
329  << " tmin= " << minKinEnergy << " max= " << maxKinEnergy
330  << G4endl;
331  }
332 
333 
334  } while( grej*G4UniformRand() > f );
335 
336  G4double deltaMomentum =
337  sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
338  G4double totalMomentum = totEnergy*sqrt(beta2);
339  G4double cost = deltaKinEnergy * (totEnergy + electron_mass_c2) /
340  (deltaMomentum * totalMomentum);
341 
342  G4double sint = sqrt(1.0 - cost*cost);
343 
344  G4double phi = twopi * G4UniformRand() ;
345 
346  G4ThreeVector deltaDirection(sint*cos(phi),sint*sin(phi), cost) ;
347  G4ThreeVector direction = dp->GetMomentumDirection();
348  deltaDirection.rotateUz(direction);
349 
350  // primary change
351  kineticEnergy -= deltaKinEnergy;
352  G4ThreeVector dir = totalMomentum*direction - deltaMomentum*deltaDirection;
353  direction = dir.unit();
354  fParticleChange->SetProposedKineticEnergy(kineticEnergy);
355  fParticleChange->SetProposedMomentumDirection(direction);
356 
357  // create G4DynamicParticle object for delta ray
358  G4DynamicParticle* delta = new G4DynamicParticle(theElectron,
359  deltaDirection,deltaKinEnergy);
360  vdp->push_back(delta);
361 }
362 
363 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4double MaxSecondaryKinEnergy(const G4DynamicParticle *dynParticle)
Definition: G4VEmModel.hh:448
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
G4double GetKineticEnergy() const
const char * p
Definition: xmltok.h:285
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *)
int G4int
Definition: G4Types.hh:78
string material
Definition: eplot.py:19
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
G4double GetElectronDensity() const
Definition: G4Material.hh:215
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy)
G4EmCorrections * EmCorrections()
const G4ThreeVector & GetMomentumDirection() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
void SetProposedKineticEnergy(G4double proposedKinEnergy)
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
float electron_mass_c2
Definition: hepunit.py:274
function dcross(V1, V2)
Definition: leptonew.f:3567
void SetProposedMomentumDirection(const G4ThreeVector &dir)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
G4double G4Log(G4double x)
Definition: G4Log.hh:227
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
G4double ShellCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double DensityCorrection(G4double x)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
Hep3Vector unit() const
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4double GetMeanExcitationEnergy() const
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
G4MuBetheBlochModel(const G4ParticleDefinition *p=0, const G4String &nam="MuBetheBloch")
double G4double
Definition: G4Types.hh:76
const G4Material * GetMaterial() const