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

#include <G4EnergyLossForExtrapolator.hh>

Public Member Functions

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

Private Member Functions

G4double ComputeValue (G4double x, const G4PhysicsTable *table, size_t idxMat)
 
const G4ParticleDefinitionFindParticle (const G4String &name)
 
const G4PhysicsTableGetPhysicsTable (ExtTableType type) const
 
G4bool SetupKinematics (const G4ParticleDefinition *, const G4Material *, G4double kinEnergy)
 

Private Attributes

G4double beta2 = 0.0
 
G4double bg2 = 0.0
 
G4double charge2 = 0.0
 
const G4MaterialcurrentMaterial = nullptr
 
const G4ParticleDefinitioncurrentParticle = nullptr
 
const G4ParticleDefinitionelectron = nullptr
 
G4double electronDensity = 0.0
 
G4double emax = 0.0
 
G4double emin = 0.0
 
G4double gam = 1.0
 
size_t index = 0
 
G4bool isMaster = false
 
G4double kineticEnergy = 0.0
 
G4double linLossLimit = 0.01
 
G4double maxEnergyTransfer = 0.0
 
const G4ParticleDefinitionmuonMinus = nullptr
 
const G4ParticleDefinitionmuonPlus = nullptr
 
G4int nbins = 80
 
size_t nmat = 0
 
const G4ParticleDefinitionpositron = nullptr
 
const G4ParticleDefinitionproton = nullptr
 
G4double radLength = 0.0
 
G4double tmax = 0.0
 
G4int verbose = 0
 

Static Private Attributes

static G4TablesForExtrapolatortables = nullptr
 

Detailed Description

Definition at line 67 of file G4EnergyLossForExtrapolator.hh.

Constructor & Destructor Documentation

◆ G4EnergyLossForExtrapolator() [1/2]

G4EnergyLossForExtrapolator::G4EnergyLossForExtrapolator ( G4int  verb = 1)
explicit

Definition at line 70 of file G4EnergyLossForExtrapolator.cc.

References emax, emin, CLHEP::MeV, and CLHEP::TeV.

◆ ~G4EnergyLossForExtrapolator()

G4EnergyLossForExtrapolator::~G4EnergyLossForExtrapolator ( )

Definition at line 79 of file G4EnergyLossForExtrapolator.cc.

80{
81 if(isMaster) {
82 delete tables;
83 tables = nullptr;
84 }
85}
static G4TablesForExtrapolator * tables

References isMaster, and tables.

◆ G4EnergyLossForExtrapolator() [2/2]

G4EnergyLossForExtrapolator::G4EnergyLossForExtrapolator ( const G4EnergyLossForExtrapolator )
delete

Member Function Documentation

◆ AverageScatteringAngle() [1/2]

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

Definition at line 312 of file G4EnergyLossForExtrapolator.cc.

317{
318 G4double theta = 0.0;
319 if(SetupKinematics(part, mat, kinEnergy)) {
320 G4double t = stepLength/radLength;
321 G4double y = std::max(0.001, t);
322 theta = 19.23*CLHEP::MeV*std::sqrt(charge2*t)*(1.0 + 0.038*G4Log(y))
323 /(beta2*gam*part->GetPDGMass());
324 }
325 return theta;
326}
G4double G4Log(G4double x)
Definition: G4Log.hh:226
double G4double
Definition: G4Types.hh:83
G4bool SetupKinematics(const G4ParticleDefinition *, const G4Material *, G4double kinEnergy)
T max(const T t1, const T t2)
brief Return the largest of the two arguments

References beta2, charge2, G4Log(), gam, G4ParticleDefinition::GetPDGMass(), G4INCL::Math::max(), CLHEP::MeV, radLength, and SetupKinematics().

Referenced by AverageScatteringAngle(), and ComputeTrueStep().

◆ AverageScatteringAngle() [2/2]

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

Definition at line 216 of file G4EnergyLossForExtrapolator.hh.

220{
221 return AverageScatteringAngle(kinEnergy,step,mat,FindParticle(name));
222}
G4double AverageScatteringAngle(G4double kinEnergy, G4double step, const G4Material *, const G4ParticleDefinition *part)
const G4ParticleDefinition * FindParticle(const G4String &name)
const char * name(G4int ptype)

References AverageScatteringAngle(), FindParticle(), and G4InuclParticleNames::name().

◆ ComputeDEDX()

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

