Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4hRDEnergyLoss Class Referenceabstract

#include <G4hRDEnergyLoss.hh>

+ Inheritance diagram for G4hRDEnergyLoss:

Public Member Functions

 G4hRDEnergyLoss (const G4String &)
 
 ~G4hRDEnergyLoss ()
 
virtual G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, enum G4ForceCondition *condition)=0
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &Step)=0
 
- Public Member Functions inherited from G4VContinuousDiscreteProcess
 G4VContinuousDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VContinuousDiscreteProcess (G4VContinuousDiscreteProcess &)
 
virtual ~G4VContinuousDiscreteProcess ()
 
G4VContinuousDiscreteProcessoperator= (const G4VContinuousDiscreteProcess &)=delete
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4VProcessoperator= (const G4VProcess &)=delete
 
G4bool operator== (const G4VProcess &right) const
 
G4bool operator!= (const G4VProcess &right) const
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual const G4VProcessGetCreatorProcess () const
 
virtual void StartTracking (G4Track *)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
virtual void ProcessDescription (std::ostream &outfile) const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Static Public Member Functions

static G4int GetNumberOfProcesses ()
 
static void SetNumberOfProcesses (G4int number)
 
static void PlusNumberOfProcesses ()
 
static void MinusNumberOfProcesses ()
 
static void SetdRoverRange (G4double value)
 
static void SetRndmStep (G4bool value)
 
static void SetEnlossFluc (G4bool value)
 
static void SetStepFunction (G4double c1, G4double c2)
 
- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 

Protected Member Functions

G4bool CutsWhereModified ()
 
- Protected Member Functions inherited from G4VContinuousDiscreteProcess
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
virtual G4double GetContinuousStepLimit (const G4Track &aTrack, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)=0
 
void SetGPILSelection (G4GPILSelection selection)
 
G4GPILSelection GetGPILSelection () const
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Static Protected Member Functions

static void BuildDEDXTable (const G4ParticleDefinition &aParticleType)
 

Protected Attributes

const G4double MaxExcitationNumber
 
const G4double probLimFluct
 
const long nmaxDirectFluct
 
const long nmaxCont1
 
const long nmaxCont2
 
G4PhysicsTabletheLossTable
 
G4double linLossLimit
 
G4double MinKineticEnergy
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft = -1.0
 
G4double currentInteractionLength = -1.0
 
G4double theInitialNumberOfInteractionLength = -1.0
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType = fNotDefined
 
G4int theProcessSubType = -1
 
G4double thePILfactor = 1.0
 
G4int verboseLevel = 0
 
G4bool enableAtRestDoIt = true
 
G4bool enableAlongStepDoIt = true
 
G4bool enablePostStepDoIt = true
 

Static Protected Attributes

static G4ThreadLocal G4PhysicsTabletheDEDXpTable = 0
 
static G4ThreadLocal G4PhysicsTabletheDEDXpbarTable = 0
 
static G4ThreadLocal G4PhysicsTabletheRangepTable = 0
 
static G4ThreadLocal G4PhysicsTabletheRangepbarTable = 0
 
static G4ThreadLocal G4PhysicsTabletheInverseRangepTable = 0
 
static G4ThreadLocal G4PhysicsTabletheInverseRangepbarTable = 0
 
static G4ThreadLocal G4PhysicsTabletheLabTimepTable = 0
 
static G4ThreadLocal G4PhysicsTabletheLabTimepbarTable = 0
 
static G4ThreadLocal G4PhysicsTabletheProperTimepTable = 0
 
static G4ThreadLocal G4PhysicsTabletheProperTimepbarTable = 0
 
static G4ThreadLocal G4PhysicsTable ** RecorderOfpProcess = 0
 
static G4ThreadLocal G4PhysicsTable ** RecorderOfpbarProcess = 0
 
static G4ThreadLocal G4int CounterOfpProcess = 0
 
static G4ThreadLocal G4int CounterOfpbarProcess = 0
 
static G4ThreadLocal G4double ParticleMass
 
static G4ThreadLocal G4double ptableElectronCutInRange = 0.0
 
