Geant4 9.6.0
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 ()
 
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 ()
 
G4int operator== (const G4VProcess &right) const
 
G4int 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 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
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 

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 previousStepSize)
 
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
 
G4VParticleChangepParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 

Static Protected Attributes

static G4PhysicsTabletheDEDXpTable = 0
 
static G4PhysicsTabletheDEDXpbarTable = 0
 
static G4PhysicsTabletheRangepTable = 0
 
static G4PhysicsTabletheRangepbarTable = 0
 
static G4PhysicsTabletheInverseRangepTable = 0
 
static G4PhysicsTabletheInverseRangepbarTable = 0
 
static G4PhysicsTabletheLabTimepTable = 0
 
static G4PhysicsTabletheLabTimepbarTable = 0
 
static G4PhysicsTabletheProperTimepTable = 0
 
static G4PhysicsTabletheProperTimepbarTable = 0
 
static G4PhysicsTable ** RecorderOfpProcess
 
static G4PhysicsTable ** RecorderOfpbarProcess
 
static G4int CounterOfpProcess = 0
 
static G4int CounterOfpbarProcess = 0
 
static G4double ParticleMass
 
static G4double ptableElectronCutInRange = 0.0*mm
 
static G4double pbartableElectronCutInRange = 0.0*mm
 
static G4double Charge
 
static G4double LowestKineticEnergy = 10.*eV
 
static G4double HighestKineticEnergy = 100.*GeV
 
static G4int TotBin = 360
 
static G4double RTable
 
static G4double LOGRTable
 
static G4double dRoverRange = 0.20
 
static G4double finalRange = 200.*micrometer
 
static G4double c1lim = dRoverRange
 
static G4double c2lim = 2.*(1.-dRoverRange)*finalRange
 
static G4double c3lim = -(1.-dRoverRange)*finalRange*finalRange
 
static G4bool rndmStepFlag = false
 
static G4bool EnlossFlucFlag = true
 

Detailed Description

Definition at line 93 of file G4hRDEnergyLoss.hh.

Constructor & Destructor Documentation

◆ G4hRDEnergyLoss()

G4hRDEnergyLoss::G4hRDEnergyLoss ( const G4String processName)

Definition at line 163 of file G4hRDEnergyLoss.cc.

164 : G4VContinuousDiscreteProcess (processName),
165 MaxExcitationNumber (1.e6),
166 probLimFluct (0.01),
167 nmaxDirectFluct (100),
168 nmaxCont1(4),
169 nmaxCont2(16),
170 theLossTable(0),
171 linLossLimit(0.05),
172 MinKineticEnergy(0.0)
173{;}
const G4double MaxExcitationNumber
G4double MinKineticEnergy
const long nmaxDirectFluct
G4PhysicsTable * theLossTable
const long nmaxCont1
const G4double probLimFluct
const long nmaxCont2

◆ ~G4hRDEnergyLoss()

G4hRDEnergyLoss::~G4hRDEnergyLoss ( )

Definition at line 177 of file G4hRDEnergyLoss.cc.

178{
179 if(theLossTable) {
181 delete theLossTable;
182 }
183}
void clearAndDestroy()

Member Function Documentation

◆ BuildDEDXTable()

void G4hRDEnergyLoss::BuildDEDXTable ( const G4ParticleDefinition aParticleType)
staticprotected

Definition at line 247 of file G4hRDEnergyLoss.cc.