Definition at line 224 of file G4EnergyLossForExtrapolator.cc.

227{
228 if(mat->GetNumberOfMaterials() != nmat) { Initialisation(); }
229 G4double x = 0.0;
230 if(part == electron) {
232 } else if(part == positron) {
234 } else if(part == muonPlus || part == muonMinus) {
236 } else {
240 }
241 return x;
242}
const G4ParticleDefinition * muonPlus
G4double ComputeValue(G4double x, const G4PhysicsTable *table, size_t idxMat)
const G4ParticleDefinition * muonMinus
const G4ParticleDefinition * electron
const G4PhysicsTable * GetPhysicsTable(ExtTableType type) const
const G4ParticleDefinition * positron
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:679
size_t GetIndex() const
Definition: G4Material.hh:256
G4double GetPDGCharge() const
static constexpr double eplus
static constexpr double proton_mass_c2

References ComputeValue(), electron, CLHEP::eplus, fDedxElectron, fDedxMuon, fDedxPositron, fDedxProton, G4Material::GetIndex(), G4Material::GetNumberOfMaterials(), G4ParticleDefinition::GetPDGCharge(), G4ParticleDefinition::GetPDGMass(), GetPhysicsTable(), Initialisation(), muonMinus, muonPlus, nmat, positron, and CLHEP::proton_mass_c2.

Referenced by EnergyAfterStep(), and EnergyBeforeStep().

◆ ComputeEnergy()

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

Definition at line 272 of file G4EnergyLossForExtrapolator.cc.

275{
276 if(mat->GetNumberOfMaterials() != nmat) { Initialisation(); }
277 G4double x = 0.0;
278 if(part == electron) {
280 } else if(part == positron) {
282 } else if(part == muonPlus || part == muonMinus) {
284 } else {
285 G4double massratio = CLHEP::proton_mass_c2/part->GetPDGMass();
287 G4double r = range*massratio*q*q;
288 x = ComputeValue(r, GetPhysicsTable(fInvRangeProton), mat->GetIndex())/massratio;
289 }
290 return x;
291}
@ fInvRangeElectron
@ fInvRangePositron

References ComputeValue(), electron, CLHEP::eplus, fInvRangeElectron, fInvRangeMuon, fInvRangePositron, fInvRangeProton, G4Material::GetIndex(), G4Material::GetNumberOfMaterials(), G4ParticleDefinition::GetPDGCharge(), G4ParticleDefinition::GetPDGMass(), GetPhysicsTable(), Initialisation(), muonMinus, muonPlus, nmat, positron, and CLHEP::proton_mass_c2.

Referenced by EnergyAfterStep(), and EnergyBeforeStep().

◆ ComputeRange()

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

Definition at line 247 of file G4EnergyLossForExtrapolator.cc.

250{
251 if(mat->GetNumberOfMaterials() != nmat) { Initialisation(); }
252 G4double x = 0.0;
253 if(part == electron) {
255 } else if(part == positron) {
257 } else if(part == muonPlus || part == muonMinus) {
259 } else {
260 G4double massratio = CLHEP::proton_mass_c2/part->GetPDGMass();
261 G4double e = ekin*massratio;
264 /(q*q*massratio);
265 }
266 return x;
267}

References ComputeValue(), electron, CLHEP::eplus, fRangeElectron, fRangeMuon, fRangePositron, fRangeProton, G4Material::GetIndex(), G4Material::GetNumberOfMaterials(), G4ParticleDefinition::GetPDGCharge(), G4ParticleDefinition::GetPDGMass(), GetPhysicsTable(), Initialisation(), muonMinus, muonPlus, nmat, positron, and CLHEP::proton_mass_c2.

Referenced by EnergyAfterStep(), EnergyBeforeStep(), G4ErrorEnergyLoss::GetContinuousStepLimit(), and TrueStepLength().

◆ ComputeTrueStep()

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

Definition at line 238 of file G4EnergyLossForExtrapolator.hh.

242{
243 G4double theta = AverageScatteringAngle(kinEnergy,stepLength,mat,part);
244 return stepLength*std::sqrt(1.0 + 0.625*theta*theta);
245}

References AverageScatteringAngle().

Referenced by EnergyDispersion(), and TrueStepLength().

◆ ComputeValue()

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

Definition at line 250 of file G4EnergyLossForExtrapolator.hh.

253{
254 return (nullptr != table) ? ((*table)[idxMat])->Value(x, index) : 0.0;
255}

