Geant4-11
Public Member Functions | Private Member Functions | Private Attributes | Static Private Attributes
G4QAOLowEnergyLoss Class Reference

#include <G4QAOLowEnergyLoss.hh>

Inheritance diagram for G4QAOLowEnergyLoss:
G4VLowEnergyModel

Public Member Functions

 G4QAOLowEnergyLoss (const G4String &name)
 
G4double HighEnergyLimit (const G4ParticleDefinition *aParticle) const override
 
G4double HighEnergyLimit (const G4ParticleDefinition *aParticle, const G4Material *material) const override
 
G4bool IsInCharge (const G4DynamicParticle *particle, const G4Material *material) const override
 
G4bool IsInCharge (const G4ParticleDefinition *aParticle, const G4Material *material) const override
 
G4double LowEnergyLimit (const G4ParticleDefinition *aParticle) const override
 
G4double LowEnergyLimit (const G4ParticleDefinition *aParticle, const G4Material *material) const override
 
G4double TheValue (const G4DynamicParticle *particle, const G4Material *material) override
 
G4double TheValue (const G4ParticleDefinition *aParticle, const G4Material *material, G4double kineticEnergy) override
 
 ~G4QAOLowEnergyLoss ()
 

Private Member Functions

G4double EnergyLoss (const G4Material *material, G4double kineticEnergy, G4double zParticle) const
 
G4double GetL0 (G4double normEnergy) const
 
G4double GetL1 (G4double normEnergy) const
 
G4double GetL2 (G4double normEnergy) const
 
G4int GetNumberOfShell (const G4Material *material) const
 
G4double GetOccupationNumber (G4int Z, G4int ShellNb) const
 
G4double GetOscillatorEnergy (const G4Material *material, G4int nbOfTheShell) const
 
G4double GetShellEnergy (const G4Material *material, G4int nbOfTheShell) const
 
G4double GetShellStrength (const G4Material *material, G4int nbOfTheShell) const
 

Private Attributes

G4int numberOfMaterials
 
G4int sizeL0
 
G4int sizeL1
 
G4int sizeL2
 

Static Private Attributes

static const G4double alShellEnergy [3] ={ 2795e-6, 202e-6, 16.9e-6}
 
static const G4double alShellStrength [3] ={ 0.1349, 0.6387, 0.2264}
 
static const G4double auShellEnergy [6] ={ 96235e-6, 25918e-6, 4116e-6, 599e-6, 87.3e-6, 36.9e-6}
 
static const G4double auShellStrength [6] ={ 0.0139, 0.0803, 0.2473, 0.423, 0.1124, 0.1231}
 
static const G4double cuShellEnergy [4] ={ 16931e-6, 1930e-6, 199e-6, 39.6e-6}
 
static const G4double cuShellStrength [4] ={ 0.0505, 0.2561, 0.4913, 0.2021}
 
static const G4int fNumberOfShells [101]
 
static const G4double L0 [67][2]
 
static const G4double L1 [22][2]
 
static const G4double L2 [14][2]
 
static const G4int materialAvailable [6] = {13,14,29,73,79,78}
 
static const G4int nbOfElectronPerSubShell [1540]
 
static const G4int nbofShellForMaterial [6] = {3,3,4,6,6,6 }
 
static const G4double ptShellEnergy [6] ={ 95017e-6, 25590e-6, 4063e-6, 576e-6, 81.9e-6, 31.4e-6}
 
static const G4double ptShellStrength [6] ={ 0.0129, 0.0745, 0.2295, 0.4627, 0.1324, 0.0879}
 
static const G4double siShellEnergy [3] ={ 3179e-6, 249e-6, 20.3e-6 }
 
static const G4double siShellStrength [3] ={ 0.1222, 0.5972, 0.2806}
 
static const G4double taShellEnergy [6] ={ 88926e-6, 18012e-6, 3210e-6, 575e-6, 108.7e-6, 30.8e-6}
 
