Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4hRDEnergyLoss Class Referenceabstract

#include <G4hRDEnergyLoss.hh>

Inheritance diagram for G4hRDEnergyLoss:
Collaboration 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 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
 
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
 
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 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 G4ThreadLocal
G4PhysicsTable
theDEDXpTable = 0
 
static G4ThreadLocal
G4PhysicsTable
theDEDXpbarTable = 0
 
static G4ThreadLocal
G4PhysicsTable
theRangepTable = 0
 
static G4ThreadLocal
G4PhysicsTable
theRangepbarTable = 0
 
static G4ThreadLocal
G4PhysicsTable
theInverseRangepTable = 0
 
static G4ThreadLocal
G4PhysicsTable
theInverseRangepbarTable = 0
 
static G4ThreadLocal
G4PhysicsTable
theLabTimepTable = 0
 
static G4ThreadLocal
G4PhysicsTable
theLabTimepbarTable = 0
 
static G4ThreadLocal
G4PhysicsTable
theProperTimepTable = 0
 
static G4ThreadLocal
G4PhysicsTable
theProperTimepbarTable = 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 93 of file G4hRDEnergyLoss.hh.

Constructor & Destructor Documentation

G4hRDEnergyLoss::G4hRDEnergyLoss ( const G4String processName)

Definition at line 162 of file G4hRDEnergyLoss.cc.

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

Definition at line 180 of file G4hRDEnergyLoss.cc.

181 {
182  if(theLossTable)
183  {
185  delete theLossTable;
186  }
187 }
G4PhysicsTable * theLossTable
void clearAndDestroy()

Here is the call graph for this function:

Member Function Documentation

void G4hRDEnergyLoss::BuildDEDXTable ( const G4ParticleDefinition aParticleType)
staticprotected

Definition at line 251 of file G4hRDEnergyLoss.cc.

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

Here is the call graph for this function:

Here is the caller graph for this function:

G4bool G4hRDEnergyLoss::CutsWhereModified ( )
protected

Definition at line 1155 of file G4hRDEnergyLoss.cc.

1156 {
1157  G4bool wasModified = false;
1158  const G4ProductionCutsTable* theCoupleTable=
1160  size_t numOfCouples = theCoupleTable->GetTableSize();
1161 
1162  for (size_t j=0; j<numOfCouples; j++){
1163  if (theCoupleTable->GetMaterialCutsCouple(j)->IsRecalcNeeded()) {
1164  wasModified = true;
1165  break;
1166  }
1167  }
1168  return wasModified;
1169 }
G4bool IsRecalcNeeded() const
bool G4bool
Definition: G4Types.hh:79
static G4ProductionCutsTable * GetProductionCutsTable()
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 191 of file G4hRDEnergyLoss.cc.

192 {
193  return NumberOfProcesses;
194 }
void G4hRDEnergyLoss::MinusNumberOfProcesses ( )
static

Definition at line 212 of file G4hRDEnergyLoss.cc.

213 {
214  NumberOfProcesses--;
215 }
void G4hRDEnergyLoss::PlusNumberOfProcesses ( )
static

Definition at line 205 of file G4hRDEnergyLoss.cc.

206 {
207  NumberOfProcesses++;
208 }
virtual G4VParticleChange* G4hRDEnergyLoss::PostStepDoIt ( const G4Track track,
const G4Step Step 
)
pure virtual

Reimplemented from G4VContinuousDiscreteProcess.

Implemented in G4hImpactIonisation.

void G4hRDEnergyLoss::SetdRoverRange ( G4double  value)
static

Definition at line 219 of file G4hRDEnergyLoss.cc.

220 {
221  dRoverRange = value;
222 }
const XML_Char int const XML_Char * value
Definition: expat.h:331
static G4ThreadLocal G4double dRoverRange
void G4hRDEnergyLoss::SetEnlossFluc ( G4bool  value)
static

Definition at line 233 of file G4hRDEnergyLoss.cc.