static G4ThreadLocal G4double pbartableElectronCutInRange = 0.0
 
static G4ThreadLocal G4double Charge
 
static G4ThreadLocal G4double LowestKineticEnergy = 1e-05
 
static G4ThreadLocal G4double HighestKineticEnergy = 1.e5
 
static G4ThreadLocal G4int TotBin = 360
 
static G4ThreadLocal G4double RTable
 
static G4ThreadLocal G4double LOGRTable
 
static G4ThreadLocal G4double dRoverRange = 0.20
 
static G4ThreadLocal G4double finalRange = 0.2
 
static G4ThreadLocal G4double c1lim = 0.20
 
static G4ThreadLocal G4double c2lim = 0.32
 
static G4ThreadLocal G4double c3lim = -0.032
 
static G4ThreadLocal G4bool rndmStepFlag = false
 
static G4ThreadLocal G4bool EnlossFlucFlag = true
 

Detailed Description

Definition at line 91 of file G4hRDEnergyLoss.hh.

Constructor & Destructor Documentation

◆ G4hRDEnergyLoss()

G4hRDEnergyLoss::G4hRDEnergyLoss ( const G4String processName)

Definition at line 161 of file G4hRDEnergyLoss.cc.

162 : G4VContinuousDiscreteProcess (processName),
163 MaxExcitationNumber (1.e6),
164 probLimFluct (0.01),
165 nmaxDirectFluct (100),
166 nmaxCont1(4),
167 nmaxCont2(16),
168 theLossTable(0),
169 linLossLimit(0.05),
170 MinKineticEnergy(0.0)
171{
174 if (!RecorderOfProcess) RecorderOfProcess = new G4PhysicsTable*[100];
175}
const G4double MaxExcitationNumber
G4double MinKineticEnergy
static G4ThreadLocal G4PhysicsTable ** RecorderOfpProcess
const long nmaxDirectFluct
G4PhysicsTable * theLossTable
const long nmaxCont1
static G4ThreadLocal G4PhysicsTable ** RecorderOfpbarProcess
const G4double probLimFluct
const long nmaxCont2

◆ ~G4hRDEnergyLoss()

G4hRDEnergyLoss::~G4hRDEnergyLoss ( )

Definition at line 179 of file G4hRDEnergyLoss.cc.

180{
181 if(theLossTable)
182 {
184 delete theLossTable;
185 }
186}
void clearAndDestroy()

Member Function Documentation

◆ BuildDEDXTable()

void G4hRDEnergyLoss::BuildDEDXTable ( const G4ParticleDefinition aParticleType)
staticprotected

Definition at line 250 of file G4hRDEnergyLoss.cc.

