Geant4  10.02.p03
G4eLowEnergyLoss Class Referenceabstract

#include <G4eLowEnergyLoss.hh>

Inheritance diagram for G4eLowEnergyLoss:
Collaboration diagram for G4eLowEnergyLoss:

Public Member Functions

 G4eLowEnergyLoss (const G4String &)
 
 ~G4eLowEnergyLoss ()
 
G4bool IsApplicable (const G4ParticleDefinition &)
 
void BuildDEDXTable (const G4ParticleDefinition &aParticleType)
 
G4double GetContinuousStepLimit (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)
 
G4VParticleChange * AlongStepDoIt (const G4Track &track, const G4Step &Step)
 
virtual G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
virtual G4VParticleChange * PostStepDoIt (const G4Track &track, const G4Step &step)=0
 
void ActivateFluorescence (G4bool val)
 
G4bool Fluorescence () const
 
- Public Member Functions inherited from G4RDVeLowEnergyLoss
 G4RDVeLowEnergyLoss (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4RDVeLowEnergyLoss (G4RDVeLowEnergyLoss &)
 
virtual ~G4RDVeLowEnergyLoss ()
 
- 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 G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChange * AtRestDoIt (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 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 void SetNbOfProcesses (G4int nb)
 
static void PlusNbOfProcesses ()
 
static void MinusNbOfProcesses ()
 
static G4int GetNbOfProcesses ()
 
static void SetLowerBoundEloss (G4double val)
 
static void SetUpperBoundEloss (G4double val)
 
static void SetNbinEloss (G4int nb)
 
static G4double GetLowerBoundEloss ()
 
static G4double GetUpperBoundEloss ()
 
static G4int GetNbinEloss ()
 
- Static Public Member Functions inherited from G4RDVeLowEnergyLoss
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

virtual std::vector< G4DynamicParticle * > * DeexciteAtom (const G4MaterialCutsCouple *, G4double, G4double)
 
- Protected Member Functions inherited from G4RDVeLowEnergyLoss
G4double GetLossWithFluct (const G4DynamicParticle *aParticle, const G4MaterialCutsCouple *couple, G4double MeanLoss, G4double step)
 
- Protected Member Functions inherited from G4VContinuousDiscreteProcess
void SetGPILSelection (G4GPILSelection selection)
 
G4GPILSelection GetGPILSelection () const
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Protected Attributes

G4PhysicsTabletheLossTable
 
G4double MinKineticEnergy
 
G4double Charge
 
G4double lastCharge
 
- Protected Attributes inherited from G4RDVeLowEnergyLoss
const G4MateriallastMaterial
 
G4int imat
 
G4double f1Fluct
 
G4double f2Fluct
 
G4double e1Fluct
 
G4double e2Fluct
 
G4double rateFluct
 
G4double ipotFluct
 
G4double e1LogFluct
 
G4double e2LogFluct
 
G4double ipotLogFluct
 
const G4int nmaxCont1
 
const G4int nmaxCont2
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChange * pParticleChange
 
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 G4PhysicsTabletheDEDXElectronTable = 0
 
static G4PhysicsTabletheDEDXPositronTable = 0
 
static G4PhysicsTabletheRangeElectronTable = 0
 
static G4PhysicsTabletheRangePositronTable = 0
 
static G4PhysicsTabletheInverseRangeElectronTable = 0
 
static G4PhysicsTabletheInverseRangePositronTable = 0
 
static G4PhysicsTabletheLabTimeElectronTable = 0
 
static G4PhysicsTabletheLabTimePositronTable = 0
 
static G4PhysicsTabletheProperTimeElectronTable = 0
 
static G4PhysicsTabletheProperTimePositronTable = 0
 
static G4int NbOfProcesses = 2
 
static G4int CounterOfElectronProcess = 0
 
static G4int CounterOfPositronProcess = 0
 
static G4PhysicsTable ** RecorderOfElectronProcess
 
static G4PhysicsTable ** RecorderOfPositronProcess
 
- Static Protected Attributes inherited from G4RDVeLowEnergyLoss
static G4double ParticleMass
 
static G4double taulow
 
static G4double tauhigh
 
static G4double ltaulow
 
static G4double ltauhigh
 
static G4double dRoverRange = 20*perCent
 
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
 

Private Member Functions

G4double GetConstraints (const G4DynamicParticle *aParticle, const G4MaterialCutsCouple *couple)
 
 G4eLowEnergyLoss (G4eLowEnergyLoss &)
 
G4eLowEnergyLossoperator= (const G4eLowEnergyLoss &right)
 

Private Attributes

G4PhysicsTabletheDEDXTable
 
G4int CounterOfProcess
 
G4PhysicsTable ** RecorderOfProcess
 
G4double fdEdx
 
G4double fRangeNow
 
G4double linLossLimit
 
G4ParticleChangeForLoss fParticleChange
 
G4bool theFluo
 

Static Private Attributes

static G4int NbinEloss = 360
 
static G4double LowerBoundEloss = 10.*eV
 
static G4double UpperBoundEloss = 100.*GeV
 
static G4double RTable
 
static G4double LOGRTable
 
static G4PhysicsTabletheeRangeCoeffATable = 0
 
static G4PhysicsTabletheeRangeCoeffBTable = 0
 
static G4PhysicsTabletheeRangeCoeffCTable = 0
 
static G4PhysicsTablethepRangeCoeffATable = 0
 
static G4PhysicsTablethepRangeCoeffBTable = 0
 
static G4PhysicsTablethepRangeCoeffCTable = 0
 
static G4EnergyLossMessengereLossMessenger = 0
 

Additional Inherited Members

- Static Protected Member Functions inherited from G4RDVeLowEnergyLoss
static G4PhysicsTableBuildRangeTable (G4PhysicsTable *theDEDXTable, G4PhysicsTable *theRangeTable, G4double Tmin, G4double Tmax, G4int nbin)
 
static G4PhysicsTableBuildLabTimeTable (G4PhysicsTable *theDEDXTable, G4PhysicsTable *theLabTimeTable, G4double Tmin, G4double Tmax, G4int nbin)
 
static G4PhysicsTableBuildProperTimeTable (G4PhysicsTable *theDEDXTable, G4PhysicsTable *ProperTimeTable, G4double Tmin, G4double Tmax, G4int nbin)
 
static G4PhysicsTableBuildRangeCoeffATable (G4PhysicsTable *theRangeTable, G4PhysicsTable *theCoeffATable, G4double Tmin, G4double Tmax, G4int nbin)
 
static G4PhysicsTableBuildRangeCoeffBTable (G4PhysicsTable *theRangeTable, G4PhysicsTable *theCoeffBTable, G4double Tmin, G4double Tmax, G4int nbin)
 
static G4PhysicsTableBuildRangeCoeffCTable (G4PhysicsTable *theRangeTable, G4PhysicsTable *theCoeffCTable, G4double Tmin, G4double Tmax, G4int nbin)
 
static G4PhysicsTableBuildInverseRangeTable (G4PhysicsTable *theRangeTable, G4PhysicsTable *theRangeCoeffATable, G4PhysicsTable *theRangeCoeffBTable, G4PhysicsTable *theRangeCoeffCTable, G4PhysicsTable *theInverseRangeTable, G4double Tmin, G4double Tmax, G4int nbin)
 

Detailed Description

Definition at line 87 of file G4eLowEnergyLoss.hh.

Constructor & Destructor Documentation

◆ G4eLowEnergyLoss() [1/2]

G4eLowEnergyLoss::G4eLowEnergyLoss ( const G4String processName)

Definition at line 127 of file G4eLowEnergyLoss.cc.

128  : G4RDVeLowEnergyLoss (processName),
129  theLossTable(0),
130  MinKineticEnergy(1.*eV),
131  Charge(-1.),
132  lastCharge(0.),
133  theDEDXTable(0),
134  CounterOfProcess(0),
136  fdEdx(0),
137  fRangeNow(0),
138  linLossLimit(0.05),
139  theFluo(false)
140 {
141 
142  //create (only once) EnergyLoss messenger
144 }
G4PhysicsTable * theDEDXTable
static G4EnergyLossMessenger * eLossMessenger
static const double eV
Definition: G4SIunits.hh:212
G4PhysicsTable * theLossTable
G4PhysicsTable ** RecorderOfProcess

◆ ~G4eLowEnergyLoss()

G4eLowEnergyLoss::~G4eLowEnergyLoss ( )

Definition at line 148 of file G4eLowEnergyLoss.cc.

149 {
150  if (theLossTable)
151  {
153  delete theLossTable;
154  }
155 }
G4PhysicsTable * theLossTable
void clearAndDestroy()
Here is the call graph for this function:

◆ G4eLowEnergyLoss() [2/2]

G4eLowEnergyLoss::G4eLowEnergyLoss ( G4eLowEnergyLoss )
private

Member Function Documentation

◆ ActivateFluorescence()

void G4eLowEnergyLoss::ActivateFluorescence ( G4bool  val)

◆ AlongStepDoIt()

G4VParticleChange * G4eLowEnergyLoss::AlongStepDoIt ( const G4Track &  track,
const G4Step &  Step 
)
virtual

Implements G4RDVeLowEnergyLoss.

Definition at line 399 of file G4eLowEnergyLoss.cc.

401 {
402  // compute the energy loss after a Step
403 
404  static const G4double faclow = 1.5 ;
405 
406  // get particle and material pointers from trackData
407  const G4DynamicParticle* aParticle = trackData.GetDynamicParticle();
408  G4double E = aParticle->GetKineticEnergy();
409 
410  const G4MaterialCutsCouple* couple = trackData.GetMaterialCutsCouple();
411 
412  G4double Step = stepData.GetStepLength();
413 
414  aParticleChange.Initialize(trackData);
415  //fParticleChange.Initialize(trackData);
416 
417  G4double MeanLoss, finalT;
418 
419  if (E < MinKineticEnergy) finalT = 0.;
420 
421  else if ( E< faclow*LowerBoundEloss)
422  {
423  if (Step >= fRangeNow) finalT = 0.;
424  // else finalT = E*(1.-Step/fRangeNow) ;
425  else finalT = E*(1.-std::sqrt(Step/fRangeNow)) ;
426  }
427 
428  else if (E>=UpperBoundEloss) finalT = E - Step*fdEdx;
429 
430  else if (Step >= fRangeNow) finalT = 0.;
431 
432  else
433  {
434  if(Step/fRangeNow < linLossLimit) finalT = E-Step*fdEdx ;
435  else
436  {
438  (G4Electron::Electron(),fRangeNow-Step,couple);
440  (G4Positron::Positron(),fRangeNow-Step,couple);
441  }
442  }
443 
444  if(finalT < MinKineticEnergy) finalT = 0. ;
445 
446  MeanLoss = E-finalT ;
447 
448  //now the loss with fluctuation
449  if ((EnlossFlucFlag) && (finalT > 0.) && (finalT < E)&&(E > LowerBoundEloss))
450  {
451  finalT = E-GetLossWithFluct(aParticle,couple,MeanLoss,Step);
452  if (finalT < 0.) finalT = 0.;
453  }
454 
455  // kill the particle if the kinetic energy <= 0
456  if (finalT <= 0. )
457  {
458  finalT = 0.;
459  if(Charge > 0.0) aParticleChange.ProposeTrackStatus(fStopButAlive);
460  else aParticleChange.ProposeTrackStatus(fStopAndKill);
461  }
462 
463  G4double edep = E - finalT;
464 
465  aParticleChange.ProposeEnergy(finalT);
466 
467  // Deexcitation of ionised atoms
468  std::vector<G4DynamicParticle*>* deexcitationProducts = 0;
469  if (theFluo) deexcitationProducts = DeexciteAtom(couple,E,edep);
470 
471  size_t nSecondaries = 0;
472  if (deexcitationProducts != 0) nSecondaries = deexcitationProducts->size();
473  aParticleChange.SetNumberOfSecondaries(nSecondaries);
474 
475  if (nSecondaries > 0) {
476 
477  const G4StepPoint* preStep = stepData.GetPreStepPoint();
478  const G4StepPoint* postStep = stepData.GetPostStepPoint();
479  G4ThreeVector r = preStep->GetPosition();
480  G4ThreeVector deltaR = postStep->GetPosition();
481  deltaR -= r;
482  G4double t = preStep->GetGlobalTime();
483  G4double deltaT = postStep->GetGlobalTime();
484  deltaT -= t;
485  G4double time, q;
487 
488  for (size_t i=0; i<nSecondaries; i++) {
489 
490  G4DynamicParticle* part = (*deexcitationProducts)[i];
491  if (part != 0) {
492  G4double eSecondary = part->GetKineticEnergy();
493  edep -= eSecondary;
494  if (edep > 0.)
495  {
496  q = G4UniformRand();
497  time = deltaT*q + t;
498  position = deltaR*q;
499  position += r;
500  G4Track* newTrack = new G4Track(part, time, position);
501  aParticleChange.AddSecondary(newTrack);
502  }
503  else
504  {
505  edep += eSecondary;
506  delete part;
507  part = 0;
508  }
509  }
510  }
511  }
512  delete deexcitationProducts;
513 
514  aParticleChange.ProposeLocalEnergyDeposit(edep);
515 
516  return &aParticleChange;
517 }
static G4double LowerBoundEloss
Double_t edep
G4double GetKineticEnergy() const
#define position
Definition: xmlparse.cc:622
#define G4UniformRand()
Definition: Randomize.hh:97
TString part[npart]
static G4Positron * Positron()
Definition: G4Positron.cc:94
static G4double UpperBoundEloss
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
Definition: Step.hh:41
static G4Electron * Electron()
Definition: G4Electron.cc:94
double G4double
Definition: G4Types.hh:76
virtual std::vector< G4DynamicParticle * > * DeexciteAtom(const G4MaterialCutsCouple *, G4double, G4double)
G4double GetLossWithFluct(const G4DynamicParticle *aParticle, const G4MaterialCutsCouple *couple, G4double MeanLoss, G4double step)
static G4double GetPreciseEnergyFromRange(const G4ParticleDefinition *aParticle, G4double range, const G4Material *aMaterial)
Here is the call graph for this function:

◆ BuildDEDXTable()

void G4eLowEnergyLoss::BuildDEDXTable ( const G4ParticleDefinition aParticleType)

Definition at line 208 of file G4eLowEnergyLoss.cc.

210 {
211  ParticleMass = aParticleType.GetPDGMass();
212  Charge = aParticleType.GetPDGCharge()/eplus;
213 
214  // calculate data members LOGRTable,RTable first
215 
216  G4double lrate = std::log(UpperBoundEloss/LowerBoundEloss);
217  LOGRTable=lrate/NbinEloss;
218  RTable =std::exp(LOGRTable);
219  // Build energy loss table as a sum of the energy loss due to the
220  // different processes.
221  //
222 
223  const G4ProductionCutsTable* theCoupleTable=
225  size_t numOfCouples = theCoupleTable->GetTableSize();
226 
227  // create table for the total energy loss
228 
229  if (&aParticleType==G4Electron::Electron())
230  {
234  {
236  {
238  delete theDEDXElectronTable;
239  }
240  theDEDXElectronTable = new G4PhysicsTable(numOfCouples);
242  }
243  }
244  if (&aParticleType==G4Positron::Positron())
245  {
249  {
251  {
253  delete theDEDXPositronTable;
254  }
255  theDEDXPositronTable = new G4PhysicsTable(numOfCouples);
257  }
258  }
259 
261  {
262  // fill the tables
263  // loop for materials
264  G4double LowEdgeEnergy , Value;
265  G4bool isOutRange;
266  G4PhysicsTable* pointer;
267 
268  for (size_t J=0; J<numOfCouples; J++)
269  {
270  // create physics vector and fill it
271 
272  G4PhysicsLogVector* aVector = new G4PhysicsLogVector(
274 
275  // loop for the kinetic energy
276  for (G4int i=0; i<NbinEloss; i++)
277  {
278  LowEdgeEnergy = aVector->GetLowEdgeEnergy(i) ;
279  //here comes the sum of the different tables created by the
280  //processes (ionisation,bremsstrahlung,etc...)
281  Value = 0.;
282  for (G4int process=0; process < NbOfProcesses; process++)
283  {
284  pointer= RecorderOfProcess[process];
285  Value += (*pointer)[J]->GetValue(LowEdgeEnergy,isOutRange);
286  }
287 
288  aVector->PutValue(i,Value) ;
289  }
290 
291  theDEDXTable->insert(aVector) ;
292 
293  }
294 
295 
296  //reset counter to zero
297  if (&aParticleType==G4Electron::Electron()) CounterOfElectronProcess=0;
298  if (&aParticleType==G4Positron::Positron()) CounterOfPositronProcess=0;
299 
300  ParticleMass = aParticleType.GetPDGMass();
301 
302  if (&aParticleType==G4Electron::Electron())
303  {
304  // Build range table
307  LowerBoundEloss,UpperBoundEloss,NbinEloss);
308 
309  // Build lab/proper time tables
312  LowerBoundEloss,UpperBoundEloss,NbinEloss);
315  LowerBoundEloss,UpperBoundEloss,NbinEloss);
316 
317  // Build coeff tables for the energy loss calculation
320  LowerBoundEloss,UpperBoundEloss,NbinEloss);
321 
324  LowerBoundEloss,UpperBoundEloss,NbinEloss);
325 
328  LowerBoundEloss,UpperBoundEloss,NbinEloss);
329 
330  // invert the range table
336  LowerBoundEloss,UpperBoundEloss,NbinEloss);
337  }
338  if (&aParticleType==G4Positron::Positron())
339  {
340  // Build range table
343  LowerBoundEloss,UpperBoundEloss,NbinEloss);
344 
345 
346  // Build lab/proper time tables
349  LowerBoundEloss,UpperBoundEloss,NbinEloss);
352  LowerBoundEloss,UpperBoundEloss,NbinEloss);
353 
354  // Build coeff tables for the energy loss calculation
357  LowerBoundEloss,UpperBoundEloss,NbinEloss);
358 
361  LowerBoundEloss,UpperBoundEloss,NbinEloss);
362 
365  LowerBoundEloss,UpperBoundEloss,NbinEloss);
366 
367  // invert the range table
373  LowerBoundEloss,UpperBoundEloss,NbinEloss);
374  }
375 
376  if(verboseLevel > 1) {
377  G4cout << (*theDEDXElectronTable) << G4endl;
378  }
379 
380 
381  // make the energy loss and the range table available
382  G4EnergyLossTables::Register(&aParticleType,
383  (&aParticleType==G4Electron::Electron())?
385  (&aParticleType==G4Electron::Electron())?
387  (&aParticleType==G4Electron::Electron())?
389  (&aParticleType==G4Electron::Electron())?
391  (&aParticleType==G4Electron::Electron())?
393  LowerBoundEloss, UpperBoundEloss, 1.,NbinEloss);
394  }
395 }
G4int verboseLevel
Definition: G4VProcess.hh:368
void insert(G4PhysicsVector *)
static G4PhysicsTable * theInverseRangePositronTable
static G4PhysicsTable * theProperTimePositronTable
static G4int CounterOfElectronProcess
static G4int NbOfProcesses
static G4PhysicsTable * thepRangeCoeffATable
static G4PhysicsTable * theRangeElectronTable
static G4double LowerBoundEloss
G4PhysicsTable * theDEDXTable
static G4double ParticleMass
int G4int
Definition: G4Types.hh:78
static G4PhysicsTable * BuildRangeTable(G4PhysicsTable *theDEDXTable, G4PhysicsTable *theRangeTable, G4double Tmin, G4double Tmax, G4int nbin)
static G4PhysicsTable * BuildProperTimeTable(G4PhysicsTable *theDEDXTable, G4PhysicsTable *ProperTimeTable, G4double Tmin, G4double Tmax, G4int nbin)
static G4PhysicsTable * BuildLabTimeTable(G4PhysicsTable *theDEDXTable, G4PhysicsTable *theLabTimeTable, G4double Tmin, G4double Tmax, G4int nbin)
G4double GetLowEdgeEnergy(size_t binNumber) const
static G4PhysicsTable * theLabTimePositronTable
G4GLOB_DLL std::ostream G4cout
static G4PhysicsTable * theeRangeCoeffATable
bool G4bool
Definition: G4Types.hh:79
static G4PhysicsTable * theDEDXElectronTable
static G4PhysicsTable * BuildRangeCoeffBTable(G4PhysicsTable *theRangeTable, G4PhysicsTable *theCoeffBTable, G4double Tmin, G4double Tmax, G4int nbin)
void PutValue(size_t index, G4double theValue)
static G4PhysicsTable * thepRangeCoeffBTable
static G4PhysicsTable * theRangePositronTable
static G4double LOGRTable
static G4PhysicsTable ** RecorderOfPositronProcess
static G4PhysicsTable * theInverseRangeElectronTable
static G4ProductionCutsTable * GetProductionCutsTable()
static G4Positron * Positron()
Definition: G4Positron.cc:94
static G4double RTable
static G4PhysicsTable * theeRangeCoeffCTable
static G4PhysicsTable ** RecorderOfElectronProcess
static G4PhysicsTable * theLabTimeElectronTable
static G4PhysicsTable * BuildInverseRangeTable(G4PhysicsTable *theRangeTable, G4PhysicsTable *theRangeCoeffATable, G4PhysicsTable *theRangeCoeffBTable, G4PhysicsTable *theRangeCoeffCTable, G4PhysicsTable *theInverseRangeTable, G4double Tmin, G4double Tmax, G4int nbin)
static G4double UpperBoundEloss
static G4PhysicsTable * theProperTimeElectronTable
static G4int NbinEloss
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
static G4int CounterOfPositronProcess
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)
static G4PhysicsTable * theDEDXPositronTable
static const double eplus
Definition: G4SIunits.hh:196
static G4PhysicsTable * BuildRangeCoeffATable(G4PhysicsTable *theRangeTable, G4PhysicsTable *theCoeffATable, G4double Tmin, G4double Tmax, G4int nbin)
G4PhysicsTable ** RecorderOfProcess
G4double GetPDGCharge() const
static G4PhysicsTable * theeRangeCoeffBTable
void clearAndDestroy()
static G4PhysicsTable * BuildRangeCoeffCTable(G4PhysicsTable *theRangeTable, G4PhysicsTable *theCoeffCTable, G4double Tmin, G4double Tmax, G4int nbin)
static G4PhysicsTable * thepRangeCoeffCTable
Here is the call graph for this function:
Here is the caller graph for this function:

◆ DeexciteAtom()

virtual std::vector<G4DynamicParticle*>* G4eLowEnergyLoss::DeexciteAtom ( const G4MaterialCutsCouple ,
G4double  ,
G4double   
)
inlineprotectedvirtual

Reimplemented in G4LowEnergyIonisation.

Definition at line 154 of file G4eLowEnergyLoss.hh.

156  { return 0; };
Here is the caller graph for this function:

◆ Fluorescence()

G4bool G4eLowEnergyLoss::Fluorescence ( ) const
Here is the caller graph for this function:

◆ GetConstraints()

G4double G4eLowEnergyLoss::GetConstraints ( const G4DynamicParticle aParticle,
const G4MaterialCutsCouple couple 
)
private

◆ GetContinuousStepLimit()

G4double G4eLowEnergyLoss::GetContinuousStepLimit ( const G4Track &  track,
G4double  previousStepSize,
G4double  currentMinimumStep,
G4double currentSafety 
)
virtual

Implements G4RDVeLowEnergyLoss.

◆ GetLowerBoundEloss()

G4double G4eLowEnergyLoss::GetLowerBoundEloss ( )
static

Definition at line 192 of file G4eLowEnergyLoss.cc.

193 {
194  return LowerBoundEloss;
195 }
static G4double LowerBoundEloss
Here is the caller graph for this function:

◆ GetMeanFreePath()

virtual G4double G4eLowEnergyLoss::GetMeanFreePath ( const G4Track &  track,
G4double  previousStepSize,
G4ForceCondition *  condition 
)
pure virtual

◆ GetNbinEloss()

G4int G4eLowEnergyLoss::GetNbinEloss ( )
static

Definition at line 202 of file G4eLowEnergyLoss.cc.

203 {
204  return NbinEloss;
205 }
static G4int NbinEloss
Here is the caller graph for this function:

◆ GetNbOfProcesses()

G4int G4eLowEnergyLoss::GetNbOfProcesses ( )
static

Definition at line 172 of file G4eLowEnergyLoss.cc.

173 {
174  return NbOfProcesses;
175 }
static G4int NbOfProcesses

◆ GetUpperBoundEloss()

G4double G4eLowEnergyLoss::GetUpperBoundEloss ( )
static

Definition at line 197 of file G4eLowEnergyLoss.cc.

198 {
199  return UpperBoundEloss;
200 }
static G4double UpperBoundEloss
Here is the caller graph for this function:

◆ IsApplicable()

G4bool G4eLowEnergyLoss::IsApplicable ( const G4ParticleDefinition )
virtual

Reimplemented from G4VProcess.

Reimplemented in G4LowEnergyBremsstrahlung, and G4LowEnergyIonisation.

