Geant4  10.02.p03
G4TablesForExtrapolator Class Reference

#include <G4TablesForExtrapolator.hh>

Collaboration diagram for G4TablesForExtrapolator:

Public Member Functions

 G4TablesForExtrapolator (G4int verb, G4int bins, G4double e1, G4double e2)
 
 ~G4TablesForExtrapolator ()
 
const G4PhysicsTableGetPhysicsTable (ExtTableType type) const
 

Private Member Functions

void Initialisation ()
 
G4PhysicsTablePrepareTable ()
 
void ComputeElectronDEDX (const G4ParticleDefinition *part, G4PhysicsTable *table)
 
void ComputeMuonDEDX (const G4ParticleDefinition *part, G4PhysicsTable *table)
 
void ComputeProtonDEDX (const G4ParticleDefinition *part, G4PhysicsTable *table)
 
void ComputeTrasportXS (const G4ParticleDefinition *part, G4PhysicsTable *table)
 
G4TablesForExtrapolatoroperator= (const G4TablesForExtrapolator &right)
 
 G4TablesForExtrapolator (const G4TablesForExtrapolator &)
 

Private Attributes

const G4ParticleDefinitioncurrentParticle
 
const G4ParticleDefinitionelectron
 
const G4ParticleDefinitionpositron
 
const G4ParticleDefinitionmuonPlus
 
const G4ParticleDefinitionmuonMinus
 
const G4ParticleDefinitionproton
 
G4DataVector cuts
 
G4ProductionCutspcuts
 
std::vector< const G4MaterialCutsCouple * > couples
 
G4PhysicsTablededxElectron
 
G4PhysicsTablededxPositron
 
G4PhysicsTablededxMuon
 
G4PhysicsTablededxProton
 
G4PhysicsTablerangeElectron
 
G4PhysicsTablerangePositron
 
G4PhysicsTablerangeMuon
 
G4PhysicsTablerangeProton
 
G4PhysicsTableinvRangeElectron
 
G4PhysicsTableinvRangePositron
 
G4PhysicsTableinvRangeMuon
 
G4PhysicsTableinvRangeProton
 
G4PhysicsTablemscElectron
 
G4int verbose
 
G4int nbins
 
G4int nmat
 
G4double emin
 
G4double emax
 
G4double mass
 
G4double charge2
 
G4bool splineFlag
 

Detailed Description

Definition at line 75 of file G4TablesForExtrapolator.hh.

Constructor & Destructor Documentation

◆ G4TablesForExtrapolator() [1/2]

G4TablesForExtrapolator::G4TablesForExtrapolator ( G4int  verb,
G4int  bins,
G4double  e1,
G4double  e2 
)

Definition at line 74 of file G4TablesForExtrapolator.cc.

Here is the call graph for this function:

◆ ~G4TablesForExtrapolator()

G4TablesForExtrapolator::~G4TablesForExtrapolator ( )

Definition at line 83 of file G4TablesForExtrapolator.cc.

84 {
85  for(G4int i=0; i<nmat; i++) {delete couples[i];}
86 
100 
101  delete dedxElectron;
102  delete dedxPositron;
103  delete dedxProton;
104  delete dedxMuon;
105  delete rangeElectron;
106  delete rangePositron;
107  delete rangeProton;
108  delete rangeMuon;
109  delete invRangeElectron;
110  delete invRangePositron;
111  delete invRangeProton;
112  delete invRangeMuon;
113  delete mscElectron;
114  delete pcuts;
115 }
std::vector< const G4MaterialCutsCouple * > couples
int G4int
Definition: G4Types.hh:78
void clearAndDestroy()
Here is the call graph for this function:

◆ G4TablesForExtrapolator() [2/2]

G4TablesForExtrapolator::G4TablesForExtrapolator ( const G4TablesForExtrapolator )
private

Member Function Documentation

◆ ComputeElectronDEDX()

void G4TablesForExtrapolator::ComputeElectronDEDX ( const G4ParticleDefinition part,
G4PhysicsTable table 
)
private

Definition at line 271 of file G4TablesForExtrapolator.cc.

