Geant4  10.02.p03
G4EnergyLossForExtrapolator Class Reference

#include <G4EnergyLossForExtrapolator.hh>

Collaboration diagram for G4EnergyLossForExtrapolator:

Public Member Functions

 G4EnergyLossForExtrapolator (G4int verb=1)
 
 ~G4EnergyLossForExtrapolator ()
 
G4double ComputeDEDX (G4double kinEnergy, const G4ParticleDefinition *)
 
G4double ComputeRange (G4double kinEnergy, const G4ParticleDefinition *)
 
G4double ComputeEnergy (G4double range, const G4ParticleDefinition *)
 
G4double EnergyAfterStep (G4double kinEnergy, G4double step, const G4Material *, const G4ParticleDefinition *)
 
G4double EnergyBeforeStep (G4double kinEnergy, G4double step, const G4Material *, const G4ParticleDefinition *)
 
G4double TrueStepLength (G4double kinEnergy, G4double step, const G4Material *, const G4ParticleDefinition *part)
 
G4double EnergyAfterStep (G4double kinEnergy, G4double step, const G4Material *, const G4String &particleName)
 
G4double EnergyBeforeStep (G4double kinEnergy, G4double step, const G4Material *, const G4String &particleName)
 
G4double AverageScatteringAngle (G4double kinEnergy, G4double step, const G4Material *, const G4ParticleDefinition *part)
 
G4double AverageScatteringAngle (G4double kinEnergy, G4double step, const G4Material *, const G4String &particleName)
 
G4double ComputeTrueStep (const G4Material *, const G4ParticleDefinition *part, G4double kinEnergy, G4double stepLength)
 
G4double EnergyDispersion (G4double kinEnergy, G4double step, const G4Material *, const G4ParticleDefinition *)
 
G4double EnergyDispersion (G4double kinEnergy, G4double step, const G4Material *, const G4String &particleName)
 
void SetVerbose (G4int val)
 
void SetMinKinEnergy (G4double)
 
void SetMaxKinEnergy (G4double)
 
void SetMaxEnergyTransfer (G4double)
 

Private Member Functions

void Initialisation ()
 
void BuildTables ()
 
G4bool SetupKinematics (const G4ParticleDefinition *, const G4Material *, G4double kinEnergy)
 
const G4ParticleDefinitionFindParticle (const G4String &name)
 
G4double ComputeValue (G4double x, const G4PhysicsTable *table, size_t idx)
 
const G4PhysicsTableGetPhysicsTable (ExtTableType type) const
 
G4EnergyLossForExtrapolatoroperator= (const G4EnergyLossForExtrapolator &right)
 
 G4EnergyLossForExtrapolator (const G4EnergyLossForExtrapolator &)
 

Private Attributes

const G4ParticleDefinitioncurrentParticle
 
const G4ParticleDefinitionelectron
 
const G4ParticleDefinitionpositron
 
const G4ParticleDefinitionmuonPlus
 
const G4ParticleDefinitionmuonMinus
 
const G4ParticleDefinitionproton
 
G4String currentParticleName
 
size_t idxDedxElectron
 
size_t idxDedxPositron
 
size_t idxDedxMuon
 
size_t idxDedxProton
 
size_t idxRangeElectron
 
size_t idxRangePositron
 
size_t idxRangeMuon
 
size_t idxRangeProton
 
size_t idxInvRangeElectron
 
size_t idxInvRangePositron
 
size_t idxInvRangeMuon
 
size_t idxInvRangeProton
 
size_t idxMscElectron
 
const G4MaterialcurrentMaterial
 
G4int index
 
G4double electronDensity
 
G4double radLength
 
G4double mass
 
G4double charge2
 
G4double kineticEnergy
 
G4double gam
 
G4double bg2
 
G4double beta2
 
G4double tmax
 
G4double linLossLimit
 
G4double emin
 
G4double emax
 
G4double maxEnergyTransfer
 
G4int nbins
 
G4int nmat
 
G4int verbose
 

Static Private Attributes

static G4TablesForExtrapolatortables = 0
 

Detailed Description

