Geant4  10.02.p03
G4QAOLowEnergyLoss Class Reference

#include <G4QAOLowEnergyLoss.hh>

Inheritance diagram for G4QAOLowEnergyLoss:
Collaboration diagram for G4QAOLowEnergyLoss:

Public Member Functions

 G4QAOLowEnergyLoss (const G4String &name)
 
 ~G4QAOLowEnergyLoss ()
 
G4double HighEnergyLimit (const G4ParticleDefinition *aParticle, const G4Material *material) const
 
G4double LowEnergyLimit (const G4ParticleDefinition *aParticle, const G4Material *material) const
 
G4double HighEnergyLimit (const G4ParticleDefinition *aParticle) const
 
G4double LowEnergyLimit (const G4ParticleDefinition *aParticle) const
 
G4bool IsInCharge (const G4DynamicParticle *particle, const G4Material *material) const
 
G4bool IsInCharge (const G4ParticleDefinition *aParticle, const G4Material *material) const
 
G4double TheValue (const G4DynamicParticle *particle, const G4Material *material)
 
G4double TheValue (const G4ParticleDefinition *aParticle, const G4Material *material, G4double kineticEnergy)
 
- Public Member Functions inherited from G4VLowEnergyModel
 G4VLowEnergyModel (const G4String &name)
 
virtual ~G4VLowEnergyModel ()
 

Private Member Functions

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

Private Attributes

G4int numberOfMaterials
 
G4int sizeL0
 
G4int sizeL1
 
G4int sizeL2
 

Static Private Attributes

static const G4int materialAvailable [6] = {13,14,29,73,79,78}
 
static const G4int nbofShellForMaterial [6] = {3,3,4,6,6,6 }
 
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 siShellEnergy [3] ={ 3179e-6, 249e-6, 20.3e-6 }
 
static const G4double siShellStrength [3] ={ 0.1222, 0.5972, 0.2806}
 
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 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}
 
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 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 L0 [67][2]
 
static const G4double L1 [22][2]
 
static const G4double L2 [14][2]
 
static const G4int nbOfElectronPerSubShell [1540]
 
static const G4int fNumberOfShells [101]
 

Detailed Description

Definition at line 53 of file G4QAOLowEnergyLoss.hh.

Constructor & Destructor Documentation

◆ G4QAOLowEnergyLoss()

G4QAOLowEnergyLoss::G4QAOLowEnergyLoss ( const G4String name)

Definition at line 70 of file G4QAOLowEnergyLoss.cc.

71  : G4VLowEnergyModel(name)
72 {
74  sizeL0 = 67;
75  sizeL1 = 22;
76  sizeL2 = 14;
77 }
G4VLowEnergyModel(const G4String &name)

◆ ~G4QAOLowEnergyLoss()

G4QAOLowEnergyLoss::~G4QAOLowEnergyLoss ( )

Definition at line 80 of file G4QAOLowEnergyLoss.cc.

81 {;}

Member Function Documentation

◆ EnergyLoss()

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

Definition at line 171 of file G4QAOLowEnergyLoss.cc.