static const G4double taShellStrength [6] ={ 0.0126, 0.0896, 0.2599, 0.3413, 0.2057, 0.0908}
 

Detailed Description

Definition at line 53 of file G4QAOLowEnergyLoss.hh.

Constructor & Destructor Documentation

◆ G4QAOLowEnergyLoss()

G4QAOLowEnergyLoss::G4QAOLowEnergyLoss ( const G4String name)
explicit

Definition at line 70 of file G4QAOLowEnergyLoss.cc.

72{
74 sizeL0 = 67;
75 sizeL1 = 22;
76 sizeL2 = 14;
77}
G4VLowEnergyModel(const G4String &name)
const char * name(G4int ptype)

References numberOfMaterials, sizeL0, sizeL1, and sizeL2.

◆ ~G4QAOLowEnergyLoss()

G4QAOLowEnergyLoss::~G4QAOLowEnergyLoss ( )

Definition at line 81 of file G4QAOLowEnergyLoss.cc.

82{;}

Member Function Documentation

◆ EnergyLoss()

G4double G4QAOLowEnergyLoss::EnergyLoss ( const G4Material material,
G4double  kineticEnergy,
G4double  zParticle 
) const
private

Definition at line 174 of file G4QAOLowEnergyLoss.cc.

177{
178 G4int nbOfShell = GetNumberOfShell(material);
179 if(nbOfShell < 1) nbOfShell = 1;
180 G4double dedx=0;
181
182 G4double v = c_light * std::sqrt( 2.0 * kineticEnergy / proton_mass_c2 );
183 G4double coeff = twopi * proton_mass_c2 *
184 (material-> GetTotNbOfElectPerVolume()) /
186 G4double fBetheVelocity = fine_structure_const * c_light / v;
188 kineticEnergy ;
189
190 G4double l0Term = 0, l1Term = 0, l2Term = 0;
191
192 for (G4int nos = 0 ; nos < nbOfShell ; nos++){
193
194 G4double l0 = 0, l1 = 0, l2 = 0;
195 G4double NormalizedEnergy = ( 2.0 * electron_mass_c2 * v * v ) /
197
198 G4double shStrength = GetShellStrength(material,nos);
199
200 l0 = GetL0(NormalizedEnergy);
201 l0Term += shStrength * l0;
202
203 l1 = GetL1(NormalizedEnergy);
204 l1Term += shStrength * l1;
205
206 l2 = GetL2(NormalizedEnergy);
207 l2Term += shStrength * l2;
208 }
209
210 dedx = coeff * zParticle * zParticle * (l0Term
211 + zParticle * fBetheVelocity * l1Term
212 + zParticle * zParticle * fBetheVelocity * fBetheVelocity * l2Term);
213
214 return dedx ;
215
216}
static constexpr double twopi
Definition: G4SIunits.hh:56
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
G4double GetShellEnergy(const G4Material *material, G4int nbOfTheShell) const
G4int GetNumberOfShell(const G4Material *material) const
G4double GetL0(G4double normEnergy) const
G4double GetShellStrength(const G4Material *material, G4int nbOfTheShell) const
G4double GetL2(G4double normEnergy) const
G4double GetL1(G4double normEnergy) const
string material
Definition: eplot.py:19
float electron_mass_c2
Definition: hepunit.py:273
float hbarc_squared
Definition: hepunit.py:265
float c_light
Definition: hepunit.py:256
float proton_mass_c2
Definition: hepunit.py:274
int fine_structure_const
Definition: hepunit.py:286
float c_squared
Definition: hepunit.py:257

References source.hepunit::c_light, source.hepunit::c_squared, source.hepunit::electron_mass_c2, source.hepunit::fine_structure_const, GetL0(), GetL1(), GetL2(), GetNumberOfShell(), GetShellEnergy(), GetShellStrength(), source.hepunit::hbarc_squared, eplot::material, source.hepunit::proton_mass_c2, and twopi.