◆ MinusNbOfProcesses()

void G4eLowEnergyLoss::MinusNbOfProcesses ( )
static

Definition at line 167 of file G4eLowEnergyLoss.cc.

168 {
169  NbOfProcesses--;
170 }
static G4int NbOfProcesses

◆ operator=()

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

◆ PlusNbOfProcesses()

void G4eLowEnergyLoss::PlusNbOfProcesses ( )
static

Definition at line 162 of file G4eLowEnergyLoss.cc.

163 {
164  NbOfProcesses++;
165 }
static G4int NbOfProcesses

◆ PostStepDoIt()

virtual G4VParticleChange* G4eLowEnergyLoss::PostStepDoIt ( const G4Track &  track,
const G4Step &  step 
)
pure virtual

◆ SetLowerBoundEloss()

void G4eLowEnergyLoss::SetLowerBoundEloss ( G4double  val)
static

Definition at line 177 of file G4eLowEnergyLoss.cc.

178 {
179  LowerBoundEloss=val;
180 }
static G4double LowerBoundEloss

◆ SetNbinEloss()

void G4eLowEnergyLoss::SetNbinEloss ( G4int  nb)
static

Definition at line 187 of file G4eLowEnergyLoss.cc.

188 {
189  NbinEloss=nb;
190 }
static G4int NbinEloss