174 {
175  G4int nbOfShell = GetNumberOfShell(material);
176  if(nbOfShell < 1) nbOfShell = 1;
177 
178  //G4int iz = (G4int)(material->GetZ());
179  //if(iz < 1) iz = 1;
180  //else if(iz > 100) iz = 100;
181  //G4int nbOfShell = fNumberOfShells[iz];
182 
183  /*
184  if(material->GetName()=="Graphite"){
185  G4cout << " E(MeV)= " << kineticEnergy/MeV << " n= "
186  << nbOfShell
187  << " for " << material->GetZ() << G4endl;
188  }
189  */
190  G4double dedx=0;
191 
192  G4double v = c_light * std::sqrt( 2.0 * kineticEnergy / proton_mass_c2 );
193  G4double coeff = twopi * proton_mass_c2 *
194  (material-> GetTotNbOfElectPerVolume()) /
196  G4double fBetheVelocity = fine_structure_const * c_light / v;
198  kineticEnergy ;
199 
200  //G4double beta = std::sqrt( 2.0 * kineticEnergy / proton_mass_c2 );
201  //G4double fBetheVelocity = std::sqrt( 2.0 * 25.0 * keV / proton_mass_c2 )/beta;
202  //G4double coeff= twopi_mc2_rcl2*(material->GetElectronDensity())/(beta*beta);
203 
204 
205  G4double l0Term = 0, l1Term = 0, l2Term = 0;
206 
207  for (G4int nos = 0 ; nos < nbOfShell ; nos++){
208 
209  G4double l0 = 0, l1 = 0, l2 = 0;
210  G4double NormalizedEnergy = ( 2.0 * electron_mass_c2 * v * v ) /
211  ( c_squared * GetShellEnergy(material,nos) );
212  //G4double NormalizedEnergy = ( 2.0 * electron_mass_c2 * beta * beta ) /
213  // GetShellEnergy(material,nos);
214  G4double shStrength = GetShellStrength(material,nos);
215 
216  l0 = GetL0(NormalizedEnergy);
217  l0Term += shStrength * l0;
218 
219  l1 = GetL1(NormalizedEnergy);
220  l1Term += shStrength * l1;
221 
222  l2 = GetL2(NormalizedEnergy);
223  l2Term += shStrength * l2;
224 
225  /*
226  if(material->GetName()=="Graphite"){
227  G4cout << nos << ". "
228  << " E(MeV)= " << kineticEnergy/MeV
229  << " normE= " << NormalizedEnergy
230  << " sh en= " << GetShellEnergy(material,nos)
231  << " str= " << shStrength
232  << " v0/v= " << fBetheVelocity
233  << " l0= " << l0Term
234  << " l1= " << l1Term
235  << " l2= " << l2Term
236  << G4endl;
237  }
238  */
239  }
240 
241  dedx = coeff * zParticle * zParticle * (l0Term
242  + zParticle * fBetheVelocity * l1Term
243  + zParticle * zParticle * fBetheVelocity * fBetheVelocity * l2Term);
244 
245  // G4cout << " E(MeV)= " << kineticEnergy/MeV
246  // << " dedx(Mev/mm)= " << dedx*mm/MeV << G4endl;
247 
248  return dedx ;
249 
250 }
float hbarc_squared
Definition: hepunit.py:266
float c_squared
Definition: hepunit.py:258
G4double GetShellStrength(const G4Material *material, G4int nbOfTheShell) const
G4double GetL0(G4double normEnergy) const
G4double GetShellEnergy(const G4Material *material, G4int nbOfTheShell) const
G4double GetL2(G4double normEnergy) const
int G4int
Definition: G4Types.hh:78
int fine_structure_const
Definition: hepunit.py:287
static const double twopi
Definition: G4SIunits.hh:75
float proton_mass_c2
Definition: hepunit.py:275
G4double GetL1(G4double normEnergy) const
float electron_mass_c2
Definition: hepunit.py:274
double G4double
Definition: G4Types.hh:76
G4int GetNumberOfShell(const G4Material *material) const
float c_light
Definition: hepunit.py:257
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetL0()

G4double G4QAOLowEnergyLoss::GetL0 ( G4double  normEnergy) const
private

Definition at line 366 of file G4QAOLowEnergyLoss.cc.

367 {
368  G4int n;
369 
370  for(n = 0; n < sizeL0; n++) {
371  if( normEnergy < L0[n][0] ) break;
372  }
373  if(0 == n) n = 1 ;
374  if(n >= sizeL0) n = sizeL0 - 1 ;
375 
376  G4double l0 = L0[n][1];
377  G4double l0p = L0[n-1][1];
378  G4double bethe = l0p + (l0 - l0p) * ( normEnergy - L0[n-1][0]) /
379  (L0[n][0] - L0[n-1][0]);
380  return bethe ;
381 
382 }
int G4int
Definition: G4Types.hh:78
Char_t n[5]
static const G4double L0[67][2]
double G4double
Definition: G4Types.hh:76
Here is the caller graph for this function:

◆ GetL1()

G4double G4QAOLowEnergyLoss::GetL1 ( G4double  normEnergy) const
private

Definition at line 384 of file G4QAOLowEnergyLoss.cc.

385 {
386  G4int n;
387 
388  for(n = 0; n < sizeL1; n++) {
389  if( normEnergy < L1[n][0] ) break;
390  }
391  if(0 == n) n = 1 ;
392  if(n >= sizeL1) n = sizeL1 - 1 ;
393 
394  G4double l1 = L1[n][1];
395  G4double l1p = L1[n-1][1];
396  G4double barkas= l1p + (l1 - l1p) * ( normEnergy - L1[n-1][0]) /
397  (L1[n][0] - L1[n-1][0]);
398 
399  return barkas;
400 
401 }
static const G4double L1[22][2]
int G4int
Definition: G4Types.hh:78
Char_t n[5]
double G4double
Definition: G4Types.hh:76
Here is the caller graph for this function:

◆ GetL2()

G4double G4QAOLowEnergyLoss::GetL2 ( G4double  normEnergy) const
private

Definition at line 404 of file G4QAOLowEnergyLoss.cc.

405 {
406  G4int n;
407  for(n = 0; n < sizeL2; n++) {
408  if( normEnergy < L2[n][0] ) break;
409  }
410  if(0 == n) n = 1 ;
411  if(n >= sizeL2) n = sizeL2 - 1 ;
412 
413  G4double l2 = L2[n][1];
414  G4double l2p = L2[n-1][1];
415  G4double bloch = l2p + (l2 - l2p) * ( normEnergy - L2[n-1][0]) /
416  (L2[n][0] - L2[n-1][0]);
417 
418  return bloch;
419 }
int G4int
Definition: G4Types.hh:78
Char_t n[5]
double G4double
Definition: G4Types.hh:76
static const G4double L2[14][2]
Here is the caller graph for this function:

◆ GetNumberOfShell()

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

Definition at line 253 of file G4QAOLowEnergyLoss.cc.

254 {
255  // Set return value from table
256  G4int Z = (G4int)(material->GetZ());
257  G4int nShell = 0;
258 
259  // Set return value if in material available from Aahrus
260  for(G4int i=0; i<numberOfMaterials; i++) {
261 
262  if(materialAvailable[i] == Z){
263  nShell = nbofShellForMaterial[i];
264  break;
265  }
266  else nShell = fNumberOfShells[Z];
267  }
268 
269  return nShell;
270 }
G4double GetZ() const
Definition: G4Material.cc:625
static const G4int nbofShellForMaterial[6]
static const G4int fNumberOfShells[101]
int G4int
Definition: G4Types.hh:78
Float_t Z
static const G4int materialAvailable[6]
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetOccupationNumber()

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

Definition at line 356 of file G4QAOLowEnergyLoss.cc.

357 {
358 
359  G4int indice = ShellNb ;
360  for (G4int z = 1 ; z < Z ; z++) {indice += fNumberOfShells[z];}
361 
362  return nbOfElectronPerSubShell[indice+1];
363 }
static const G4int fNumberOfShells[101]
int G4int
Definition: G4Types.hh:78
static const G4int nbOfElectronPerSubShell[1540]
Float_t Z
Here is the caller graph for this function:

◆ GetOscillatorEnergy()

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

Definition at line 297 of file G4QAOLowEnergyLoss.cc.

