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

#include <G4VRestContinuousDiscreteProcess.hh>

+ Inheritance diagram for G4VRestContinuousDiscreteProcess:

Public Member Functions

 G4VRestContinuousDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VRestContinuousDiscreteProcess (G4VRestContinuousDiscreteProcess &)
 
virtual ~G4VRestContinuousDiscreteProcess ()
 
G4VRestContinuousDiscreteProcessoperator= (const G4VRestContinuousDiscreteProcess &)=delete
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (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 AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (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 GetMeanLifeTime (const G4Track &aTrack, G4ForceCondition *condition)=0
 
virtual G4double GetContinuousStepLimit (const G4Track &aTrack, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)=0
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, 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 G4VRestContinuousDiscreteProcess.hh.

Constructor & Destructor Documentation

◆ G4VRestContinuousDiscreteProcess() [1/2]

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

Definition at line 50 of file G4VRestContinuousDiscreteProcess.cc.

52 : G4VProcess(aName, aType)
53{
54}

◆ G4VRestContinuousDiscreteProcess() [2/2]

G4VRestContinuousDiscreteProcess::G4VRestContinuousDiscreteProcess ( G4VRestContinuousDiscreteProcess & right)

Definition at line 62 of file G4VRestContinuousDiscreteProcess.cc.

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

◆ ~G4VRestContinuousDiscreteProcess()

G4VRestContinuousDiscreteProcess::~G4VRestContinuousDiscreteProcess ( )
virtual

Definition at line 57 of file G4VRestContinuousDiscreteProcess.cc.

58{
59}

Member Function Documentation

◆ AlongStepDoIt()

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

Implements G4VProcess.

Definition at line 146 of file G4VRestContinuousDiscreteProcess.cc.

148{
149 return pParticleChange;
150}
G4VParticleChange * pParticleChange

◆ AlongStepGetPhysicalInteractionLength()

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

Implements G4VProcess.

Definition at line 114 of file G4VRestContinuousDiscreteProcess.cc.

120{
121 // GPILSelection is set to defaule value of CandidateForSelection
122 valueGPILSelection = CandidateForSelection;
123
124 // get Step limit proposed by the process
125 G4double steplength = GetContinuousStepLimit(track, previousStepSize,
126 currentMinimumStep, currentSafety);
127
128 // set return value for G4GPILSelection
129 *selection = valueGPILSelection;
130
131#ifdef G4VERBOSE
132 if (verboseLevel>1)
133 {
134 G4cout << "G4VRestContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength() - ";
135 G4cout << "[ " << GetProcessName() << "]" << G4endl;
136 track.GetDynamicParticle()->DumpInfo();
137 G4cout << " in Material " << track.GetMaterial()->GetName() << G4endl;
138 G4cout << "IntractionLength= " << steplength/cm <<"[cm] " << G4endl;
139 }
140#endif
141 return steplength;
142}
@ CandidateForSelection
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 * G4VRestContinuousDiscreteProcess::AtRestDoIt ( const G4Track & ,
const G4Step &  )
virtual

Implements G4VProcess.

Definition at line 105 of file G4VRestContinuousDiscreteProcess.cc.

107{
109
110 return pParticleChange;
111}
void ClearNumberOfInteractionLengthLeft()

◆ AtRestGetPhysicalInteractionLength()

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

Implements G4VProcess.

Definition at line 70 of file G4VRestContinuousDiscreteProcess.cc.

73{
74 // beginning of tracking
76
77 // condition is set to "Not Forced"
79
80 // get mean life time
82
85
86#ifdef G4VERBOSE
87 if ((currentInteractionLength <0.0) || (verboseLevel>2))
88 {
91 G4cout << "G4VRestContinuousDiscreteProcess::AtRestGetPhysicalInteractionLength() - ";
92 G4cout << "[ " << GetProcessName() << "]" << G4endl;
94 G4cout << " in Material " << track.GetMaterial()->GetName() << G4endl;
95 G4cout << "MeanLifeTime = " << t << " [ns]"
96 << G4endl;
97 }
98#endif
99
100 return time;
101}
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
#define DBL_MAX
Definition templates.hh:62

◆ GetContinuousStepLimit()

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

◆ GetGPILSelection()

G4GPILSelection G4VRestContinuousDiscreteProcess::GetGPILSelection ( ) const
inlineprotected

Definition at line 115 of file G4VRestContinuousDiscreteProcess.hh.

116 { return valueGPILSelection; }

◆ GetMeanFreePath()

virtual G4double G4VRestContinuousDiscreteProcess::GetMeanFreePath ( const G4Track & aTrack,
G4double previousStepSize,
G4ForceCondition * condition )
protectedpure virtual

◆ GetMeanLifeTime()

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

◆ operator=()

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

◆ PostStepDoIt()

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

Implements G4VProcess.

Definition at line 205 of file G4VRestContinuousDiscreteProcess.cc.

207{
209
210 return pParticleChange;
211}

◆ PostStepGetPhysicalInteractionLength()

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

Implements G4VProcess.

Definition at line 153 of file G4VRestContinuousDiscreteProcess.cc.

157{
158 if ( (previousStepSize < 0.0) || (theNumberOfInteractionLengthLeft<=0.0))
159 {
160 // beginning of tracking (or just after DoIt() of this process)
162 }
163 else if ( previousStepSize > 0.0)
164 {
165 // subtract NumberOfInteractionLengthLeft
167 }
168 else
169 {
170 // zero step
171 // DO NOTHING
172 }
173
174 // condition is set to "Not Forced"
176
177 // get mean free path
178 currentInteractionLength = GetMeanFreePath(track,previousStepSize,condition);
179
180
181 G4double value;
183 {
185 }
186 else
187 {
188 value = DBL_MAX;
189 }
190#ifdef G4VERBOSE
191 if (verboseLevel>1)
192 {
193 G4cout << "G4VRestContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength() - ";
194 G4cout << "[ " << GetProcessName() << "]" << G4endl;
195 track.GetDynamicParticle()->DumpInfo();
196 G4cout << " in Material " << track.GetMaterial()->GetName() << G4endl;
197 G4cout << "InteractionLength= " << value/cm <<"[cm] " << G4endl;
198 }
199#endif
200 return value;
201}
void SubtractNumberOfInteractionLengthLeft(G4double prevStepSize)
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0

◆ SetGPILSelection()

void G4VRestContinuousDiscreteProcess::SetGPILSelection ( G4GPILSelection selection)
inlineprotected

Definition at line 112 of file G4VRestContinuousDiscreteProcess.hh.

113 { valueGPILSelection = selection; }

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