Geant4_10
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 76600 2013-11-13 08:30:02Z 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 #include "G4ProductionCutsTable.hh"
65 #include "G4MaterialCutsCouple.hh"
66 #include "G4Log.hh"
67 
68 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
69 
70 using namespace std;
71 
72 std::vector<G4double>* G4mplIonisationWithDeltaModel::dedx0 = 0;
73 
75  const G4String& nam)
77  magCharge(mCharge),
78  twoln10(log(100.0)),
79  betalow(0.01),
80  betalim(0.1),
81  beta2lim(betalim*betalim),
82  bg2lim(beta2lim*(1.0 + beta2lim))
83 {
84  nmpl = G4lrint(std::fabs(magCharge) * 2 * fine_structure_const);
85  if(nmpl > 6) { nmpl = 6; }
86  else if(nmpl < 1) { nmpl = 1; }
87  pi_hbarc2_over_mc2 = pi * hbarc * hbarc / electron_mass_c2;
88  chargeSquare = magCharge * magCharge;
89  dedxlim = 45.*nmpl*nmpl*GeV*cm2/g;
90  fParticleChange = 0;
91  theElectron = G4Electron::Electron();
92  G4cout << "### Monopole ionisation model with d-electron production, Gmag= "
93  << magCharge/eplus << G4endl;
94  monopole = 0;
95  mass = 0.0;
96 }
97 
98 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
99 
101 {
102  if(IsMaster()) { delete dedx0; }
103 }
104 
105 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
106 
108 {
109  monopole = p;
110  mass = monopole->GetPDGMass();
111  G4double emin =
112  std::min(LowEnergyLimit(),0.1*mass*(1/sqrt(1 - betalow*betalow) - 1));
113  G4double emax =
114  std::max(HighEnergyLimit(),10*mass*(1/sqrt(1 - beta2lim) - 1));
115  SetLowEnergyLimit(emin);
116  SetHighEnergyLimit(emax);
117 }
118 
119 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
120 
121 void
123  const G4DataVector&)
124 {
125  if(!monopole) { SetParticle(p); }
126  if(!fParticleChange) { fParticleChange = GetParticleChangeForLoss(); }
127  if(IsMaster()) {
128  if(!dedx0) { dedx0 = new std::vector<G4double>; }
129  G4ProductionCutsTable* theCoupleTable=
131  G4int numOfCouples = theCoupleTable->GetTableSize();
132  G4int n = dedx0->size();
133  if(n < numOfCouples) { dedx0->resize(numOfCouples); }
134 
135  // initialise vector
136  for(G4int i=0; i<numOfCouples; ++i) {
137 
138  const G4Material* material =
139  theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
140  G4double eDensity = material->GetElectronDensity();
141  G4double vF = electron_Compton_length*pow(3*pi*pi*eDensity,0.3333333333);
142  (*dedx0)[i] = pi_hbarc2_over_mc2*eDensity*nmpl*nmpl*
143  (G4Log(2*vF/fine_structure_const) - 0.5)/vF;
144  }
145  }
146 }
147 
148 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
149 
150 G4double
152  const G4ParticleDefinition* p,
153  G4double kineticEnergy,
154  G4double maxEnergy)
155 {
156  if(!monopole) { SetParticle(p); }
157  G4double tmax = MaxSecondaryEnergy(p,kineticEnergy);
158  G4double cutEnergy = std::min(tmax, maxEnergy);
159  cutEnergy = std::max(LowEnergyLimit(), cutEnergy);
160  G4double tau = kineticEnergy / mass;
161  G4double gam = tau + 1.0;
162  G4double bg2 = tau * (tau + 2.0);
163  G4double beta2 = bg2 / (gam * gam);
164  G4double beta = sqrt(beta2);
165 
166  // low-energy asymptotic formula
167  //G4double dedx = dedxlim*beta*material->GetDensity();
168  G4double dedx = (*dedx0)[CurrentCouple()->GetIndex()]*beta;
169 
170  // above asymptotic
171  if(beta > betalow) {
172 
173  // high energy
174  if(beta >= betalim) {
175  dedx = ComputeDEDXAhlen(material, bg2, cutEnergy);
176 
177  } else {
178 
179  //G4double dedx1 = dedxlim*betalow*material->GetDensity();
180  G4double dedx1 = (*dedx0)[CurrentCouple()->GetIndex()]*betalow;
181  G4double dedx2 = ComputeDEDXAhlen(material, bg2lim, cutEnergy);
182 
183  // extrapolation between two formula
184  G4double kapa2 = beta - betalow;
185  G4double kapa1 = betalim - beta;
186  dedx = (kapa1*dedx1 + kapa2*dedx2)/(kapa1 + kapa2);
187  }
188  }
189  return dedx;
190 }
191 
192 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
193 
194 G4double
195 G4mplIonisationWithDeltaModel::ComputeDEDXAhlen(const G4Material* material,
196  G4double bg2,
197  G4double cutEnergy)
198 {
199  G4double eDensity = material->GetElectronDensity();
200  G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy();
201 
202  // Ahlen's formula for nonconductors, [1]p157, f(5.7)
203  G4double dedx =
204  0.5*(log(2.0 * electron_mass_c2 * bg2*cutEnergy / (eexc*eexc)) - 1.0);
205 
206  // Kazama et al. cross-section correction
207  G4double k = 0.406;
208  if(nmpl > 1) { k = 0.346; }
209 
210  // Bloch correction
211  const G4double B[7] = { 0.0, 0.248, 0.672, 1.022, 1.243, 1.464, 1.685};
212 
213  dedx += 0.5 * k - B[nmpl];
214 
215  // density effect correction
216  G4double x = G4Log(bg2)/twoln10;
217  dedx -= material->GetIonisation()->DensityCorrection(x);
218 
219  // now compute the total ionization loss
220  dedx *= pi_hbarc2_over_mc2 * eDensity * nmpl * nmpl;
221 
222  if (dedx < 0.0) { dedx = 0; }
223  return dedx;
224 }
225 
226 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
227 
228 G4double
230  const G4ParticleDefinition* p,
231  G4double kineticEnergy,
232  G4double cut,
233  G4double maxKinEnergy)
234 {
235  if(!monopole) { SetParticle(p); }
236  G4double cross = 0.0;
237  G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
238  G4double maxEnergy = std::min(tmax,maxKinEnergy);
239  G4double cutEnergy = std::max(LowEnergyLimit(), cut);
240  if(cutEnergy < maxEnergy) {
241  cross = (0.5/cutEnergy - 0.5/maxEnergy)*pi_hbarc2_over_mc2 * nmpl * nmpl;
242  }
243  return cross;
244 }
245 
246 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
247 
248 G4double
250  const G4ParticleDefinition* p,
251  G4double kineticEnergy,
253  G4double cutEnergy,
254  G4double maxEnergy)
255 {
256  G4double cross =
257  Z*ComputeCrossSectionPerElectron(p,kineticEnergy,cutEnergy,maxEnergy);
258  return cross;
259 }
260 
261 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
262 
263 void
264 G4mplIonisationWithDeltaModel::SampleSecondaries(vector<G4DynamicParticle*>* vdp,
265  const G4MaterialCutsCouple*,
266  const G4DynamicParticle* dp,
267  G4double minKinEnergy,
268  G4double maxEnergy)
269 {
270  G4double kineticEnergy = dp->GetKineticEnergy();
271  G4double tmax = MaxSecondaryEnergy(dp->GetDefinition(),kineticEnergy);
272 
273  G4double maxKinEnergy = std::min(maxEnergy,tmax);
274  if(minKinEnergy >= maxKinEnergy) { return; }
275 
276  //G4cout << "G4mplIonisationWithDeltaModel::SampleSecondaries: E(GeV)= "
277  // << kineticEnergy/GeV << " M(GeV)= " << mass/GeV
278  // << " tmin(MeV)= " << minKinEnergy/MeV << G4endl;
279 
280  G4double totEnergy = kineticEnergy + mass;
281  G4double etot2 = totEnergy*totEnergy;
282  G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/etot2;
283 
284  // sampling without nuclear size effect
285  G4double q = G4UniformRand();
286  G4double deltaKinEnergy = minKinEnergy*maxKinEnergy
287  /(minKinEnergy*(1.0 - q) + maxKinEnergy*q);
288 
289  // delta-electron is produced
290  G4double totMomentum = totEnergy*sqrt(beta2);
291  G4double deltaMomentum =
292  sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
293  G4double cost = deltaKinEnergy * (totEnergy + electron_mass_c2) /
294  (deltaMomentum * totMomentum);
295  if(cost > 1.0) { cost = 1.0; }
296 
297  G4double sint = sqrt((1.0 - cost)*(1.0 + cost));
298 
299  G4double phi = twopi * G4UniformRand() ;
300 
301  G4ThreeVector deltaDirection(sint*cos(phi),sint*sin(phi), cost);
302  G4ThreeVector direction = dp->GetMomentumDirection();
303  deltaDirection.rotateUz(direction);
304 
305  // create G4DynamicParticle object for delta ray
306  G4DynamicParticle* delta =
307  new G4DynamicParticle(theElectron,deltaDirection,deltaKinEnergy);
308 
309  vdp->push_back(delta);
310 
311  // Change kinematics of primary particle
312  kineticEnergy -= deltaKinEnergy;
313  G4ThreeVector finalP = direction*totMomentum - deltaDirection*deltaMomentum;
314  finalP = finalP.unit();
315 
316  fParticleChange->SetProposedKineticEnergy(kineticEnergy);
317  fParticleChange->SetProposedMomentumDirection(finalP);
318 }
319 
320 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
321 
323  const G4MaterialCutsCouple* couple,
324  const G4DynamicParticle* dp,
325  G4double tmax,
326  G4double length,
327  G4double meanLoss)
328 {
329  G4double siga = Dispersion(couple->GetMaterial(),dp,tmax,length);
330  G4double loss = meanLoss;
331  siga = sqrt(siga);
332  G4double twomeanLoss = meanLoss + meanLoss;
333 
334  if(twomeanLoss < siga) {
335  G4double x;
336  do {
337  loss = twomeanLoss*G4UniformRand();
338  x = (loss - meanLoss)/siga;
339  } while (1.0 - 0.5*x*x < G4UniformRand());
340  } else {
341  do {
342  loss = G4RandGauss::shoot(meanLoss,siga);
343  } while (0.0 > loss || loss > twomeanLoss);
344  }
345  return loss;
346 }
347 
348 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
349 
350 G4double
352  const G4DynamicParticle* dp,
353  G4double tmax,
354  G4double length)
355 {
356  G4double siga = 0.0;
357  G4double tau = dp->GetKineticEnergy()/mass;
358  if(tau > 0.0) {
359  G4double electronDensity = material->GetElectronDensity();
360  G4double gam = tau + 1.0;
361  G4double invbeta2 = (gam*gam)/(tau * (tau+2.0));
362  siga = (invbeta2 - 0.5) * twopi_mc2_rcl2 * tmax * length
363  * electronDensity * chargeSquare;
364  }
365  return siga;
366 }
367 
368 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
369 
370 G4double
372  G4double kinEnergy)
373 {
374  G4double tau = kinEnergy/mass;
375  return 2.0*electron_mass_c2*tau*(tau + 2.);
376 }
377 
378 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
ThreeVector shoot(const G4int Ap, const G4int Af)
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:224
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:599
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:107
G4mplIonisationWithDeltaModel(G4double mCharge, const G4String &nam="mplIonisationWithDelta")
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
G4double GetKineticEnergy() const
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:592
const char * p
Definition: xmltok.h:285
virtual G4double Dispersion(const G4Material *, const G4DynamicParticle *, G4double tmax, G4double length)
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
virtual G4double SampleFluctuations(const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmax, G4double length, G4double meanLoss)
void SetParticle(const G4ParticleDefinition *p)
G4bool IsMaster() const
Definition: G4VEmModel.hh:676
G4ParticleDefinition * GetDefinition() const
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy)
tuple x
Definition: test.py:50
int G4int
Definition: G4Types.hh:78
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:683
string material
Definition: eplot.py:19
function g(Y1, Y2, PT2)
Definition: hijing1.383.f:5205
Char_t n[5]
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:426
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
Float_t Z
Definition: plot.C:39
G4double GetElectronDensity() const
Definition: G4Material.hh:215
const G4ThreeVector & GetMomentumDirection() const
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
void SetProposedKineticEnergy(G4double proposedKinEnergy)
float electron_mass_c2
Definition: hepunit.py:274
void SetProposedMomentumDirection(const G4ThreeVector &dir)
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
G4double G4Log(G4double x)
Definition: G4Log.hh:227
electron_Compton_length
Definition: hepunit.py:289
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy)
static G4ProductionCutsTable * GetProductionCutsTable()
G4double DensityCorrection(G4double x)
G4double GetPDGMass() const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
int G4lrint(double ad)
Definition: templates.hh:163
T max(const T t1, const T t2)
brief Return the largest of the two arguments
Hep3Vector unit() const
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4double GetMeanExcitationEnergy() const
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:690
const G4Material * GetMaterial() const