249{
250 // calculate data members TotBin,LOGRTable,RTable first
251
252 const G4ProductionCutsTable* theCoupleTable=
254 size_t numOfCouples = theCoupleTable->GetTableSize();
255
256 // create/fill proton or antiproton tables depending on the charge
257 Charge = aParticleType.GetPDGCharge()/eplus;
258 ParticleMass = aParticleType.GetPDGMass() ;
259
260 if (Charge>0.) {theDEDXTable= theDEDXpTable;}
261 else {theDEDXTable= theDEDXpbarTable;}
262
263 if( ((Charge>0.) && (theDEDXTable==0)) ||
264 ((Charge<0.) && (theDEDXTable==0))
265 )
266 {
267
268 // Build energy loss table as a sum of the energy loss due to the
269 // different processes.
270 if( Charge >0.)
271 {
272 RecorderOfProcess=RecorderOfpProcess;
273 CounterOfProcess=CounterOfpProcess;
274
275 if(CounterOfProcess == NumberOfProcesses)
276 {
277 if(theDEDXpTable)
279 delete theDEDXpTable; }
280 theDEDXpTable = new G4PhysicsTable(numOfCouples);
281 theDEDXTable = theDEDXpTable;
282 }
283 }
284 else
285 {
286 RecorderOfProcess=RecorderOfpbarProcess;
287 CounterOfProcess=CounterOfpbarProcess;
288
289 if(CounterOfProcess == NumberOfProcesses)
290 {
293 delete theDEDXpbarTable; }
294 theDEDXpbarTable = new G4PhysicsTable(numOfCouples);
295 theDEDXTable = theDEDXpbarTable;
296 }
297 }
298
299 if(CounterOfProcess == NumberOfProcesses)
300 {
301 // loop for materials
302 G4double LowEdgeEnergy , Value ;
303 G4bool isOutRange ;
304 G4PhysicsTable* pointer ;
305
306 for (size_t J=0; J<numOfCouples; J++)
307 {
308 // create physics vector and fill it
311
312 // loop for the kinetic energy
313 for (G4int i=0; i<TotBin; i++)
314 {
315 LowEdgeEnergy = aVector->GetLowEdgeEnergy(i) ;
316 Value = 0. ;
317
318 // loop for the contributing processes
319 for (G4int process=0; process < NumberOfProcesses; process++)
320 {
321 pointer= RecorderOfProcess[process];
322 Value += (*pointer)[J]->
323 GetValue(LowEdgeEnergy,isOutRange) ;
324 }
325
326 aVector->PutValue(i,Value) ;
327 }
328
329 theDEDXTable->insert(aVector) ;
330 }
331
332 // reset counter to zero ..................
333 if( Charge >0.)
335 else
337
338 // Build range table
339 BuildRangeTable( aParticleType);
340
341 // Build lab/proper time tables
342 BuildTimeTables( aParticleType) ;
343
344 // Build coeff tables for the energy loss calculation
345 BuildRangeCoeffATable( aParticleType);
346 BuildRangeCoeffBTable( aParticleType);
347 BuildRangeCoeffCTable( aParticleType);
348
349 // invert the range table
350
351 BuildInverseRangeTable(aParticleType);
352 }
353 }
354 // make the energy loss and the range table available
355
356 G4EnergyLossTables::Register(&aParticleType,
357 (Charge>0)?
359 (Charge>0)?
361 (Charge>0)?
363 (Charge>0)?
365 (Charge>0)?
368 proton_mass_c2/aParticleType.GetPDGMass(),
369 TotBin);
370
371}
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
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 *)
virtual G4double GetLowEdgeEnergy(size_t binNumber) const
void PutValue(size_t index, G4double theValue)
static G4ProductionCutsTable * GetProductionCutsTable()
static G4PhysicsTable * theDEDXpTable
static G4PhysicsTable * theInverseRangepbarTable
static G4PhysicsTable * theInverseRangepTable
static G4PhysicsTable ** RecorderOfpbarProcess
static G4int CounterOfpbarProcess
static G4PhysicsTable * theRangepTable
static G4PhysicsTable * theDEDXpbarTable
static G4double LowestKineticEnergy
static G4PhysicsTable * theLabTimepTable
static G4double Charge
static G4PhysicsTable * theProperTimepbarTable
static G4PhysicsTable * theProperTimepTable
static G4PhysicsTable ** RecorderOfpProcess
static G4PhysicsTable * theLabTimepbarTable
static G4int TotBin
static G4double HighestKineticEnergy
static G4int CounterOfpProcess
static G4double ParticleMass
static G4PhysicsTable * theRangepbarTable

Referenced by G4hImpactIonisation::BuildPhysicsTable().

◆ CutsWhereModified()

G4bool G4hRDEnergyLoss::CutsWhereModified ( )
protected

Definition at line 1141 of file G4hRDEnergyLoss.cc.

1142{
1143 G4bool wasModified = false;
1144 const G4ProductionCutsTable* theCoupleTable=
1146 size_t numOfCouples = theCoupleTable->GetTableSize();
1147
1148 for (size_t j=0; j<numOfCouples; j++){
1149 if (theCoupleTable->GetMaterialCutsCouple(j)->IsRecalcNeeded()) {
1150 wasModified = true;
1151 break;
1152 }
1153 }
1154 return wasModified;
1155}
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 187 of file G4hRDEnergyLoss.cc.

188{
189 return NumberOfProcesses;
190}

◆ MinusNumberOfProcesses()

void G4hRDEnergyLoss::MinusNumberOfProcesses ( )
static

Definition at line 208 of file G4hRDEnergyLoss.cc.

209{
210 NumberOfProcesses--;
211}

