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