References index.

Referenced by ComputeDEDX(), ComputeEnergy(), ComputeRange(), and TrueStepLength().

◆ EnergyAfterStep() [1/2]

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

Definition at line 90 of file G4EnergyLossForExtrapolator.cc.

94{
95 G4double kinEnergyFinal = kinEnergy;
96 if(SetupKinematics(part, mat, kinEnergy)) {
97 G4double step = TrueStepLength(kinEnergy,stepLength,mat,part);
98 G4double r = ComputeRange(kinEnergy,part,mat);
99 if(r <= step) {
100 kinEnergyFinal = 0.0;
101 } else if(step < linLossLimit*r) {
102 kinEnergyFinal -= step*ComputeDEDX(kinEnergy,part,mat);
103 } else {
104 G4double r1 = r - step;
105 kinEnergyFinal = ComputeEnergy(r1,part,mat);
106 }
107 }
108 return kinEnergyFinal;
109}
G4double TrueStepLength(G4double kinEnergy, G4double step, const G4Material *, const G4ParticleDefinition *part)
G4double ComputeRange(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *)
G4double ComputeDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *)
G4double ComputeEnergy(G4double range, const G4ParticleDefinition *, const G4Material *)

References ComputeDEDX(), ComputeEnergy(), ComputeRange(), linLossLimit, SetupKinematics(), and TrueStepLength().

Referenced by G4ErrorEnergyLoss::AlongStepDoIt(), and EnergyAfterStep().

◆ EnergyAfterStep() [2/2]

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

Definition at line 194 of file G4EnergyLossForExtrapolator.hh.

198{
199 return EnergyAfterStep(kinEnergy,step,mat,FindParticle(name));
200}
G4double EnergyAfterStep(G4double kinEnergy, G4double step, const G4Material *, const G4ParticleDefinition *)

References EnergyAfterStep(), FindParticle(), and G4InuclParticleNames::name().

◆ EnergyBeforeStep() [1/2]

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

Definition at line 114 of file G4EnergyLossForExtrapolator.cc.

118{
119 //G4cout << "G4EnergyLossForExtrapolator::EnergyBeforeStep" << G4endl;
120 G4double kinEnergyFinal = kinEnergy;
121
122 if(SetupKinematics(part, mat, kinEnergy)) {
123 G4double step = TrueStepLength(kinEnergy,stepLength,mat,part);
124 G4double r = ComputeRange(kinEnergy,part,mat);
125
126 if(step < linLossLimit*r) {
127 kinEnergyFinal += step*ComputeDEDX(kinEnergy,part,mat);
128 } else {
129 G4double r1 = r + step;
130 kinEnergyFinal = ComputeEnergy(r1,part,mat);
131 }
132 }
133 return kinEnergyFinal;
134}

References ComputeDEDX(), ComputeEnergy(), ComputeRange(), linLossLimit, SetupKinematics(), and TrueStepLength().

Referenced by G4ErrorEnergyLoss::AlongStepDoIt(), and EnergyBeforeStep().

◆ EnergyBeforeStep() [2/2]

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

Definition at line 205 of file G4EnergyLossForExtrapolator.hh.

209{
210 return EnergyBeforeStep(kinEnergy,step,mat,FindParticle(name));
211}
G4double EnergyBeforeStep(G4double kinEnergy, G4double step, const G4Material *, const G4ParticleDefinition *)

References EnergyBeforeStep(), FindParticle(), and G4InuclParticleNames::name().

◆ EnergyDispersion() [1/2]

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

Definition at line 296 of file G4EnergyLossForExtrapolator.cc.

300{
301 G4double sig2 = 0.0;
302 if(SetupKinematics(part, mat, kinEnergy)) {
303 G4double step = ComputeTrueStep(mat,part,kinEnergy,stepLength);
304 sig2 = (1.0/beta2 - 0.5)
306 }
307 return sig2;
308}
G4double ComputeTrueStep(const G4Material *, const G4ParticleDefinition *part, G4double kinEnergy, G4double stepLength)
static constexpr double twopi_mc2_rcl2

References beta2, charge2, ComputeTrueStep(), electronDensity, SetupKinematics(), tmax, and CLHEP::twopi_mc2_rcl2.

Referenced by EnergyDispersion().

◆ EnergyDispersion() [2/2]

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

Definition at line 227 of file G4EnergyLossForExtrapolator.hh.

