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

#include <G4ErrorStepLengthLimitProcess.hh>

Inheritance diagram for G4ErrorStepLengthLimitProcess:

Public Member Functions

 G4ErrorStepLengthLimitProcess (const G4String &processName="G4ErrorStepLengthLimit")
 ~G4ErrorStepLengthLimitProcess ()
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
Public Member Functions inherited from G4VErrorLimitProcess
 G4VErrorLimitProcess (const G4String &processName)
 ~G4VErrorLimitProcess ()
virtual G4double GetMeanFreePath (const class G4Track &, G4double, enum G4ForceCondition *)
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
G4double GetStepLimit () const
void SetStepLimit (G4double val)
Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &aName, G4ProcessType aType=fNotDefined)
 G4VDiscreteProcess (G4VDiscreteProcess &)
virtual ~G4VDiscreteProcess ()
G4VDiscreteProcessoperator= (const G4VDiscreteProcess &)=delete
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
virtual G4double GetCrossSection (const G4double, const G4MaterialCutsCouple *)
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *)
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
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 const G4VProcessGetCreatorProcess () const
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 G4VDiscreteProcess
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
void ClearNumberOfInteractionLengthLeft ()
Protected Attributes inherited from G4VErrorLimitProcess
G4double theStepLimit
G4double theStepLength
G4VParticleChange theParticleChange
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 54 of file G4ErrorStepLengthLimitProcess.hh.

Constructor & Destructor Documentation

◆ G4ErrorStepLengthLimitProcess()

G4ErrorStepLengthLimitProcess::G4ErrorStepLengthLimitProcess ( const G4String & processName = "G4ErrorStepLengthLimit")

Definition at line 41 of file G4ErrorStepLengthLimitProcess.cc.

43 : G4VErrorLimitProcess(processName)
44{
45 theStepLimit = 1000. * mm; // kInfinity;
46}
G4VErrorLimitProcess(const G4String &processName)

◆ ~G4ErrorStepLengthLimitProcess()

G4ErrorStepLengthLimitProcess::~G4ErrorStepLengthLimitProcess ( )

Definition at line 49 of file G4ErrorStepLengthLimitProcess.cc.

49{}

Member Function Documentation

◆ PostStepGetPhysicalInteractionLength()

G4double G4ErrorStepLengthLimitProcess::PostStepGetPhysicalInteractionLength ( const G4Track & track,
G4double previousStepSize,
G4ForceCondition * condition )
virtual

Implements G4VErrorLimitProcess.

Definition at line 52 of file G4ErrorStepLengthLimitProcess.cc.

54{
56
57#ifdef G4VERBOSE
59 {
60 G4cout
61 << "G4ErrorStepLengthLimitProcess::PostStepGetPhysicalInteractionLength "
62 << theStepLimit << G4endl;
63 }
64#endif
65
66 return theStepLimit;
67}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout

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