◆ SetNbOfProcesses()

void G4eLowEnergyLoss::SetNbOfProcesses ( G4int  nb)
static

Definition at line 157 of file G4eLowEnergyLoss.cc.

158 {
159  NbOfProcesses=nb;
160 }
static G4int NbOfProcesses

◆ SetUpperBoundEloss()

void G4eLowEnergyLoss::SetUpperBoundEloss ( G4double  val)
static

Definition at line 182 of file G4eLowEnergyLoss.cc.

183 {
184  UpperBoundEloss=val;
185 }
static G4double UpperBoundEloss

Member Data Documentation

◆ Charge

G4double G4eLowEnergyLoss::Charge
protected

Definition at line 164 of file G4eLowEnergyLoss.hh.

◆ CounterOfElectronProcess

G4int G4eLowEnergyLoss::CounterOfElectronProcess = 0
staticprotected

Definition at line 187 of file G4eLowEnergyLoss.hh.

◆ CounterOfPositronProcess

G4int G4eLowEnergyLoss::CounterOfPositronProcess = 0
staticprotected

Definition at line 188 of file G4eLowEnergyLoss.hh.

◆ CounterOfProcess

G4int G4eLowEnergyLoss::CounterOfProcess
private

Definition at line 205 of file G4eLowEnergyLoss.hh.