Referenced by TheValue().

◆ GetL0()

G4double G4QAOLowEnergyLoss::GetL0 ( G4double  normEnergy) const
private

Definition at line 333 of file G4QAOLowEnergyLoss.cc.

334{
335 G4int n;
336
337 for(n = 0; n < sizeL0; n++) {
338 if( normEnergy < L0[n][0] ) break;
339 }
340 if(0 == n) n = 1 ;
341 if(n >= sizeL0) n = sizeL0 - 1 ;
342
343 G4double l0 = L0[n][1];
344 G4double l0p = L0[n-1][1];
345 G4double bethe = l0p + (l0 - l0p) * ( normEnergy - L0[n-1][0]) /
346 (L0[n][0] - L0[n-1][0]);
347 return bethe ;
348}
static const G4double L0[67][2]

References L0, CLHEP::detail::n, and sizeL0.

Referenced by EnergyLoss().

◆ GetL1()

G4double G4QAOLowEnergyLoss::GetL1 ( G4double  normEnergy) const
private

Definition at line 352 of file G4QAOLowEnergyLoss.cc.

353{
354 G4int n;
355
356 for(n = 0; n < sizeL1; n++) {
357 if( normEnergy < L1[n][0] ) break;
358 }
359 if(0 == n) n = 1 ;
360 if(n >= sizeL1) n = sizeL1 - 1 ;
361
362 G4double l1 = L1[n][1];
363 G4double l1p = L1[n-1][1];
364 G4double barkas= l1p + (l1 - l1p) * ( normEnergy - L1[n-1][0]) /
365 (L1[n][0] - L1[n-1][0]);
366 return barkas;
367}
static const G4double L1[22][2]

References L1, CLHEP::detail::n, and sizeL1.

Referenced by EnergyLoss().

◆ GetL2()

G4double G4QAOLowEnergyLoss::GetL2 ( G4double  normEnergy) const
private

Definition at line 371 of file G4QAOLowEnergyLoss.cc.

372{
373 G4int n;
374 for(n = 0; n < sizeL2; n++) {
375 if( normEnergy < L2[n][0] ) break;
376 }
377 if(0 == n) n = 1 ;
378 if(n >= sizeL2) n = sizeL2 - 1 ;
379
380 G4double l2 = L2[n][1];
381 G4double l2p = L2[n-1][1];
382 G4double bloch = l2p + (l2 - l2p) * ( normEnergy - L2[n-1][0]) /
383 (L2[n][0] - L2[n-1][0]);
384 return bloch;
385}
static const G4double L2[14][2]

References L2, CLHEP::detail::n, and sizeL2.

Referenced by EnergyLoss().

◆ GetNumberOfShell()

G4int G4QAOLowEnergyLoss::GetNumberOfShell ( const G4Material material) const
private

Definition at line 220 of file G4QAOLowEnergyLoss.cc.

221{
222 // Set return value from table
223 G4int Z = (G4int)(material->GetZ());
224 G4int nShell = 0;
225
226 // Set return value if in material available from Aahrus
227 for(G4int i=0; i<numberOfMaterials; i++) {
228
229 if(materialAvailable[i] == Z){
230 nShell = nbofShellForMaterial[i];
231 break;
232 }
233 else
234 nShell = fNumberOfShells[Z];
235 }
236 return nShell;
237}
const G4int Z[17]
static const G4int nbofShellForMaterial[6]
static const G4int materialAvailable[6]
static const G4int fNumberOfShells[101]

References fNumberOfShells, eplot::material, materialAvailable, nbofShellForMaterial, numberOfMaterials, and Z.

Referenced by EnergyLoss().

◆ GetOccupationNumber()

G4double G4QAOLowEnergyLoss::GetOccupationNumber ( G4int  Z,
G4int  ShellNb 
) const
private