299 {
300 
301  const G4Element* element = material->GetElement(0);
302 
303  G4int Z = (G4int)(element->GetZ());
304 
305 G4double squaredPlasmonEnergy = 28.816 * 28.816 * 1e-6
306  * material->GetDensity()/g/cm3
307  * (Z/element->GetN()) ;
308 //G4double squaredPlasmonEnergy = 28.16 * 28.16 * 1e-6
309 // * (material->GetDensity()) * (cm3/g)
310 // * Z / (element->GetN()) ;
311 
312  G4double plasmonTerm = 0.66667 * GetOccupationNumber(Z,nbOfTheShell)
313  * squaredPlasmonEnergy / (Z*Z) ;
314 
315  G4double ionTerm = G4Exp(0.5) * (element->GetAtomicShell(nbOfTheShell)) ;
316 
317  ionTerm = ionTerm*ionTerm ;
318 
319  G4double oscShellEnergy = std::sqrt( ionTerm + plasmonTerm );
320 
321 /* if(material->GetName()=="Graphite"){
322  G4cout << "\t" << Z
323  << "\t" << nbOfTheShell
324  << "\t" << squaredPlasmonEnergy
325  << "\t" << plasmonTerm
326  << "\t" << ionTerm
327  << "\t" << GetOccupationNumber(Z,nbOfTheShell)
328  << "\t" << element->GetAtomicShell(nbOfTheShell)
329  << "\t" << oscShellEnergy << G4endl;}
330 */
331  return oscShellEnergy;
332 }
G4double GetAtomicShell(G4int index) const
Definition: G4Element.cc:365
G4double GetDensity() const
Definition: G4Material.hh:180
int G4int
Definition: G4Types.hh:78
G4double GetN() const
Definition: G4Element.hh:134
function g(Y1, Y2, PT2)
Definition: hijing1.383.f:5206
Float_t Z
static const double cm3
Definition: G4SIunits.hh:120
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:202
G4double GetOccupationNumber(G4int Z, G4int ShellNb) const
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
double G4double
Definition: G4Types.hh:76
G4double GetZ() const
Definition: G4Element.hh:131
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetShellEnergy()

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

Definition at line 274 of file G4QAOLowEnergyLoss.cc.

276 {
277  //
278  G4double shellEnergy = alShellEnergy[0];
279 
280  if(material->GetZ() == 13) shellEnergy = alShellEnergy[nbOfTheShell];
281  else if(material->GetZ() == 14)shellEnergy = siShellEnergy[nbOfTheShell];
282  else if(material->GetZ() == 29)shellEnergy = cuShellEnergy[nbOfTheShell];
283  else if(material->GetZ() == 73)shellEnergy = taShellEnergy[nbOfTheShell];
284  else if(material->GetZ() == 79)shellEnergy = auShellEnergy[nbOfTheShell];
285  else if(material->GetZ() == 78)shellEnergy = ptShellEnergy[nbOfTheShell];
286  else if (material->GetNumberOfElements() == 1)
287  shellEnergy = GetOscillatorEnergy(material, nbOfTheShell);
288  else G4cout << "WARNING - G4QAOLowEnergyLoss::GetShellEnergy - "
289  << "The model is not available for "
290  << material->GetName()
291  << G4endl;
292 
293  return shellEnergy;
294 }
G4double GetZ() const
Definition: G4Material.cc:625
static const G4double cuShellEnergy[4]
G4GLOB_DLL std::ostream G4cout
G4double GetOscillatorEnergy(const G4Material *material, G4int nbOfTheShell) const
static const G4double taShellEnergy[6]
static const G4double siShellEnergy[3]
static const G4double auShellEnergy[6]
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
#define G4endl
Definition: G4ios.hh:61
static const G4double alShellEnergy[3]
double G4double
Definition: G4Types.hh:76
const G4String & GetName() const
Definition: G4Material.hh:178
static const G4double ptShellEnergy[6]
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetShellStrength()

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

Definition at line 334 of file G4QAOLowEnergyLoss.cc.