◆ eLossMessenger

G4EnergyLossMessenger * G4eLowEnergyLoss::eLossMessenger = 0
staticprivate

Definition at line 237 of file G4eLowEnergyLoss.hh.

◆ fdEdx

G4double G4eLowEnergyLoss::fdEdx
private

Definition at line 208 of file G4eLowEnergyLoss.hh.

◆ fParticleChange

G4ParticleChangeForLoss G4eLowEnergyLoss::fParticleChange
private

Definition at line 215 of file G4eLowEnergyLoss.hh.

◆ fRangeNow

G4double G4eLowEnergyLoss::fRangeNow
private

Definition at line 209 of file G4eLowEnergyLoss.hh.

◆ lastCharge

G4double G4eLowEnergyLoss::lastCharge
protected

Definition at line 164 of file G4eLowEnergyLoss.hh.

◆ linLossLimit

G4double G4eLowEnergyLoss::linLossLimit
private

Definition at line 211 of file G4eLowEnergyLoss.hh.

◆ LOGRTable

G4double G4eLowEnergyLoss::LOGRTable
staticprivate

Definition at line 225 of file G4eLowEnergyLoss.hh.

◆ LowerBoundEloss

G4double G4eLowEnergyLoss::LowerBoundEloss = 10.*eV
staticprivate