251{
252 // calculate data members TotBin,LOGRTable,RTable first
253
256 if (!RecorderOfProcess) RecorderOfProcess = new G4PhysicsTable*[100];
257
258 const G4ProductionCutsTable* theCoupleTable=
260 std::size_t numOfCouples = theCoupleTable->GetTableSize();
261
262 // create/fill proton or antiproton tables depending on the charge
263 Charge = aParticleType.GetPDGCharge()/eplus;
264 ParticleMass = aParticleType.GetPDGMass() ;
265
266 if (Charge>0.) {theDEDXTable= theDEDXpTable;}
267 else {theDEDXTable= theDEDXpbarTable;}
268
269 if( ((Charge>0.) && (theDEDXTable==0)) ||
270 ((Charge<0.) && (theDEDXTable==0))
271 )
272 {
273
274 // Build energy loss table as a sum of the energy loss due to the
275 // different processes.
276 if( Charge >0.)
277 {
278 RecorderOfProcess=RecorderOfpProcess;
279 CounterOfProcess=CounterOfpProcess;
280
281 if(CounterOfProcess == NumberOfProcesses)
282 {
283 if(theDEDXpTable)
285 delete theDEDXpTable; }
286 theDEDXpTable = new G4PhysicsTable(numOfCouples);
287 theDEDXTable = theDEDXpTable;
288 }
289 }
290 else
291 {
292 RecorderOfProcess=RecorderOfpbarProcess;
293 CounterOfProcess=CounterOfpbarProcess;
294
295 if(CounterOfProcess == NumberOfProcesses)
296 {
299 delete theDEDXpbarTable; }
300 theDEDXpbarTable = new G4PhysicsTable(numOfCouples);
301 theDEDXTable = theDEDXpbarTable;
302 }
303 }
304
305 if(CounterOfProcess == NumberOfProcesses)
306 {
307 // loop for materials
308 G4double LowEdgeEnergy , Value ;
309 G4bool isOutRange ;
310 G4PhysicsTable* pointer ;
311
312 for (std::size_t J=0; J<numOfCouples; ++J)
313 {
314 // create physics vector and fill it
315 G4PhysicsLogVector* aVector =
318
319 // loop for the kinetic energy
320 for (G4int i=0; i<TotBin; i++)
321 {
322 LowEdgeEnergy = aVector->GetLowEdgeEnergy(i) ;
323 Value = 0. ;
324
325 // loop for the contributing processes
326 for (G4int process=0; process < NumberOfProcesses; process++)
327 {
328 pointer= RecorderOfProcess[process];
329 Value += (*pointer)[J]->
330 GetValue(LowEdgeEnergy,isOutRange) ;
331 }
332
333 aVector->PutValue(i,Value) ;
334 }
335
336 theDEDXTable->insert(aVector) ;
337 }
338
339 // reset counter to zero ..................
340 if( Charge >0.)
342 else
344
345 // Build range table
346 BuildRangeTable( aParticleType);
347
348 // Build lab/proper time tables
349 BuildTimeTables( aParticleType) ;
350
351 // Build coeff tables for the energy loss calculation
352 BuildRangeCoeffATable( aParticleType);
353 BuildRangeCoeffBTable( aParticleType);
354 BuildRangeCoeffCTable( aParticleType);
355
356 // invert the range table
357
358 BuildInverseRangeTable(aParticleType);
359 }
360 }
361 // make the energy loss and the range table available
362
363 G4EnergyLossTables::Register(&aParticleType,
364 (Charge>0)?
366 (Charge>0)?
368 (Charge>0)?
370 (Charge>0)?
372 (Charge>0)?
375 proton_mass_c2/aParticleType.GetPDGMass(),
376 TotBin);
377
378}
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
static void Register(const G4ParticleDefinition *p, const G4PhysicsTable *tDEDX, const G4PhysicsTable *tRange, const G4PhysicsTable *tInverseRange, const G4PhysicsTable *tLabTime, const G4PhysicsTable *tProperTime, G4double lowestKineticEnergy, G4double highestKineticEnergy, G4double massRatio, G4int NumberOfBins)
G4double GetPDGCharge() const
void insert(G4PhysicsVector *)
G4double GetLowEdgeEnergy(const std::size_t index) const
void PutValue(const std::size_t index, const G4double value)
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
static G4ThreadLocal G4PhysicsTable * theDEDXpTable
static G4ThreadLocal G4PhysicsTable * theInverseRangepTable
static G4ThreadLocal G4double Charge
static G4ThreadLocal G4PhysicsTable * theDEDXpbarTable
static G4ThreadLocal G4int TotBin
static G4ThreadLocal G4PhysicsTable * theProperTimepTable
static G4ThreadLocal G4PhysicsTable * theLabTimepTable
static G4ThreadLocal G4PhysicsTable * theRangepTable
static G4ThreadLocal G4PhysicsTable * theLabTimepbarTable
static G4ThreadLocal G4PhysicsTable * theInverseRangepbarTable
static G4ThreadLocal G4int CounterOfpProcess
static G4ThreadLocal G4double ParticleMass
static G4ThreadLocal G4double HighestKineticEnergy
static G4ThreadLocal G4int CounterOfpbarProcess
static G4ThreadLocal G4PhysicsTable * theProperTimepbarTable
static G4ThreadLocal G4double LowestKineticEnergy
static G4ThreadLocal G4PhysicsTable * theRangepbarTable

Referenced by G4hImpactIonisation::BuildPhysicsTable().

◆ CutsWhereModified()