231{
232 return EnergyDispersion(kinEnergy,step,mat,FindParticle(name));
233}
G4double EnergyDispersion(G4double kinEnergy, G4double step, const G4Material *, const G4ParticleDefinition *)

References EnergyDispersion(), FindParticle(), and G4InuclParticleNames::name().

◆ FindParticle()

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

Definition at line 211 of file G4EnergyLossForExtrapolator.cc.

212{
214 if(nullptr == currentParticle) {
215 G4cout << "### G4EnergyLossForExtrapolator WARNING: "
216 << "FindParticle() fails to find " << name << G4endl;
217 }
218 return currentParticle;
219}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
const G4ParticleDefinition * currentParticle
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()

References currentParticle, G4ParticleTable::FindParticle(), G4cout, G4endl, G4ParticleTable::GetParticleTable(), and G4InuclParticleNames::name().

Referenced by AverageScatteringAngle(), EnergyAfterStep(), EnergyBeforeStep(), and EnergyDispersion().

◆ GetPhysicsTable()

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

Definition at line 186 of file G4EnergyLossForExtrapolator.hh.

187{
188 return tables->GetPhysicsTable(type);
189}
const G4PhysicsTable * GetPhysicsTable(ExtTableType type) const

References G4TablesForExtrapolator::GetPhysicsTable(), and tables.

Referenced by ComputeDEDX(), ComputeEnergy(), ComputeRange(), and TrueStepLength().

◆ Initialisation()

void G4EnergyLossForExtrapolator::Initialisation ( )

Definition at line 330 of file G4EnergyLossForExtrapolator.cc.

331{
332 if(verbose>0) {
333 G4cout << "### G4EnergyLossForExtrapolator::Initialisation tables= "
334 << tables << G4endl;
335 }
341
342 // initialisation for the 1st run
343 if(nullptr == tables) {
344#ifdef G4MULTITHREADED
345 G4MUTEXLOCK(&extrMutex);
346 if(nullptr == tables) {
347#endif
348 isMaster = true;
352 if(verbose > 0) {
353 G4cout << "### G4EnergyLossForExtrapolator::BuildTables for "
354 << nmat << " materials Nbins= "
355 << nbins << " Emin(MeV)= " << emin << " Emax(MeV)= " << emax
356 << G4endl;
357 }
358#ifdef G4MULTITHREADED
359 }
360 G4MUTEXUNLOCK(&extrMutex);
361#endif
362 }
363
364 // initialisation for the next run
366#ifdef G4MULTITHREADED
367 G4MUTEXLOCK(&extrMutex);
368#endif
370#ifdef G4MULTITHREADED
371 G4MUTEXUNLOCK(&extrMutex);
372#endif
373 }
375}
#define G4MUTEXLOCK(mutex)
Definition: G4Threading.hh:251
#define G4MUTEXUNLOCK(mutex)
Definition: G4Threading.hh:254
static G4Electron * Electron()
Definition: G4Electron.cc:93
const G4ParticleDefinition * proton
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:99
static G4MuonPlus * MuonPlus()
Definition: G4MuonPlus.cc:98
static G4Positron * Positron()
Definition: G4Positron.cc:93
static G4Proton * Proton()
Definition: G4Proton.cc:92

References G4Electron::Electron(), electron, emax, emin, G4cout, G4endl, G4MUTEXLOCK, G4MUTEXUNLOCK, G4Material::GetNumberOfMaterials(), G4TablesForExtrapolator::Initialisation(), isMaster, G4MuonMinus::MuonMinus(), muonMinus, G4MuonPlus::MuonPlus(), muonPlus, nbins, nmat, G4Positron::Positron(), positron, G4Proton::Proton(), proton, tables, and verbose.

Referenced by ComputeDEDX(), ComputeEnergy(), ComputeRange(), and SetupKinematics().

◆ operator=()

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

◆ SetMaxEnergyTransfer()

void G4EnergyLossForExtrapolator::SetMaxEnergyTransfer ( G4double  val)
inline

Definition at line 280 of file G4EnergyLossForExtrapolator.hh.

281{
282 maxEnergyTransfer = val;
283}

References maxEnergyTransfer.

◆ SetMaxKinEnergy()

void G4EnergyLossForExtrapolator::SetMaxKinEnergy ( G4double  val)
inline

Definition at line 273 of file G4EnergyLossForExtrapolator.hh.

274{
275 emax = val;
276}

References emax.

