Geant4 9.6.0
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 ()
 
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 ()
 
G4int operator== (const G4VProcess &right) const
 
G4int 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
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 

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 previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChangepParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 

Detailed Description

Definition at line 55 of file G4VContinuousDiscreteProcess.hh.

Constructor & Destructor Documentation

◆ G4VContinuousDiscreteProcess() [1/2]

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

Definition at line 55 of file G4VContinuousDiscreteProcess.cc.

56 : G4VProcess(aName, aType),
57 valueGPILSelection(CandidateForSelection)
58{
59 enableAtRestDoIt = false;
60}
@ CandidateForSelection
G4bool enableAtRestDoIt
Definition: G4VProcess.hh:350

Referenced by G4VContinuousDiscreteProcess().

◆ G4VContinuousDiscreteProcess() [2/2]

G4VContinuousDiscreteProcess::G4VContinuousDiscreteProcess ( G4VContinuousDiscreteProcess right)

Definition at line 66 of file G4VContinuousDiscreteProcess.cc.

67 : G4VProcess(right),
68 valueGPILSelection(right.valueGPILSelection)
69{
70}

◆ ~G4VContinuousDiscreteProcess()

G4VContinuousDiscreteProcess::~G4VContinuousDiscreteProcess ( )
virtual

Definition at line 62 of file G4VContinuousDiscreteProcess.cc.

63{
64}

Member Function Documentation

◆ AlongStepDoIt()

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

Implements G4VProcess.

Reimplemented in G4VEnergyLossProcess, G4VMultipleScattering, G4VeLowEnergyLoss, G4VEnergyLoss, and G4hImpactIonisation.

Definition at line 124 of file G4VContinuousDiscreteProcess.cc.

128{
129// clear NumberOfInteractionLengthLeft
131 return pParticleChange;
132}
void ClearNumberOfInteractionLengthLeft()
Definition: G4VProcess.hh:418
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283

◆ AlongStepGetPhysicalInteractionLength()

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

Implements G4VProcess.

Reimplemented in G4VMultipleScattering, and G4VEnergyLossProcess.

Definition at line 134 of file G4VContinuousDiscreteProcess.cc.

141{
142 // GPILSelection is set to defaule value of CandidateForSelection
143 valueGPILSelection = CandidateForSelection;
144
145 // get Step limit proposed by the process
146 G4double steplength = GetContinuousStepLimit(track,previousStepSize,currentMinimumStep, currentSafety);
147
148 // set return value for G4GPILSelection
149 *selection = valueGPILSelection;
150
151#ifdef G4VERBOSE
152 if (verboseLevel>1){
153 G4cout << "G4VContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength ";
154 G4cout << "[ " << GetProcessName() << "]" <<G4endl;
155 track.GetDynamicParticle()->DumpInfo();
156 G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
157 G4cout << "IntractionLength= " << steplength/cm <<"[cm] " <<G4endl;
158 }
159#endif
160 return steplength ;
161}
double G4double
Definition: G4Types.hh:64
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
void DumpInfo(G4int mode=0) const
const G4String & GetName() const
Definition: G4Material.hh:177
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
virtual G4double GetContinuousStepLimit(const G4Track &aTrack, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)=0
G4int verboseLevel
Definition: G4VProcess.hh:368
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379

◆ AtRestDoIt()

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

Implements G4VProcess.

Definition at line 99 of file G4VContinuousDiscreteProcess.hh.

102 {return 0;};

◆ AtRestGetPhysicalInteractionLength()

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

Implements G4VProcess.

Definition at line 93 of file G4VContinuousDiscreteProcess.hh.

96 { 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 128 of file G4VContinuousDiscreteProcess.hh.

128{return valueGPILSelection;};

◆ GetMeanFreePath()

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

◆ PostStepDoIt()

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

Implements G4VProcess.

Reimplemented in G4VEnergyLossProcess, G4VMultipleScattering, G4hImpactIonisation, G4VeLowEnergyLoss, G4hRDEnergyLoss, and G4VEnergyLoss.

Definition at line 114 of file G4VContinuousDiscreteProcess.cc.

118{
119// clear NumberOfInteractionLengthLeft
121 return pParticleChange;
122}

Referenced by G4hImpactIonisation::PostStepDoIt().

◆ PostStepGetPhysicalInteractionLength()

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

Implements G4VProcess.

Reimplemented in G4VMultipleScattering, G4ePolarizedIonisation, and G4VEnergyLossProcess.

Definition at line 72 of file G4VContinuousDiscreteProcess.cc.

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

◆ SetGPILSelection()

void G4VContinuousDiscreteProcess::SetGPILSelection ( G4GPILSelection  selection)
inlineprotected

Definition at line 125 of file G4VContinuousDiscreteProcess.hh.

126 { valueGPILSelection = selection;};

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