G4bool G4hRDEnergyLoss::CutsWhereModified ( )
protected

Definition at line 1154 of file G4hRDEnergyLoss.cc.

1155{
1156 G4bool wasModified = false;
1157 const G4ProductionCutsTable* theCoupleTable=
1159 G4int numOfCouples = (G4int)theCoupleTable->GetTableSize();
1160
1161 for (G4int j=0; j<numOfCouples; ++j){
1162 if (theCoupleTable->GetMaterialCutsCouple(j)->IsRecalcNeeded()) {
1163 wasModified = true;
1164 break;
1165 }
1166 }
1167 return wasModified;
1168}
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const

Referenced by G4hImpactIonisation::BuildPhysicsTable().

◆ GetMeanFreePath()

virtual G4double G4hRDEnergyLoss::GetMeanFreePath ( const G4Track track,
G4double  previousStepSize,
enum G4ForceCondition condition 
)
pure virtual

◆ GetNumberOfProcesses()

G4int G4hRDEnergyLoss::GetNumberOfProcesses ( )
static

Definition at line 190 of file G4hRDEnergyLoss.cc.

191{
192 return NumberOfProcesses;
193}

◆ MinusNumberOfProcesses()

void G4hRDEnergyLoss::MinusNumberOfProcesses ( )
static

Definition at line 211 of file G4hRDEnergyLoss.cc.

212{
213 NumberOfProcesses--;
214}

◆ PlusNumberOfProcesses()

void G4hRDEnergyLoss::PlusNumberOfProcesses ( )
static

Definition at line 204 of file G4hRDEnergyLoss.cc.

205{
206 NumberOfProcesses++;
207}

◆ PostStepDoIt()

virtual G4VParticleChange * G4hRDEnergyLoss::PostStepDoIt ( const G4Track track,
const G4Step Step 
)
pure virtual

Reimplemented from G4VContinuousDiscreteProcess.

Implemented in G4hImpactIonisation.

◆ SetdRoverRange()

void G4hRDEnergyLoss::SetdRoverRange ( G4double  value)
static

Definition at line 218 of file G4hRDEnergyLoss.cc.

219{
220 dRoverRange = value;
221}
static G4ThreadLocal G4double dRoverRange

◆ SetEnlossFluc()

void G4hRDEnergyLoss::SetEnlossFluc ( G4bool  value)
static

Definition at line 232 of file G4hRDEnergyLoss.cc.

233{
234 EnlossFlucFlag = value;
235}
static G4ThreadLocal G4bool EnlossFlucFlag

◆ SetNumberOfProcesses()

void G4hRDEnergyLoss::SetNumberOfProcesses ( G4int  number)
static

Definition at line 197 of file G4hRDEnergyLoss.cc.

198{
199 NumberOfProcesses=number;
200}

◆ SetRndmStep()

void G4hRDEnergyLoss::SetRndmStep ( G4bool  value)
static

Definition at line 225 of file G4hRDEnergyLoss.cc.

226{
227 rndmStepFlag = value;
228}
static G4ThreadLocal G4bool rndmStepFlag

◆ SetStepFunction()

void G4hRDEnergyLoss::SetStepFunction ( G4double  c1,
G4double  c2 
)
static

Definition at line 239 of file G4hRDEnergyLoss.cc.

240{
241 dRoverRange = c1;
242 finalRange = c2;
246}
static G4ThreadLocal G4double finalRange
static G4ThreadLocal G4double c3lim
static G4ThreadLocal G4double c1lim
static G4ThreadLocal G4double c2lim

Member Data Documentation

◆ c1lim

G4ThreadLocal G4double G4hRDEnergyLoss::c1lim = 0.20
staticprotected

Definition at line 192 of file G4hRDEnergyLoss.hh.

Referenced by SetStepFunction().

◆ c2lim

G4ThreadLocal G4double G4hRDEnergyLoss::c2lim = 0.32
staticprotected

Definition at line 192 of file G4hRDEnergyLoss.hh.

Referenced by SetStepFunction().

◆ c3lim

G4ThreadLocal G4double G4hRDEnergyLoss::c3lim = -0.032
staticprotected

