Geant4-11
Public Member Functions | Private Member Functions | Private Attributes
G4IonFluctuations Class Reference

#include <G4IonFluctuations.hh>

Inheritance diagram for G4IonFluctuations:
G4VEmFluctuationModel

Public Member Functions

G4double Dispersion (const G4Material *, const G4DynamicParticle *, const G4double tcut, const G4double tmax, const G4double length) override
 
 G4IonFluctuations (const G4IonFluctuations &)=delete
 
 G4IonFluctuations (const G4String &nam="IonFluc")
 
const G4StringGetName () const
 
void InitialiseMe (const G4ParticleDefinition *) override
 
G4IonFluctuationsoperator= (const G4IonFluctuations &right)=delete
 
G4double SampleFluctuations (const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double tcut, const G4double tmax, const G4double length, const G4double meanLoss) override
 
void SetParticleAndCharge (const G4ParticleDefinition *, G4double q2) override
 
 ~G4IonFluctuations () override
 

Private Member Functions

G4double Factor (const G4Material *, G4double Zeff)
 
G4double RelativisticFactor (const G4Material *, G4double Zeff)
 

Private Attributes

G4double beta2 = 0.0
 
G4double charge = 1.0
 
G4double chargeSquare = 1.0
 
G4double effChargeSquare = 1.0
 
G4LossTableManagerfManager
 
G4Powg4calc
 
G4double kineticEnergy = 0.0
 
G4double minFraction = 0.2
 
G4double minLoss
 
const G4String name
 
G4double parameter
 
const G4ParticleDefinitionparticle = nullptr
 
G4double particleMass
 
G4double theBohrBeta2
 
G4UniversalFluctuationuniFluct
 
G4double xmin = 0.2
 

Detailed Description

Definition at line 59 of file G4IonFluctuations.hh.

Constructor & Destructor Documentation

◆ G4IonFluctuations() [1/2]

G4IonFluctuations::G4IonFluctuations ( const G4String nam = "IonFluc")
explicit

Definition at line 73 of file G4IonFluctuations.cc.

78 minLoss(0.001*CLHEP::eV)
79{
82}
static constexpr double keV
Definition: G4SIunits.hh:202
G4UniversalFluctuation * uniFluct
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4VEmFluctuationModel(const G4String &nam)
static constexpr double proton_mass_c2
static constexpr double MeV
static constexpr double eV

References g4calc, G4Pow::GetInstance(), and uniFluct.

◆ ~G4IonFluctuations()

G4IonFluctuations::~G4IonFluctuations ( )
override

Definition at line 86 of file G4IonFluctuations.cc.

87{}

◆ G4IonFluctuations() [2/2]

G4IonFluctuations::G4IonFluctuations ( const G4IonFluctuations )
delete

Member Function Documentation

◆ Dispersion()

G4double G4IonFluctuations::Dispersion ( const G4Material material,
const G4DynamicParticle dp,
const G4double  tcut,
const G4double  tmax,
const G4double  length 
)
overridevirtual

Implements G4VEmFluctuationModel.

Definition at line 171 of file G4IonFluctuations.cc.

176{
177 if(dp->GetDefinition() != particle) { InitialiseMe(dp->GetDefinition()); }
178
179 const G4double beta = dp->GetBeta();
181 beta2 = beta*beta;
182
183 G4double siga = (tmax/beta2 - 0.5*tcut)*CLHEP::twopi_mc2_rcl2*length*
184 material->GetElectronDensity()*effChargeSquare;
185
186 // Low velocity - additional ion charge fluctuations according to
187 // Q.Yang et al., NIM B61(1991)149-155.
188 //G4cout << "sigE= " << sqrt(siga) << " charge= " << charge <<G4endl;
189
190 G4double Z = material->GetIonisation()->GetZeffective();
192
193 // heavy ion correction
194 // G4double f1 = 1.065e-4*chargeSquare;
195 // if(beta2 > theBohrBeta2) f1/= beta2;
196 // else f1/= theBohrBeta2;
197 // if(f1 > 2.5) f1 = 2.5;
198 // fac *= (1.0 + f1);
199
200 // taking into account the cutg
201 G4double fac_cut = 1.0 + (fac - 1.0)*2.0*CLHEP::electron_mass_c2*beta2
202 /(tmax*(1.0 - beta2));
203 if(fac_cut > 0.01 && fac > 0.01) {
204 siga *= fac_cut;
205 }
206 /*
207 G4cout << "siga(keV)= " << sqrt(siga)/keV << " fac= " << fac
208 << " f1= " << fac_cut << G4endl;
209 */
210 return siga;
211}
static const G4double fac
double G4double
Definition: G4Types.hh:83
const G4int Z[17]
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double GetBeta() const
const G4ParticleDefinition * particle
G4double Factor(const G4Material *, G4double Zeff)
void InitialiseMe(const G4ParticleDefinition *) override
static constexpr double electron_mass_c2
static constexpr double twopi_mc2_rcl2
string material
Definition: eplot.py:19

