Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 76600 2013-11-13 08:30: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"
60 #include "G4LossTableManager.hh"
62 #include "G4ProductionCutsTable.hh"
63 #include "G4MaterialCutsCouple.hh"
64 #include "G4Log.hh"
65 
66 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
67 
68 using namespace std;
69 
70 std::vector<G4double>* G4mplIonisationModel::dedx0 = 0;
71 
74  magCharge(mCharge),
75  twoln10(log(100.0)),
76  betalow(0.01),
77  betalim(0.1),
78  beta2lim(betalim*betalim),
79  bg2lim(beta2lim*(1.0 + beta2lim))
80 {
81  nmpl = G4int(abs(magCharge) * 2 * fine_structure_const + 0.5);
82  if(nmpl > 6) { nmpl = 6; }
83  else if(nmpl < 1) { nmpl = 1; }
84  pi_hbarc2_over_mc2 = pi * hbarc * hbarc / electron_mass_c2;
85  chargeSquare = magCharge * magCharge;
86  dedxlim = 45.*nmpl*nmpl*GeV*cm2/g;
87  fParticleChange = 0;
88  monopole = 0;
89  mass = 0.0;
90 }
91 
92 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
93 
95 {
96  if(IsMaster()) { delete dedx0; }
97 }
98 
99 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
100 
102 {
103  monopole = p;
104  mass = monopole->GetPDGMass();
105  G4double emin =
106  std::min(LowEnergyLimit(),0.1*mass*(1/sqrt(1 - betalow*betalow) - 1));
107  G4double emax =
108  std::max(HighEnergyLimit(),10*mass*(1/sqrt(1 - beta2lim) - 1));
109  SetLowEnergyLimit(emin);
110  SetHighEnergyLimit(emax);
111 }
112 
113 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
114 
116  const G4DataVector&)
117 {
118  if(!monopole) { SetParticle(p); }
119  if(!fParticleChange) { fParticleChange = GetParticleChangeForLoss(); }
120  if(IsMaster()) {
121  if(!dedx0) { dedx0 = new std::vector<G4double>; }
122  G4ProductionCutsTable* theCoupleTable=
124  G4int numOfCouples = theCoupleTable->GetTableSize();
125  G4int n = dedx0->size();
126  if(n < numOfCouples) { dedx0->resize(numOfCouples); }
127 
128  // initialise vector
129  for(G4int i=0; i<numOfCouples; ++i) {
130 
131  const G4Material* material =
132  theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
133  G4double eDensity = material->GetElectronDensity();
134  G4double vF = electron_Compton_length*pow(3*pi*pi*eDensity,0.3333333333);
135  (*dedx0)[i] = pi_hbarc2_over_mc2*eDensity*nmpl*nmpl*
136  (G4Log(2*vF/fine_structure_const) - 0.5)/vF;
137  }
138  }
139 }
140 
141 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
142 
144  const G4ParticleDefinition* p,
145  G4double kineticEnergy,
146  G4double)
147 {
148  if(!monopole) { SetParticle(p); }
149  G4double tau = kineticEnergy / mass;
150  G4double gam = tau + 1.0;
151  G4double bg2 = tau * (tau + 2.0);
152  G4double beta2 = bg2 / (gam * gam);
153  G4double beta = sqrt(beta2);
154 
155  // low-energy asymptotic formula
156  //G4double dedx = dedxlim*beta*material->GetDensity();
157  G4double dedx = (*dedx0)[CurrentCouple()->GetIndex()]*beta;
158 
159  // above asymptotic
160  if(beta > betalow) {
161 
162  // high energy
163  if(beta >= betalim) {
164  dedx = ComputeDEDXAhlen(material, bg2);
165 
166  } else {
167 
168  //G4double dedx1 = dedxlim*betalow*material->GetDensity();
169  G4double dedx1 = (*dedx0)[CurrentCouple()->GetIndex()]*betalow;
170  G4double dedx2 = ComputeDEDXAhlen(material, bg2lim);
171 
172  // extrapolation between two formula
173  G4double kapa2 = beta - betalow;
174  G4double kapa1 = betalim - beta;
175  dedx = (kapa1*dedx1 + kapa2*dedx2)/(kapa1 + kapa2);
176  }
177  }
178  return dedx;
179 }
180 
181 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
182 
183 G4double G4mplIonisationModel::ComputeDEDXAhlen(const G4Material* material,
184  G4double bg2)
185 {
186  G4double eDensity = material->GetElectronDensity();
187  G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy();
188  G4double cden = material->GetIonisation()->GetCdensity();
189  G4double mden = material->GetIonisation()->GetMdensity();
190  G4double aden = material->GetIonisation()->GetAdensity();
191  G4double x0den = material->GetIonisation()->GetX0density();
192  G4double x1den = material->GetIonisation()->GetX1density();
193 
194  // Ahlen's formula for nonconductors, [1]p157, f(5.7)
195  G4double dedx = log(2.0 * electron_mass_c2 * bg2 / eexc) - 0.5;
196 
197  // Kazama et al. cross-section correction
198  G4double k = 0.406;
199  if(nmpl > 1) k = 0.346;
200 
201  // Bloch correction
202  const G4double B[7] = { 0.0, 0.248, 0.672, 1.022, 1.243, 1.464, 1.685};
203 
204  dedx += 0.5 * k - B[nmpl];
205 
206  // density effect correction
207  G4double deltam;
208  G4double x = log(bg2) / twoln10;
209  if ( x >= x0den ) {
210  deltam = twoln10 * x - cden;
211  if ( x < x1den ) deltam += aden * pow((x1den-x), mden);
212  dedx -= 0.5 * deltam;
213  }
214 
215  // now compute the total ionization loss
216  dedx *= pi_hbarc2_over_mc2 * eDensity * nmpl * nmpl;
217 
218  if (dedx < 0.0) dedx = 0;
219  return dedx;
220 }
221 
222 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
223 
224 void G4mplIonisationModel::SampleSecondaries(std::vector<G4DynamicParticle*>*,
225  const G4MaterialCutsCouple*,
226  const G4DynamicParticle*,
227  G4double,
228  G4double)
229 {}
230 
231 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
232 
234  const G4MaterialCutsCouple* couple,
235  const G4DynamicParticle* dp,
236  G4double tmax,
237  G4double length,
238  G4double meanLoss)
239 {
240  G4double siga = Dispersion(couple->GetMaterial(),dp,tmax,length);
241  G4double loss = meanLoss;
242  siga = sqrt(siga);
243  G4double twomeanLoss = meanLoss + meanLoss;
244 
245  if(twomeanLoss < siga) {
246  G4double x;
247  do {
248  loss = twomeanLoss*G4UniformRand();
249  x = (loss - meanLoss)/siga;
250  } while (1.0 - 0.5*x*x < G4UniformRand());
251  } else {
252  do {
253  loss = G4RandGauss::shoot(meanLoss,siga);
254  } while (0.0 > loss || loss > twomeanLoss);
255  }
256  return loss;
257 }
258 
259 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
260 
262  const G4DynamicParticle* dp,
263  G4double tmax,
264  G4double length)
265 {
266  G4double siga = 0.0;
267  G4double tau = dp->GetKineticEnergy()/mass;
268  if(tau > 0.0) {
269  G4double electronDensity = material->GetElectronDensity();
270  G4double gam = tau + 1.0;
271  G4double invbeta2 = (gam*gam)/(tau * (tau+2.0));
272  siga = (invbeta2 - 0.5) * twopi_mc2_rcl2 * tmax * length
273  * electronDensity * chargeSquare;
274  }
275  return siga;
276 }
277 
278 //....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
virtual G4double Dispersion(const G4Material *, const G4DynamicParticle *, G4double tmax, G4double length)
G4double GetAdensity() const
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:107
G4double GetKineticEnergy() const
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:592
const char * p
Definition: xmltok.h:285
G4double GetX1density() const
G4bool IsMaster() const
Definition: G4VEmModel.hh:676
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
int G4int
Definition: G4Types.hh:78
G4mplIonisationModel(G4double mCharge, const G4String &nam="mplIonisation")
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:683
string material
Definition: eplot.py:19
function g(Y1, Y2, PT2)
Definition: hijing1.383.f:5205
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:426
#define G4UniformRand()
Definition: Randomize.hh:87
G4double GetElectronDensity() const
Definition: G4Material.hh:215
float electron_mass_c2
Definition: hepunit.py:274
const G4int n
void SetParticle(const G4ParticleDefinition *p)
G4double G4Log(G4double x)
Definition: G4Log.hh:227
electron_Compton_length
Definition: hepunit.py:289
static G4ProductionCutsTable * GetProductionCutsTable()
G4double GetX0density() const
G4double GetPDGMass() const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
G4double GetCdensity() const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4double GetMeanExcitationEnergy() const
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:690
G4double GetMdensity() const
virtual G4double SampleFluctuations(const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmax, G4double length, G4double meanLoss)
const G4Material * GetMaterial() const