274 {
277  ioni->Initialise(part, cuts);
278  brem->Initialise(part, cuts);
279 
281  charge2 = 1.0;
283 
285  if(0<verbose) {
286  G4cout << "G4TablesForExtrapolator::ComputeElectronDEDX for "
287  << part->GetParticleName()
288  << G4endl;
289  }
290  for(G4int i=0; i<nmat; ++i) {
291 
292  const G4Material* mat = (*mtable)[i];
293  if(1<verbose) {
294  G4cout << "i= " << i << " mat= " << mat->GetName() << G4endl;
295  }
296  const G4MaterialCutsCouple* couple = couples[i];
297  G4PhysicsVector* aVector = (*table)[i];
298 
299  for(G4int j=0; j<=nbins; ++j) {
300 
301  G4double e = aVector->Energy(j);
302  G4double dedx = ioni->ComputeDEDX(couple,part,e,e)
303  + brem->ComputeDEDX(couple,part,e,e);
304  if(1<verbose) {
305  G4cout << "j= " << j
306  << " e(MeV)= " << e/MeV
307  << " dedx(Mev/cm)= " << dedx*cm/MeV
308  << " dedx(Mev.cm2/g)= "
309  << dedx/((MeV*mat->GetDensity())/(g/cm2))
310  << G4endl;
311  }
312  aVector->PutValue(j,dedx);
313  }
314  if(splineFlag) { aVector->FillSecondDerivatives(); }
315  }
316  delete ioni;
317  delete brem;
318 }
std::vector< const G4MaterialCutsCouple * > couples
static const double cm
Definition: G4SIunits.hh:118
static const double MeV
Definition: G4SIunits.hh:211
static const double cm2
Definition: G4SIunits.hh:119
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:589
std::vector< G4Material * > G4MaterialTable
G4double GetDensity() const
Definition: G4Material.hh:180
int G4int
Definition: G4Types.hh:78
void FillSecondDerivatives()
Float_t mat
function g(Y1, Y2, PT2)
Definition: hijing1.383.f:5206
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
TString part[npart]
void PutValue(size_t index, G4double theValue)
float electron_mass_c2
Definition: hepunit.py:274
G4double ComputeDEDX(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
Definition: G4VEmModel.hh:490
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
#define G4endl
Definition: G4ios.hh:61
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
G4double Energy(size_t index) const
double G4double
Definition: G4Types.hh:76
const G4String & GetName() const
Definition: G4Material.hh:178
const G4ParticleDefinition * currentParticle
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ComputeMuonDEDX()

void G4TablesForExtrapolator::ComputeMuonDEDX ( const G4ParticleDefinition part,
G4PhysicsTable table 
)
private

Definition at line 323 of file G4TablesForExtrapolator.cc.

325 {
326  G4BetheBlochModel* ioni = new G4BetheBlochModel();
329  ioni->Initialise(part, cuts);
330  pair->Initialise(part, cuts);
331  brem->Initialise(part, cuts);
332 
333  mass = part->GetPDGMass();
334  charge2 = 1.0;
336 
338 
339  if(0<verbose) {
340  G4cout << "G4TablesForExtrapolator::ComputeMuonDEDX for "
341  << part->GetParticleName()
342  << G4endl;
343  }
344 
345  for(G4int i=0; i<nmat; ++i) {
346 
347  const G4Material* mat = (*mtable)[i];
348  if(1<verbose) {
349  G4cout << "i= " << i << " mat= " << mat->GetName() << G4endl;
350  }
351  const G4MaterialCutsCouple* couple = couples[i];
352  G4PhysicsVector* aVector = (*table)[i];
353  for(G4int j=0; j<=nbins; j++) {
354 
355  G4double e = aVector->Energy(j);
356  G4double dedx = ioni->ComputeDEDX(couple,part,e,e) +
357  pair->ComputeDEDX(couple,part,e,e) +
358  brem->ComputeDEDX(couple,part,e,e);
359  aVector->PutValue(j,dedx);
360  if(1<verbose) {
361  G4cout << "j= " << j
362  << " e(MeV)= " << e/MeV
363  << " dedx(Mev/cm)= " << dedx*cm/MeV
364  << " dedx(Mev/(g/cm2)= "
365  << dedx/((MeV*mat->GetDensity())/(g/cm2))
366  << G4endl;
367  }
368  }
369  if(splineFlag) { aVector->FillSecondDerivatives(); }
370  }
371  delete ioni;
372 }
std::vector< const G4MaterialCutsCouple * > couples
static const double cm
Definition: G4SIunits.hh:118
static const double MeV
Definition: G4SIunits.hh:211
static const double cm2
Definition: G4SIunits.hh:119
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:589
std::vector< G4Material * > G4MaterialTable
G4double GetDensity() const
Definition: G4Material.hh:180
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
int G4int
Definition: G4Types.hh:78
void FillSecondDerivatives()
Float_t mat
function g(Y1, Y2, PT2)
Definition: hijing1.383.f:5206
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
TString part[npart]
void PutValue(size_t index, G4double theValue)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
G4double ComputeDEDX(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
Definition: G4VEmModel.hh:490
#define G4endl
Definition: G4ios.hh:61
G4double Energy(size_t index) const
double G4double
Definition: G4Types.hh:76
const G4String & GetName() const
Definition: G4Material.hh:178
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
const G4ParticleDefinition * currentParticle
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ComputeProtonDEDX()

void G4TablesForExtrapolator::ComputeProtonDEDX ( const G4ParticleDefinition part,
G4PhysicsTable table 
)
private

Definition at line 377 of file G4TablesForExtrapolator.cc.

379 {
380  G4BetheBlochModel* ioni = new G4BetheBlochModel();
381  ioni->Initialise(part, cuts);
382 
383  mass = part->GetPDGMass();
384  charge2 = 1.0;
386 
388 
389  if(0<verbose) {
390  G4cout << "G4TablesForExtrapolator::ComputeProtonDEDX for "
391  << part->GetParticleName()
392  << G4endl;
393  }
394 
395  for(G4int i=0; i<nmat; ++i) {
396 
397  const G4Material* mat = (*mtable)[i];
398  if(1<verbose)
399  G4cout << "i= " << i << " mat= " << mat->GetName() << G4endl;
400  const G4MaterialCutsCouple* couple = couples[i];
401  G4PhysicsVector* aVector = (*table)[i];
402  for(G4int j=0; j<=nbins; j++) {
403 
404  G4double e = aVector->Energy(j);
405  G4double dedx = ioni->ComputeDEDX(couple,part,e,e);
406  aVector->PutValue(j,dedx);
407  if(1<verbose) {
408  G4cout << "j= " << j
409  << " e(MeV)= " << e/MeV
410  << " dedx(Mev/cm)= " << dedx*cm/MeV
411  << " dedx(Mev.cm2/g)= " << dedx/((mat->GetDensity())/(g/cm2))
412  << G4endl;
413  }
414  }
415  if(splineFlag) { aVector->FillSecondDerivatives(); }
416  }
417  delete ioni;
418 }
std::vector< const G4MaterialCutsCouple * > couples
static const double cm
Definition: G4SIunits.hh:118
static const double MeV
Definition: G4SIunits.hh:211
static const double cm2
Definition: G4SIunits.hh:119
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:589
std::vector< G4Material * > G4MaterialTable
G4double GetDensity() const
Definition: G4Material.hh:180
int G4int
Definition: G4Types.hh:78
void FillSecondDerivatives()
Float_t mat
function g(Y1, Y2, PT2)
Definition: hijing1.383.f:5206
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
TString part[npart]
void PutValue(size_t index, G4double theValue)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
G4double ComputeDEDX(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
Definition: G4VEmModel.hh:490
#define G4endl
Definition: G4ios.hh:61
G4double Energy(size_t index) const
double G4double
Definition: G4Types.hh:76
const G4String & GetName() const
Definition: G4Material.hh:178
const G4ParticleDefinition * currentParticle
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ComputeTrasportXS()

void G4TablesForExtrapolator::ComputeTrasportXS ( const G4ParticleDefinition part,
G4PhysicsTable table 
)
private

Definition at line 423 of file G4TablesForExtrapolator.cc.

425 {
426  G4WentzelVIModel* msc = new G4WentzelVIModel();
428  msc->Initialise(part, cuts);
429 
430  mass = part->GetPDGMass();
431  charge2 = 1.0;
433 
435 
436  if(0<verbose) {
437  G4cout << "G4TablesForExtrapolator::ComputeProtonDEDX for "
438  << part->GetParticleName()
439  << G4endl;
440  }
441 
442  for(G4int i=0; i<nmat; ++i) {
443 
444  const G4Material* mat = (*mtable)[i];
445  msc->SetCurrentCouple(couples[i]);
446  if(1<verbose)
447  G4cout << "i= " << i << " mat= " << mat->GetName() << G4endl;
448  G4PhysicsVector* aVector = (*table)[i];
449  for(G4int j=0; j<=nbins; j++) {
450 
451  G4double e = aVector->Energy(j);
452  G4double xs = msc->CrossSectionPerVolume(mat,part,e);
453  aVector->PutValue(j,xs);
454  if(1<verbose) {
455  G4cout << "j= " << j << " e(MeV)= " << e/MeV
456  << " xs(1/mm)= " << xs*mm << G4endl;
457  }
458  }
459  if(splineFlag) { aVector->FillSecondDerivatives(); }
460  }
461  delete msc;
462 }
std::vector< const G4MaterialCutsCouple * > couples
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:258
static const double MeV
Definition: G4SIunits.hh:211
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:589
std::vector< G4Material * > G4MaterialTable
int G4int
Definition: G4Types.hh:78
void FillSecondDerivatives()
Float_t mat
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
static const double pi
Definition: SystemOfUnits.h:53
TString part[npart]
void PutValue(size_t index, G4double theValue)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:445
#define G4endl
Definition: G4ios.hh:61
G4double Energy(size_t index) const
double G4double
Definition: G4Types.hh:76
const G4String & GetName() const
Definition: G4Material.hh:178
static const double mm
Definition: G4SIunits.hh:114
void SetPolarAngleLimit(G4double)
Definition: G4VEmModel.hh:760
const G4ParticleDefinition * currentParticle
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetPhysicsTable()

const G4PhysicsTable * G4TablesForExtrapolator::GetPhysicsTable ( ExtTableType  type) const

Definition at line 120 of file G4TablesForExtrapolator.cc.

121 {
122  const G4PhysicsTable* table = 0;
123  switch (type)
124  {
125  case fDedxElectron:
126  table = dedxElectron;
127  break;
128  case fDedxPositron:
129  table = dedxPositron;
130  break;
131  case fDedxProton:
132  table = dedxProton;
133  break;
134  case fDedxMuon:
135  table = dedxMuon;
136  break;
137  case fRangeElectron:
138  table = rangeElectron;
139  break;
140  case fRangePositron:
141  table = rangePositron;
142  break;
143  case fRangeProton:
144  table = rangeProton;
145  break;
146  case fRangeMuon:
147  table = rangeMuon;
148  break;
149  case fInvRangeElectron:
150  table = invRangeElectron;
151  break;
152  case fInvRangePositron:
153  table = invRangePositron;
154  break;
155  case fInvRangeProton:
156  table = invRangeProton;
157  break;
158  case fInvRangeMuon:
159  table = invRangeMuon;
160  break;
161  case fMscElectron:
162  table = mscElectron;
163  }
164  return table;
165 }
Here is the caller graph for this function:

◆ Initialisation()

void G4TablesForExtrapolator::Initialisation ( )
private

Definition at line 169 of file G4TablesForExtrapolator.cc.

170 {
171  if(verbose>1) {
172  G4cout << "### G4TablesForExtrapolator::Initialisation" << G4endl;
173  }
174  currentParticle = 0;
175 
181 
182  mass = charge2 = 0.0;
183 
186  pcuts = new G4ProductionCuts();
187 
188  couples.resize(nmat,0);
189  cuts.resize(nmat,DBL_MAX);
190 
191  for(G4int i=0; i<nmat; ++i) {
192  couples[i] = new G4MaterialCutsCouple((*mtable)[i],pcuts);
193  }
194 
196 
210 
211  G4LossTableBuilder builder;
212 
213  if(verbose>1) {
214  G4cout << "### G4TablesForExtrapolator Builds electron tables"
215  << G4endl;
216  }
220 
221  if(verbose>1) {
222  G4cout << "### G4TablesForExtrapolator Builds positron tables"
223  << G4endl;
224  }
228 
229  if(verbose>1) {
230  G4cout << "### G4TablesForExtrapolator Builds muon tables" << G4endl;
231  }
235  /*
236  G4cout << "DEDX MUON" << G4endl
237  G4cout << *dedxMuon << G4endl;
238  G4cout << "RANGE MUON" << G4endl
239  G4cout << *rangeMuon << G4endl;
240  G4cout << "INVRANGE MUON" << G4endl
241  G4cout << *invRangeMuon << G4endl;
242  */
243  if(verbose>1) {
244  G4cout << "### G4TablesForExtrapolator Builds proton tables"
245  << G4endl;
246  }
250 
252 }
std::vector< const G4MaterialCutsCouple * > couples
void BuildRangeTable(const G4PhysicsTable *dedxTable, G4PhysicsTable *rangeTable, G4bool isIonisation=false)
static G4MuonPlus * MuonPlus()
Definition: G4MuonPlus.cc:99
const G4ParticleDefinition * electron
void ComputeMuonDEDX(const G4ParticleDefinition *part, G4PhysicsTable *table)
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:589
std::vector< G4Material * > G4MaterialTable
const G4ParticleDefinition * positron
int G4int
Definition: G4Types.hh:78
const G4ParticleDefinition * muonPlus
G4bool Spline() const
G4GLOB_DLL std::ostream G4cout
void ComputeProtonDEDX(const G4ParticleDefinition *part, G4PhysicsTable *table)
void ComputeTrasportXS(const G4ParticleDefinition *part, G4PhysicsTable *table)
static G4Proton * Proton()
Definition: G4Proton.cc:93
void BuildInverseRangeTable(const G4PhysicsTable *rangeTable, G4PhysicsTable *invRangeTable, G4bool isIonisation=false)
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:596
static G4Positron * Positron()
Definition: G4Positron.cc:94
static G4EmParameters * Instance()
void ComputeElectronDEDX(const G4ParticleDefinition *part, G4PhysicsTable *table)
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:100
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
#define DBL_MAX
Definition: templates.hh:83
const G4ParticleDefinition * proton
const G4ParticleDefinition * muonMinus
const G4ParticleDefinition * currentParticle
Here is the call graph for this function:
Here is the caller graph for this function:

◆ operator=()

G4TablesForExtrapolator& G4TablesForExtrapolator::operator= ( const G4TablesForExtrapolator right)
private

◆ PrepareTable()

G4PhysicsTable * G4TablesForExtrapolator::PrepareTable ( )
private

Definition at line 256 of file G4TablesForExtrapolator.cc.

257 {
258  G4PhysicsTable* table = new G4PhysicsTable();
259 
260  for(G4int i=0; i<nmat; ++i) {
261 
263  v->SetSpline(splineFlag);
264  table->push_back(v);
265  }
266  return table;
267 }
void push_back(G4PhysicsVector *)
int G4int
Definition: G4Types.hh:78
void SetSpline(G4bool)
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ charge2

G4double G4TablesForExtrapolator::charge2
private

Definition at line 140 of file G4TablesForExtrapolator.hh.

◆ couples

std::vector<const G4MaterialCutsCouple*> G4TablesForExtrapolator::couples
private

Definition at line 117 of file G4TablesForExtrapolator.hh.

◆ currentParticle

const G4ParticleDefinition* G4TablesForExtrapolator::currentParticle
private

Definition at line 107 of file G4TablesForExtrapolator.hh.

◆ cuts

G4DataVector G4TablesForExtrapolator::cuts
private

Definition at line 114 of file G4TablesForExtrapolator.hh.

◆ dedxElectron

G4PhysicsTable* G4TablesForExtrapolator::dedxElectron
private

Definition at line 119 of file G4TablesForExtrapolator.hh.

◆ dedxMuon

G4PhysicsTable* G4TablesForExtrapolator::dedxMuon
private

Definition at line 121 of file G4TablesForExtrapolator.hh.

◆ dedxPositron

G4PhysicsTable* G4TablesForExtrapolator::dedxPositron
private

Definition at line 120 of file G4TablesForExtrapolator.hh.

◆ dedxProton

G4PhysicsTable* G4TablesForExtrapolator::dedxProton
private

Definition at line 122 of file G4TablesForExtrapolator.hh.

◆ electron

const G4ParticleDefinition* G4TablesForExtrapolator::electron
private

Definition at line 108 of file G4TablesForExtrapolator.hh.

◆ emax

G4double G4TablesForExtrapolator::emax
private

Definition at line 138 of file G4TablesForExtrapolator.hh.

◆ emin

G4double G4TablesForExtrapolator::emin
private

Definition at line 137 of file G4TablesForExtrapolator.hh.

◆ invRangeElectron

G4PhysicsTable* G4TablesForExtrapolator::invRangeElectron
private

Definition at line 127 of file G4TablesForExtrapolator.hh.

◆ invRangeMuon

G4PhysicsTable* G4TablesForExtrapolator::invRangeMuon
private

Definition at line 129 of file G4TablesForExtrapolator.hh.

◆ invRangePositron

G4PhysicsTable* G4TablesForExtrapolator::invRangePositron
private

Definition at line 128 of file G4TablesForExtrapolator.hh.

◆ invRangeProton

G4PhysicsTable* G4TablesForExtrapolator::invRangeProton
private

Definition at line 130 of file G4TablesForExtrapolator.hh.

◆ mass

G4double G4TablesForExtrapolator::mass
private

Definition at line 139 of file G4TablesForExtrapolator.hh.

◆ mscElectron

G4PhysicsTable* G4TablesForExtrapolator::mscElectron
private

Definition at line 131 of file G4TablesForExtrapolator.hh.

◆ muonMinus

const G4ParticleDefinition* G4TablesForExtrapolator::muonMinus
private

Definition at line 111 of file G4TablesForExtrapolator.hh.

◆ muonPlus

const G4ParticleDefinition* G4TablesForExtrapolator::muonPlus
private

Definition at line 110 of file G4TablesForExtrapolator.hh.

◆ nbins

G4int G4TablesForExtrapolator::nbins
private

Definition at line 134 of file G4TablesForExtrapolator.hh.

◆ nmat

G4int G4TablesForExtrapolator::nmat
private

Definition at line 135 of file G4TablesForExtrapolator.hh.

◆ pcuts

G4ProductionCuts* G4TablesForExtrapolator::pcuts
private

Definition at line 116 of file G4TablesForExtrapolator.hh.

◆ positron

const G4ParticleDefinition* G4TablesForExtrapolator::positron
private

Definition at line 109 of file G4TablesForExtrapolator.hh.

◆ proton

const G4ParticleDefinition* G4TablesForExtrapolator::proton
private

Definition at line 112 of file G4TablesForExtrapolator.hh.

◆ rangeElectron

G4PhysicsTable* G4TablesForExtrapolator::rangeElectron
private

Definition at line 123 of file G4TablesForExtrapolator.hh.

◆ rangeMuon

G4PhysicsTable* G4TablesForExtrapolator::rangeMuon
private

Definition at line 125 of file G4TablesForExtrapolator.hh.

◆ rangePositron

G4PhysicsTable* G4TablesForExtrapolator::rangePositron
private

Definition at line 124 of file G4TablesForExtrapolator.hh.

◆ rangeProton

G4PhysicsTable* G4TablesForExtrapolator::rangeProton
private

Definition at line 126 of file G4TablesForExtrapolator.hh.

◆ splineFlag

G4bool G4TablesForExtrapolator::splineFlag
private

Definition at line 142 of file G4TablesForExtrapolator.hh.

◆ verbose

G4int G4TablesForExtrapolator::verbose
private

Definition at line 133 of file G4TablesForExtrapolator.hh.


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