Definition at line 223 of file G4eLowEnergyLoss.hh.

◆ MinKineticEnergy

G4double G4eLowEnergyLoss::MinKineticEnergy
protected

Definition at line 160 of file G4eLowEnergyLoss.hh.

◆ NbinEloss

G4int G4eLowEnergyLoss::NbinEloss = 360
staticprivate

Definition at line 221 of file G4eLowEnergyLoss.hh.

◆ NbOfProcesses

G4int G4eLowEnergyLoss::NbOfProcesses = 2
staticprotected

Definition at line 186 of file G4eLowEnergyLoss.hh.

◆ RecorderOfElectronProcess

G4PhysicsTable ** G4eLowEnergyLoss::RecorderOfElectronProcess
staticprotected
Initial value:
=
new G4PhysicsTable*[10]

Definition at line 189 of file G4eLowEnergyLoss.hh.

◆ RecorderOfPositronProcess

G4PhysicsTable ** G4eLowEnergyLoss::RecorderOfPositronProcess
staticprotected
Initial value:
=
new G4PhysicsTable*[10]

Definition at line 190 of file G4eLowEnergyLoss.hh.

◆ RecorderOfProcess

G4PhysicsTable** G4eLowEnergyLoss::RecorderOfProcess
private

Definition at line 206 of file G4eLowEnergyLoss.hh.