234 {
236 }
const XML_Char int const XML_Char * value
Definition: expat.h:331
static G4ThreadLocal G4bool EnlossFlucFlag
void G4hRDEnergyLoss::SetNumberOfProcesses ( G4int  number)
static

Definition at line 198 of file G4hRDEnergyLoss.cc.

199 {
200  NumberOfProcesses=number;
201 }
void G4hRDEnergyLoss::SetRndmStep ( G4bool  value)
static

Definition at line 226 of file G4hRDEnergyLoss.cc.

227 {
229 }
static G4ThreadLocal G4bool rndmStepFlag
const XML_Char int const XML_Char * value
Definition: expat.h:331
void G4hRDEnergyLoss::SetStepFunction ( G4double  c1,
G4double  c2 
)
static

Definition at line 240 of file G4hRDEnergyLoss.cc.

241 {
242  dRoverRange = c1;
243  finalRange = c2;
247 }
static c2_factory< G4double > c2
static G4ThreadLocal G4double c3lim
static G4ThreadLocal G4double c1lim
static G4ThreadLocal G4double finalRange
static G4ThreadLocal G4double dRoverRange
static G4ThreadLocal G4double c2lim
tuple c1
Definition: plottest35.py:14

Member Data Documentation

G4ThreadLocal G4double G4hRDEnergyLoss::c1lim = 0.20
staticprotected

Definition at line 194 of file G4hRDEnergyLoss.hh.

G4ThreadLocal G4double G4hRDEnergyLoss::c2lim = 0.32
staticprotected

Definition at line 194 of file G4hRDEnergyLoss.hh.

G4ThreadLocal G4double G4hRDEnergyLoss::c3lim = -0.032
staticprotected

Definition at line 194 of file G4hRDEnergyLoss.hh.

G4ThreadLocal G4double G4hRDEnergyLoss::Charge
staticprotected

Definition at line 174 of file G4hRDEnergyLoss.hh.

G4ThreadLocal G4int G4hRDEnergyLoss::CounterOfpbarProcess = 0
staticprotected

Definition at line 165 of file G4hRDEnergyLoss.hh.

G4ThreadLocal G4int G4hRDEnergyLoss::CounterOfpProcess = 0
staticprotected

Definition at line 164 of file G4hRDEnergyLoss.hh.

G4ThreadLocal G4double G4hRDEnergyLoss::dRoverRange = 0.20
staticprotected

Definition at line 191 of file G4hRDEnergyLoss.hh.

G4ThreadLocal G4bool G4hRDEnergyLoss::EnlossFlucFlag = true
staticprotected

Definition at line 197 of file G4hRDEnergyLoss.hh.

G4ThreadLocal G4double G4hRDEnergyLoss::finalRange = 0.2
staticprotected

Definition at line 193 of file G4hRDEnergyLoss.hh.

G4ThreadLocal G4double G4hRDEnergyLoss::HighestKineticEnergy = 1.e5
staticprotected

Definition at line 177 of file G4hRDEnergyLoss.hh.

G4double G4hRDEnergyLoss::linLossLimit
protected

Definition at line 187 of file G4hRDEnergyLoss.hh.

G4ThreadLocal G4double G4hRDEnergyLoss::LOGRTable
staticprotected

Definition at line 181 of file G4hRDEnergyLoss.hh.

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

Definition at line 176 of file G4hRDEnergyLoss.hh.

const G4double G4hRDEnergyLoss::MaxExcitationNumber
protected

Definition at line 136 of file G4hRDEnergyLoss.hh.

G4double G4hRDEnergyLoss::MinKineticEnergy
protected

Definition at line 189 of file G4hRDEnergyLoss.hh.

const long G4hRDEnergyLoss::nmaxCont1
protected

Definition at line 138 of file G4hRDEnergyLoss.hh.

const long G4hRDEnergyLoss::nmaxCont2
protected