Definition at line 192 of file G4hRDEnergyLoss.hh.

Referenced by SetStepFunction().

◆ Charge

G4ThreadLocal G4double G4hRDEnergyLoss::Charge
staticprotected

Definition at line 172 of file G4hRDEnergyLoss.hh.

Referenced by BuildDEDXTable().

◆ CounterOfpbarProcess

G4ThreadLocal G4int G4hRDEnergyLoss::CounterOfpbarProcess = 0
staticprotected

Definition at line 163 of file G4hRDEnergyLoss.hh.

Referenced by BuildDEDXTable(), and G4hImpactIonisation::BuildPhysicsTable().

◆ CounterOfpProcess

G4ThreadLocal G4int G4hRDEnergyLoss::CounterOfpProcess = 0
staticprotected

Definition at line 162 of file G4hRDEnergyLoss.hh.

Referenced by BuildDEDXTable(), and G4hImpactIonisation::BuildPhysicsTable().

◆ dRoverRange

G4ThreadLocal G4double G4hRDEnergyLoss::dRoverRange = 0.20
staticprotected

Definition at line 189 of file G4hRDEnergyLoss.hh.

Referenced by SetdRoverRange(), and SetStepFunction().

◆ EnlossFlucFlag

G4ThreadLocal G4bool G4hRDEnergyLoss::EnlossFlucFlag = true
staticprotected

Definition at line 195 of file G4hRDEnergyLoss.hh.

Referenced by G4hImpactIonisation::AlongStepDoIt(), and SetEnlossFluc().

◆ finalRange

G4ThreadLocal G4double G4hRDEnergyLoss::finalRange = 0.2
staticprotected

Definition at line 191 of file G4hRDEnergyLoss.hh.

Referenced by SetStepFunction().

◆ HighestKineticEnergy

◆ linLossLimit

G4double G4hRDEnergyLoss::linLossLimit
protected

Definition at line 185 of file G4hRDEnergyLoss.hh.

Referenced by G4hImpactIonisation::AlongStepDoIt().

◆ LOGRTable

G4ThreadLocal G4double G4hRDEnergyLoss::LOGRTable
staticprotected

Definition at line 179 of file G4hRDEnergyLoss.hh.

◆ LowestKineticEnergy

G4ThreadLocal G4double G4hRDEnergyLoss::LowestKineticEnergy = 1e-05
staticprotected

◆ MaxExcitationNumber

const G4double G4hRDEnergyLoss::MaxExcitationNumber
protected

Definition at line 134 of file G4hRDEnergyLoss.hh.

◆ MinKineticEnergy

G4double G4hRDEnergyLoss::MinKineticEnergy
protected

◆ nmaxCont1

const long G4hRDEnergyLoss::nmaxCont1
protected

Definition at line 136 of file G4hRDEnergyLoss.hh.

◆ nmaxCont2

const long G4hRDEnergyLoss::nmaxCont2
protected

Definition at line 136 of file G4hRDEnergyLoss.hh.

◆ nmaxDirectFluct

const long G4hRDEnergyLoss::nmaxDirectFluct
protected

Definition at line 136 of file G4hRDEnergyLoss.hh.

◆ ParticleMass

G4ThreadLocal G4double G4hRDEnergyLoss::ParticleMass
staticprotected

Definition at line 166 of file G4hRDEnergyLoss.hh.

Referenced by BuildDEDXTable().

◆ pbartableElectronCutInRange

G4ThreadLocal G4double G4hRDEnergyLoss::pbartableElectronCutInRange = 0.0
staticprotected

Definition at line 170 of file G4hRDEnergyLoss.hh.

◆ probLimFluct

const G4double G4hRDEnergyLoss::probLimFluct
protected

Definition at line 135 of file G4hRDEnergyLoss.hh.

◆ ptableElectronCutInRange

G4ThreadLocal G4double G4hRDEnergyLoss::ptableElectronCutInRange = 0.0
staticprotected

Definition at line 169 of file G4hRDEnergyLoss.hh.

◆ RecorderOfpbarProcess

