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

#include <G4VRestContinuousProcess.hh>

+ Inheritance diagram for G4VRestContinuousProcess:

Public Member Functions

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

Protected Member Functions

virtual G4double GetContinuousStepLimit (const G4Track &aTrack, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)=0
 
virtual G4double GetMeanLifeTime (const G4Track &aTrack, G4ForceCondition *condition)=0
 
void SetGPILSelection (G4GPILSelection selection)
 
G4GPILSelection GetGPILSelection () const
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- 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 45 of file G4VRestContinuousProcess.hh.

Constructor & Destructor Documentation

◆ G4VRestContinuousProcess() [1/2]

G4VRestContinuousProcess::G4VRestContinuousProcess ( const G4String & aName,
G4ProcessType aType = fNotDefined )

Definition at line 50 of file G4VRestContinuousProcess.cc.

52 : G4VProcess(aName, aType), valueGPILSelection(CandidateForSelection)
53{
54 enablePostStepDoIt = false;
55}
@ CandidateForSelection
G4bool enablePostStepDoIt

◆ G4VRestContinuousProcess() [2/2]

G4VRestContinuousProcess::G4VRestContinuousProcess ( G4VRestContinuousProcess & right)

Definition at line 63 of file G4VRestContinuousProcess.cc.

65 : G4VProcess(right), valueGPILSelection(right.valueGPILSelection)
66{
67}

◆ ~G4VRestContinuousProcess()

G4VRestContinuousProcess::~G4VRestContinuousProcess ( )
virtual

Definition at line 58 of file G4VRestContinuousProcess.cc.

59{
60}

Member Function Documentation

◆ AlongStepDoIt()

G4VParticleChange * G4VRestContinuousProcess::AlongStepDoIt ( const G4Track & ,
const G4Step &  )
virtual

Implements G4VProcess.

Definition at line 140 of file G4VRestContinuousProcess.cc.

142{
143 return pParticleChange;
144}
G4VParticleChange * pParticleChange

◆ AlongStepGetPhysicalInteractionLength()

G4double G4VRestContinuousProcess::AlongStepGetPhysicalInteractionLength ( const G4Track & track,
G4double previousStepSize,
G4double currentMinimumStep,
G4double & currentSafety,
G4GPILSelection * selection )
virtual

Implements G4VProcess.

Definition at line 70 of file G4VRestContinuousProcess.cc.

76{
77 // GPILSelection is set to defaule value of CandidateForSelection
78 valueGPILSelection = CandidateForSelection;
79
80 // get Step limit proposed by the process
81 G4double steplength = GetContinuousStepLimit(track, previousStepSize,
82 currentMinimumStep, currentSafety);
83
84 // set return value for G4GPILSelection
85 *selection = valueGPILSelection;
86#ifdef G4VERBOSE
87 if (verboseLevel>1)
88 {
89 G4cout << "G4VRestContinuousProcess::AlongStepGetPhysicalInteractionLength() - ";
90 G4cout << "[ " << GetProcessName() << "]" << G4endl;
92 G4cout << " in Material " << track.GetMaterial()->GetName() << G4endl;
93 G4cout << "IntractionLength= " << steplength/cm <<"[cm] " << G4endl;
94 }
95#endif
96 return steplength ;
97}
double G4double
Definition G4Types.hh:83
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
void DumpInfo(G4int mode=0) const
const G4String & GetName() const
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
G4int verboseLevel
const G4String & GetProcessName() const
virtual G4double GetContinuousStepLimit(const G4Track &aTrack, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)=0

◆ AtRestDoIt()

G4VParticleChange * G4VRestContinuousProcess::AtRestDoIt ( const G4Track & ,
const G4Step &  )
virtual

Implements G4VProcess.

Definition at line 130 of file G4VRestContinuousProcess.cc.

132{
134
135 return pParticleChange;
136}
void ClearNumberOfInteractionLengthLeft()

◆ AtRestGetPhysicalInteractionLength()

G4double G4VRestContinuousProcess::AtRestGetPhysicalInteractionLength ( const G4Track & track,
G4ForceCondition * condition )
virtual

Implements G4VProcess.

Definition at line 100 of file G4VRestContinuousProcess.cc.

103{
104 // beginning of tracking
106
107 // condition is set to "Not Forced"
109
110 // get mean life time
112
113#ifdef G4VERBOSE
114 if ((currentInteractionLength <0.0) || (verboseLevel>2))
115 {
116 G4cout << "G4VRestContinuousProcess::AtRestGetPhysicalInteractionLength() - ";
117 G4cout << "[ " << GetProcessName() << "]" << G4endl;
118 track.GetDynamicParticle()->DumpInfo();
119 G4cout << " in Material " << track.GetMaterial()->GetName() << G4endl;
120 G4cout << "MeanLifeTime = " << currentInteractionLength/ns << "[ns]"
121 << G4endl;
122 }
123#endif
124
126}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced
G4double currentInteractionLength
virtual void ResetNumberOfInteractionLengthLeft()
Definition G4VProcess.cc:80
G4double theNumberOfInteractionLengthLeft
virtual G4double GetMeanLifeTime(const G4Track &aTrack, G4ForceCondition *condition)=0

◆ GetContinuousStepLimit()

virtual G4double G4VRestContinuousProcess::GetContinuousStepLimit ( const G4Track & aTrack,
G4double previousStepSize,
G4double currentMinimumStep,
G4double & currentSafety )
protectedpure virtual

◆ GetGPILSelection()

G4GPILSelection G4VRestContinuousProcess::GetGPILSelection ( ) const
inlineprotected

Definition at line 112 of file G4VRestContinuousProcess.hh.

113 { return valueGPILSelection; }

◆ GetMeanLifeTime()

virtual G4double G4VRestContinuousProcess::GetMeanLifeTime ( const G4Track & aTrack,
G4ForceCondition * condition )
protectedpure virtual

◆ operator=()

G4VRestContinuousProcess & G4VRestContinuousProcess::operator= ( const G4VRestContinuousProcess & )
delete

◆ PostStepDoIt()

virtual G4VParticleChange * G4VRestContinuousProcess::PostStepDoIt ( const G4Track & ,
const G4Step &  )
inlinevirtual

Implements G4VProcess.

Definition at line 87 of file G4VRestContinuousProcess.hh.

90 { return 0; }

◆ PostStepGetPhysicalInteractionLength()

virtual G4double G4VRestContinuousProcess::PostStepGetPhysicalInteractionLength ( const G4Track & ,
G4double ,
G4ForceCondition *  )
inlinevirtual

Implements G4VProcess.

Definition at line 80 of file G4VRestContinuousProcess.hh.

84 { return -1.0; }

◆ SetGPILSelection()

void G4VRestContinuousProcess::SetGPILSelection ( G4GPILSelection selection)
inlineprotected

Definition at line 109 of file G4VRestContinuousProcess.hh.

110 { valueGPILSelection = selection; }

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