336 {
337  G4double shellStrength = alShellStrength[0];
338 
339  if(material->GetZ() == 13) shellStrength = alShellStrength[nbOfTheShell];
340  else if(material->GetZ() == 14)shellStrength =siShellStrength[nbOfTheShell];
341  else if(material->GetZ() == 29)shellStrength =cuShellStrength[nbOfTheShell];
342  else if(material->GetZ() == 73)shellStrength =taShellStrength[nbOfTheShell];
343  else if(material->GetZ() == 79)shellStrength =auShellStrength[nbOfTheShell];
344  else if(material->GetZ() == 78)shellStrength =ptShellStrength[nbOfTheShell];
345  else if (material->GetNumberOfElements() == 1){
346  G4int Z = (G4int)(material->GetZ());
347  shellStrength = GetOccupationNumber(Z,nbOfTheShell) / Z ;}
348  else G4cout << "WARNING - G4QAOLowEnergyLoss::GetShellEnergy - "
349  << "The model is not available for "
350  << material->GetName()
351  << G4endl;
352 
353  return shellStrength;
354 }
G4double GetZ() const
Definition: G4Material.cc:625
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
Float_t Z
G4double GetOccupationNumber(G4int Z, G4int ShellNb) const
static const G4double siShellStrength[3]
static const G4double cuShellStrength[4]
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
static const G4double ptShellStrength[6]
#define G4endl
Definition: G4ios.hh:61
static const G4double taShellStrength[6]
double G4double
Definition: G4Types.hh:76
const G4String & GetName() const
Definition: G4Material.hh:178
static const G4double auShellStrength[6]
static const G4double alShellStrength[3]
Here is the call graph for this function:
Here is the caller graph for this function:

◆ HighEnergyLimit() [1/2]

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

Implements G4VLowEnergyModel.

Definition at line 84 of file G4QAOLowEnergyLoss.cc.

86 {
87  return 2.0*MeV ;
88 }
static const double MeV
Definition: G4SIunits.hh:211

◆ HighEnergyLimit() [2/2]

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

Implements G4VLowEnergyModel.

Definition at line 99 of file G4QAOLowEnergyLoss.cc.

100 {
101  return 2.0*MeV ;
102 }
static const double MeV
Definition: G4SIunits.hh:211

◆ IsInCharge() [1/2]

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

Implements G4VLowEnergyModel.

Definition at line 111 of file G4QAOLowEnergyLoss.cc.

113 {
114  G4bool isInCharge = false;
115 
116  G4bool hasMaterial = false;
117 
118  if (material->GetNumberOfElements() == 1) hasMaterial = true;
119 
120  if ((particle->GetDefinition()) == (G4AntiProton::AntiProtonDefinition())
121  && hasMaterial) isInCharge = true;
122 
123  return isInCharge;
124 
125 }
static G4AntiProton * AntiProtonDefinition()
Definition: G4AntiProton.cc:88
bool G4bool
Definition: G4Types.hh:79
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
G4ParticleDefinition * GetDefinition() const
Here is the call graph for this function:

◆ IsInCharge() [2/2]

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

Implements G4VLowEnergyModel.

Definition at line 127 of file G4QAOLowEnergyLoss.cc.

129 {
130 
131  G4bool isInCharge = false;
132 
133  G4bool hasMaterial = false;
134 
135  if (material->GetNumberOfElements() == 1) hasMaterial = true;
136 
137 
138  if (aParticle == (G4AntiProton::AntiProtonDefinition())
139  && hasMaterial) isInCharge = true;
140 
141  return isInCharge;
142 
143 }
static G4AntiProton * AntiProtonDefinition()
Definition: G4AntiProton.cc:88
bool G4bool
Definition: G4Types.hh:79
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
Here is the call graph for this function:

◆ LowEnergyLimit() [1/2]

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

Implements G4VLowEnergyModel.

Definition at line 91 of file G4QAOLowEnergyLoss.cc.

93 {
94  // return 50.0*keV ;
95  return 5.0*keV ;
96 }
static const double keV
Definition: G4SIunits.hh:213

◆ LowEnergyLimit() [2/2]

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

Implements G4VLowEnergyModel.

Definition at line 104 of file G4QAOLowEnergyLoss.cc.

105 {
106  // return 50.0*keV ;
107  return 5.0*keV ;
108 }
static const double keV
Definition: G4SIunits.hh:213

◆ TheValue() [1/2]

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

Implements G4VLowEnergyModel.

Definition at line 146 of file G4QAOLowEnergyLoss.cc.