Definition at line 67 of file G4EnergyLossForExtrapolator.hh.

Constructor & Destructor Documentation

◆ G4EnergyLossForExtrapolator() [1/2]

G4EnergyLossForExtrapolator::G4EnergyLossForExtrapolator ( G4int  verb = 1)

Definition at line 67 of file G4EnergyLossForExtrapolator.cc.

69 {
70  currentParticle = 0;
71  currentMaterial = 0;
72 
73  linLossLimit = 0.01;
74  emin = 1.*MeV;
75  emax = 10.*TeV;
76  nbins = 70;
77 
78  nmat = index = 0;
79 
81  = kineticEnergy = tmax = 0.0;
82  gam = 1.0;
83 
88 
90 }
static const double MeV
Definition: G4SIunits.hh:211
const G4ParticleDefinition * muonMinus
const G4ParticleDefinition * proton
const G4ParticleDefinition * muonPlus
const G4ParticleDefinition * electron
static const double TeV
Definition: G4SIunits.hh:215
const G4ParticleDefinition * positron
#define DBL_MAX
Definition: templates.hh:83
const G4ParticleDefinition * currentParticle

◆ ~G4EnergyLossForExtrapolator()

G4EnergyLossForExtrapolator::~G4EnergyLossForExtrapolator ( )

Definition at line 346 of file G4EnergyLossForExtrapolator.cc.

347 {
348  G4AutoLock l(&G4EnergyLossForExtrapolatorMutex);
349  delete tables;
350  tables = 0;
351 }
static G4TablesForExtrapolator * tables

◆ G4EnergyLossForExtrapolator() [2/2]

G4EnergyLossForExtrapolator::G4EnergyLossForExtrapolator ( const G4EnergyLossForExtrapolator )
private

Member Function Documentation

◆ AverageScatteringAngle() [1/2]

G4double G4EnergyLossForExtrapolator::AverageScatteringAngle ( G4double  kinEnergy,
G4double  step,
const G4Material mat,
const G4ParticleDefinition part 
)
inline

Definition at line 248 of file G4EnergyLossForExtrapolator.hh.

253 {
254  G4double theta = 0.0;
255  if(SetupKinematics(part, mat, kinEnergy)) {
256  G4double t = stepLength/radLength;
257  G4double y = std::max(0.001, t);
258  theta = 19.23*CLHEP::MeV*std::sqrt(charge2*t)*(1.0 + 0.038*G4Log(y))
259  /(beta2*gam*mass);
260  }
261  return theta;
262 }
G4bool SetupKinematics(const G4ParticleDefinition *, const G4Material *, G4double kinEnergy)
Double_t y
G4double G4Log(G4double x)
Definition: G4Log.hh:230
double G4double
Definition: G4Types.hh:76
static const double MeV
Here is the call graph for this function:
Here is the caller graph for this function:

◆ AverageScatteringAngle() [2/2]

G4double G4EnergyLossForExtrapolator::AverageScatteringAngle ( G4double  kinEnergy,
G4double  step,
const G4Material mat,
const G4String particleName 
)
inline

Definition at line 227 of file G4EnergyLossForExtrapolator.hh.

231 {
232  return AverageScatteringAngle(kinEnergy,step,mat,FindParticle(name));
233 }
G4double AverageScatteringAngle(G4double kinEnergy, G4double step, const G4Material *, const G4ParticleDefinition *part)
G4String name
Definition: TRTMaterials.hh:40
const G4ParticleDefinition * FindParticle(const G4String &name)
Here is the call graph for this function:

◆ BuildTables()

void G4EnergyLossForExtrapolator::BuildTables ( )
private

Definition at line 330 of file G4EnergyLossForExtrapolator.cc.

331 {
332  G4AutoLock l(&G4EnergyLossForExtrapolatorMutex);
333 
334  if(!tables) {
335  if(verbose > 0) {
336  G4cout << "### G4EnergyLossForExtrapolator::BuildTables for "
337  << nmat << " materials Nbins= " << nbins
338  << " Emin(MeV)= " << emin << " Emax(MeV)= " << emax << G4endl;
339  }
341  }
342 }
G4GLOB_DLL std::ostream G4cout
static G4TablesForExtrapolator * tables
#define G4endl
Definition: G4ios.hh:61
Here is the caller graph for this function:

◆ ComputeDEDX()

G4double G4EnergyLossForExtrapolator::ComputeDEDX ( G4double  kinEnergy,
const G4ParticleDefinition part 
)

Definition at line 240 of file G4EnergyLossForExtrapolator.cc.

242 {
243  G4double x = 0.0;
244  if(part == electron) {
246  } else if(part == positron) {
248  } else if(part == muonPlus || part == muonMinus) {
250  } else {
251  G4double e = ekin*proton_mass_c2/mass;
253  }
254  return x;
255 }
const G4ParticleDefinition * muonMinus
const G4PhysicsTable * GetPhysicsTable(ExtTableType type) const
G4double ComputeValue(G4double x, const G4PhysicsTable *table, size_t idx)
float proton_mass_c2
Definition: hepunit.py:275
const G4ParticleDefinition * muonPlus
const G4ParticleDefinition * electron
const G4ParticleDefinition * positron
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ComputeEnergy()

G4double G4EnergyLossForExtrapolator::ComputeEnergy ( G4double  range,
const G4ParticleDefinition part 
)

Definition at line 282 of file G4EnergyLossForExtrapolator.cc.

284 {
285  G4double x = 0.0;
286  if(part == electron) {
289  } else if(part == positron) {
292  } else if(part == muonPlus || part == muonMinus) {
294  } else {
295  G4double massratio = proton_mass_c2/mass;
296  G4double r = range*massratio*charge2;
298  idxInvRangeProton)/massratio;
299  }
300  return x;
301 }
const G4ParticleDefinition * muonMinus
const G4PhysicsTable * GetPhysicsTable(ExtTableType type) const
G4double ComputeValue(G4double x, const G4PhysicsTable *table, size_t idx)
float proton_mass_c2
Definition: hepunit.py:275
const G4ParticleDefinition * muonPlus
const G4ParticleDefinition * electron
const G4ParticleDefinition * positron
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ComputeRange()

G4double G4EnergyLossForExtrapolator::ComputeRange ( G4double  kinEnergy,
const G4ParticleDefinition part 
)

Definition at line 260 of file G4EnergyLossForExtrapolator.cc.

262 {
263  G4double x = 0.0;
264  if(part == electron) {
266  } else if(part == positron) {
268  } else if(part == muonPlus || part == muonMinus) {
270  } else {
271  G4double massratio = proton_mass_c2/mass;
272  G4double e = ekin*massratio;
274  /(charge2*massratio);
275  }
276  return x;
277 }
const G4ParticleDefinition * muonMinus
const G4PhysicsTable * GetPhysicsTable(ExtTableType type) const
G4double ComputeValue(G4double x, const G4PhysicsTable *table, size_t idx)
float proton_mass_c2
Definition: hepunit.py:275
const G4ParticleDefinition * muonPlus
const G4ParticleDefinition * electron
const G4ParticleDefinition * positron
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ComputeTrueStep()

G4double G4EnergyLossForExtrapolator::ComputeTrueStep ( const G4Material mat,
const G4ParticleDefinition part,
G4double  kinEnergy,
G4double  stepLength 
)
inline

Definition at line 267 of file G4EnergyLossForExtrapolator.hh.

271 {
272  G4double theta = AverageScatteringAngle(kinEnergy,stepLength,mat,part);
273  return stepLength*std::sqrt(1.0 + 0.625*theta*theta);
274 }
G4double AverageScatteringAngle(G4double kinEnergy, G4double step, const G4Material *, const G4ParticleDefinition *part)
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ComputeValue()

G4double G4EnergyLossForExtrapolator::ComputeValue ( G4double  x,
const G4PhysicsTable table,
size_t  idx 
)
inlineprivate

Definition at line 296 of file G4EnergyLossForExtrapolator.hh.

299 {
300  G4double res = 0.0;
301  if(table) { res = ((*table)[index])->Value(x, idx); }
302  return res;
303 }
double G4double
Definition: G4Types.hh:76
Here is the caller graph for this function:

◆ EnergyAfterStep() [1/2]

G4double G4EnergyLossForExtrapolator::EnergyAfterStep ( G4double  kinEnergy,
G4double  step,
const G4Material mat,
const G4ParticleDefinition part 
)

Definition at line 95 of file G4EnergyLossForExtrapolator.cc.

99 {
100  if(0 == nmat) { Initialisation(); }
101  G4double kinEnergyFinal = kinEnergy;
102  if(SetupKinematics(part, mat, kinEnergy)) {
103  G4double step = TrueStepLength(kinEnergy,stepLength,mat,part);
104  G4double r = ComputeRange(kinEnergy,part);
105  if(r <= step) {
106  kinEnergyFinal = 0.0;
107  } else if(step < linLossLimit*r) {
108  kinEnergyFinal -= step*ComputeDEDX(kinEnergy,part);
109  } else {
110  G4double r1 = r - step;
111  kinEnergyFinal = ComputeEnergy(r1,part);
112  }
113  }
114  return kinEnergyFinal;
115 }
G4bool SetupKinematics(const G4ParticleDefinition *, const G4Material *, G4double kinEnergy)
G4double TrueStepLength(G4double kinEnergy, G4double step, const G4Material *, const G4ParticleDefinition *part)
G4double ComputeDEDX(G4double kinEnergy, const G4ParticleDefinition *)
G4double ComputeRange(G4double kinEnergy, const G4ParticleDefinition *)
G4double ComputeEnergy(G4double range, const G4ParticleDefinition *)
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ EnergyAfterStep() [2/2]

G4double G4EnergyLossForExtrapolator::EnergyAfterStep ( G4double  kinEnergy,
G4double  step,
const G4Material mat,
const G4String particleName 
)
inline

Definition at line 205 of file G4EnergyLossForExtrapolator.hh.

209 {
210  return EnergyAfterStep(kinEnergy,step,mat,FindParticle(name));
211 }
G4String name
Definition: TRTMaterials.hh:40
const G4ParticleDefinition * FindParticle(const G4String &name)
G4double EnergyAfterStep(G4double kinEnergy, G4double step, const G4Material *, const G4ParticleDefinition *)
Here is the call graph for this function:

◆ EnergyBeforeStep() [1/2]

G4double G4EnergyLossForExtrapolator::EnergyBeforeStep ( G4double  kinEnergy,
G4double  step,
const G4Material mat,
const G4ParticleDefinition part 
)

Definition at line 120 of file G4EnergyLossForExtrapolator.cc.

124 {
125  //G4cout << "G4EnergyLossForExtrapolator::EnergyBeforeStep" << G4endl;
126  if(0 == nmat) { Initialisation(); }
127  G4double kinEnergyFinal = kinEnergy;
128 
129  if(SetupKinematics(part, mat, kinEnergy)) {
130  G4double step = TrueStepLength(kinEnergy,stepLength,mat,part);
131  G4double r = ComputeRange(kinEnergy,part);
132 
133  if(step < linLossLimit*r) {
134  kinEnergyFinal += step*ComputeDEDX(kinEnergy,part);
135  } else {
136  G4double r1 = r + step;
137  kinEnergyFinal = ComputeEnergy(r1,part);
138  }
139  }
140  return kinEnergyFinal;
141 }
G4bool SetupKinematics(const G4ParticleDefinition *, const G4Material *, G4double kinEnergy)
G4double TrueStepLength(G4double kinEnergy, G4double step, const G4Material *, const G4ParticleDefinition *part)
G4double ComputeDEDX(G4double kinEnergy, const G4ParticleDefinition *)
G4double ComputeRange(G4double kinEnergy, const G4ParticleDefinition *)
G4double ComputeEnergy(G4double range, const G4ParticleDefinition *)
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ EnergyBeforeStep() [2/2]

G4double G4EnergyLossForExtrapolator::EnergyBeforeStep ( G4double  kinEnergy,
G4double  step,
const G4Material mat,
const G4String particleName 
)
inline

Definition at line 216 of file G4EnergyLossForExtrapolator.hh.

220 {
221  return EnergyBeforeStep(kinEnergy,step,mat,FindParticle(name));
222 }
G4String name
Definition: TRTMaterials.hh:40
const G4ParticleDefinition * FindParticle(const G4String &name)
G4double EnergyBeforeStep(G4double kinEnergy, G4double step, const G4Material *, const G4ParticleDefinition *)
Here is the call graph for this function:

◆ EnergyDispersion() [1/2]

G4double G4EnergyLossForExtrapolator::EnergyDispersion ( G4double  kinEnergy,
G4double  step,
const G4Material mat,
const G4ParticleDefinition part 
)
inline

Definition at line 279 of file G4EnergyLossForExtrapolator.hh.

283 {
284  G4double sig2 = 0.0;
285  if(SetupKinematics(part, mat, kinEnergy)) {
286  G4double step = ComputeTrueStep(mat,part,kinEnergy,stepLength);
287  sig2 = (1.0/beta2 - 0.5)
289  }
290  return sig2;
291 }
G4bool SetupKinematics(const G4ParticleDefinition *, const G4Material *, G4double kinEnergy)
G4double ComputeTrueStep(const G4Material *, const G4ParticleDefinition *part, G4double kinEnergy, G4double stepLength)
static const double twopi_mc2_rcl2
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ EnergyDispersion() [2/2]

G4double G4EnergyLossForExtrapolator::EnergyDispersion ( G4double  kinEnergy,
G4double  step,
const G4Material mat,
const G4String particleName 
)
inline

Definition at line 238 of file G4EnergyLossForExtrapolator.hh.

242 {
243  return EnergyDispersion(kinEnergy,step,mat,FindParticle(name));
244 }
G4double EnergyDispersion(G4double kinEnergy, G4double step, const G4Material *, const G4ParticleDefinition *)
G4String name
Definition: TRTMaterials.hh:40
const G4ParticleDefinition * FindParticle(const G4String &name)
Here is the call graph for this function:

◆ FindParticle()

const G4ParticleDefinition * G4EnergyLossForExtrapolator::FindParticle ( const G4String name)
private

Definition at line 221 of file G4EnergyLossForExtrapolator.cc.

222 {
223  const G4ParticleDefinition* p = 0;
224  if(name != currentParticleName) {
226  if(!p) {
227  G4cout << "### G4EnergyLossForExtrapolator WARNING: "
228  << "FindParticle() fails to find "
229  << name << G4endl;
230  }
231  } else {
232  p = currentParticle;
233  }
234  return p;
235 }
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
G4GLOB_DLL std::ostream G4cout
static G4ParticleTable * GetParticleTable()
#define G4endl
Definition: G4ios.hh:61
const G4ParticleDefinition * currentParticle
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetPhysicsTable()

const G4PhysicsTable * G4EnergyLossForExtrapolator::GetPhysicsTable ( ExtTableType  type) const
inlineprivate

Definition at line 197 of file G4EnergyLossForExtrapolator.hh.

198 {
199  return tables->GetPhysicsTable(type);
200 }
static G4TablesForExtrapolator * tables
const G4PhysicsTable * GetPhysicsTable(ExtTableType type) const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ Initialisation()

void G4EnergyLossForExtrapolator::Initialisation ( )
private

Definition at line 305 of file G4EnergyLossForExtrapolator.cc.

306 {
307  if(verbose>1) {
308  G4cout << "### G4EnergyLossForExtrapolator::Initialisation" << G4endl;
309  }
310  currentParticle = 0;
311  currentMaterial = 0;
312  kineticEnergy = 0.0;
313 
319 
321  currentParticleName = "";
322  BuildTables();
323 }
static G4MuonPlus * MuonPlus()
Definition: G4MuonPlus.cc:99
const G4ParticleDefinition * muonMinus
G4GLOB_DLL std::ostream G4cout
static G4Proton * Proton()
Definition: G4Proton.cc:93
const G4ParticleDefinition * proton
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:596
static G4Positron * Positron()
Definition: G4Positron.cc:94
const G4ParticleDefinition * muonPlus
const G4ParticleDefinition * electron
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:100
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
const G4ParticleDefinition * positron
const G4ParticleDefinition * currentParticle
Here is the call graph for this function:
Here is the caller graph for this function:

◆ operator=()

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

◆ SetMaxEnergyTransfer()

void G4EnergyLossForExtrapolator::SetMaxEnergyTransfer ( G4double  val)
inline

Definition at line 328 of file G4EnergyLossForExtrapolator.hh.

329 {
330  maxEnergyTransfer = val;
331 }

◆ SetMaxKinEnergy()

void G4EnergyLossForExtrapolator::SetMaxKinEnergy ( G4double  val)
inline

Definition at line 321 of file G4EnergyLossForExtrapolator.hh.

322 {
323  emax = val;
324 }

◆ SetMinKinEnergy()

void G4EnergyLossForExtrapolator::SetMinKinEnergy ( G4double  val)
inline

Definition at line 314 of file G4EnergyLossForExtrapolator.hh.

315 {
316  emin = val;
317 }

◆ SetupKinematics()

G4bool G4EnergyLossForExtrapolator::SetupKinematics ( const G4ParticleDefinition part,
const G4Material mat,
G4double  kinEnergy 
)
private

Definition at line 173 of file G4EnergyLossForExtrapolator.cc.

176 {
177  if(0 == nmat) { Initialisation(); }
178  if(!part || !mat || kinEnergy < keV) { return false; }
179  G4bool flag = false;
180  if(part != currentParticle) {
181  flag = true;
183  mass = part->GetPDGMass();
184  G4double q = part->GetPDGCharge()/eplus;
185  charge2 = q*q;
186  }
187  if(mat != currentMaterial) {
188  G4int i = mat->GetIndex();
189  if(i >= nmat) {
190  G4cout << "### G4EnergyLossForExtrapolator WARNING:index i= "
191  << i << " is out of table - NO extrapolation" << G4endl;
192  } else {
193  flag = true;
196  radLength = mat->GetRadlen();
197  index = i;
198  }
199  }
200  if(flag || kinEnergy != kineticEnergy) {
201  kineticEnergy = kinEnergy;
202  G4double tau = kinEnergy/mass;
203 
204  gam = tau + 1.0;
205  bg2 = tau * (tau + 2.0);
206  beta2 = bg2/(gam*gam);
207  tmax = kinEnergy;
208  if(part == electron) tmax *= 0.5;
209  else if(part != positron) {
211  tmax = 2.0*bg2*electron_mass_c2/(1.0 + 2.0*gam*r + r*r);
212  }
214  }
215  return true;
216 }
size_t GetIndex() const
Definition: G4Material.hh:262
int G4int
Definition: G4Types.hh:78
Float_t mat
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
TString part[npart]
float electron_mass_c2
Definition: hepunit.py:274
G4double GetElectronDensity() const
Definition: G4Material.hh:217
const G4ParticleDefinition * electron
#define G4endl
Definition: G4ios.hh:61
static const double keV
Definition: G4SIunits.hh:213
const G4ParticleDefinition * positron
double G4double
Definition: G4Types.hh:76
static const double eplus
Definition: G4SIunits.hh:196
G4double GetRadlen() const
Definition: G4Material.hh:220
G4double GetPDGCharge() const
const G4ParticleDefinition * currentParticle
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetVerbose()

void G4EnergyLossForExtrapolator::SetVerbose ( G4int  val)
inline

Definition at line 307 of file G4EnergyLossForExtrapolator.hh.

308 {
309  verbose = val;
310 }

◆ TrueStepLength()

G4double G4EnergyLossForExtrapolator::TrueStepLength ( G4double  kinEnergy,
G4double  step,
const G4Material mat,
const G4ParticleDefinition part 
)

Definition at line 146 of file G4EnergyLossForExtrapolator.cc.

150 {
151  G4double res = stepLength;
152  if(0 == nmat) { Initialisation(); }
153  //G4cout << "## G4EnergyLossForExtrapolator::TrueStepLength L= " << res
154  // << " " << part->GetParticleName() << G4endl;
155  if(SetupKinematics(part, mat, kinEnergy)) {
156  if(part == electron || part == positron) {
157  G4double x = stepLength*
159  //G4cout << " x= " << x << G4endl;
160  if(x < 0.2) { res *= (1.0 + 0.5*x + x*x/3.0); }
161  else if(x < 0.9999) { res = -G4Log(1.0 - x)*stepLength/x; }
162  else { res = ComputeRange(kinEnergy,part); }
163  } else {
164  res = ComputeTrueStep(mat,part,kinEnergy,stepLength);
165  }
166  }
167  return res;
168 }
const G4PhysicsTable * GetPhysicsTable(ExtTableType type) const
G4bool SetupKinematics(const G4ParticleDefinition *, const G4Material *, G4double kinEnergy)
G4double ComputeValue(G4double x, const G4PhysicsTable *table, size_t idx)
G4double ComputeTrueStep(const G4Material *, const G4ParticleDefinition *part, G4double kinEnergy, G4double stepLength)
G4double ComputeRange(G4double kinEnergy, const G4ParticleDefinition *)
G4double G4Log(G4double x)
Definition: G4Log.hh:230
const G4ParticleDefinition * electron
const G4ParticleDefinition * positron
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ beta2

G4double G4EnergyLossForExtrapolator::beta2
private

Definition at line 181 of file G4EnergyLossForExtrapolator.hh.

◆ bg2

G4double G4EnergyLossForExtrapolator::bg2
private

Definition at line 180 of file G4EnergyLossForExtrapolator.hh.

◆ charge2

G4double G4EnergyLossForExtrapolator::charge2
private

Definition at line 177 of file G4EnergyLossForExtrapolator.hh.

◆ currentMaterial

const G4Material* G4EnergyLossForExtrapolator::currentMaterial
private

Definition at line 171 of file G4EnergyLossForExtrapolator.hh.

◆ currentParticle

const G4ParticleDefinition* G4EnergyLossForExtrapolator::currentParticle
private

Definition at line 148 of file G4EnergyLossForExtrapolator.hh.

◆ currentParticleName

G4String G4EnergyLossForExtrapolator::currentParticleName
private

Definition at line 155 of file G4EnergyLossForExtrapolator.hh.

◆ electron

const G4ParticleDefinition* G4EnergyLossForExtrapolator::electron
private

Definition at line 149 of file G4EnergyLossForExtrapolator.hh.

◆ electronDensity

G4double G4EnergyLossForExtrapolator::electronDensity
private

Definition at line 174 of file G4EnergyLossForExtrapolator.hh.

◆ emax

G4double G4EnergyLossForExtrapolator::emax
private

Definition at line 186 of file G4EnergyLossForExtrapolator.hh.

◆ emin

G4double G4EnergyLossForExtrapolator::emin
private

Definition at line 185 of file G4EnergyLossForExtrapolator.hh.

◆ gam

G4double G4EnergyLossForExtrapolator::gam
private

Definition at line 179 of file G4EnergyLossForExtrapolator.hh.

◆ idxDedxElectron

size_t G4EnergyLossForExtrapolator::idxDedxElectron
private

Definition at line 157 of file G4EnergyLossForExtrapolator.hh.

◆ idxDedxMuon

size_t G4EnergyLossForExtrapolator::idxDedxMuon
private

Definition at line 159 of file G4EnergyLossForExtrapolator.hh.

◆ idxDedxPositron

size_t G4EnergyLossForExtrapolator::idxDedxPositron
private

Definition at line 158 of file G4EnergyLossForExtrapolator.hh.

◆ idxDedxProton

size_t G4EnergyLossForExtrapolator::idxDedxProton
private

Definition at line 160 of file G4EnergyLossForExtrapolator.hh.

◆ idxInvRangeElectron

size_t G4EnergyLossForExtrapolator::idxInvRangeElectron
private

Definition at line 165 of file G4EnergyLossForExtrapolator.hh.

◆ idxInvRangeMuon

size_t G4EnergyLossForExtrapolator::idxInvRangeMuon
private

Definition at line 167 of file G4EnergyLossForExtrapolator.hh.

◆ idxInvRangePositron

size_t G4EnergyLossForExtrapolator::idxInvRangePositron
private

Definition at line 166 of file G4EnergyLossForExtrapolator.hh.

◆ idxInvRangeProton

size_t G4EnergyLossForExtrapolator::idxInvRangeProton
private

Definition at line 168 of file G4EnergyLossForExtrapolator.hh.

◆ idxMscElectron

size_t G4EnergyLossForExtrapolator::idxMscElectron
private

Definition at line 169 of file G4EnergyLossForExtrapolator.hh.

◆ idxRangeElectron

size_t G4EnergyLossForExtrapolator::idxRangeElectron
private

Definition at line 161 of file G4EnergyLossForExtrapolator.hh.

◆ idxRangeMuon

size_t G4EnergyLossForExtrapolator::idxRangeMuon
private

Definition at line 163 of file G4EnergyLossForExtrapolator.hh.

◆ idxRangePositron

size_t G4EnergyLossForExtrapolator::idxRangePositron
private

Definition at line 162 of file G4EnergyLossForExtrapolator.hh.

◆ idxRangeProton

size_t G4EnergyLossForExtrapolator::idxRangeProton
private

Definition at line 164 of file G4EnergyLossForExtrapolator.hh.

◆ index

G4int G4EnergyLossForExtrapolator::index
private

Definition at line 172 of file G4EnergyLossForExtrapolator.hh.

◆ kineticEnergy

G4double G4EnergyLossForExtrapolator::kineticEnergy
private

Definition at line 178 of file G4EnergyLossForExtrapolator.hh.

◆ linLossLimit

G4double G4EnergyLossForExtrapolator::linLossLimit
private

Definition at line 184 of file G4EnergyLossForExtrapolator.hh.

◆ mass

G4double G4EnergyLossForExtrapolator::mass
private

Definition at line 176 of file G4EnergyLossForExtrapolator.hh.

◆ maxEnergyTransfer

G4double G4EnergyLossForExtrapolator::maxEnergyTransfer
private

Definition at line 187 of file G4EnergyLossForExtrapolator.hh.

◆ muonMinus

const G4ParticleDefinition* G4EnergyLossForExtrapolator::muonMinus
private

Definition at line 152 of file G4EnergyLossForExtrapolator.hh.

◆ muonPlus

const G4ParticleDefinition* G4EnergyLossForExtrapolator::muonPlus
private

Definition at line 151 of file G4EnergyLossForExtrapolator.hh.

◆ nbins

G4int G4EnergyLossForExtrapolator::nbins
private

Definition at line 189 of file G4EnergyLossForExtrapolator.hh.

◆ nmat

G4int G4EnergyLossForExtrapolator::nmat
private

Definition at line 190 of file G4EnergyLossForExtrapolator.hh.

◆ positron

const G4ParticleDefinition* G4EnergyLossForExtrapolator::positron
private

Definition at line 150 of file G4EnergyLossForExtrapolator.hh.

◆ proton

const G4ParticleDefinition* G4EnergyLossForExtrapolator::proton
private

Definition at line 153 of file G4EnergyLossForExtrapolator.hh.

◆ radLength

G4double G4EnergyLossForExtrapolator::radLength
private

Definition at line 175 of file G4EnergyLossForExtrapolator.hh.

◆ tables

G4TablesForExtrapolator * G4EnergyLossForExtrapolator::tables = 0
staticprivate

Definition at line 146 of file G4EnergyLossForExtrapolator.hh.

◆ tmax

G4double G4EnergyLossForExtrapolator::tmax
private

Definition at line 182 of file G4EnergyLossForExtrapolator.hh.

◆ verbose

G4int G4EnergyLossForExtrapolator::verbose
private

Definition at line 191 of file G4EnergyLossForExtrapolator.hh.


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