References anonymous_namespace{G4PionRadiativeDecayChannel.cc}::beta, beta2, effChargeSquare, CLHEP::electron_mass_c2, fac, Factor(), G4DynamicParticle::GetBeta(), G4DynamicParticle::GetDefinition(), G4DynamicParticle::GetKineticEnergy(), InitialiseMe(), kineticEnergy, eplot::material, particle, CLHEP::twopi_mc2_rcl2, and Z.

Referenced by SampleFluctuations().

◆ Factor()

G4double G4IonFluctuations::Factor ( const G4Material material,
G4double  Zeff 
)
private

Definition at line 215 of file G4IonFluctuations.cc.

216{
217 // The aproximation of energy loss fluctuations
218 // Q.Yang et al., NIM B61(1991)149-155.
219
220 // Reduced energy in MeV/AMU
222
223 // simple approximation for higher beta2
225
226 // tabulation for lower beta2
227 if( beta2 < 3.0*theBohrBeta2*Z ) {
228
229 static const G4double a[96][4] = {
230 {-0.3291, -0.8312, 0.2460, -1.0220},
231 {-0.5615, -0.5898, 0.5205, -0.7258},
232 {-0.5280, -0.4981, 0.5519, -0.5865},
233 {-0.5125, -0.4625, 0.5660, -0.5190},
234 {-0.5127, -0.8595, 0.5626, -0.8721},
235 {-0.5174, -1.1930, 0.5565, -1.1980},
236 {-0.5179, -1.1850, 0.5560, -1.2070},
237 {-0.5209, -0.9355, 0.5590, -1.0250},
238 {-0.5255, -0.7766, 0.5720, -0.9412},
239
240 {-0.5776, -0.6665, 0.6598, -0.8484},
241 {-0.6013, -0.6045, 0.7321, -0.7671},
242 {-0.5781, -0.5518, 0.7605, -0.6919},
243 {-0.5587, -0.4981, 0.7835, -0.6195},
244 {-0.5466, -0.4656, 0.7978, -0.5771},
245 {-0.5406, -0.4690, 0.8031, -0.5718},
246 {-0.5391, -0.5061, 0.8024, -0.5974},
247 {-0.5380, -0.6483, 0.7962, -0.6970},
248 {-0.5355, -0.7722, 0.7962, -0.7839},
249 {-0.5329, -0.7720, 0.7988, -0.7846},
250
251 {-0.5335, -0.7671, 0.7984, -0.7933},
252 {-0.5324, -0.7612, 0.7998, -0.8031},
253 {-0.5305, -0.7300, 0.8031, -0.7990},
254 {-0.5307, -0.7178, 0.8049, -0.8216},
255 {-0.5248, -0.6621, 0.8165, -0.7919},
256 {-0.5180, -0.6502, 0.8266, -0.7986},
257 {-0.5084, -0.6408, 0.8396, -0.8048},
258 {-0.4967, -0.6331, 0.8549, -0.8093},
259 {-0.4861, -0.6508, 0.8712, -0.8432},
260 {-0.4700, -0.6186, 0.8961, -0.8132},
261
262 {-0.4545, -0.5720, 0.9227, -0.7710},
263 {-0.4404, -0.5226, 0.9481, -0.7254},
264 {-0.4288, -0.4778, 0.9701, -0.6850},
265 {-0.4199, -0.4425, 0.9874, -0.6539},
266 {-0.4131, -0.4188, 0.9998, -0.6332},
267 {-0.4089, -0.4057, 1.0070, -0.6218},
268 {-0.4039, -0.3913, 1.0150, -0.6107},
269 {-0.3987, -0.3698, 1.0240, -0.5938},
270 {-0.3977, -0.3608, 1.0260, -0.5852},
271 {-0.3972, -0.3600, 1.0260, -0.5842},
272
273 {-0.3985, -0.3803, 1.0200, -0.6013},
274 {-0.3985, -0.3979, 1.0150, -0.6168},
275 {-0.3968, -0.3990, 1.0160, -0.6195},
276 {-0.3971, -0.4432, 1.0050, -0.6591},
277 {-0.3944, -0.4665, 1.0010, -0.6825},
278 {-0.3924, -0.5109, 0.9921, -0.7235},
279 {-0.3882, -0.5158, 0.9947, -0.7343},
280 {-0.3838, -0.5125, 0.9999, -0.7370},
281 {-0.3786, -0.4976, 1.0090, -0.7310},
282 {-0.3741, -0.4738, 1.0200, -0.7155},
283
284 {-0.3969, -0.4496, 1.0320, -0.6982},
285 {-0.3663, -0.4297, 1.0430, -0.6828},
286 {-0.3630, -0.4120, 1.0530, -0.6689},
287 {-0.3597, -0.3964, 1.0620, -0.6564},
288 {-0.3555, -0.3809, 1.0720, -0.6454},
289 {-0.3525, -0.3607, 1.0820, -0.6289},
290 {-0.3505, -0.3465, 1.0900, -0.6171},
291 {-0.3397, -0.3570, 1.1020, -0.6384},
292 {-0.3314, -0.3552, 1.1130, -0.6441},
293 {-0.3235, -0.3531, 1.1230, -0.6498},
294
295 {-0.3150, -0.3483, 1.1360, -0.6539},
296 {-0.3060, -0.3441, 1.1490, -0.6593},
297 {-0.2968, -0.3396, 1.1630, -0.6649},
298 {-0.2935, -0.3225, 1.1760, -0.6527},
299 {-0.2797, -0.3262, 1.1940, -0.6722},
300 {-0.2704, -0.3202, 1.2100, -0.6770},
301 {-0.2815, -0.3227, 1.2480, -0.6775},
302 {-0.2880, -0.3245, 1.2810, -0.6801},
303 {-0.3034, -0.3263, 1.3270, -0.6778},
304 {-0.2936, -0.3215, 1.3430, -0.6835},
305
306 {-0.3282, -0.3200, 1.3980, -0.6650},
307 {-0.3260, -0.3070, 1.4090, -0.6552},
308 {-0.3511, -0.3074, 1.4470, -0.6442},
309 {-0.3501, -0.3064, 1.4500, -0.6442},
310 {-0.3490, -0.3027, 1.4550, -0.6418},
311 {-0.3487, -0.3048, 1.4570, -0.6447},
312 {-0.3478, -0.3074, 1.4600, -0.6483},
313 {-0.3501, -0.3283, 1.4540, -0.6669},
314 {-0.3494, -0.3373, 1.4550, -0.6765},
315 {-0.3485, -0.3373, 1.4570, -0.6774},
316
317 {-0.3462, -0.3300, 1.4630, -0.6728},
318 {-0.3462, -0.3225, 1.4690, -0.6662},
319 {-0.3453, -0.3094, 1.4790, -0.6553},
320 {-0.3844, -0.3134, 1.5240, -0.6412},
321 {-0.3848, -0.3018, 1.5310, -0.6303},
322 {-0.3862, -0.2955, 1.5360, -0.6237},
323 {-0.4262, -0.2991, 1.5860, -0.6115},
324 {-0.4278, -0.2910, 1.5900, -0.6029},
325 {-0.4303, -0.2817, 1.5940, -0.5927},
326 {-0.4315, -0.2719, 1.6010, -0.5829},
327
328 {-0.4359, -0.2914, 1.6050, -0.6010},
329 {-0.4365, -0.2982, 1.6080, -0.6080},
330 {-0.4253, -0.3037, 1.6120, -0.6150},
331 {-0.4335, -0.3245, 1.6160, -0.6377},
332 {-0.4307, -0.3292, 1.6210, -0.6447},
333 {-0.4284, -0.3204, 1.6290, -0.6380},
334 {-0.4227, -0.3217, 1.6360, -0.6438}
335 } ;
336
337 G4int iz = G4lrint(Z) - 2;
338 if( 0 > iz ) { iz = 0; }
339 else if(95 < iz ) { iz = 95; }
340
341 const G4double ss = 1.0 + a[iz][0]*g4calc->powA(energy,a[iz][1])+
342 + a[iz][2]*g4calc->powA(energy,a[iz][3]);
343
344 // protection for the validity range for low beta
345 static const G4double slim = 0.001;
346 if(ss < slim) { s1 = 1.0/slim; }
347 // for high value of beta
348 else if(s1*ss < 1.0) { s1 = 1.0/ss; }
349 }
350 G4int i = 0 ;
351 G4double factor = 1.0 ;
352
353 // The index of set of parameters i = 0 for protons(hadrons) in gases
354 // 1 for protons(hadrons) in solids
355 // 2 for ions in atomic gases
356 // 3 for ions in molecular gases
357 // 4 for ions in solids
358 static const G4double b[5][4] = {
359 {0.1014, 0.3700, 0.9642, 3.987},
360 {0.1955, 0.6941, 2.522, 1.040},
361 {0.05058, 0.08975, 0.1419, 10.80},
362 {0.05009, 0.08660, 0.2751, 3.787},
363 {0.01273, 0.03458, 0.3951, 3.812}
364 } ;
365
366 // protons (hadrons)
367 if(1.5 > charge) {
368 if( kStateGas != material->GetState() ) { i = 1; }
369
370 // ions
371 } else {
372
373 factor = charge * g4calc->A13(charge/Z);
374
375 if( kStateGas == material->GetState() ) {
376 energy /= (charge * std::sqrt(charge)) ;
377
378 if(1 == (material->GetNumberOfElements())) {
379 i = 2 ;
380 } else {
381 i = 3 ;
382 }
383
384 } else {
385 energy /= (charge * std::sqrt(charge*Z)) ;
386 i = 4 ;
387 }
388 }
389
390 G4double x = b[i][2];
391 G4double y = energy * b[i][3];
392 if(y <= 0.2) x *= (y*(1.0 - 0.5*y));
393 else x *= (1.0 - g4calc->expA(-y));
394
395 y = energy - b[i][1];
396
397 const G4double s2 = factor * x * b[i][0] / (y*y + x*x);
398 /*
399 G4cout << "s1= " << s1 << " s2= " << s2 << " q^2= " << effChargeSquare
400 << " e= " << energy << G4endl;
401 */
402 return s1*effChargeSquare/chargeSquare + s2;
403}
@ kStateGas
Definition: G4Material.hh:111
int G4int
Definition: G4Types.hh:85
G4double RelativisticFactor(const G4Material *, G4double Zeff)
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
static constexpr double amu_c2
G4double energy(const ThreeVector &p, const G4double m)
int G4lrint(double ad)
Definition: templates.hh:134