Definition at line 322 of file G4QAOLowEnergyLoss.cc.

323{
324
325 G4int indice = ShellNb ;
326 for (G4int z = 1 ; z < Z ; z++) {indice += fNumberOfShells[z];}
327
328 return nbOfElectronPerSubShell[indice+1];
329}
static const G4int nbOfElectronPerSubShell[1540]

References fNumberOfShells, nbOfElectronPerSubShell, and Z.

Referenced by GetOscillatorEnergy(), and GetShellStrength().

◆ GetOscillatorEnergy()

G4double G4QAOLowEnergyLoss::GetOscillatorEnergy ( const G4Material material,
G4int  nbOfTheShell 
) const
private

Definition at line 270 of file G4QAOLowEnergyLoss.cc.

272{
273
274 const G4Element* element = material->GetElement(0);
275
276 G4int Z = (G4int)(element->GetZ());
277
278 G4double squaredPlasmonEnergy = 28.816 * 28.816 * 1e-6
279 * material->GetDensity()/g/cm3
280 * (Z/element->GetN()) ;
281
282 G4double plasmonTerm = 0.66667 * GetOccupationNumber(Z,nbOfTheShell)
283 * squaredPlasmonEnergy / (Z*Z) ;
284
285 G4double ionTerm = G4Exp(0.5) * (element->GetAtomicShell(nbOfTheShell)) ;
286 ionTerm = ionTerm*ionTerm ;
287
288 G4double oscShellEnergy = std::sqrt( ionTerm + plasmonTerm );
289 return oscShellEnergy;
290}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
static constexpr double cm3
Definition: G4SIunits.hh:101
static constexpr double g
Definition: G4SIunits.hh:168
G4double GetZ() const
Definition: G4Element.hh:131
G4double GetN() const
Definition: G4Element.hh:135
G4double GetAtomicShell(G4int index) const
Definition: G4Element.cc:366
G4double GetOccupationNumber(G4int Z, G4int ShellNb) const

References cm3, g, G4Exp(), G4Element::GetAtomicShell(), G4Element::GetN(), GetOccupationNumber(), G4Element::GetZ(), eplot::material, and Z.

Referenced by GetShellEnergy().

◆ GetShellEnergy()

G4double G4QAOLowEnergyLoss::GetShellEnergy ( const G4Material material,
G4int  nbOfTheShell 
) const
private

Definition at line 241 of file G4QAOLowEnergyLoss.cc.

243{
244 //
245 G4double shellEnergy = alShellEnergy[0];
246
247 if(material->GetZ() == 13) shellEnergy = alShellEnergy[nbOfTheShell];
248 else if(material->GetZ() == 14)shellEnergy = siShellEnergy[nbOfTheShell];
249 else if(material->GetZ() == 29)shellEnergy = cuShellEnergy[nbOfTheShell];
250 else if(material->GetZ() == 73)shellEnergy = taShellEnergy[nbOfTheShell];
251 else if(material->GetZ() == 79)shellEnergy = auShellEnergy[nbOfTheShell];
252 else if(material->GetZ() == 78)shellEnergy = ptShellEnergy[nbOfTheShell];
253 else if (material->GetNumberOfElements() == 1)
254 shellEnergy = GetOscillatorEnergy(material, nbOfTheShell);
255 else
256 {
258 ed << "The model is not available for "
259 << material->GetName()
260 << G4endl;
261 G4Exception("G4QAOLowEnergyLoss::GetShellEnergy()",
262 "em2638",JustWarning,ed);
263 }
264
265 return shellEnergy;
266}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
#define G4endl
Definition: G4ios.hh:57
G4double GetOscillatorEnergy(const G4Material *material, G4int nbOfTheShell) const
static const G4double alShellEnergy[3]
static const G4double siShellEnergy[3]
static const G4double cuShellEnergy[4]
static const G4double ptShellEnergy[6]
static const G4double auShellEnergy[6]
static const G4double taShellEnergy[6]