◆ SetMinKinEnergy()

void G4EnergyLossForExtrapolator::SetMinKinEnergy ( G4double  val)
inline

Definition at line 266 of file G4EnergyLossForExtrapolator.hh.

267{
268 emin = val;
269}

References emin.

◆ SetupKinematics()

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

Definition at line 165 of file G4EnergyLossForExtrapolator.cc.

168{
169 if(mat->GetNumberOfMaterials() != nmat) { Initialisation(); }
170 if(nullptr == part || nullptr == mat || kinEnergy < CLHEP::keV)
171 { return false; }
172 if(part != currentParticle) {
173 currentParticle = part;
174 G4double q = part->GetPDGCharge()/eplus;
175 charge2 = q*q;
176 }
177 if(mat != currentMaterial) {
178 size_t i = mat->GetIndex();
179 if(i >= nmat) {
180 G4cout << "### G4EnergyLossForExtrapolator WARNING: material index i= "
181 << i << " above number of materials " << nmat << G4endl;
182 return false;
183 } else {
184 currentMaterial = mat;
186 radLength = mat->GetRadlen();
187 }
188 }
189 if(kinEnergy != kineticEnergy) {
190 kineticEnergy = kinEnergy;
191 G4double mass = part->GetPDGMass();
192 G4double tau = kinEnergy/mass;
193
194 gam = tau + 1.0;
195 bg2 = tau * (tau + 2.0);
196 beta2 = bg2/(gam*gam);
197 tmax = kinEnergy;
198 if(part == electron) tmax *= 0.5;
199 else if(part != positron) {
201 tmax = 2.0*bg2*CLHEP::electron_mass_c2/(1.0 + 2.0*gam*r + r*r);
202 }
204 }
205 return true;
206}
static constexpr double eplus
Definition: G4SIunits.hh:184
G4double GetElectronDensity() const
Definition: G4Material.hh:213
G4double GetRadlen() const
Definition: G4Material.hh:216
static constexpr double electron_mass_c2
static constexpr double keV
T min(const T t1, const T t2)
brief Return the smallest of the two arguments

References beta2, bg2, charge2, currentMaterial, currentParticle, electron, CLHEP::electron_mass_c2, electronDensity, eplus, G4cout, G4endl, gam, G4Material::GetElectronDensity(), G4Material::GetIndex(), G4Material::GetNumberOfMaterials(), G4ParticleDefinition::GetPDGCharge(), G4ParticleDefinition::GetPDGMass(), G4Material::GetRadlen(), Initialisation(), CLHEP::keV, kineticEnergy, maxEnergyTransfer, G4INCL::Math::min(), nmat, positron, radLength, and tmax.

Referenced by AverageScatteringAngle(), EnergyAfterStep(), EnergyBeforeStep(), EnergyDispersion(), and TrueStepLength().

◆ SetVerbose()

void G4EnergyLossForExtrapolator::SetVerbose ( G4int  val)
inline

Definition at line 259 of file G4EnergyLossForExtrapolator.hh.

260{
261 verbose = val;
262}

References verbose.

◆ TrueStepLength()

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

Definition at line 139 of file G4EnergyLossForExtrapolator.cc.

143{
144 G4double res = stepLength;
145 //G4cout << "## G4EnergyLossForExtrapolator::TrueStepLength L= " << res
146 // << " " << part->GetParticleName() << G4endl;
147 if(SetupKinematics(part, mat, kinEnergy)) {
148 if(part == electron || part == positron) {
149 const G4double x = stepLength*
151 //G4cout << " x= " << x << G4endl;
152 if(x < 0.2) { res *= (1.0 + 0.5*x + x*x/3.0); }
153 else if(x < 0.9999) { res = -G4Log(1.0 - x)*stepLength/x; }
154 else { res = ComputeRange(kinEnergy, part, mat); }
155 } else {
156 res = ComputeTrueStep(mat,part,kinEnergy,stepLength);
157 }
158 }
159 return res;
160}

References ComputeRange(), ComputeTrueStep(), ComputeValue(), electron, fMscElectron, G4Log(), G4Material::GetIndex(), GetPhysicsTable(), positron, and SetupKinematics().

Referenced by EnergyAfterStep(), EnergyBeforeStep(), and G4EnergySplitter::SplitEnergyInVolumes().

Field Documentation

◆ beta2

G4double G4EnergyLossForExtrapolator::beta2 = 0.0
private

