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

#include <G4VContinuousDiscreteProcess.hh>

+ Inheritance diagram for G4VContinuousDiscreteProcess:

Public Member Functions

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

Constructor & Destructor Documentation

◆ G4VContinuousDiscreteProcess() [1/2]

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

Definition at line 50 of file G4VContinuousDiscreteProcess.cc.

52 : G4VProcess(aName, aType)
53{
54 enableAtRestDoIt = false;
55}
G4bool enableAtRestDoIt

◆ G4VContinuousDiscreteProcess() [2/2]

G4VContinuousDiscreteProcess::G4VContinuousDiscreteProcess ( G4VContinuousDiscreteProcess & right)

Definition at line 63 of file G4VContinuousDiscreteProcess.cc.

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

◆ ~G4VContinuousDiscreteProcess()

G4VContinuousDiscreteProcess::~G4VContinuousDiscreteProcess ( )
virtual

Definition at line 58 of file G4VContinuousDiscreteProcess.cc.

59{
60}

Member Function Documentation

◆ AlongStepDoIt()

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

Implements G4VProcess.

Reimplemented in G4AdjointForcedInteractionForGamma, G4hImpactIonisation, G4VEnergyLossProcess, and G4VMultipleScattering.

Definition at line 133 of file G4VContinuousDiscreteProcess.cc.

135{
136 // clear NumberOfInteractionLengthLeft
138
139 return pParticleChange;
140}
void ClearNumberOfInteractionLengthLeft()
G4VParticleChange * pParticleChange

◆ AlongStepGetPhysicalInteractionLength()

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

Implements G4VProcess.

Reimplemented in G4VEnergyLossProcess, and G4VMultipleScattering.

Definition at line 143 of file G4VContinuousDiscreteProcess.cc.

149{
150 // GPILSelection is set to defaule value of CandidateForSelection
151 valueGPILSelection = CandidateForSelection;
152
153 // get Step limit proposed by the process
154 G4double steplength = GetContinuousStepLimit(track,previousStepSize,currentMinimumStep,currentSafety);
155
156 // set return value for G4GPILSelection
157 *selection = valueGPILSelection;
158
159#ifdef G4VERBOSE
160 if (verboseLevel>1)
161 {
162 G4cout << "G4VContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength() - ";
163 G4cout << "[ " << GetProcessName() << "]" << G4endl;
164 track.GetDynamicParticle()->DumpInfo();
165 G4cout << " in Material " << track.GetMaterial()->GetName() << G4endl;
166 G4cout << "IntractionLength= " << steplength/cm <<"[cm] " << G4endl;
167 }
168#endif
169 return steplength;
170}
@ 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
virtual G4double GetContinuousStepLimit(const G4Track &aTrack, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)=0
G4int verboseLevel
const G4String & GetProcessName() const

◆ AtRestDoIt()

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

Implements G4VProcess.

Definition at line 87 of file G4VContinuousDiscreteProcess.hh.

90 { return 0; }

◆ AtRestGetPhysicalInteractionLength()

virtual G4double G4VContinuousDiscreteProcess::AtRestGetPhysicalInteractionLength ( const G4Track & ,
G4ForceCondition *  )
inlinevirtual

Implements G4VProcess.

Definition at line 81 of file G4VContinuousDiscreteProcess.hh.

84 { return -1.0; }

◆ GetContinuousStepLimit()

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

◆ GetGPILSelection()

G4GPILSelection G4VContinuousDiscreteProcess::GetGPILSelection ( ) const
inlineprotected

Definition at line 109 of file G4VContinuousDiscreteProcess.hh.

110 { return valueGPILSelection; }

◆ GetMeanFreePath()

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

◆ operator=()

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

◆ PostStepDoIt()

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

Implements G4VProcess.

Reimplemented in G4AdjointForcedInteractionForGamma, G4hImpactIonisation, G4hRDEnergyLoss, and G4VEnergyLossProcess.

Definition at line 122 of file G4VContinuousDiscreteProcess.cc.

124{
125 // clear NumberOfInteractionLengthLeft
127
128 return pParticleChange;
129}

Referenced by G4hImpactIonisation::PostStepDoIt().

◆ PostStepGetPhysicalInteractionLength()

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

Implements G4VProcess.

Reimplemented in G4AdjointForcedInteractionForGamma, G4PolarizedIonisation, G4VEnergyLossProcess, and G4VMultipleScattering.

Definition at line 71 of file G4VContinuousDiscreteProcess.cc.

75{
76 if ( (previousStepSize <=0.0) || (theNumberOfInteractionLengthLeft<=0.0))
77 {
78 // beginning of tracking (or just after DoIt() of this process)
80 }
81 else if ( previousStepSize > 0.0)
82 {
83 // subtract NumberOfInteractionLengthLeft
85 }
86 else
87 {
88 // zero step
89 // DO NOTHING
90 }
91
92 // condition is set to "Not Forced"
94
95 // get mean free path
96 currentInteractionLength = GetMeanFreePath(track,previousStepSize,condition);
97
98 G4double value;
100 {
102 }
103 else
104 {
105 value = DBL_MAX;
106 }
107#ifdef G4VERBOSE
108 if (verboseLevel>1)
109 {
110 G4cout << "G4VContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength() - ";
111 G4cout << "[ " << GetProcessName() << "]" << G4endl;
112 track.GetDynamicParticle()->DumpInfo();
113 G4cout << " in Material " << track.GetMaterial()->GetName() << G4endl;
114 G4cout << "InteractionLength= " << value/cm <<"[cm] " << G4endl;
115 }
116#endif
117 return value;
118}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
G4double currentInteractionLength
void SubtractNumberOfInteractionLengthLeft(G4double prevStepSize)
virtual void ResetNumberOfInteractionLengthLeft()
Definition G4VProcess.cc:80
G4double theNumberOfInteractionLengthLeft
#define DBL_MAX
Definition templates.hh:62

◆ SetGPILSelection()

void G4VContinuousDiscreteProcess::SetGPILSelection ( G4GPILSelection selection)
inlineprotected

Definition at line 106 of file G4VContinuousDiscreteProcess.hh.

107 { valueGPILSelection = selection; }

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