◆ PlusNumberOfProcesses()

void G4hRDEnergyLoss::PlusNumberOfProcesses ( )
static

Definition at line 201 of file G4hRDEnergyLoss.cc.

202{
203 NumberOfProcesses++;
204}

◆ 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 215 of file G4hRDEnergyLoss.cc.

216{
217 dRoverRange = value;
218}
static G4double dRoverRange

◆ SetEnlossFluc()

void G4hRDEnergyLoss::SetEnlossFluc ( G4bool  value)
static

Definition at line 229 of file G4hRDEnergyLoss.cc.

230{
231 EnlossFlucFlag = value;
232}
static G4bool EnlossFlucFlag

◆ SetNumberOfProcesses()

void G4hRDEnergyLoss::SetNumberOfProcesses ( G4int  number)
static

Definition at line 194 of file G4hRDEnergyLoss.cc.

195{
196 NumberOfProcesses=number;
197}

◆ SetRndmStep()

void G4hRDEnergyLoss::SetRndmStep ( G4bool  value)
static

Definition at line 222 of file G4hRDEnergyLoss.cc.

223{
224 rndmStepFlag = value;
225}
static G4bool rndmStepFlag

◆ SetStepFunction()

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

Definition at line 236 of file G4hRDEnergyLoss.cc.

237{
238 dRoverRange = c1;
239 finalRange = c2;
243}
static G4double finalRange
static G4double c2lim
static G4double c3lim
static G4double c1lim

Member Data Documentation

◆ c1lim

G4double G4hRDEnergyLoss::c1lim = dRoverRange
staticprotected

Definition at line 194 of file G4hRDEnergyLoss.hh.

Referenced by SetStepFunction().

◆ c2lim

G4double G4hRDEnergyLoss::c2lim = 2.*(1.-dRoverRange)*finalRange
staticprotected

Definition at line 194 of file G4hRDEnergyLoss.hh.

Referenced by SetStepFunction().

◆ c3lim

G4double G4hRDEnergyLoss::c3lim = -(1.-dRoverRange)*finalRange*finalRange
staticprotected

Definition at line 194 of file G4hRDEnergyLoss.hh.

Referenced by SetStepFunction().

◆ Charge

G4double G4hRDEnergyLoss::Charge
staticprotected

Definition at line 174 of file G4hRDEnergyLoss.hh.

Referenced by BuildDEDXTable().

◆ CounterOfpbarProcess

G4int G4hRDEnergyLoss::CounterOfpbarProcess = 0
staticprotected

Definition at line 165 of file G4hRDEnergyLoss.hh.

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

◆ CounterOfpProcess

G4int G4hRDEnergyLoss::CounterOfpProcess = 0
staticprotected

Definition at line 164 of file G4hRDEnergyLoss.hh.

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

◆ dRoverRange

G4double G4hRDEnergyLoss::dRoverRange = 0.20
staticprotected

Definition at line 191 of file G4hRDEnergyLoss.hh.

Referenced by SetdRoverRange(), and SetStepFunction().

◆ EnlossFlucFlag

G4bool G4hRDEnergyLoss::EnlossFlucFlag = true
staticprotected

Definition at line 197 of file G4hRDEnergyLoss.hh.

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

◆ finalRange

G4double G4hRDEnergyLoss::finalRange = 200.*micrometer
staticprotected

Definition at line 193 of file G4hRDEnergyLoss.hh.

Referenced by SetStepFunction().

◆ HighestKineticEnergy

◆ linLossLimit

G4double G4hRDEnergyLoss::linLossLimit
protected

Definition at line 187 of file G4hRDEnergyLoss.hh.

Referenced by G4hImpactIonisation::AlongStepDoIt().

◆ LOGRTable

G4double G4hRDEnergyLoss::LOGRTable
staticprotected

Definition at line 181 of file G4hRDEnergyLoss.hh.

◆ LowestKineticEnergy

G4double G4hRDEnergyLoss::LowestKineticEnergy = 10.*eV
staticprotected

◆ MaxExcitationNumber

const G4double G4hRDEnergyLoss::MaxExcitationNumber
protected

Definition at line 136 of file G4hRDEnergyLoss.hh.

◆ MinKineticEnergy

G4double G4hRDEnergyLoss::MinKineticEnergy
protected

◆ nmaxCont1

const long G4hRDEnergyLoss::nmaxCont1
protected

Definition at line 138 of file G4hRDEnergyLoss.hh.

◆ nmaxCont2

const long G4hRDEnergyLoss::nmaxCont2
protected

