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
G4IonFluctuations.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 file
30//
31//
32// File name: G4IonFluctuation
33//
34// Author: Vladimir Ivanchenko
35//
36// Creation date: 03.01.2002
37//
38// Modifications:
39//
40// 28-12-02 add method Dispersion (V.Ivanchenko)
41// 07-02-03 change signature (V.Ivanchenko)
42// 13-02-03 Add name (V.Ivanchenko)
43// 23-05-03 Add control on parthalogical cases (V.Ivanchenko)
44// 16-10-03 Changed interface to Initialisation (V.Ivanchenko)
45// 27-09-07 Use FermiEnergy from material, add cut dependence (V.Ivanchenko)
46// 01-02-08 Add protection for small energies and optimise the code (V.Ivanchenko)
47// 01-06-08 Added initialisation of effective charge prestep (V.Ivanchenko)
48//
49// Class Description:
50//
51// -------------------------------------------------------------------
52//
53
54//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
55//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
56
57#include "G4IonFluctuations.hh"
60#include "G4SystemOfUnits.hh"
61#include "Randomize.hh"
62#include "G4Poisson.hh"
64#include "G4Material.hh"
65#include "G4DynamicParticle.hh"
66#include "G4Pow.hh"
67#include "G4Log.hh"
68
69//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
70
71using namespace std;
72
75 particleMass(CLHEP::proton_mass_c2),
76 parameter(10.0*CLHEP::MeV/CLHEP::proton_mass_c2),
77 theBohrBeta2(50.0*keV/CLHEP::proton_mass_c2),
78 minLoss(0.001*CLHEP::eV)
79{
80 uniFluct = new G4UniversalFluctuation();
81 g4calc = G4Pow::GetInstance();
82}
83
84//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
85
87
88//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
89
91{
92 particle = part;
93 particleMass = part->GetPDGMass();
94 charge = part->GetPDGCharge()/eplus;
95 effChargeSquare = chargeSquare = charge*charge;
96 uniFluct->InitialiseMe(part);
97}
98
99//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
100
103 const G4DynamicParticle* dp,
104 const G4double tcut,
105 const G4double tmax,
106 const G4double length,
107 const G4double meanLoss)
108{
109 // G4cout << "### meanLoss= " << meanLoss << G4endl;
110 if(meanLoss <= minLoss) return meanLoss;
111
112 //G4cout << "G4IonFluctuations::SampleFluctuations E(MeV)= "
113 // << dp->GetKineticEnergy()
114 // << " Elim(MeV)= " << parameter*charge*particleMass << G4endl;
115
116 // Vavilov fluctuations above energy threshold
117 if(dp->GetKineticEnergy() > parameter*charge*particleMass) {
118 return uniFluct->SampleFluctuations(couple,dp,tcut,tmax,length,meanLoss);
119 }
120
121 const G4Material* material = couple->GetMaterial();
122 G4double siga = Dispersion(material,dp,tcut,tmax,length);
123 G4double loss = meanLoss;
124
125 //G4cout << "### siga= " << sqrt(siga) << " navr= " << navr << G4endl;
126
127 // Gaussian fluctuation
128
129 // Increase fluctuations for big fractional energy loss
130 //G4cout << "siga= " << siga << G4endl;
131 if ( meanLoss > minFraction*kineticEnergy ) {
132 G4double gam = (kineticEnergy - meanLoss)/particleMass + 1.0;
133 G4double b2 = 1.0 - 1.0/(gam*gam);
134 if(b2 < xmin*beta2) b2 = xmin*beta2;
135 G4double x = b2/beta2;
136 G4double x3 = 1.0/(x*x*x);
137 siga *= 0.25*(1.0 + x)*(x3 + (1.0/b2 - 0.5)/(1.0/beta2 - 0.5) );
138 }
139 siga = std::sqrt(siga);
140 G4double sn = meanLoss/siga;
141 G4double twomeanLoss = meanLoss + meanLoss;
142 // G4cout << "siga= " << siga << " sn= " << sn << G4endl;
143
144 CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
145 // thick target case
146 if (sn >= 2.0) {
147
148 do {
149 loss = G4RandGauss::shoot(rndmEngine,meanLoss,siga);
150 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
151 } while (0.0 > loss || twomeanLoss < loss);
152
153 // Gamma distribution
154 } else if(sn > 0.1) {
155
156 G4double neff = sn*sn;
157 loss = meanLoss*G4RandGamma::shoot(rndmEngine,neff,1.0)/neff;
158
159 // uniform distribution for very small steps
160 } else {
161 loss = twomeanLoss*rndmEngine->flat();
162 }
163
164 //G4cout << "meanLoss= " << meanLoss << " loss= " << loss << G4endl;
165 return loss;
166}
167
168//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
169
171 const G4DynamicParticle* dp,
172 const G4double tcut,
173 const G4double tmax,
174 const G4double length)
175{
176 if(dp->GetDefinition() != particle) { InitialiseMe(dp->GetDefinition()); }
177
178 const G4double beta = dp->GetBeta();
179 kineticEnergy = dp->GetKineticEnergy();
180 beta2 = beta*beta;
181
182 G4double siga = (tmax/beta2 - 0.5*tcut)*CLHEP::twopi_mc2_rcl2*length*
183 material->GetElectronDensity()*effChargeSquare;
184
185 // Low velocity - additional ion charge fluctuations according to
186 // Q.Yang et al., NIM B61(1991)149-155.
187 //G4cout << "sigE= " << sqrt(siga) << " charge= " << charge <<G4endl;
188
189 G4double Z = material->GetIonisation()->GetZeffective();
190 G4double fac = Factor(material, Z);
191
192 // heavy ion correction
193 // G4double f1 = 1.065e-4*chargeSquare;
194 // if(beta2 > theBohrBeta2) f1/= beta2;
195 // else f1/= theBohrBeta2;
196 // if(f1 > 2.5) f1 = 2.5;
197 // fac *= (1.0 + f1);
198
199 // taking into account the cutg
200 G4double fac_cut = 1.0 + (fac - 1.0)*2.0*CLHEP::electron_mass_c2*beta2
201 /(tmax*(1.0 - beta2));
202 if(fac_cut > 0.01 && fac > 0.01) {
203 siga *= fac_cut;
204 }
205 /*
206 G4cout << "siga(keV)= " << sqrt(siga)/keV << " fac= " << fac
207 << " f1= " << fac_cut << G4endl;
208 */
209 return siga;
210}
211
212//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
213
214G4double G4IonFluctuations::Factor(const G4Material* material, G4double Z)
215{
216 // The aproximation of energy loss fluctuations
217 // Q.Yang et al., NIM B61(1991)149-155.
218
219 // Reduced energy in MeV/AMU
220 G4double energy = kineticEnergy*CLHEP::amu_c2/(particleMass*CLHEP::MeV);
221
222 // simple approximation for higher beta2
223 G4double s1 = RelativisticFactor(material, Z);
224
225 // tabulation for lower beta2
226 if( beta2 < 3.0*theBohrBeta2*Z ) {
227
228 static const G4double a[96][4] = {
229 {-0.3291, -0.8312, 0.2460, -1.0220},
230 {-0.5615, -0.5898, 0.5205, -0.7258},
231 {-0.5280, -0.4981, 0.5519, -0.5865},
232 {-0.5125, -0.4625, 0.5660, -0.5190},
233 {-0.5127, -0.8595, 0.5626, -0.8721},
234 {-0.5174, -1.1930, 0.5565, -1.1980},
235 {-0.5179, -1.1850, 0.5560, -1.2070},
236 {-0.5209, -0.9355, 0.5590, -1.0250},
237 {-0.5255, -0.7766, 0.5720, -0.9412},
238
239 {-0.5776, -0.6665, 0.6598, -0.8484},
240 {-0.6013, -0.6045, 0.7321, -0.7671},
241 {-0.5781, -0.5518, 0.7605, -0.6919},
242 {-0.5587, -0.4981, 0.7835, -0.6195},
243 {-0.5466, -0.4656, 0.7978, -0.5771},
244 {-0.5406, -0.4690, 0.8031, -0.5718},
245 {-0.5391, -0.5061, 0.8024, -0.5974},
246 {-0.5380, -0.6483, 0.7962, -0.6970},
247 {-0.5355, -0.7722, 0.7962, -0.7839},
248 {-0.5329, -0.7720, 0.7988, -0.7846},
249
250 {-0.5335, -0.7671, 0.7984, -0.7933},
251 {-0.5324, -0.7612, 0.7998, -0.8031},
252 {-0.5305, -0.7300, 0.8031, -0.7990},
253 {-0.5307, -0.7178, 0.8049, -0.8216},
254 {-0.5248, -0.6621, 0.8165, -0.7919},
255 {-0.5180, -0.6502, 0.8266, -0.7986},
256 {-0.5084, -0.6408, 0.8396, -0.8048},
257 {-0.4967, -0.6331, 0.8549, -0.8093},
258 {-0.4861, -0.6508, 0.8712, -0.8432},
259 {-0.4700, -0.6186, 0.8961, -0.8132},
260
261 {-0.4545, -0.5720, 0.9227, -0.7710},
262 {-0.4404, -0.5226, 0.9481, -0.7254},
263 {-0.4288, -0.4778, 0.9701, -0.6850},
264 {-0.4199, -0.4425, 0.9874, -0.6539},
265 {-0.4131, -0.4188, 0.9998, -0.6332},
266 {-0.4089, -0.4057, 1.0070, -0.6218},
267 {-0.4039, -0.3913, 1.0150, -0.6107},
268 {-0.3987, -0.3698, 1.0240, -0.5938},
269 {-0.3977, -0.3608, 1.0260, -0.5852},
270 {-0.3972, -0.3600, 1.0260, -0.5842},
271
272 {-0.3985, -0.3803, 1.0200, -0.6013},
273 {-0.3985, -0.3979, 1.0150, -0.6168},
274 {-0.3968, -0.3990, 1.0160, -0.6195},
275 {-0.3971, -0.4432, 1.0050, -0.6591},
276 {-0.3944, -0.4665, 1.0010, -0.6825},
277 {-0.3924, -0.5109, 0.9921, -0.7235},
278 {-0.3882, -0.5158, 0.9947, -0.7343},
279 {-0.3838, -0.5125, 0.9999, -0.7370},
280 {-0.3786, -0.4976, 1.0090, -0.7310},
281 {-0.3741, -0.4738, 1.0200, -0.7155},
282
283 {-0.3969, -0.4496, 1.0320, -0.6982},
284 {-0.3663, -0.4297, 1.0430, -0.6828},
285 {-0.3630, -0.4120, 1.0530, -0.6689},
286 {-0.3597, -0.3964, 1.0620, -0.6564},
287 {-0.3555, -0.3809, 1.0720, -0.6454},
288 {-0.3525, -0.3607, 1.0820, -0.6289},
289 {-0.3505, -0.3465, 1.0900, -0.6171},
290 {-0.3397, -0.3570, 1.1020, -0.6384},
291 {-0.3314, -0.3552, 1.1130, -0.6441},
292 {-0.3235, -0.3531, 1.1230, -0.6498},
293
294 {-0.3150, -0.3483, 1.1360, -0.6539},
295 {-0.3060, -0.3441, 1.1490, -0.6593},
296 {-0.2968, -0.3396, 1.1630, -0.6649},
297 {-0.2935, -0.3225, 1.1760, -0.6527},
298 {-0.2797, -0.3262, 1.1940, -0.6722},
299 {-0.2704, -0.3202, 1.2100, -0.6770},
300 {-0.2815, -0.3227, 1.2480, -0.6775},
301 {-0.2880, -0.3245, 1.2810, -0.6801},
302 {-0.3034, -0.3263, 1.3270, -0.6778},
303 {-0.2936, -0.3215, 1.3430, -0.6835},
304
305 {-0.3282, -0.3200, 1.3980, -0.6650},
306 {-0.3260, -0.3070, 1.4090, -0.6552},
307 {-0.3511, -0.3074, 1.4470, -0.6442},
308 {-0.3501, -0.3064, 1.4500, -0.6442},
309 {-0.3490, -0.3027, 1.4550, -0.6418},
310 {-0.3487, -0.3048, 1.4570, -0.6447},
311 {-0.3478, -0.3074, 1.4600, -0.6483},
312 {-0.3501, -0.3283, 1.4540, -0.6669},
313 {-0.3494, -0.3373, 1.4550, -0.6765},
314 {-0.3485, -0.3373, 1.4570, -0.6774},
315
316 {-0.3462, -0.3300, 1.4630, -0.6728},
317 {-0.3462, -0.3225, 1.4690, -0.6662},
318 {-0.3453, -0.3094, 1.4790, -0.6553},
319 {-0.3844, -0.3134, 1.5240, -0.6412},
320 {-0.3848, -0.3018, 1.5310, -0.6303},
321 {-0.3862, -0.2955, 1.5360, -0.6237},
322 {-0.4262, -0.2991, 1.5860, -0.6115},
323 {-0.4278, -0.2910, 1.5900, -0.6029},
324 {-0.4303, -0.2817, 1.5940, -0.5927},
325 {-0.4315, -0.2719, 1.6010, -0.5829},
326
327 {-0.4359, -0.2914, 1.6050, -0.6010},
328 {-0.4365, -0.2982, 1.6080, -0.6080},
329 {-0.4253, -0.3037, 1.6120, -0.6150},
330 {-0.4335, -0.3245, 1.6160, -0.6377},
331 {-0.4307, -0.3292, 1.6210, -0.6447},
332 {-0.4284, -0.3204, 1.6290, -0.6380},
333 {-0.4227, -0.3217, 1.6360, -0.6438}
334 } ;
335
336 G4int iz = G4lrint(Z) - 2;
337 if( 0 > iz ) { iz = 0; }
338 else if(95 < iz ) { iz = 95; }
339
340 const G4double ss = 1.0 + a[iz][0]*g4calc->powA(energy,a[iz][1])+
341 + a[iz][2]*g4calc->powA(energy,a[iz][3]);
342
343 // protection for the validity range for low beta
344 static const G4double slim = 0.001;
345 if(ss < slim) { s1 = 1.0/slim; }
346 // for high value of beta
347 else if(s1*ss < 1.0) { s1 = 1.0/ss; }
348 }
349 G4int i = 0 ;
350 G4double factor = 1.0 ;
351
352 // The index of set of parameters i = 0 for protons(hadrons) in gases
353 // 1 for protons(hadrons) in solids
354 // 2 for ions in atomic gases
355 // 3 for ions in molecular gases
356 // 4 for ions in solids
357 static const G4double b[5][4] = {
358 {0.1014, 0.3700, 0.9642, 3.987},
359 {0.1955, 0.6941, 2.522, 1.040},
360 {0.05058, 0.08975, 0.1419, 10.80},
361 {0.05009, 0.08660, 0.2751, 3.787},
362 {0.01273, 0.03458, 0.3951, 3.812}
363 } ;
364
365 // protons (hadrons)
366 if(1.5 > charge) {
367 if( kStateGas != material->GetState() ) { i = 1; }
368
369 // ions
370 } else {
371
372 factor = charge * g4calc->A13(charge/Z);
373
374 if( kStateGas == material->GetState() ) {
375 energy /= (charge * std::sqrt(charge)) ;
376
377 if(1 == (material->GetNumberOfElements())) {
378 i = 2 ;
379 } else {
380 i = 3 ;
381 }
382
383 } else {
384 energy /= (charge * std::sqrt(charge*Z)) ;
385 i = 4 ;
386 }
387 }
388
389 G4double x = b[i][2];
390 G4double y = energy * b[i][3];
391 if(y <= 0.2) x *= (y*(1.0 - 0.5*y));
392 else x *= (1.0 - g4calc->expA(-y));
393
394 y = energy - b[i][1];
395
396 const G4double s2 = factor * x * b[i][0] / (y*y + x*x);
397 /*
398 G4cout << "s1= " << s1 << " s2= " << s2 << " q^2= " << effChargeSquare
399 << " e= " << energy << G4endl;
400 */
401 return s1*effChargeSquare/chargeSquare + s2;
402}
403
404//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
405
406G4double G4IonFluctuations::RelativisticFactor(const G4Material* mat,
407 G4double Z)
408{
411
412 // H.Geissel et al. NIM B, 195 (2002) 3.
413 G4double bF2= 2.0*eF/CLHEP::electron_mass_c2;
414 G4double f = 0.4*(1.0 - beta2)/((1.0 - 0.5*beta2)*Z);
415 if(beta2 > bF2) f *= G4Log(2.0*CLHEP::electron_mass_c2*beta2/I)*bF2/beta2;
416 else f *= G4Log(4.0*eF/I);
417
418 // G4cout << "f= " << f << " beta2= " << beta2
419 // << " bf2= " << bF2 << " q^2= " << chargeSquare << G4endl;
420
421 return 1.0 + f;
422}
423
424//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
425
427 G4double q2)
428{
429 if(part != particle) {
430 particle = part;
431 particleMass = part->GetPDGMass();
432 charge = part->GetPDGCharge()/eplus;
433 chargeSquare = charge*charge;
434 }
435 effChargeSquare = q2;
436 uniFluct->SetParticleAndCharge(part, q2);
437}
438
439//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
G4double G4Log(G4double x)
Definition: G4Log.hh:227
@ kStateGas
Definition: G4Material.hh:110
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
virtual double flat()=0
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double GetBeta() const
~G4IonFluctuations() 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
G4IonFluctuations(const G4String &nam="IonFluc")
void SetParticleAndCharge(const G4ParticleDefinition *, G4double q2) override
void InitialiseMe(const G4ParticleDefinition *) override
G4double GetMeanExcitationEnergy() const
G4double GetZeffective() const
G4double GetFermiEnergy() const
const G4Material * GetMaterial() const
G4State GetState() const
Definition: G4Material.hh:176
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:221
size_t GetNumberOfElements() const
Definition: G4Material.hh:181
G4double GetElectronDensity() const
Definition: G4Material.hh:212
G4double GetPDGCharge() const
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4double A13(G4double A) const
Definition: G4Pow.cc:116
G4double expA(G4double A) const
Definition: G4Pow.hh:203
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:230
G4double SampleFluctuations(const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double, const G4double, const G4double, const G4double) override
void InitialiseMe(const G4ParticleDefinition *) override
void SetParticleAndCharge(const G4ParticleDefinition *, G4double q2) override
Definition: DoubConv.h:17
G4double energy(const ThreeVector &p, const G4double m)
int G4lrint(double ad)
Definition: templates.hh:134