◆ RTable

G4double G4eLowEnergyLoss::RTable
staticprivate

Definition at line 225 of file G4eLowEnergyLoss.hh.

◆ theDEDXElectronTable

G4PhysicsTable * G4eLowEnergyLoss::theDEDXElectronTable = 0
staticprotected

Definition at line 167 of file G4eLowEnergyLoss.hh.

◆ theDEDXPositronTable

G4PhysicsTable * G4eLowEnergyLoss::theDEDXPositronTable = 0
staticprotected

Definition at line 168 of file G4eLowEnergyLoss.hh.

◆ theDEDXTable

G4PhysicsTable* G4eLowEnergyLoss::theDEDXTable
private

Definition at line 203 of file G4eLowEnergyLoss.hh.

◆ theeRangeCoeffATable

G4PhysicsTable * G4eLowEnergyLoss::theeRangeCoeffATable = 0
staticprivate

Definition at line 230 of file G4eLowEnergyLoss.hh.

◆ theeRangeCoeffBTable

G4PhysicsTable * G4eLowEnergyLoss::theeRangeCoeffBTable = 0
staticprivate

Definition at line 231 of file G4eLowEnergyLoss.hh.

◆ theeRangeCoeffCTable

G4PhysicsTable * G4eLowEnergyLoss::theeRangeCoeffCTable = 0
staticprivate

Definition at line 232 of file G4eLowEnergyLoss.hh.

◆ theFluo

G4bool G4eLowEnergyLoss::theFluo
private

Definition at line 239 of file G4eLowEnergyLoss.hh.

◆ theInverseRangeElectronTable

G4PhysicsTable * G4eLowEnergyLoss::theInverseRangeElectronTable = 0
staticprotected

Definition at line 173 of file G4eLowEnergyLoss.hh.

◆ theInverseRangePositronTable

G4PhysicsTable * G4eLowEnergyLoss::theInverseRangePositronTable = 0
staticprotected

Definition at line 174 of file G4eLowEnergyLoss.hh.

◆ theLabTimeElectronTable

G4PhysicsTable * G4eLowEnergyLoss::theLabTimeElectronTable = 0
staticprotected

Definition at line 177 of file G4eLowEnergyLoss.hh.

◆ theLabTimePositronTable

G4PhysicsTable * G4eLowEnergyLoss::theLabTimePositronTable = 0
staticprotected

Definition at line 178 of file G4eLowEnergyLoss.hh.

◆ theLossTable

G4PhysicsTable* G4eLowEnergyLoss::theLossTable
protected

Definition at line 156 of file G4eLowEnergyLoss.hh.

◆ thepRangeCoeffATable

G4PhysicsTable * G4eLowEnergyLoss::thepRangeCoeffATable = 0
staticprivate

Definition at line 233 of file G4eLowEnergyLoss.hh.

◆ thepRangeCoeffBTable

G4PhysicsTable * G4eLowEnergyLoss::thepRangeCoeffBTable = 0
staticprivate

Definition at line 234 of file G4eLowEnergyLoss.hh.

◆ thepRangeCoeffCTable

G4PhysicsTable * G4eLowEnergyLoss::thepRangeCoeffCTable = 0
staticprivate

Definition at line 235 of file G4eLowEnergyLoss.hh.

◆ theProperTimeElectronTable

G4PhysicsTable * G4eLowEnergyLoss::theProperTimeElectronTable = 0
staticprotected

Definition at line 179 of file G4eLowEnergyLoss.hh.

◆ theProperTimePositronTable

G4PhysicsTable * G4eLowEnergyLoss::theProperTimePositronTable = 0
staticprotected

Definition at line 180 of file G4eLowEnergyLoss.hh.

◆ theRangeElectronTable

G4PhysicsTable * G4eLowEnergyLoss::theRangeElectronTable = 0
staticprotected

Definition at line 169 of file G4eLowEnergyLoss.hh.

◆ theRangePositronTable

G4PhysicsTable * G4eLowEnergyLoss::theRangePositronTable = 0
staticprotected

Definition at line 170 of file G4eLowEnergyLoss.hh.

◆ UpperBoundEloss

G4double G4eLowEnergyLoss::UpperBoundEloss = 100.*GeV
staticprivate

Definition at line 224 of file G4eLowEnergyLoss.hh.


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