Definition at line 138 of file G4hRDEnergyLoss.hh.

◆ nmaxDirectFluct

const long G4hRDEnergyLoss::nmaxDirectFluct
protected

Definition at line 138 of file G4hRDEnergyLoss.hh.

◆ ParticleMass

G4double G4hRDEnergyLoss::ParticleMass
staticprotected

Definition at line 168 of file G4hRDEnergyLoss.hh.

Referenced by BuildDEDXTable().

◆ pbartableElectronCutInRange

G4double G4hRDEnergyLoss::pbartableElectronCutInRange = 0.0*mm
staticprotected

Definition at line 172 of file G4hRDEnergyLoss.hh.

◆ probLimFluct

const G4double G4hRDEnergyLoss::probLimFluct
protected

Definition at line 137 of file G4hRDEnergyLoss.hh.

◆ ptableElectronCutInRange

G4double G4hRDEnergyLoss::ptableElectronCutInRange = 0.0*mm
staticprotected

Definition at line 171 of file G4hRDEnergyLoss.hh.

◆ RecorderOfpbarProcess

G4PhysicsTable ** G4hRDEnergyLoss::RecorderOfpbarProcess
staticprotected
Initial value:
=
new G4PhysicsTable*[100]

Definition at line 163 of file G4hRDEnergyLoss.hh.

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

◆ RecorderOfpProcess

G4PhysicsTable ** G4hRDEnergyLoss::RecorderOfpProcess
staticprotected
Initial value:
=
new G4PhysicsTable*[100]

Definition at line 162 of file G4hRDEnergyLoss.hh.

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

◆ rndmStepFlag

G4bool G4hRDEnergyLoss::rndmStepFlag = false
staticprotected

Definition at line 196 of file G4hRDEnergyLoss.hh.

Referenced by SetRndmStep().

◆ RTable

G4double G4hRDEnergyLoss::RTable
staticprotected

Definition at line 181 of file G4hRDEnergyLoss.hh.

◆ theDEDXpbarTable

G4PhysicsTable * G4hRDEnergyLoss::theDEDXpbarTable = 0
staticprotected

Definition at line 145 of file G4hRDEnergyLoss.hh.

Referenced by BuildDEDXTable().

◆ theDEDXpTable

G4PhysicsTable * G4hRDEnergyLoss::theDEDXpTable = 0
staticprotected

Definition at line 144 of file G4hRDEnergyLoss.hh.

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

◆ theInverseRangepbarTable

G4PhysicsTable * G4hRDEnergyLoss::theInverseRangepbarTable = 0
staticprotected

Definition at line 151 of file G4hRDEnergyLoss.hh.

Referenced by BuildDEDXTable().

◆ theInverseRangepTable

G4PhysicsTable * G4hRDEnergyLoss::theInverseRangepTable = 0
staticprotected

Definition at line 150 of file G4hRDEnergyLoss.hh.

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

◆ theLabTimepbarTable

G4PhysicsTable * G4hRDEnergyLoss::theLabTimepbarTable = 0
staticprotected

Definition at line 155 of file G4hRDEnergyLoss.hh.

Referenced by BuildDEDXTable().

◆ theLabTimepTable

G4PhysicsTable * G4hRDEnergyLoss::theLabTimepTable = 0
staticprotected

Definition at line 154 of file G4hRDEnergyLoss.hh.

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

◆ theLossTable

G4PhysicsTable* G4hRDEnergyLoss::theLossTable
protected

Definition at line 185 of file G4hRDEnergyLoss.hh.

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

◆ theProperTimepbarTable

G4PhysicsTable * G4hRDEnergyLoss::theProperTimepbarTable = 0
staticprotected

Definition at line 158 of file G4hRDEnergyLoss.hh.

Referenced by BuildDEDXTable().

◆ theProperTimepTable

G4PhysicsTable * G4hRDEnergyLoss::theProperTimepTable = 0
staticprotected

Definition at line 157 of file G4hRDEnergyLoss.hh.

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

◆ theRangepbarTable

G4PhysicsTable * G4hRDEnergyLoss::theRangepbarTable = 0
staticprotected

Definition at line 147 of file G4hRDEnergyLoss.hh.

Referenced by BuildDEDXTable().

◆ theRangepTable

G4PhysicsTable * G4hRDEnergyLoss::theRangepTable = 0
staticprotected

Definition at line 146 of file G4hRDEnergyLoss.hh.

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

◆ TotBin

G4int G4hRDEnergyLoss::TotBin = 360
staticprotected

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