Geant4-11
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//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class header file
30//
31//
32// File name: G4mplIonisationWithDeltaModel
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"
60#include "G4Electron.hh"
61#include "G4DynamicParticle.hh"
64#include "G4Log.hh"
65#include "G4Pow.hh"
66
67//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
68
69using namespace std;
70
71std::vector<G4double>* G4mplIonisationWithDeltaModel::dedx0 = nullptr;
72
74 const G4String& nam)
76 magCharge(mCharge),
77 twoln10(std::log(100.0)),
78 betalow(0.01),
79 betalim(0.1),
80 beta2lim(betalim*betalim),
81 bg2lim(beta2lim*(1.0 + beta2lim))
82{
83 nmpl = G4lrint(std::abs(magCharge) * 2 * fine_structure_const);
84 if(nmpl > 6) { nmpl = 6; }
85 else if(nmpl < 1) { nmpl = 1; }
88 dedxlim = 45.*nmpl*nmpl*GeV*cm2/g;
89 fParticleChange = nullptr;
91 G4cout << "### Monopole ionisation model with d-electron production, Gmag= "
93 monopole = nullptr;
94 mass = 0.0;
95}
96
97//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
98
100{
101 if(IsMaster()) { delete dedx0; }
102}
103
104//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
105
107{
108 monopole = p;
110 G4double emin =
111 std::min(LowEnergyLimit(),0.1*mass*(1./sqrt(1. - betalow*betalow) - 1.));
112 G4double emax =
113 std::max(HighEnergyLimit(),10*mass*(1./sqrt(1. - beta2lim) - 1.));
114 SetLowEnergyLimit(emin);
116}
117
118//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
119
120void
122 const G4DataVector&)
123{
124 if(!monopole) { SetParticle(p); }
126 if(IsMaster()) {
127 if(!dedx0) { dedx0 = new std::vector<G4double>; }
128 G4ProductionCutsTable* theCoupleTable=
130 G4int numOfCouples = theCoupleTable->GetTableSize();
131 G4int n = dedx0->size();
132 if(n < numOfCouples) { dedx0->resize(numOfCouples); }
133 G4Pow* g4calc = G4Pow::GetInstance();
134
135 // initialise vector assuming low conductivity
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 vF2 = 2*electron_Compton_length*g4calc->A13(3.*pi*pi*eDensity);
142 (*dedx0)[i] = pi_hbarc2_over_mc2*eDensity*nmpl*nmpl*
143 (G4Log(vF2/fine_structure_const) - 0.5)/vF2;
144 }
145 }
146}
147
148//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
149
152 const G4MaterialCutsCouple* couple)
153{
155}
156
157//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
158
161 const G4ParticleDefinition* p,
162 G4double kineticEnergy,
163 G4double maxEnergy)
164{
165 if(!monopole) { SetParticle(p); }
166 G4double tmax = MaxSecondaryEnergy(p,kineticEnergy);
167 G4double cutEnergy = std::min(tmax, maxEnergy);
168 cutEnergy = std::max(LowEnergyLimit(), cutEnergy);
169 G4double tau = kineticEnergy / mass;
170 G4double gam = tau + 1.0;
171 G4double bg2 = tau * (tau + 2.0);
172 G4double beta2 = bg2 / (gam * gam);
173 G4double beta = sqrt(beta2);
174
175 // low-energy asymptotic formula
176 G4double dedx = (*dedx0)[CurrentCouple()->GetIndex()]*beta;
177
178 // above asymptotic
179 if(beta > betalow) {
180
181 // high energy
182 if(beta >= betalim) {
183 dedx = ComputeDEDXAhlen(material, bg2, cutEnergy);
184
185 } else {
186 G4double dedx1 = (*dedx0)[CurrentCouple()->GetIndex()]*betalow;
187 G4double dedx2 = ComputeDEDXAhlen(material, bg2lim, cutEnergy);
188
189 // extrapolation between two formula
190 G4double kapa2 = beta - betalow;
191 G4double kapa1 = betalim - beta;
192 dedx = (kapa1*dedx1 + kapa2*dedx2)/(kapa1 + kapa2);
193 }
194 }
195 return dedx;
196}
197
198//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
199
202 G4double bg2,
203 G4double cutEnergy)
204{
205 G4double eDensity = material->GetElectronDensity();
206 G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy();
207
208 // Ahlen's formula for nonconductors, [1]p157, f(5.7)
209 G4double dedx =
210 0.5*(G4Log(2.0*electron_mass_c2*bg2*cutEnergy/(eexc*eexc)) -1.0);
211
212 // Kazama et al. cross-section correction
213 G4double k = 0.406;
214 if(nmpl > 1) { k = 0.346; }
215
216 // Bloch correction
217 const G4double B[7] = { 0.0, 0.248, 0.672, 1.022, 1.243, 1.464, 1.685};
218
219 dedx += 0.5 * k - B[nmpl];
220
221 // density effect correction
222 G4double x = G4Log(bg2)/twoln10;
223 dedx -= material->GetIonisation()->DensityCorrection(x);
224
225 // now compute the total ionization loss
226 dedx *= pi_hbarc2_over_mc2 * eDensity * nmpl * nmpl;
227
228 dedx = std::max(dedx, 0.0);
229 return dedx;
230}
231
232//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
233
236 const G4ParticleDefinition* p,
237 G4double kineticEnergy,
238 G4double cut,
239 G4double maxKinEnergy)
240{
241 if(!monopole) { SetParticle(p); }
242 G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
243 G4double maxEnergy = std::min(tmax, maxKinEnergy);
244 G4double cutEnergy = std::max(LowEnergyLimit(), cut);
245 G4double cross = (cutEnergy < maxEnergy)
246 ? (0.5/cutEnergy - 0.5/maxEnergy)*pi_hbarc2_over_mc2 * nmpl * nmpl : 0.0;
247 return cross;
248}
249
250//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
251
254 const G4ParticleDefinition* p,
255 G4double kineticEnergy,
257 G4double cutEnergy,
258 G4double maxEnergy)
259{
260 G4double cross =
261 Z*ComputeCrossSectionPerElectron(p,kineticEnergy,cutEnergy,maxEnergy);
262 return cross;
263}
264
265//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
266
267void
270 const G4DynamicParticle* dp,
271 G4double minKinEnergy,
272 G4double maxEnergy)
273{
274 G4double kineticEnergy = dp->GetKineticEnergy();
275 G4double tmax = MaxSecondaryEnergy(dp->GetDefinition(),kineticEnergy);
276
277 G4double maxKinEnergy = std::min(maxEnergy,tmax);
278 if(minKinEnergy >= maxKinEnergy) { return; }
279
280 //G4cout << "G4mplIonisationWithDeltaModel::SampleSecondaries: E(GeV)= "
281 // << kineticEnergy/GeV << " M(GeV)= " << mass/GeV
282 // << " tmin(MeV)= " << minKinEnergy/MeV << G4endl;
283
284 G4double totEnergy = kineticEnergy + mass;
285 G4double etot2 = totEnergy*totEnergy;
286 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/etot2;
287
288 // sampling without nuclear size effect
290 G4double deltaKinEnergy = minKinEnergy*maxKinEnergy
291 /(minKinEnergy*(1.0 - q) + maxKinEnergy*q);
292
293 // delta-electron is produced
294 G4double totMomentum = totEnergy*sqrt(beta2);
295 G4double deltaMomentum =
296 sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
297 G4double cost = deltaKinEnergy * (totEnergy + electron_mass_c2) /
298 (deltaMomentum * totMomentum);
299 cost = std::min(cost, 1.0);
300
301 G4double sint = sqrt((1.0 - cost)*(1.0 + cost));
302
303 G4double phi = twopi * G4UniformRand() ;
304
305 G4ThreeVector deltaDirection(sint*cos(phi),sint*sin(phi), cost);
306 G4ThreeVector direction = dp->GetMomentumDirection();
307 deltaDirection.rotateUz(direction);
308
309 // create G4DynamicParticle object for delta ray
310 G4DynamicParticle* delta =
311 new G4DynamicParticle(theElectron,deltaDirection,deltaKinEnergy);
312
313 vdp->push_back(delta);
314
315 // Change kinematics of primary particle
316 kineticEnergy -= deltaKinEnergy;
317 G4ThreeVector finalP = direction*totMomentum - deltaDirection*deltaMomentum;
318 finalP = finalP.unit();
319
322}
323
324//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
325
327 const G4MaterialCutsCouple* couple,
328 const G4DynamicParticle* dp,
329 const G4double tcut,
330 const G4double tmax,
331 const G4double length,
332 const G4double meanLoss)
333{
334 G4double siga = Dispersion(couple->GetMaterial(),dp,tcut,tmax,length);
335 G4double loss = meanLoss;
336 siga = std::sqrt(siga);
337 G4double twomeanLoss = meanLoss + meanLoss;
338
339 if(twomeanLoss < siga) {
340 G4double x;
341 do {
342 loss = twomeanLoss*G4UniformRand();
343 x = (loss - meanLoss)/siga;
344 // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
345 } while (1.0 - 0.5*x*x < G4UniformRand());
346 } else {
347 do {
348 loss = G4RandGauss::shoot(meanLoss,siga);
349 // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
350 } while (0.0 > loss || loss > twomeanLoss);
351 }
352 return loss;
353}
354
355//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
356
359 const G4DynamicParticle* dp,
360 const G4double tcut,
361 const G4double tmax,
362 const G4double length)
363{
364 G4double siga = 0.0;
365 G4double tau = dp->GetKineticEnergy()/mass;
366 if(tau > 0.0) {
367 const G4double beta = dp->GetBeta();
368 siga = (tmax/(beta*beta) - 0.5*tcut) * twopi_mc2_rcl2 * length
369 * material->GetElectronDensity() * chargeSquare;
370 }
371 return siga;
372}
373
374//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
375
378 G4double kinEnergy)
379{
380 G4double tau = kinEnergy/mass;
381 return 2.0*electron_mass_c2*tau*(tau + 2.);
382}
383
384//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
static const G4double emax
G4double B(G4double temperature)
G4double G4Log(G4double x)
Definition: G4Log.hh:226
static constexpr double twopi
Definition: G4SIunits.hh:56
static constexpr double cm2
Definition: G4SIunits.hh:100
static constexpr double eplus
Definition: G4SIunits.hh:184
static constexpr double g
Definition: G4SIunits.hh:168
static constexpr double GeV
Definition: G4SIunits.hh:203
static constexpr double pi
Definition: G4SIunits.hh:55
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector unit() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double GetBeta() const
static G4Electron * Electron()
Definition: G4Electron.cc:93
G4double GetMeanExcitationEnergy() const
const G4Material * GetMaterial() const
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:222
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void SetProposedMomentumDirection(const G4ThreeVector &dir)
Definition: G4Pow.hh:49
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4double A13(G4double A) const
Definition: G4Pow.cc:120
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:767
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:662
G4bool IsMaster() const
Definition: G4VEmModel.hh:746
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:655
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:774
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:490
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:108
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *couple) override
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4mplIonisationWithDeltaModel(G4double mCharge, const G4String &nam="mplIonisationWithDelta")
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
static std::vector< G4double > * dedx0
virtual G4double SampleFluctuations(const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double tcut, const G4double tmax, const G4double length, const G4double meanLoss) override
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) override
void SetParticle(const G4ParticleDefinition *p)
virtual G4double Dispersion(const G4Material *, const G4DynamicParticle *, const G4double tcut, const G4double tmax, const G4double length) override
G4double ComputeDEDXAhlen(const G4Material *material, G4double bg2, G4double cut)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
ThreeVector shoot(const G4int Ap, const G4int Af)
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
string material
Definition: eplot.py:19
float electron_mass_c2
Definition: hepunit.py:273
int twopi_mc2_rcl2
Definition: hepunit.py:293
float hbarc
Definition: hepunit.py:264
float electron_Compton_length
Definition: hepunit.py:288
int fine_structure_const
Definition: hepunit.py:286
int G4lrint(double ad)
Definition: templates.hh:134