148 {
149  G4double zParticle = (G4int)(particle->GetCharge())/eplus;
150 
151  G4double energy = particle->GetKineticEnergy() ;
152  G4double eloss = EnergyLoss(material,energy,zParticle) ;
153 
154  return eloss ;
155 }
int G4int
Definition: G4Types.hh:78
G4double GetKineticEnergy() const
double energy
Definition: plottest35.C:25
G4double GetCharge() const
G4double EnergyLoss(const G4Material *material, G4double kineticEnergy, G4double zParticle) const
double G4double
Definition: G4Types.hh:76
static const double eplus
Definition: G4SIunits.hh:196
Here is the call graph for this function:

◆ TheValue() [2/2]

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

Implements G4VLowEnergyModel.

Definition at line 158 of file G4QAOLowEnergyLoss.cc.

161 {
162  G4double zParticle = (aParticle->GetPDGCharge())/eplus;
163 
164  G4double eloss = EnergyLoss(material,kineticEnergy,zParticle) ;
165 
166  return eloss ;
167 }
G4double EnergyLoss(const G4Material *material, G4double kineticEnergy, G4double zParticle) const
double G4double
Definition: G4Types.hh:76
static const double eplus
Definition: G4SIunits.hh:196
G4double GetPDGCharge() const
Here is the call graph for this function:

Member Data Documentation

◆ alShellEnergy

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

Definition at line 123 of file G4QAOLowEnergyLoss.hh.

◆ alShellStrength

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

Definition at line 124 of file G4QAOLowEnergyLoss.hh.

◆ auShellEnergy

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

Definition at line 131 of file G4QAOLowEnergyLoss.hh.

◆ auShellStrength

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

Definition at line 132 of file G4QAOLowEnergyLoss.hh.

◆ cuShellEnergy

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

Definition at line 127 of file G4QAOLowEnergyLoss.hh.

◆ cuShellStrength

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

Definition at line 128 of file G4QAOLowEnergyLoss.hh.

◆ 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 143 of file G4QAOLowEnergyLoss.hh.

◆ L0

const G4double G4QAOLowEnergyLoss::L0
staticprivate

Definition at line 139 of file G4QAOLowEnergyLoss.hh.

◆ 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 140 of file G4QAOLowEnergyLoss.hh.

◆ 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 141 of file G4QAOLowEnergyLoss.hh.

◆ materialAvailable

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

Definition at line 118 of file G4QAOLowEnergyLoss.hh.

◆ nbOfElectronPerSubShell

const G4int G4QAOLowEnergyLoss::nbOfElectronPerSubShell
staticprivate

Definition at line 142 of file G4QAOLowEnergyLoss.hh.

◆ nbofShellForMaterial

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

Definition at line 122 of file G4QAOLowEnergyLoss.hh.

◆ numberOfMaterials

G4int G4QAOLowEnergyLoss::numberOfMaterials
private

Definition at line 136 of file G4QAOLowEnergyLoss.hh.

◆ ptShellEnergy

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

Definition at line 133 of file G4QAOLowEnergyLoss.hh.

◆ ptShellStrength

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

Definition at line 134 of file G4QAOLowEnergyLoss.hh.

◆ siShellEnergy

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

Definition at line 125 of file G4QAOLowEnergyLoss.hh.

◆ siShellStrength

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

Definition at line 126 of file G4QAOLowEnergyLoss.hh.

◆ sizeL0

G4int G4QAOLowEnergyLoss::sizeL0
private

Definition at line 145 of file G4QAOLowEnergyLoss.hh.

◆ sizeL1

G4int G4QAOLowEnergyLoss::sizeL1
private

Definition at line 146 of file G4QAOLowEnergyLoss.hh.

◆ sizeL2

G4int G4QAOLowEnergyLoss::sizeL2
private

Definition at line 147 of file G4QAOLowEnergyLoss.hh.

◆ taShellEnergy

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

Definition at line 129 of file G4QAOLowEnergyLoss.hh.

◆ taShellStrength

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

Definition at line 130 of file G4QAOLowEnergyLoss.hh.


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