Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4StepLimiter Class Reference

#include <G4StepLimiter.hh>

+ Inheritance diagram for G4StepLimiter:

Public Member Functions

 G4StepLimiter (const G4String &processName="StepLimiter")
 
virtual ~G4StepLimiter ()
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual 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 ()
 
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 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 &)
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- 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
 

Detailed Description

Definition at line 44 of file G4StepLimiter.hh.

Constructor & Destructor Documentation

◆ G4StepLimiter()

G4StepLimiter::G4StepLimiter ( const G4String processName = "StepLimiter")

Definition at line 43 of file G4StepLimiter.cc.

44 : G4VProcess(aName, fGeneral )
45{
46 // set Process Sub Type
47 SetProcessSubType(static_cast<int>(STEP_LIMITER));
48
49 if (verboseLevel>1) {
50 G4cout << GetProcessName() << " is created "<< G4endl;
51 }
52}
@ fGeneral
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4int verboseLevel
Definition: G4VProcess.hh:356
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:406
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382

◆ ~G4StepLimiter()

G4StepLimiter::~G4StepLimiter ( )
virtual

Definition at line 56 of file G4StepLimiter.cc.

57{
58}

Member Function Documentation

◆ AlongStepDoIt()

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

Implements G4VProcess.

Definition at line 87 of file G4StepLimiter.hh.

90 {return 0;};

◆ AlongStepGetPhysicalInteractionLength()

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

Implements G4VProcess.

Definition at line 78 of file G4StepLimiter.hh.

84 { return -1.0; };

◆ AtRestDoIt()

virtual G4VParticleChange * G4StepLimiter::AtRestDoIt ( const G4Track ,
const G4Step  
)
inlinevirtual

Implements G4VProcess.

Definition at line 72 of file G4StepLimiter.hh.

75 {return 0;};

◆ AtRestGetPhysicalInteractionLength()

virtual G4double G4StepLimiter::AtRestGetPhysicalInteractionLength ( const G4Track ,
G4ForceCondition  
)
inlinevirtual

Implements G4VProcess.

Definition at line 66 of file G4StepLimiter.hh.

69 { return -1.0; };

◆ PostStepDoIt()

G4VParticleChange * G4StepLimiter::PostStepDoIt ( const G4Track aTrack,
const G4Step  
)
virtual

Implements G4VProcess.

Definition at line 92 of file G4StepLimiter.cc.

96{
98 return &aParticleChange;
99}
virtual void Initialize(const G4Track &)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:327

◆ PostStepGetPhysicalInteractionLength()

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

Implements G4VProcess.

Definition at line 70 of file G4StepLimiter.cc.

74{
75 // condition is set to "Not Forced"
77
78 G4double proposedStep = DBL_MAX;
79 G4UserLimits* pUserLimits =
80 aTrack.GetVolume()->GetLogicalVolume()->GetUserLimits();
81 if (pUserLimits) {
82 // max allowed step length
83 //
84 proposedStep = pUserLimits->GetMaxAllowedStep(aTrack);
85 if (proposedStep < 0.) proposedStep =0.;
86 }
87 return proposedStep;
88}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced
double G4double
Definition: G4Types.hh:83
virtual G4double GetMaxAllowedStep(const G4Track &)
#define DBL_MAX
Definition: templates.hh:62

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