Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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$
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"
59#include "G4SystemOfUnits.hh"
60#include "G4LossTableManager.hh"
62#include "G4Electron.hh"
63#include "G4DynamicParticle.hh"
64
65//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
66
67using 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 = G4int(abs(magCharge) * 2 * fine_structure_const + 0.5);
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
114void
116 const G4DataVector&)
117{
118 if(!monopole) { SetParticle(p); }
119 if(!fParticleChange) { fParticleChange = GetParticleChangeForLoss(); }
120}
121
122//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
123
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 G4double tau = kineticEnergy / mass;
134 G4double gam = tau + 1.0;
135 G4double bg2 = tau * (tau + 2.0);
136 G4double beta2 = bg2 / (gam * gam);
137 G4double beta = sqrt(beta2);
138
139 // low-energy asymptotic formula
140 G4double dedx = dedxlim*beta*material->GetDensity();
141
142 // above asymptotic
143 if(beta > betalow) {
144
145 // high energy
146 if(beta >= betalim) {
147 dedx = ComputeDEDXAhlen(material, bg2, cutEnergy);
148
149 } else {
150
151 G4double dedx1 = dedxlim*betalow*material->GetDensity();
152 G4double dedx2 = ComputeDEDXAhlen(material, bg2lim, cutEnergy);
153
154 // extrapolation between two formula
155 G4double kapa2 = beta - betalow;
156 G4double kapa1 = betalim - beta;
157 dedx = (kapa1*dedx1 + kapa2*dedx2)/(kapa1 + kapa2);
158 }
159 }
160 return dedx;
161}
162
163//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
164
166G4mplIonisationWithDeltaModel::ComputeDEDXAhlen(const G4Material* material,
167 G4double bg2,
168 G4double cutEnergy)
169{
170 G4double eDensity = material->GetElectronDensity();
171 G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy();
172
173 // Ahlen's formula for nonconductors, [1]p157, f(5.7)
174 G4double dedx =
175 0.5*(log(2.0 * electron_mass_c2 * bg2*cutEnergy / (eexc*eexc)) - 1.0);
176
177 // Kazama et al. cross-section correction
178 G4double k = 0.406;
179 if(nmpl > 1) { k = 0.346; }
180
181 // Bloch correction
182 const G4double B[7] = { 0.0, 0.248, 0.672, 1.022, 1.243, 1.464, 1.685};
183
184 dedx += 0.5 * k - B[nmpl];
185
186 // density effect correction
187 G4double x = log(bg2)/twoln10;
188 dedx -= material->GetIonisation()->DensityCorrection(x);
189
190 // now compute the total ionization loss
191 dedx *= pi_hbarc2_over_mc2 * eDensity * nmpl * nmpl;
192
193 if (dedx < 0.0) { dedx = 0; }
194 return dedx;
195}
196
197//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
198
201 const G4ParticleDefinition* p,
202 G4double kineticEnergy,
203 G4double cutEnergy,
204 G4double maxKinEnergy)
205{
206 if(!monopole) { SetParticle(p); }
207 G4double cross = 0.0;
208 G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
209 G4double maxEnergy = min(tmax,maxKinEnergy);
210 if(cutEnergy < maxEnergy) {
211 cross = (1.0/cutEnergy - 1.0/maxEnergy)*twopi_mc2_rcl2*chargeSquare;
212 }
213 return cross;
214}
215
216//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
217
220 const G4ParticleDefinition* p,
221 G4double kineticEnergy,
223 G4double cutEnergy,
224 G4double maxEnergy)
225{
226 G4double cross =
227 Z*ComputeCrossSectionPerElectron(p,kineticEnergy,cutEnergy,maxEnergy);
228 return cross;
229}
230
231//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
232
233void
236 const G4DynamicParticle* dp,
237 G4double minKinEnergy,
238 G4double maxEnergy)
239{
240 G4double kineticEnergy = dp->GetKineticEnergy();
241 G4double tmax = MaxSecondaryEnergy(dp->GetDefinition(),kineticEnergy);
242
243 G4double maxKinEnergy = std::min(maxEnergy,tmax);
244 if(minKinEnergy >= maxKinEnergy) { return; }
245
246 //G4cout << "G4mplIonisationWithDeltaModel::SampleSecondaries: E(GeV)= "
247 // << kineticEnergy/GeV << " M(GeV)= " << mass/GeV
248 // << " tmin(MeV)= " << minKinEnergy/MeV << G4endl;
249
250 G4double totEnergy = kineticEnergy + mass;
251 G4double etot2 = totEnergy*totEnergy;
252 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/etot2;
253
254 // sampling without nuclear size effect
256 G4double deltaKinEnergy = minKinEnergy*maxKinEnergy
257 /(minKinEnergy*(1.0 - q) + maxKinEnergy*q);
258
259 // delta-electron is produced
260 G4double totMomentum = totEnergy*sqrt(beta2);
261 G4double deltaMomentum =
262 sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
263 G4double cost = deltaKinEnergy * (totEnergy + electron_mass_c2) /
264 (deltaMomentum * totMomentum);
265 if(cost > 1.0) { cost = 1.0; }
266
267 G4double sint = sqrt((1.0 - cost)*(1.0 + cost));
268
269 G4double phi = twopi * G4UniformRand() ;
270
271 G4ThreeVector deltaDirection(sint*cos(phi),sint*sin(phi), cost);
272 G4ThreeVector direction = dp->GetMomentumDirection();
273 deltaDirection.rotateUz(direction);
274
275 // create G4DynamicParticle object for delta ray
276 G4DynamicParticle* delta =
277 new G4DynamicParticle(theElectron,deltaDirection,deltaKinEnergy);
278
279 vdp->push_back(delta);
280
281 // Change kinematics of primary particle
282 kineticEnergy -= deltaKinEnergy;
283 G4ThreeVector finalP = direction*totMomentum - deltaDirection*deltaMomentum;
284 finalP = finalP.unit();
285
286 fParticleChange->SetProposedKineticEnergy(kineticEnergy);
287 fParticleChange->SetProposedMomentumDirection(finalP);
288}
289
290//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
291
293 const G4Material* material,
294 const G4DynamicParticle* dp,
295 G4double& tmax,
296 G4double& length,
297 G4double& meanLoss)
298{
299 G4double siga = Dispersion(material,dp,tmax,length);
300 G4double loss = meanLoss;
301 siga = sqrt(siga);
302 G4double twomeanLoss = meanLoss + meanLoss;
303
304 if(twomeanLoss < siga) {
305 G4double x;
306 do {
307 loss = twomeanLoss*G4UniformRand();
308 x = (loss - meanLoss)/siga;
309 } while (1.0 - 0.5*x*x < G4UniformRand());
310 } else {
311 do {
312 loss = G4RandGauss::shoot(meanLoss,siga);
313 } while (0.0 > loss || loss > twomeanLoss);
314 }
315 return loss;
316}
317
318//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
319
322 const G4DynamicParticle* dp,
323 G4double& tmax,
324 G4double& length)
325{
326 G4double siga = 0.0;
327 G4double tau = dp->GetKineticEnergy()/mass;
328 if(tau > 0.0) {
329 G4double electronDensity = material->GetElectronDensity();
330 G4double gam = tau + 1.0;
331 G4double invbeta2 = (gam*gam)/(tau * (tau+2.0));
332 siga = (invbeta2 - 0.5) * twopi_mc2_rcl2 * tmax * length
333 * electronDensity * chargeSquare;
334 }
335 return siga;
336}
337
338//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
339
342 G4double kinEnergy)
343{
344 G4double tau = kinEnergy/mass;
345 return 2.0*electron_mass_c2*tau*(tau + 2.);
346}
347
348//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
Hep3Vector unit() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4double GetMeanExcitationEnergy() const
G4double DensityCorrection(G4double x)
G4double GetDensity() const
Definition: G4Material.hh:179
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:225
G4double GetElectronDensity() const
Definition: G4Material.hh:216
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void SetProposedMomentumDirection(const G4ThreeVector &dir)
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:585
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:529
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:522
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:592
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:95
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy)
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
virtual G4double Dispersion(const G4Material *, const G4DynamicParticle *, G4double &tmax, G4double &length)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy)
G4mplIonisationWithDeltaModel(G4double mCharge, const G4String &nam="mplIonisationWithDelta")
void SetParticle(const G4ParticleDefinition *p)
virtual G4double SampleFluctuations(const G4Material *, const G4DynamicParticle *, G4double &tmax, G4double &length, G4double &meanLoss)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)