Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4mplIonisationWithDeltaModel.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: G4mplIonisationWithDeltaModel.cc 66996 2013-01-29 14:50:52Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class header file
31 //
32 //
33 // File name: G4mplIonisationWithDeltaModel
34 //
35 // Author: Vladimir Ivanchenko
36 //
37 // Creation date: 06.09.2005
38 //
39 // Modifications:
40 // 12.08.2007 Changing low energy approximation and extrapolation.
41 // Small bug fixing and refactoring (M. Vladymyrov)
42 // 13.11.2007 Use low-energy asymptotic from [3] (V.Ivanchenko)
43 //
44 //
45 // -------------------------------------------------------------------
46 // References
47 // [1] Steven P. Ahlen: Energy loss of relativistic heavy ionizing particles,
48 // S.P. Ahlen, Rev. Mod. Phys 52(1980), p121
49 // [2] K.A. Milton arXiv:hep-ex/0602040
50 // [3] S.P. Ahlen and K. Kinoshita, Phys. Rev. D26 (1982) 2347
51 
52 
53 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
54 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
55 
57 #include "Randomize.hh"
58 #include "G4PhysicalConstants.hh"
59 #include "G4SystemOfUnits.hh"
60 #include "G4LossTableManager.hh"
62 #include "G4Electron.hh"
63 #include "G4DynamicParticle.hh"
64 
65 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
66 
67 using namespace std;
68 
70  const G4String& nam)
72  magCharge(mCharge),
73  twoln10(log(100.0)),
74  betalow(0.01),
75  betalim(0.1),
76  beta2lim(betalim*betalim),
77  bg2lim(beta2lim*(1.0 + beta2lim))
78 {
79  nmpl = G4lrint(std::fabs(magCharge) * 2 * fine_structure_const);
80  if(nmpl > 6) { nmpl = 6; }
81  else if(nmpl < 1) { nmpl = 1; }
82  pi_hbarc2_over_mc2 = pi * hbarc * hbarc / electron_mass_c2;
83  chargeSquare = magCharge * magCharge;
84  dedxlim = 45.*nmpl*nmpl*GeV*cm2/g;
85  fParticleChange = 0;
86  theElectron = G4Electron::Electron();
87  G4cout << "### Monopole ionisation model with d-electron production, Gmag= "
88  << magCharge/eplus << G4endl;
89  monopole = 0;
90  mass = 0.0;
91 }
92 
93 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
94 
96 {}
97 
98 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
99 
101 {
102  monopole = p;
103  mass = monopole->GetPDGMass();
104  G4double emin =
105  std::min(LowEnergyLimit(),0.1*mass*(1/sqrt(1 - betalow*betalow) - 1));
106  G4double emax =
107  std::max(HighEnergyLimit(),10*mass*(1/sqrt(1 - beta2lim) - 1));
108  SetLowEnergyLimit(emin);
109  SetHighEnergyLimit(emax);
110 }
111 
112 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
113 
114 void
116  const G4DataVector&)
117 {
118  if(!monopole) { SetParticle(p); }
119  if(!fParticleChange) { fParticleChange = GetParticleChangeForLoss(); }
120 }
121 
122 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
123 
124 G4double
126  const G4ParticleDefinition* p,
127  G4double kineticEnergy,
128  G4double maxEnergy)
129 {
130  if(!monopole) { SetParticle(p); }
131  G4double tmax = MaxSecondaryEnergy(p,kineticEnergy);
132  G4double cutEnergy = std::min(tmax, maxEnergy);
133  cutEnergy = std::max(LowEnergyLimit(), cutEnergy);
134  G4double tau = kineticEnergy / mass;
135  G4double gam = tau + 1.0;
136  G4double bg2 = tau * (tau + 2.0);
137  G4double beta2 = bg2 / (gam * gam);
138  G4double beta = sqrt(beta2);
139 
140  // low-energy asymptotic formula
141  G4double dedx = dedxlim*beta*material->GetDensity();
142 
143  // above asymptotic
144  if(beta > betalow) {
145 
146  // high energy
147  if(beta >= betalim) {
148  dedx = ComputeDEDXAhlen(material, bg2, cutEnergy);
149 
150  } else {
151 
152  G4double dedx1 = dedxlim*betalow*material->GetDensity();
153  G4double dedx2 = ComputeDEDXAhlen(material, bg2lim, cutEnergy);
154 
155  // extrapolation between two formula
156  G4double kapa2 = beta - betalow;
157  G4double kapa1 = betalim - beta;
158  dedx = (kapa1*dedx1 + kapa2*dedx2)/(kapa1 + kapa2);
159  }
160  }
161  return dedx;
162 }
163 
164 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
165 
166 G4double
167 G4mplIonisationWithDeltaModel::ComputeDEDXAhlen(const G4Material* material,
168  G4double bg2,
169  G4double cutEnergy)
170 {
171  G4double eDensity = material->GetElectronDensity();
172  G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy();
173 
174  // Ahlen's formula for nonconductors, [1]p157, f(5.7)
175  G4double dedx =
176  0.5*(log(2.0 * electron_mass_c2 * bg2*cutEnergy / (eexc*eexc)) - 1.0);
177 
178  // Kazama et al. cross-section correction
179  G4double k = 0.406;
180  if(nmpl > 1) { k = 0.346; }
181 
182  // Bloch correction
183  const G4double B[7] = { 0.0, 0.248, 0.672, 1.022, 1.243, 1.464, 1.685};
184 
185  dedx += 0.5 * k - B[nmpl];
186 
187  // density effect correction
188  G4double x = log(bg2)/twoln10;
189  dedx -= material->GetIonisation()->DensityCorrection(x);
190 
191  // now compute the total ionization loss
192  dedx *= pi_hbarc2_over_mc2 * eDensity * nmpl * nmpl;
193 
194  if (dedx < 0.0) { dedx = 0; }
195  return dedx;
196 }
197 
198 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
199 
200 G4double
202  const G4ParticleDefinition* p,
203  G4double kineticEnergy,
204  G4double cut,
205  G4double maxKinEnergy)
206 {
207  if(!monopole) { SetParticle(p); }
208  G4double cross = 0.0;
209  G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
210  G4double maxEnergy = std::min(tmax,maxKinEnergy);
211  G4double cutEnergy = std::max(LowEnergyLimit(), cut);
212  if(cutEnergy < maxEnergy) {
213  cross = (0.5/cutEnergy - 0.5/maxEnergy)*pi_hbarc2_over_mc2 * nmpl * nmpl;
214  }
215  return cross;
216 }
217 
218 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
219 
220 G4double
222  const G4ParticleDefinition* p,
223  G4double kineticEnergy,
225  G4double cutEnergy,
226  G4double maxEnergy)
227 {
228  G4double cross =
229  Z*ComputeCrossSectionPerElectron(p,kineticEnergy,cutEnergy,maxEnergy);
230  return cross;
231 }
232 
233 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
234 
235 void
236 G4mplIonisationWithDeltaModel::SampleSecondaries(vector<G4DynamicParticle*>* vdp,
237  const G4MaterialCutsCouple*,
238  const G4DynamicParticle* dp,
239  G4double minKinEnergy,
240  G4double maxEnergy)
241 {
242  G4double kineticEnergy = dp->GetKineticEnergy();
243  G4double tmax = MaxSecondaryEnergy(dp->GetDefinition(),kineticEnergy);
244 
245  G4double maxKinEnergy = std::min(maxEnergy,tmax);
246  if(minKinEnergy >= maxKinEnergy) { return; }
247 
248  //G4cout << "G4mplIonisationWithDeltaModel::SampleSecondaries: E(GeV)= "
249  // << kineticEnergy/GeV << " M(GeV)= " << mass/GeV
250  // << " tmin(MeV)= " << minKinEnergy/MeV << G4endl;
251 
252  G4double totEnergy = kineticEnergy + mass;
253  G4double etot2 = totEnergy*totEnergy;
254  G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/etot2;
255 
256  // sampling without nuclear size effect
257  G4double q = G4UniformRand();
258  G4double deltaKinEnergy = minKinEnergy*maxKinEnergy
259  /(minKinEnergy*(1.0 - q) + maxKinEnergy*q);
260 
261  // delta-electron is produced
262  G4double totMomentum = totEnergy*sqrt(beta2);
263  G4double deltaMomentum =
264  sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
265  G4double cost = deltaKinEnergy * (totEnergy + electron_mass_c2) /
266  (deltaMomentum * totMomentum);
267  if(cost > 1.0) { cost = 1.0; }
268 
269  G4double sint = sqrt((1.0 - cost)*(1.0 + cost));
270 
271  G4double phi = twopi * G4UniformRand() ;
272 
273  G4ThreeVector deltaDirection(sint*cos(phi),sint*sin(phi), cost);
274  G4ThreeVector direction = dp->GetMomentumDirection();
275  deltaDirection.rotateUz(direction);
276 
277  // create G4DynamicParticle object for delta ray
278  G4DynamicParticle* delta =
279  new G4DynamicParticle(theElectron,deltaDirection,deltaKinEnergy);
280 
281  vdp->push_back(delta);
282 
283  // Change kinematics of primary particle
284  kineticEnergy -= deltaKinEnergy;
285  G4ThreeVector finalP = direction*totMomentum - deltaDirection*deltaMomentum;
286  finalP = finalP.unit();
287 
288  fParticleChange->SetProposedKineticEnergy(kineticEnergy);
289  fParticleChange->SetProposedMomentumDirection(finalP);
290 }
291 
292 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
293 
295  const G4Material* material,
296  const G4DynamicParticle* dp,
297  G4double& tmax,
298  G4double& length,
299  G4double& meanLoss)
300 {
301  G4double siga = Dispersion(material,dp,tmax,length);
302  G4double loss = meanLoss;
303  siga = sqrt(siga);
304  G4double twomeanLoss = meanLoss + meanLoss;
305 
306  if(twomeanLoss < siga) {
307  G4double x;
308  do {
309  loss = twomeanLoss*G4UniformRand();
310  x = (loss - meanLoss)/siga;
311  } while (1.0 - 0.5*x*x < G4UniformRand());
312  } else {
313  do {
314  loss = G4RandGauss::shoot(meanLoss,siga);
315  } while (0.0 > loss || loss > twomeanLoss);
316  }
317  return loss;
318 }
319 
320 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
321 
322 G4double
324  const G4DynamicParticle* dp,
325  G4double& tmax,
326  G4double& length)
327 {
328  G4double siga = 0.0;
329  G4double tau = dp->GetKineticEnergy()/mass;
330  if(tau > 0.0) {
331  G4double electronDensity = material->GetElectronDensity();
332  G4double gam = tau + 1.0;
333  G4double invbeta2 = (gam*gam)/(tau * (tau+2.0));
334  siga = (invbeta2 - 0.5) * twopi_mc2_rcl2 * tmax * length
335  * electronDensity * chargeSquare;
336  }
337  return siga;
338 }
339 
340 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
341 
342 G4double
344  G4double kinEnergy)
345 {
346  G4double tau = kinEnergy/mass;
347  return 2.0*electron_mass_c2*tau*(tau + 2.);
348 }
349 
350 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....