◆ bg2

G4double G4EnergyLossForExtrapolator::bg2 = 0.0
private

Definition at line 166 of file G4EnergyLossForExtrapolator.hh.

Referenced by SetupKinematics().

◆ charge2

G4double G4EnergyLossForExtrapolator::charge2 = 0.0
private

◆ currentMaterial

const G4Material* G4EnergyLossForExtrapolator::currentMaterial = nullptr
private

Definition at line 159 of file G4EnergyLossForExtrapolator.hh.

Referenced by SetupKinematics().

◆ currentParticle

const G4ParticleDefinition* G4EnergyLossForExtrapolator::currentParticle = nullptr
private

Definition at line 153 of file G4EnergyLossForExtrapolator.hh.

Referenced by FindParticle(), and SetupKinematics().

◆ electron

const G4ParticleDefinition* G4EnergyLossForExtrapolator::electron = nullptr
private

◆ electronDensity

G4double G4EnergyLossForExtrapolator::electronDensity = 0.0
private

Definition at line 161 of file G4EnergyLossForExtrapolator.hh.

Referenced by EnergyDispersion(), and SetupKinematics().

◆ emax

G4double G4EnergyLossForExtrapolator::emax = 0.0
private

◆ emin

G4double G4EnergyLossForExtrapolator::emin = 0.0
private

◆ gam

G4double G4EnergyLossForExtrapolator::gam = 1.0
private

Definition at line 165 of file G4EnergyLossForExtrapolator.hh.

Referenced by AverageScatteringAngle(), and SetupKinematics().

◆ index

size_t G4EnergyLossForExtrapolator::index = 0
private

Definition at line 175 of file G4EnergyLossForExtrapolator.hh.

Referenced by ComputeValue().

◆ isMaster

G4bool G4EnergyLossForExtrapolator::isMaster = false
private

Definition at line 180 of file G4EnergyLossForExtrapolator.hh.

Referenced by Initialisation(), and ~G4EnergyLossForExtrapolator().

◆ kineticEnergy

G4double G4EnergyLossForExtrapolator::kineticEnergy = 0.0
private

Definition at line 164 of file G4EnergyLossForExtrapolator.hh.

Referenced by SetupKinematics().

◆ linLossLimit

G4double G4EnergyLossForExtrapolator::linLossLimit = 0.01
private

Definition at line 170 of file G4EnergyLossForExtrapolator.hh.

Referenced by EnergyAfterStep(), and EnergyBeforeStep().

◆ maxEnergyTransfer

G4double G4EnergyLossForExtrapolator::maxEnergyTransfer = 0.0
private

Definition at line 173 of file G4EnergyLossForExtrapolator.hh.

Referenced by SetMaxEnergyTransfer(), and SetupKinematics().

◆ muonMinus

const G4ParticleDefinition* G4EnergyLossForExtrapolator::muonMinus = nullptr
private

◆ muonPlus

const G4ParticleDefinition* G4EnergyLossForExtrapolator::muonPlus = nullptr
private

◆ nbins

G4int G4EnergyLossForExtrapolator::nbins = 80
private

Definition at line 177 of file G4EnergyLossForExtrapolator.hh.

Referenced by Initialisation().

◆ nmat

size_t G4EnergyLossForExtrapolator::nmat = 0
private

◆ positron

const G4ParticleDefinition* G4EnergyLossForExtrapolator::positron = nullptr
private

◆ proton

const G4ParticleDefinition* G4EnergyLossForExtrapolator::proton = nullptr
private

Definition at line 158 of file G4EnergyLossForExtrapolator.hh.

Referenced by Initialisation().

◆ radLength

G4double G4EnergyLossForExtrapolator::radLength = 0.0
private

Definition at line 162 of file G4EnergyLossForExtrapolator.hh.

Referenced by AverageScatteringAngle(), and SetupKinematics().

◆ tables

G4TablesForExtrapolator * G4EnergyLossForExtrapolator::tables = nullptr
staticprivate

◆ tmax

G4double G4EnergyLossForExtrapolator::tmax = 0.0
private

Definition at line 168 of file G4EnergyLossForExtrapolator.hh.

Referenced by EnergyDispersion(), and SetupKinematics().

◆ verbose

G4int G4EnergyLossForExtrapolator::verbose = 0
private

Definition at line 178 of file G4EnergyLossForExtrapolator.hh.

Referenced by Initialisation(), and SetVerbose().


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