References G4Pow::A13(), CLHEP::amu_c2, beta2, charge, chargeSquare, effChargeSquare, G4INCL::KinematicsUtils::energy(), G4Pow::expA(), g4calc, G4lrint(), kineticEnergy, kStateGas, eplot::material, CLHEP::MeV, particleMass, G4Pow::powA(), RelativisticFactor(), theBohrBeta2, and Z.

Referenced by Dispersion().

◆ GetName()

const G4String & G4VEmFluctuationModel::GetName ( ) const
inlineinherited

Definition at line 121 of file G4VEmFluctuationModel.hh.

122{
123 return name;
124}

References G4VEmFluctuationModel::name.

Referenced by G4LossTableManager::Register(), and G4EmConfigurator::SetModelForRegion().

◆ InitialiseMe()

void G4IonFluctuations::InitialiseMe ( const G4ParticleDefinition part)
overridevirtual

Reimplemented from G4VEmFluctuationModel.

Definition at line 91 of file G4IonFluctuations.cc.

92{
93 particle = part;
94 particleMass = part->GetPDGMass();
95 charge = part->GetPDGCharge()/eplus;
98}
static constexpr double eplus
Definition: G4SIunits.hh:184
G4double GetPDGCharge() const
void InitialiseMe(const G4ParticleDefinition *) override