G4ThreadLocal G4PhysicsTable ** G4hRDEnergyLoss::RecorderOfpbarProcess = 0
staticprotected

◆ RecorderOfpProcess

G4ThreadLocal G4PhysicsTable ** G4hRDEnergyLoss::RecorderOfpProcess = 0
staticprotected

◆ rndmStepFlag

G4ThreadLocal G4bool G4hRDEnergyLoss::rndmStepFlag = false
staticprotected

Definition at line 194 of file G4hRDEnergyLoss.hh.

Referenced by SetRndmStep().

◆ RTable

G4ThreadLocal G4double G4hRDEnergyLoss::RTable
staticprotected

Definition at line 179 of file G4hRDEnergyLoss.hh.

◆ theDEDXpbarTable

G4ThreadLocal G4PhysicsTable * G4hRDEnergyLoss::theDEDXpbarTable = 0
staticprotected

Definition at line 143 of file G4hRDEnergyLoss.hh.

Referenced by BuildDEDXTable().

◆ theDEDXpTable

G4ThreadLocal G4PhysicsTable * G4hRDEnergyLoss::theDEDXpTable = 0
staticprotected

Definition at line 142 of file G4hRDEnergyLoss.hh.

Referenced by BuildDEDXTable(), and G4hImpactIonisation::BuildPhysicsTable().

◆ theInverseRangepbarTable

G4ThreadLocal G4PhysicsTable * G4hRDEnergyLoss::theInverseRangepbarTable = 0
staticprotected

Definition at line 149 of file G4hRDEnergyLoss.hh.

Referenced by BuildDEDXTable().

◆ theInverseRangepTable

G4ThreadLocal G4PhysicsTable * G4hRDEnergyLoss::theInverseRangepTable = 0
staticprotected

Definition at line 148 of file G4hRDEnergyLoss.hh.

Referenced by BuildDEDXTable(), and G4hImpactIonisation::BuildPhysicsTable().

◆ theLabTimepbarTable

G4ThreadLocal G4PhysicsTable * G4hRDEnergyLoss::theLabTimepbarTable = 0
staticprotected

Definition at line 153 of file G4hRDEnergyLoss.hh.

Referenced by BuildDEDXTable().

◆ theLabTimepTable

G4ThreadLocal G4PhysicsTable * G4hRDEnergyLoss::theLabTimepTable = 0
staticprotected

Definition at line 152 of file G4hRDEnergyLoss.hh.

Referenced by BuildDEDXTable(), and G4hImpactIonisation::BuildPhysicsTable().

◆ theLossTable

G4PhysicsTable* G4hRDEnergyLoss::theLossTable
protected

Definition at line 183 of file G4hRDEnergyLoss.hh.

Referenced by G4hImpactIonisation::BuildPhysicsTable(), and ~G4hRDEnergyLoss().

◆ theProperTimepbarTable

G4ThreadLocal G4PhysicsTable * G4hRDEnergyLoss::theProperTimepbarTable = 0
staticprotected

Definition at line 156 of file G4hRDEnergyLoss.hh.

Referenced by BuildDEDXTable().

◆ theProperTimepTable

G4ThreadLocal G4PhysicsTable * G4hRDEnergyLoss::theProperTimepTable = 0
staticprotected

Definition at line 155 of file G4hRDEnergyLoss.hh.

Referenced by BuildDEDXTable(), and G4hImpactIonisation::BuildPhysicsTable().

◆ theRangepbarTable

G4ThreadLocal G4PhysicsTable * G4hRDEnergyLoss::theRangepbarTable = 0
staticprotected

Definition at line 145 of file G4hRDEnergyLoss.hh.

Referenced by BuildDEDXTable().

◆ theRangepTable

G4ThreadLocal G4PhysicsTable * G4hRDEnergyLoss::theRangepTable = 0
staticprotected

Definition at line 144 of file G4hRDEnergyLoss.hh.

Referenced by BuildDEDXTable(), and G4hImpactIonisation::BuildPhysicsTable().

◆ TotBin

G4ThreadLocal G4int G4hRDEnergyLoss::TotBin = 360
staticprotected

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