Geant4  10.02.p01
G4mplIonisationModel.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: G4mplIonisationModel.cc 91869 2015-08-07 15:21:02Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class header file
31 //
32 //
33 // File name: G4mplIonisationModel
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 
56 #include "G4mplIonisationModel.hh"
57 #include "Randomize.hh"
58 #include "G4PhysicalConstants.hh"
59 #include "G4SystemOfUnits.hh"
61 #include "G4ProductionCutsTable.hh"
62 #include "G4MaterialCutsCouple.hh"
63 #include "G4Log.hh"
64 
65 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
66 
67 using namespace std;
68 
69 std::vector<G4double>* G4mplIonisationModel::dedx0 = 0;
70 
73  magCharge(mCharge),
74  twoln10(log(100.0)),
75  betalow(0.01),
76  betalim(0.1),
77  beta2lim(betalim*betalim),
78  bg2lim(beta2lim*(1.0 + beta2lim))
79 {
80  nmpl = G4int(abs(magCharge) * 2 * fine_structure_const + 0.5);
81  if(nmpl > 6) { nmpl = 6; }
82  else if(nmpl < 1) { nmpl = 1; }
83  pi_hbarc2_over_mc2 = pi * hbarc * hbarc / electron_mass_c2;
85  dedxlim = 45.*nmpl*nmpl*GeV*cm2/g;
86  fParticleChange = 0;
87  monopole = 0;
88  mass = 0.0;
89 }
90 
91 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
92 
94 {
95  if(IsMaster()) { delete dedx0; }
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 
115  const G4DataVector&)
116 {
117  if(!monopole) { SetParticle(p); }
119  if(IsMaster()) {
120  if(!dedx0) { dedx0 = new std::vector<G4double>; }
121  G4ProductionCutsTable* theCoupleTable=
123  G4int numOfCouples = theCoupleTable->GetTableSize();
124  G4int n = dedx0->size();
125  if(n < numOfCouples) { dedx0->resize(numOfCouples); }
126 
127  // initialise vector
128  for(G4int i=0; i<numOfCouples; ++i) {
129 
130  const G4Material* material =
131  theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
132  G4double eDensity = material->GetElectronDensity();
133  G4double vF = electron_Compton_length*pow(3*pi*pi*eDensity,0.3333333333);
134  (*dedx0)[i] = pi_hbarc2_over_mc2*eDensity*nmpl*nmpl*
135  (G4Log(2*vF/fine_structure_const) - 0.5)/vF;
136  }
137  }
138 }
139 
140 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
141 
143  const G4ParticleDefinition* p,
144  G4double kineticEnergy,
145  G4double)
146 {
147  if(!monopole) { SetParticle(p); }
148  G4double tau = kineticEnergy / mass;
149  G4double gam = tau + 1.0;
150  G4double bg2 = tau * (tau + 2.0);
151  G4double beta2 = bg2 / (gam * gam);
152  G4double beta = sqrt(beta2);
153 
154  // low-energy asymptotic formula
155  //G4double dedx = dedxlim*beta*material->GetDensity();
156  G4double dedx = (*dedx0)[CurrentCouple()->GetIndex()]*beta;
157 
158  // above asymptotic
159  if(beta > betalow) {
160 
161  // high energy
162  if(beta >= betalim) {
163  dedx = ComputeDEDXAhlen(material, bg2);
164 
165  } else {
166 
167  //G4double dedx1 = dedxlim*betalow*material->GetDensity();
168  G4double dedx1 = (*dedx0)[CurrentCouple()->GetIndex()]*betalow;
169  G4double dedx2 = ComputeDEDXAhlen(material, bg2lim);
170 
171  // extrapolation between two formula
172  G4double kapa2 = beta - betalow;
173  G4double kapa1 = betalim - beta;
174  dedx = (kapa1*dedx1 + kapa2*dedx2)/(kapa1 + kapa2);
175  }
176  }
177  return dedx;
178 }
179 
180 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
181 
183  G4double bg2)
184 {
185  G4double eDensity = material->GetElectronDensity();
186  G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy();
187  G4double cden = material->GetIonisation()->GetCdensity();
188  G4double mden = material->GetIonisation()->GetMdensity();
189  G4double aden = material->GetIonisation()->GetAdensity();
190  G4double x0den = material->GetIonisation()->GetX0density();
191  G4double x1den = material->GetIonisation()->GetX1density();
192 
193  // Ahlen's formula for nonconductors, [1]p157, f(5.7)
194  G4double dedx = log(2.0 * electron_mass_c2 * bg2 / eexc) - 0.5;
195 
196  // Kazama et al. cross-section correction
197  G4double k = 0.406;
198  if(nmpl > 1) k = 0.346;
199 
200  // Bloch correction
201  const G4double B[7] = { 0.0, 0.248, 0.672, 1.022, 1.243, 1.464, 1.685};
202 
203  dedx += 0.5 * k - B[nmpl];
204 
205  // density effect correction
206  G4double deltam;
207  G4double x = log(bg2) / twoln10;
208  if ( x >= x0den ) {
209  deltam = twoln10 * x - cden;
210  if ( x < x1den ) deltam += aden * pow((x1den-x), mden);
211  dedx -= 0.5 * deltam;
212  }
213 
214  // now compute the total ionization loss
215  dedx *= pi_hbarc2_over_mc2 * eDensity * nmpl * nmpl;
216 
217  if (dedx < 0.0) dedx = 0;
218  return dedx;
219 }
220 
221 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
222 
223 void G4mplIonisationModel::SampleSecondaries(std::vector<G4DynamicParticle*>*,
224  const G4MaterialCutsCouple*,
225  const G4DynamicParticle*,
226  G4double,
227  G4double)
228 {}
229 
230 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
231 
233  const G4MaterialCutsCouple* couple,
234  const G4DynamicParticle* dp,
235  G4double tmax,
236  G4double length,
237  G4double meanLoss)
238 {
239  G4double siga = Dispersion(couple->GetMaterial(),dp,tmax,length);
240  G4double loss = meanLoss;
241  siga = sqrt(siga);
242  G4double twomeanLoss = meanLoss + meanLoss;
243 
244  if(twomeanLoss < siga) {
245  G4double x;
246  do {
247  loss = twomeanLoss*G4UniformRand();
248  x = (loss - meanLoss)/siga;
249  // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
250  } while (1.0 - 0.5*x*x < G4UniformRand());
251  } else {
252  do {
253  loss = G4RandGauss::shoot(meanLoss,siga);
254  // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
255  } while (0.0 > loss || loss > twomeanLoss);
256  }
257  return loss;
258 }
259 
260 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
261 
263  const G4DynamicParticle* dp,
264  G4double tmax,
265  G4double length)
266 {
267  G4double siga = 0.0;
268  G4double tau = dp->GetKineticEnergy()/mass;
269  if(tau > 0.0) {
270  G4double electronDensity = material->GetElectronDensity();
271  G4double gam = tau + 1.0;
272  G4double invbeta2 = (gam*gam)/(tau * (tau+2.0));
273  siga = (invbeta2 - 0.5) * twopi_mc2_rcl2 * tmax * length
274  * electronDensity * chargeSquare;
275  }
276  return siga;
277 }
278 
279 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
ThreeVector shoot(const G4int Ap, const G4int Af)
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:226
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
virtual G4double Dispersion(const G4Material *, const G4DynamicParticle *, G4double tmax, G4double length)
G4double GetAdensity() const
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:120
static const double cm2
Definition: G4SIunits.hh:119
G4double GetKineticEnergy() const
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
G4double GetX1density() const
G4bool IsMaster() const
Definition: G4VEmModel.hh:718
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
double B(double temperature)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
int G4int
Definition: G4Types.hh:78
static const G4double bg2lim
G4mplIonisationModel(G4double mCharge, const G4String &nam="mplIonisation")
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:725
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:452
#define G4UniformRand()
Definition: Randomize.hh:97
G4double GetElectronDensity() const
Definition: G4Material.hh:217
G4double ComputeDEDXAhlen(const G4Material *material, G4double bg2)
static const G4double beta2lim
static const double GeV
Definition: G4SIunits.hh:214
const G4int n
void SetParticle(const G4ParticleDefinition *p)
static const G4double emax
G4double G4Log(G4double x)
Definition: G4Log.hh:230
static G4ProductionCutsTable * GetProductionCutsTable()
G4double GetX0density() const
G4double GetPDGMass() const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
static const double pi
Definition: G4SIunits.hh:74
G4double GetCdensity() const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
static const double g
Definition: G4SIunits.hh:180
G4ParticleChangeForLoss * fParticleChange
const G4double x[NPOINTSGL]
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4double GetMeanExcitationEnergy() const
const G4ParticleDefinition * monopole
static const G4double twoln10
static std::vector< G4double > * dedx0
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
G4double GetMdensity() const
virtual G4double SampleFluctuations(const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmax, G4double length, G4double meanLoss)
const G4Material * GetMaterial() const