References charge, chargeSquare, effChargeSquare, eplus, G4ParticleDefinition::GetPDGCharge(), G4ParticleDefinition::GetPDGMass(), G4UniversalFluctuation::InitialiseMe(), particle, particleMass, and uniFluct.

Referenced by Dispersion().

◆ operator=()

G4IonFluctuations & G4IonFluctuations::operator= ( const G4IonFluctuations right)
delete

◆ RelativisticFactor()

G4double G4IonFluctuations::RelativisticFactor ( const G4Material mat,
G4double  Zeff 
)
private

Definition at line 407 of file G4IonFluctuations.cc.

409{
412
413 // H.Geissel et al. NIM B, 195 (2002) 3.
415 G4double f = 0.4*(1.0 - beta2)/((1.0 - 0.5*beta2)*Z);
416 if(beta2 > bF2) f *= G4Log(2.0*CLHEP::electron_mass_c2*beta2/I)*bF2/beta2;
417 else f *= G4Log(4.0*eF/I);
418
419 // G4cout << "f= " << f << " beta2= " << beta2
420 // << " bf2= " << bF2 << " q^2= " << chargeSquare << G4endl;
421
422 return 1.0 + f;
423}
G4double G4Log(G4double x)
Definition: G4Log.hh:226
G4double GetMeanExcitationEnergy() const
G4double GetFermiEnergy() const
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:222

