Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros 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//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class header file
30//
31//
32// File name: G4mplIonisationModel
33//
34// Author: Vladimir Ivanchenko
35//
36// Creation date: 06.09.2005
37//
38// Modifications:
39// 12.08.2007 Changing low energy approximation and extrapolation.
40// Small bug fixing and refactoring (M. Vladymyrov)
41// 13.11.2007 Use low-energy asymptotic from [3] (V.Ivanchenko)
42//
43//
44// -------------------------------------------------------------------
45// References
46// [1] Steven P. Ahlen: Energy loss of relativistic heavy ionizing particles,
47// S.P. Ahlen, Rev. Mod. Phys 52(1980), p121
48// [2] K.A. Milton arXiv:hep-ex/0602040
49// [3] S.P. Ahlen and K. Kinoshita, Phys. Rev. D26 (1982) 2347
50
51
52//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
53//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
54
56#include "Randomize.hh"
58#include "G4SystemOfUnits.hh"
62#include "G4Log.hh"
63#include "G4Pow.hh"
64
65//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
66
67std::vector<G4double>* G4mplIonisationModel::dedx0 = nullptr;
68
71 magCharge(mCharge),
72 twoln10(G4Log(100.0)),
73 betalow(0.01),
74 betalim(0.1),
75 beta2lim(betalim*betalim),
76 bg2lim(beta2lim*(1.0 + beta2lim))
77{
78 nmpl = G4int(std::abs(magCharge) * 2 * CLHEP::fine_structure_const + 0.5);
79 if(nmpl > 6) { nmpl = 6; }
80 else if(nmpl < 1) { nmpl = 1; }
81 pi_hbarc2_over_mc2 = CLHEP::pi*CLHEP::hbarc*CLHEP::hbarc/CLHEP::electron_mass_c2;
82 chargeSquare = magCharge * magCharge;
83 dedxlim = 45.*nmpl*nmpl*CLHEP::GeV*CLHEP::cm2/CLHEP::g;
84}
85
86//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
87
89{
90 if(IsMaster()) { delete dedx0; }
91}
92
93//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
94
96{
97 monopole = p;
98 mass = monopole->GetPDGMass();
99 G4double emin =
100 std::min(LowEnergyLimit(),0.1*mass*(1./std::sqrt(1. - betalow*betalow) - 1.));
101 G4double emax =
102 std::max(HighEnergyLimit(),10.*mass*(1./std::sqrt(1. - beta2lim) - 1.));
103 SetLowEnergyLimit(emin);
104 SetHighEnergyLimit(emax);
105}
106
107//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
108
110 const G4DataVector&)
111{
112 if(nullptr == monopole) { SetParticle(p); }
113 if(nullptr == fParticleChange) { fParticleChange = GetParticleChangeForLoss(); }
114 if(IsMaster()) {
115 if(nullptr == dedx0) { dedx0 = new std::vector<G4double>; }
116 G4ProductionCutsTable* theCoupleTable=
118 G4int numOfCouples = (G4int)theCoupleTable->GetTableSize();
119 G4int n = (G4int)dedx0->size();
120 if(n < numOfCouples) { dedx0->resize(numOfCouples); }
121
122 G4Pow* g4calc = G4Pow::GetInstance();
123
124 // initialise vector assuming low conductivity
125 for(G4int i=0; i<numOfCouples; ++i) {
126
127 const G4Material* material =
128 theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
129 G4double eDensity = material->GetElectronDensity();
130 G4double vF2 = 2*electron_Compton_length*g4calc->A13(3.*pi*pi*eDensity);
131 (*dedx0)[i] = pi_hbarc2_over_mc2*eDensity*nmpl*nmpl*
132 (G4Log(vF2/fine_structure_const) - 0.5)/vF2;
133 }
134 }
135}
136
137//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
138
140 const G4ParticleDefinition* p,
141 G4double kineticEnergy,
142 G4double)
143{
144 if(nullptr == monopole) { SetParticle(p); }
145 G4double tau = kineticEnergy / mass;
146 G4double gam = tau + 1.0;
147 G4double bg2 = tau * (tau + 2.0);
148 G4double beta2 = bg2 / (gam * gam);
149 G4double beta = std::sqrt(beta2);
150
151 // low-energy asymptotic formula
152 //G4double dedx = dedxlim*beta*material->GetDensity();
153 G4double dedx = (*dedx0)[CurrentCouple()->GetIndex()]*beta;
154
155 // above asymptotic
156 if(beta > betalow) {
157
158 // high energy
159 if(beta >= betalim) {
160 dedx = ComputeDEDXAhlen(material, bg2);
161
162 } else {
163
164 //G4double dedx1 = dedxlim*betalow*material->GetDensity();
165 G4double dedx1 = (*dedx0)[CurrentCouple()->GetIndex()]*betalow;
166 G4double dedx2 = ComputeDEDXAhlen(material, bg2lim);
167
168 // extrapolation between two formula
169 G4double kapa2 = beta - betalow;
170 G4double kapa1 = betalim - beta;
171 dedx = (kapa1*dedx1 + kapa2*dedx2)/(kapa1 + kapa2);
172 }
173 }
174 return dedx;
175}
176
177//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
178
179G4double G4mplIonisationModel::ComputeDEDXAhlen(const G4Material* material,
180 G4double bg2)
181{
182 G4double eDensity = material->GetElectronDensity();
183 G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy();
184 G4double cden = material->GetIonisation()->GetCdensity();
185 G4double mden = material->GetIonisation()->GetMdensity();
186 G4double aden = material->GetIonisation()->GetAdensity();
187 G4double x0den = material->GetIonisation()->GetX0density();
188 G4double x1den = material->GetIonisation()->GetX1density();
189
190 // Ahlen's formula for nonconductors, [1]p157, f(5.7)
191 G4double dedx = std::log(2.0 * electron_mass_c2 * bg2 / eexc) - 0.5;
192
193 // Kazama et al. cross-section correction
194 G4double k = 0.406;
195 if(nmpl > 1) k = 0.346;
196
197 // Bloch correction
198 const G4double B[7] = { 0.0, 0.248, 0.672, 1.022, 1.243, 1.464, 1.685};
199
200 dedx += 0.5 * k - B[nmpl];
201
202 // density effect correction
203 G4double deltam;
204 G4double x = std::log(bg2) / twoln10;
205 if ( x >= x0den ) {
206 deltam = twoln10 * x - cden;
207 if ( x < x1den ) deltam += aden * std::pow((x1den-x), mden);
208 dedx -= 0.5 * deltam;
209 }
210
211 // now compute the total ionization loss
212 dedx *= pi_hbarc2_over_mc2 * eDensity * nmpl * nmpl;
213
214 if (dedx < 0.0) dedx = 0.;
215 return dedx;
216}
217
218//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
219
220void G4mplIonisationModel::SampleSecondaries(std::vector<G4DynamicParticle*>*,
222 const G4DynamicParticle*,
223 G4double,
224 G4double)
225{}
226
227//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
228
230 const G4MaterialCutsCouple* couple,
231 const G4DynamicParticle* dp,
232 const G4double tcut,
233 const G4double tmax,
234 const G4double length,
235 const G4double meanLoss)
236{
237 G4double siga = Dispersion(couple->GetMaterial(),dp,tcut,tmax,length);
238 G4double loss = meanLoss;
239 siga = std::sqrt(siga);
240 G4double twomeanLoss = meanLoss + meanLoss;
241
242 if(twomeanLoss < siga) {
243 G4double x;
244 do {
245 loss = twomeanLoss*G4UniformRand();
246 x = (loss - meanLoss)/siga;
247 // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
248 } while (1.0 - 0.5*x*x < G4UniformRand());
249 } else {
250 do {
251 loss = G4RandGauss::shoot(meanLoss,siga);
252 // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
253 } while (0.0 > loss || loss > twomeanLoss);
254 }
255 return loss;
256}
257
258//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
259
261 const G4DynamicParticle* dp,
262 const G4double tcut,
263 const G4double tmax,
264 const G4double length)
265{
266 G4double siga = 0.0;
267 G4double tau = dp->GetKineticEnergy()/mass;
268 if(tau > 0.0) {
269 const G4double beta = dp->GetBeta();
270 siga = (tmax/(beta*beta) - 0.5*tcut) * twopi_mc2_rcl2 * length
271 * material->GetElectronDensity() * chargeSquare;
272 }
273 return siga;
274}
275
276//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
G4double B(G4double temperature)
G4double G4Log(G4double x)
Definition: G4Log.hh:227
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4UniformRand()
Definition: Randomize.hh:52
G4double GetKineticEnergy() const
G4double GetBeta() const
G4double GetMdensity() const
G4double GetX1density() const
G4double GetX0density() const
G4double GetCdensity() const
G4double GetMeanExcitationEnergy() const
G4double GetAdensity() const
const G4Material * GetMaterial() const
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:221
G4double GetElectronDensity() const
Definition: G4Material.hh:212
Definition: G4Pow.hh:49
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4double A13(G4double A) const
Definition: G4Pow.cc:116
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:746
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
G4bool IsMaster() const
Definition: G4VEmModel.hh:725
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:753
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:486
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:109
G4mplIonisationModel(G4double mCharge, const G4String &nam="mplIonisation")
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
void SetParticle(const G4ParticleDefinition *p)
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4double SampleFluctuations(const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double tcut, const G4double tmax, const G4double length, const G4double meanLoss) override
G4double Dispersion(const G4Material *, const G4DynamicParticle *, const G4double tcut, const G4double tmax, const G4double length) override