Definition at line 138 of file G4hRDEnergyLoss.hh.

const long G4hRDEnergyLoss::nmaxDirectFluct
protected

Definition at line 138 of file G4hRDEnergyLoss.hh.

G4ThreadLocal G4double G4hRDEnergyLoss::ParticleMass
staticprotected

Definition at line 168 of file G4hRDEnergyLoss.hh.

G4ThreadLocal G4double G4hRDEnergyLoss::pbartableElectronCutInRange = 0.0
staticprotected

Definition at line 172 of file G4hRDEnergyLoss.hh.

const G4double G4hRDEnergyLoss::probLimFluct
protected

Definition at line 137 of file G4hRDEnergyLoss.hh.

G4ThreadLocal G4double G4hRDEnergyLoss::ptableElectronCutInRange = 0.0
staticprotected

Definition at line 171 of file G4hRDEnergyLoss.hh.

G4ThreadLocal G4PhysicsTable ** G4hRDEnergyLoss::RecorderOfpbarProcess = 0
staticprotected

Definition at line 163 of file G4hRDEnergyLoss.hh.

G4ThreadLocal G4PhysicsTable ** G4hRDEnergyLoss::RecorderOfpProcess = 0
staticprotected

Definition at line 162 of file G4hRDEnergyLoss.hh.

G4ThreadLocal G4bool G4hRDEnergyLoss::rndmStepFlag = false
staticprotected

Definition at line 196 of file G4hRDEnergyLoss.hh.

G4ThreadLocal G4double G4hRDEnergyLoss::RTable
staticprotected

Definition at line 181 of file G4hRDEnergyLoss.hh.

G4ThreadLocal G4PhysicsTable * G4hRDEnergyLoss::theDEDXpbarTable = 0
staticprotected

Definition at line 145 of file G4hRDEnergyLoss.hh.

G4ThreadLocal G4PhysicsTable * G4hRDEnergyLoss::theDEDXpTable = 0
staticprotected

Definition at line 144 of file G4hRDEnergyLoss.hh.

G4ThreadLocal G4PhysicsTable * G4hRDEnergyLoss::theInverseRangepbarTable = 0
staticprotected

Definition at line 151 of file G4hRDEnergyLoss.hh.

G4ThreadLocal G4PhysicsTable * G4hRDEnergyLoss::theInverseRangepTable = 0
staticprotected

Definition at line 150 of file G4hRDEnergyLoss.hh.

G4ThreadLocal G4PhysicsTable * G4hRDEnergyLoss::theLabTimepbarTable = 0
staticprotected

Definition at line 155 of file G4hRDEnergyLoss.hh.

G4ThreadLocal G4PhysicsTable * G4hRDEnergyLoss::theLabTimepTable = 0
staticprotected

Definition at line 154 of file G4hRDEnergyLoss.hh.

G4PhysicsTable* G4hRDEnergyLoss::theLossTable
protected

Definition at line 185 of file G4hRDEnergyLoss.hh.

G4ThreadLocal G4PhysicsTable * G4hRDEnergyLoss::theProperTimepbarTable = 0
staticprotected

Definition at line 158 of file G4hRDEnergyLoss.hh.

G4ThreadLocal G4PhysicsTable * G4hRDEnergyLoss::theProperTimepTable = 0
staticprotected

Definition at line 157 of file G4hRDEnergyLoss.hh.

G4ThreadLocal G4PhysicsTable * G4hRDEnergyLoss::theRangepbarTable = 0
staticprotected

Definition at line 147 of file G4hRDEnergyLoss.hh.

G4ThreadLocal G4PhysicsTable * G4hRDEnergyLoss::theRangepTable = 0
staticprotected

Definition at line 146 of file G4hRDEnergyLoss.hh.

G4ThreadLocal G4int G4hRDEnergyLoss::TotBin = 360
staticprotected

Definition at line 178 of file G4hRDEnergyLoss.hh.


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