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

#include <G4IVRestDiscreteProcess.hh>

+ Inheritance diagram for G4IVRestDiscreteProcess:

Public Member Functions

 G4IVRestDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4IVRestDiscreteProcess (G4IVRestDiscreteProcess &)
 
virtual ~G4IVRestDiscreteProcess ()
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4VParticleChangeAlongStepDoIt (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
 

Protected Member Functions

virtual void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
virtual G4double GetMeanLifeTime (const G4Track &aTrack, G4ForceCondition *condition)=0
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Protected Attributes

G4PhysicsTabletheNlambdaTable
 
G4PhysicsTabletheInverseNlambdaTable
 
const G4double BIGSTEP
 
- 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
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 

Detailed Description

Definition at line 59 of file G4IVRestDiscreteProcess.hh.

Constructor & Destructor Documentation

◆ G4IVRestDiscreteProcess() [1/2]

G4IVRestDiscreteProcess::G4IVRestDiscreteProcess ( const G4String aName,
G4ProcessType  aType = fNotDefined 
)

Definition at line 49 of file G4IVRestDiscreteProcess.cc.

50 : G4VProcess(aName, aType),
52 BIGSTEP(1.e10)
53{
54 enableAlongStepDoIt = false;
55}
G4bool enableAlongStepDoIt
Definition: G4VProcess.hh:351

Referenced by G4IVRestDiscreteProcess().

◆ G4IVRestDiscreteProcess() [2/2]

G4IVRestDiscreteProcess::G4IVRestDiscreteProcess ( G4IVRestDiscreteProcess right)

Definition at line 61 of file G4IVRestDiscreteProcess.cc.

62 : G4VProcess(right),
64 BIGSTEP(1.e10)
65{
66}

◆ ~G4IVRestDiscreteProcess()

G4IVRestDiscreteProcess::~G4IVRestDiscreteProcess ( )
virtual

Definition at line 57 of file G4IVRestDiscreteProcess.cc.

58{
59}

Member Function Documentation

◆ AlongStepDoIt()

virtual G4VParticleChange * G4IVRestDiscreteProcess::AlongStepDoIt ( const G4Track ,
const G4Step  
)
inlinevirtual

Implements G4VProcess.

Definition at line 107 of file G4IVRestDiscreteProcess.hh.

110 {return 0;}

◆ AlongStepGetPhysicalInteractionLength()

virtual G4double G4IVRestDiscreteProcess::AlongStepGetPhysicalInteractionLength ( const G4Track ,
G4double  ,
G4double  ,
G4double ,
G4GPILSelection  
)
inlinevirtual

Implements G4VProcess.

Definition at line 98 of file G4IVRestDiscreteProcess.hh.

104 { return -1.0; }

◆ AtRestDoIt()

G4VParticleChange * G4IVRestDiscreteProcess::AtRestDoIt ( const G4Track ,
const G4Step  
)
inlinevirtual

Implements G4VProcess.

Definition at line 187 of file G4IVRestDiscreteProcess.hh.

191{
192// clear NumberOfInteractionLengthLeft
194
195 return pParticleChange;
196}
void ClearNumberOfInteractionLengthLeft()
Definition: G4VProcess.hh:418
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283

◆ AtRestGetPhysicalInteractionLength()

G4double G4IVRestDiscreteProcess::AtRestGetPhysicalInteractionLength ( const G4Track track,
G4ForceCondition condition 
)
inlinevirtual

Implements G4VProcess.

Definition at line 159 of file G4IVRestDiscreteProcess.hh.

163{
164 // beggining of tracking
166
167 // condition is set to "Not Forced"
169
170 // get mean life time
172
173#ifdef G4VERBOSE
174 if ((currentInteractionLength <0.0) || (verboseLevel>2)){
175 G4cout << "G4IVRestDiscreteProcess::AtRestGetPhysicalInteractionLength ";
176 G4cout << "[ " << GetProcessName() << "]" <<G4endl;
177 track.GetDynamicParticle()->DumpInfo();
178 G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
179 G4cout << "MeanLifeTime = " << currentInteractionLength/CLHEP::ns << "[ns]" <<G4endl;
180 }
181#endif
182
184}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
void DumpInfo(G4int mode=0) const
virtual G4double GetMeanLifeTime(const G4Track &aTrack, G4ForceCondition *condition)=0
const G4String & GetName() const
Definition: G4Material.hh:177
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
G4double currentInteractionLength
Definition: G4VProcess.hh:297
virtual void ResetNumberOfInteractionLengthLeft()
Definition: G4VProcess.cc:92
G4int verboseLevel
Definition: G4VProcess.hh:368
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:293
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379

◆ GetMeanLifeTime()

virtual G4double G4IVRestDiscreteProcess::GetMeanLifeTime ( const G4Track aTrack,
G4ForceCondition condition 
)
protectedpure virtual

◆ PostStepDoIt()

G4VParticleChange * G4IVRestDiscreteProcess::PostStepDoIt ( const G4Track ,
const G4Step  
)
inlinevirtual

Implements G4VProcess.

Definition at line 148 of file G4IVRestDiscreteProcess.hh.

152{
153// reset NumberOfInteractionLengthLeft
155
156 return pParticleChange;
157}

◆ PostStepGetPhysicalInteractionLength()

G4double G4IVRestDiscreteProcess::PostStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
virtual

Implements G4VProcess.

Definition at line 68 of file G4IVRestDiscreteProcess.cc.

73{
74 G4double value = DBL_MAX ;
75
76 return value;
77}
double G4double
Definition: G4Types.hh:64
#define DBL_MAX
Definition: templates.hh:83

◆ SubtractNumberOfInteractionLengthLeft()

void G4IVRestDiscreteProcess::SubtractNumberOfInteractionLengthLeft ( G4double  previousStepSize)
inlineprotectedvirtual

Definition at line 141 of file G4IVRestDiscreteProcess.hh.

143 {
144 // dummy routine
145 ;
146 }

Member Data Documentation

◆ BIGSTEP

const G4double G4IVRestDiscreteProcess::BIGSTEP
protected

Definition at line 133 of file G4IVRestDiscreteProcess.hh.

◆ theInverseNlambdaTable

G4PhysicsTable* G4IVRestDiscreteProcess::theInverseNlambdaTable
protected

Definition at line 131 of file G4IVRestDiscreteProcess.hh.

◆ theNlambdaTable

G4PhysicsTable* G4IVRestDiscreteProcess::theNlambdaTable
protected

Definition at line 130 of file G4IVRestDiscreteProcess.hh.


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