References alShellEnergy, auShellEnergy, cuShellEnergy, G4endl, G4Exception(), GetOscillatorEnergy(), JustWarning, eplot::material, ptShellEnergy, siShellEnergy, and taShellEnergy.

Referenced by EnergyLoss().

◆ GetShellStrength()

G4double G4QAOLowEnergyLoss::GetShellStrength ( const G4Material material,
G4int  nbOfTheShell 
) const
private

Definition at line 294 of file G4QAOLowEnergyLoss.cc.

296{
297 G4double shellStrength = alShellStrength[0];
298
299 if(material->GetZ() == 13) shellStrength = alShellStrength[nbOfTheShell];
300 else if(material->GetZ() == 14)shellStrength =siShellStrength[nbOfTheShell];
301 else if(material->GetZ() == 29)shellStrength =cuShellStrength[nbOfTheShell];
302 else if(material->GetZ() == 73)shellStrength =taShellStrength[nbOfTheShell];
303 else if(material->GetZ() == 79)shellStrength =auShellStrength[nbOfTheShell];
304 else if(material->GetZ() == 78)shellStrength =ptShellStrength[nbOfTheShell];
305 else if (material->GetNumberOfElements() == 1){
306 G4int Z = (G4int)(material->GetZ());
307 shellStrength = GetOccupationNumber(Z,nbOfTheShell) / Z ;}
308 else
309 {
311 ed << "The model is not available for "
312 << material->GetName()
313 << G4endl;
314 G4Exception("G4QAOLowEnergyLoss::GetShellStrength()",
315 "em2639",JustWarning,ed);
316 }
317 return shellStrength;
318}
static const G4double siShellStrength[3]
static const G4double ptShellStrength[6]
static const G4double cuShellStrength[4]
static const G4double taShellStrength[6]
static const G4double alShellStrength[3]
static const G4double auShellStrength[6]

References alShellStrength, auShellStrength, cuShellStrength, G4endl, G4Exception(), GetOccupationNumber(), JustWarning, eplot::material, ptShellStrength, siShellStrength, taShellStrength, and Z.

Referenced by EnergyLoss().

◆ HighEnergyLimit() [1/2]

G4double G4QAOLowEnergyLoss::HighEnergyLimit ( const G4ParticleDefinition aParticle) const
overridevirtual

Implements G4VLowEnergyModel.

Definition at line 102 of file G4QAOLowEnergyLoss.cc.

103{
104 return 2.0*MeV ;
105}
static constexpr double MeV
Definition: G4SIunits.hh:200

References MeV.

◆ HighEnergyLimit() [2/2]

G4double G4QAOLowEnergyLoss::HighEnergyLimit ( const G4ParticleDefinition aParticle,
const G4Material material 
) const
overridevirtual

Implements G4VLowEnergyModel.

Definition at line 86 of file G4QAOLowEnergyLoss.cc.

88{
89 return 2.0*MeV ;
90}

References MeV.

◆ IsInCharge() [1/2]

G4bool G4QAOLowEnergyLoss::IsInCharge ( const G4DynamicParticle particle,
const G4Material material 
) const
overridevirtual

Implements G4VLowEnergyModel.

Definition at line 116 of file G4QAOLowEnergyLoss.cc.

118{
119 G4bool isInCharge = false;
120 G4bool hasMaterial = false;
121
122 if (material->GetNumberOfElements() == 1) hasMaterial = true;
123
125 && hasMaterial) isInCharge = true;
126
127 return isInCharge;
128}
bool G4bool
Definition: G4Types.hh:86
static G4AntiProton * AntiProtonDefinition()
Definition: G4AntiProton.cc:87
G4ParticleDefinition * GetDefinition() const

References G4AntiProton::AntiProtonDefinition(), G4DynamicParticle::GetDefinition(), and eplot::material.

◆ IsInCharge() [2/2]

G4bool G4QAOLowEnergyLoss::IsInCharge ( const G4ParticleDefinition aParticle,
const G4Material material 
) const
overridevirtual

Implements G4VLowEnergyModel.

Definition at line 132 of file G4QAOLowEnergyLoss.cc.

134{
135
136 G4bool isInCharge = false;
137 G4bool hasMaterial = false;
138
139 if (material->GetNumberOfElements() == 1) hasMaterial = true;
140
141 if (aParticle == (G4AntiProton::AntiProtonDefinition())
142 && hasMaterial) isInCharge = true;
143 return isInCharge;
144}

References G4AntiProton::AntiProtonDefinition(), and eplot::material.

◆ LowEnergyLimit() [1/2]

G4double G4QAOLowEnergyLoss::LowEnergyLimit ( const G4ParticleDefinition aParticle) const
overridevirtual

Implements G4VLowEnergyModel.

Definition at line 109 of file G4QAOLowEnergyLoss.cc.

110{
111 return 5.0*keV ;
112}
static constexpr double keV
Definition: G4SIunits.hh:202

References keV.

◆ LowEnergyLimit() [2/2]

G4double G4QAOLowEnergyLoss::LowEnergyLimit ( const G4ParticleDefinition aParticle,
const G4Material material 
) const
overridevirtual

Implements G4VLowEnergyModel.

Definition at line 94 of file G4QAOLowEnergyLoss.cc.

96{
97 return 5.0*keV ;
98}

References keV.

◆ TheValue() [1/2]

G4double G4QAOLowEnergyLoss::TheValue ( const G4DynamicParticle particle,
const G4Material material 
)
overridevirtual

Implements G4VLowEnergyModel.

Definition at line 148 of file G4QAOLowEnergyLoss.cc.

150{
151 G4double zParticle = (G4int)(particle->GetCharge())/eplus;
152
153 G4double energy = particle->GetKineticEnergy() ;
154 G4double eloss = EnergyLoss(material,energy,zParticle) ;
155
156 return eloss ;
157}
static constexpr double eplus
Definition: G4SIunits.hh:184
G4double GetCharge() const
G4double GetKineticEnergy() const
G4double EnergyLoss(const G4Material *material, G4double kineticEnergy, G4double zParticle) const
G4double energy(const ThreeVector &p, const G4double m)

References G4INCL::KinematicsUtils::energy(), EnergyLoss(), eplus, G4DynamicParticle::GetCharge(), G4DynamicParticle::GetKineticEnergy(), and eplot::material.

◆ TheValue() [2/2]

G4double G4QAOLowEnergyLoss::TheValue ( const G4ParticleDefinition aParticle,
const G4Material material,
G4double  kineticEnergy 
)
overridevirtual

Implements G4VLowEnergyModel.

Definition at line 161 of file G4QAOLowEnergyLoss.cc.

164{
165 G4double zParticle = (aParticle->GetPDGCharge())/eplus;
166
167 G4double eloss = EnergyLoss(material,kineticEnergy,zParticle) ;
168
169 return eloss ;
170}
G4double GetPDGCharge() const

References EnergyLoss(), eplus, G4ParticleDefinition::GetPDGCharge(), and eplot::material.

Field Documentation

◆ alShellEnergy

const G4double G4QAOLowEnergyLoss::alShellEnergy ={ 2795e-6, 202e-6, 16.9e-6}
staticprivate

Definition at line 117 of file G4QAOLowEnergyLoss.hh.

Referenced by GetShellEnergy().

◆ alShellStrength

const G4double G4QAOLowEnergyLoss::alShellStrength ={ 0.1349, 0.6387, 0.2264}
staticprivate

Definition at line 118 of file G4QAOLowEnergyLoss.hh.

Referenced by GetShellStrength().

◆ auShellEnergy

const G4double G4QAOLowEnergyLoss::auShellEnergy ={ 96235e-6, 25918e-6, 4116e-6, 599e-6, 87.3e-6, 36.9e-6}
staticprivate

Definition at line 125 of file G4QAOLowEnergyLoss.hh.

Referenced by GetShellEnergy().

◆ auShellStrength

const G4double G4QAOLowEnergyLoss::auShellStrength ={ 0.0139, 0.0803, 0.2473, 0.423, 0.1124, 0.1231}
staticprivate

Definition at line 126 of file G4QAOLowEnergyLoss.hh.

Referenced by GetShellStrength().

◆ cuShellEnergy

const G4double G4QAOLowEnergyLoss::cuShellEnergy ={ 16931e-6, 1930e-6, 199e-6, 39.6e-6}
staticprivate

Definition at line 121 of file G4QAOLowEnergyLoss.hh.

Referenced by GetShellEnergy().

◆ cuShellStrength

const G4double G4QAOLowEnergyLoss::cuShellStrength ={ 0.0505, 0.2561, 0.4913, 0.2021}
staticprivate

Definition at line 122 of file G4QAOLowEnergyLoss.hh.

Referenced by GetShellStrength().

◆ fNumberOfShells

const G4int G4QAOLowEnergyLoss::fNumberOfShells
staticprivate
Initial value:
=
{
0 ,
1 , 1 , 2 , 2 , 3 , 3 , 4 , 4 , 3 , 4 ,
5 , 5 , 6 , 6 , 6 , 6 , 6 , 7 , 8 , 8 ,
9 , 9 , 9 , 9 , 9 , 9 , 9 , 10 , 10 , 10 ,
11 , 11 , 11 , 11 , 11 , 12 , 13 , 13 , 14 , 14 ,
14 , 14 , 14 , 14 , 14 , 15 , 15 , 15 , 16 , 16 ,
16 , 16 , 16 , 17 , 18 , 18 , 19 , 19 , 19 , 19 ,
19 , 19 , 19 , 20 , 19 , 19 , 19 , 19 , 19 , 20 ,
21 , 21 , 21 , 21 , 21 , 21 , 21 , 21 , 22 , 22 ,
23 , 23 , 23 , 23 , 24 , 24 , 25 , 25 , 26 , 26 ,
27 , 27 , 27 , 26 , 26 , 27 , 27 , 26 , 26 , 26
}

Definition at line 134 of file G4QAOLowEnergyLoss.hh.

Referenced by GetNumberOfShell(), and GetOccupationNumber().

◆ L0

const G4double G4QAOLowEnergyLoss::L0
staticprivate

Definition at line 130 of file G4QAOLowEnergyLoss.hh.

Referenced by GetL0().

◆ L1

const G4double G4QAOLowEnergyLoss::L1
staticprivate
Initial value:
=
{
{0.00, -0.000001},
{0.10, -0.00001},
{0.20, -0.00049},
{0.30, -0.00084},
{0.40, 0.00085},
{0.50, 0.00519},
{0.60, 0.01198},
{0.70, 0.02074},
{0.80, 0.03133},
{0.90, 0.04369},
{1.00, 0.06035},
{2.00, 0.24023},
{3.00, 0.44284},
{4.00, 0.62012},
{5.00, 0.77031},
{6.00, 0.90390},
{7.00, 1.02705},
{8.00, 1.10867},
{9.00, 1.17546},
{10.00, 1.21599},
{15.00, 1.24349},
{20.00, 1.16752}
}

Definition at line 131 of file G4QAOLowEnergyLoss.hh.

Referenced by GetL1().

◆ L2

const G4double G4QAOLowEnergyLoss::L2
staticprivate
Initial value:
=
{
{0.00, 0.000001},
{0.10, 0.00001},
{0.20, 0.00000},
{0.40, -0.00120},
{0.60, -0.00036},
{0.80, 0.00372},
{1.00, 0.01298},
{2.00, 0.08296},
{4.00, 0.21953},
{6.00, 0.23903},
{8.00, 0.20893},
{10.00, 0.10879},
{20.00, -0.88409},
{40.00, -1.13902}
}

Definition at line 132 of file G4QAOLowEnergyLoss.hh.

Referenced by GetL2().

◆ materialAvailable

const G4int G4QAOLowEnergyLoss::materialAvailable = {13,14,29,73,79,78}
staticprivate

Definition at line 137 of file G4QAOLowEnergyLoss.hh.

Referenced by GetNumberOfShell().

◆ nbOfElectronPerSubShell

const G4int G4QAOLowEnergyLoss::nbOfElectronPerSubShell
staticprivate

Definition at line 133 of file G4QAOLowEnergyLoss.hh.

Referenced by GetOccupationNumber().

◆ nbofShellForMaterial

const G4int G4QAOLowEnergyLoss::nbofShellForMaterial = {3,3,4,6,6,6 }
staticprivate

Definition at line 116 of file G4QAOLowEnergyLoss.hh.

Referenced by GetNumberOfShell().

◆ numberOfMaterials

G4int G4QAOLowEnergyLoss::numberOfMaterials
private

Definition at line 139 of file G4QAOLowEnergyLoss.hh.

Referenced by G4QAOLowEnergyLoss(), and GetNumberOfShell().

◆ ptShellEnergy

const G4double G4QAOLowEnergyLoss::ptShellEnergy ={ 95017e-6, 25590e-6, 4063e-6, 576e-6, 81.9e-6, 31.4e-6}
staticprivate

Definition at line 127 of file G4QAOLowEnergyLoss.hh.

Referenced by GetShellEnergy().

◆ ptShellStrength

const G4double G4QAOLowEnergyLoss::ptShellStrength ={ 0.0129, 0.0745, 0.2295, 0.4627, 0.1324, 0.0879}
staticprivate

Definition at line 128 of file G4QAOLowEnergyLoss.hh.

Referenced by GetShellStrength().

◆ siShellEnergy

const G4double G4QAOLowEnergyLoss::siShellEnergy ={ 3179e-6, 249e-6, 20.3e-6 }
staticprivate

Definition at line 119 of file G4QAOLowEnergyLoss.hh.

Referenced by GetShellEnergy().

◆ siShellStrength

const G4double G4QAOLowEnergyLoss::siShellStrength ={ 0.1222, 0.5972, 0.2806}
staticprivate

Definition at line 120 of file G4QAOLowEnergyLoss.hh.

Referenced by GetShellStrength().

◆ sizeL0

G4int G4QAOLowEnergyLoss::sizeL0
private

Definition at line 140 of file G4QAOLowEnergyLoss.hh.

Referenced by G4QAOLowEnergyLoss(), and GetL0().

◆ sizeL1

G4int G4QAOLowEnergyLoss::sizeL1
private

Definition at line 141 of file G4QAOLowEnergyLoss.hh.

Referenced by G4QAOLowEnergyLoss(), and GetL1().

◆ sizeL2

G4int G4QAOLowEnergyLoss::sizeL2
private

Definition at line 142 of file G4QAOLowEnergyLoss.hh.

Referenced by G4QAOLowEnergyLoss(), and GetL2().

◆ taShellEnergy

const G4double G4QAOLowEnergyLoss::taShellEnergy ={ 88926e-6, 18012e-6, 3210e-6, 575e-6, 108.7e-6, 30.8e-6}
staticprivate

Definition at line 123 of file G4QAOLowEnergyLoss.hh.

Referenced by GetShellEnergy().

◆ taShellStrength

const G4double G4QAOLowEnergyLoss::taShellStrength ={ 0.0126, 0.0896, 0.2599, 0.3413, 0.2057, 0.0908}
staticprivate

Definition at line 124 of file G4QAOLowEnergyLoss.hh.

Referenced by GetShellStrength().


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