References beta2, CLHEP::electron_mass_c2, G4Log(), G4IonisParamMat::GetFermiEnergy(), G4Material::GetIonisation(), G4IonisParamMat::GetMeanExcitationEnergy(), and Z.

Referenced by Factor().

◆ SampleFluctuations()

G4double G4IonFluctuations::SampleFluctuations ( const G4MaterialCutsCouple couple,
const G4DynamicParticle dp,
const G4double  tcut,
const G4double  tmax,
const G4double  length,
const G4double  meanLoss 
)
overridevirtual

Implements G4VEmFluctuationModel.

Definition at line 103 of file G4IonFluctuations.cc.

109{
110 // G4cout << "### meanLoss= " << meanLoss << G4endl;
111 if(meanLoss <= minLoss) return meanLoss;
112
113 //G4cout << "G4IonFluctuations::SampleFluctuations E(MeV)= "
114 // << dp->GetKineticEnergy()
115 // << " Elim(MeV)= " << parameter*charge*particleMass << G4endl;
116
117 // Vavilov fluctuations above energy threshold
119 return uniFluct->SampleFluctuations(couple,dp,tcut,tmax,length,meanLoss);
120 }
121
122 const G4Material* material = couple->GetMaterial();
123 G4double siga = Dispersion(material,dp,tcut,tmax,length);
124 G4double loss = meanLoss;
125
126 //G4cout << "### siga= " << sqrt(siga) << " navr= " << navr << G4endl;
127
128 // Gaussian fluctuation
129
130 // Increase fluctuations for big fractional energy loss
131 //G4cout << "siga= " << siga << G4endl;
132 if ( meanLoss > minFraction*kineticEnergy ) {
133 G4double gam = (kineticEnergy - meanLoss)/particleMass + 1.0;
134 G4double b2 = 1.0 - 1.0/(gam*gam);
135 if(b2 < xmin*beta2) b2 = xmin*beta2;
136 G4double x = b2/beta2;
137 G4double x3 = 1.0/(x*x*x);
138 siga *= 0.25*(1.0 + x)*(x3 + (1.0/b2 - 0.5)/(1.0/beta2 - 0.5) );
139 }
140 siga = std::sqrt(siga);
141 G4double sn = meanLoss/siga;
142 G4double twomeanLoss = meanLoss + meanLoss;
143 // G4cout << "siga= " << siga << " sn= " << sn << G4endl;
144
145 CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
146 // thick target case
147 if (sn >= 2.0) {
148
149 do {
150 loss = G4RandGauss::shoot(rndmEngine,meanLoss,siga);
151 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
152 } while (0.0 > loss || twomeanLoss < loss);
153
154 // Gamma distribution
155 } else if(sn > 0.1) {
156
157 G4double neff = sn*sn;
158 loss = meanLoss*G4RandGamma::shoot(rndmEngine,neff,1.0)/neff;
159
160 // uniform distribution for very small steps
161 } else {
162 loss = twomeanLoss*rndmEngine->flat();
163 }
164
165 //G4cout << "meanLoss= " << meanLoss << " loss= " << loss << G4endl;
166 return loss;
167}
virtual double flat()=0
G4double Dispersion(const G4Material *, const G4DynamicParticle *, const G4double tcut, const G4double tmax, const G4double length) override
const G4Material * GetMaterial() const
G4double SampleFluctuations(const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double, const G4double, const G4double, const G4double) override
ThreeVector shoot(const G4int Ap, const G4int Af)

References beta2, charge, Dispersion(), CLHEP::HepRandomEngine::flat(), G4InuclParticleNames::gam, G4DynamicParticle::GetKineticEnergy(), G4MaterialCutsCouple::GetMaterial(), kineticEnergy, eplot::material, minFraction, minLoss, parameter, particleMass, G4UniversalFluctuation::SampleFluctuations(), G4INCL::DeJongSpin::shoot(), uniFluct, and xmin.

◆ SetParticleAndCharge()

void G4IonFluctuations::SetParticleAndCharge ( const G4ParticleDefinition part,
G4double  q2 
)
overridevirtual

Reimplemented from G4VEmFluctuationModel.

Definition at line 427 of file G4IonFluctuations.cc.

429{
430 if(part != particle) {
431 particle = part;
432 particleMass = part->GetPDGMass();
433 charge = part->GetPDGCharge()/eplus;
435 }
436 effChargeSquare = q2;
438}
void SetParticleAndCharge(const G4ParticleDefinition *, G4double q2) override

References charge, chargeSquare, effChargeSquare, eplus, G4ParticleDefinition::GetPDGCharge(), G4ParticleDefinition::GetPDGMass(), particle, particleMass, G4UniversalFluctuation::SetParticleAndCharge(), and uniFluct.

Field Documentation

◆ beta2

G4double G4IonFluctuations::beta2 = 0.0
private

Definition at line 116 of file G4IonFluctuations.hh.

Referenced by Dispersion(), Factor(), RelativisticFactor(), and SampleFluctuations().

◆ charge

G4double G4IonFluctuations::charge = 1.0
private

◆ chargeSquare

G4double G4IonFluctuations::chargeSquare = 1.0
private

Definition at line 105 of file G4IonFluctuations.hh.

Referenced by Factor(), InitialiseMe(), and SetParticleAndCharge().

◆ effChargeSquare

G4double G4IonFluctuations::effChargeSquare = 1.0
private

Definition at line 106 of file G4IonFluctuations.hh.

Referenced by Dispersion(), Factor(), InitialiseMe(), and SetParticleAndCharge().

◆ fManager

G4LossTableManager* G4VEmFluctuationModel::fManager
privateinherited

◆ g4calc

G4Pow* G4IonFluctuations::g4calc
private

Definition at line 101 of file G4IonFluctuations.hh.

Referenced by Factor(), and G4IonFluctuations().

◆ kineticEnergy

G4double G4IonFluctuations::kineticEnergy = 0.0
private

Definition at line 115 of file G4IonFluctuations.hh.

Referenced by Dispersion(), Factor(), and SampleFluctuations().

◆ minFraction

G4double G4IonFluctuations::minFraction = 0.2
private

Definition at line 111 of file G4IonFluctuations.hh.

Referenced by SampleFluctuations().

◆ minLoss

G4double G4IonFluctuations::minLoss
private

Definition at line 113 of file G4IonFluctuations.hh.

Referenced by SampleFluctuations().

◆ name

const G4String G4VEmFluctuationModel::name
privateinherited

◆ parameter

G4double G4IonFluctuations::parameter
private

Definition at line 109 of file G4IonFluctuations.hh.

Referenced by SampleFluctuations().

◆ particle

const G4ParticleDefinition* G4IonFluctuations::particle = nullptr
private

◆ particleMass

G4double G4IonFluctuations::particleMass
private

◆ theBohrBeta2

G4double G4IonFluctuations::theBohrBeta2
private

Definition at line 110 of file G4IonFluctuations.hh.

Referenced by Factor().

◆ uniFluct

G4UniversalFluctuation* G4IonFluctuations::uniFluct
private

◆ xmin

G4double G4IonFluctuations::xmin = 0.2
private

Definition at line 112 of file G4IonFluctuations.hh.

Referenced by SampleFluctuations().


The documentation for this class was generated from the following files: