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

#include <G4IVContinuousDiscreteProcess.hh>

+ Inheritance diagram for G4IVContinuousDiscreteProcess:

Public Member Functions

 G4IVContinuousDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4IVContinuousDiscreteProcess (G4IVContinuousDiscreteProcess &)
 
virtual ~G4IVContinuousDiscreteProcess ()
 
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 GetContinuousStepLimit (const G4Track &aTrack, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)=0
 
virtual void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void SetGPILSelection (G4GPILSelection selection)
 
G4GPILSelection GetGPILSelection () const
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Protected Attributes

G4PhysicsTabletheNlambdaTable
 
G4PhysicsTabletheInverseNlambdaTable
 
const G4double BIGSTEP
 
- 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
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 

Detailed Description

Definition at line 60 of file G4IVContinuousDiscreteProcess.hh.

Constructor & Destructor Documentation

◆ G4IVContinuousDiscreteProcess() [1/2]

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

Definition at line 50 of file G4IVContinuousDiscreteProcess.cc.

51 : G4VProcess(aName, aType),
52 valueGPILSelection(CandidateForSelection),
54 BIGSTEP(1.e10)
55{
56 enableAtRestDoIt = false;
57}
@ CandidateForSelection
G4bool enableAtRestDoIt
Definition: G4VProcess.hh:350

Referenced by G4IVContinuousDiscreteProcess().

◆ G4IVContinuousDiscreteProcess() [2/2]

G4IVContinuousDiscreteProcess::G4IVContinuousDiscreteProcess ( G4IVContinuousDiscreteProcess right)

Definition at line 62 of file G4IVContinuousDiscreteProcess.cc.

63 : G4VProcess(right),
64 valueGPILSelection(right.valueGPILSelection),
66 BIGSTEP(right.BIGSTEP)
67{
68}

◆ ~G4IVContinuousDiscreteProcess()

G4IVContinuousDiscreteProcess::~G4IVContinuousDiscreteProcess ( )
virtual

Definition at line 58 of file G4IVContinuousDiscreteProcess.cc.

59{
60}

Member Function Documentation

◆ AlongStepDoIt()

G4VParticleChange * G4IVContinuousDiscreteProcess::AlongStepDoIt ( const G4Track ,
const G4Step  
)
inlinevirtual

Implements G4VProcess.

Definition at line 166 of file G4IVContinuousDiscreteProcess.hh.

170{
171// clear NumberOfInteractionLengthLeft
173 return pParticleChange;
174}
void ClearNumberOfInteractionLengthLeft()
Definition: G4VProcess.hh:418
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283

◆ AlongStepGetPhysicalInteractionLength()

G4double G4IVContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4double  currentMinimumStep,
G4double currentSafety,
G4GPILSelection selection 
)
inlinevirtual

Implements G4VProcess.

Definition at line 176 of file G4IVContinuousDiscreteProcess.hh.

183{
184 // GPILSelection is set to defaule value of CandidateForSelection
185 valueGPILSelection = CandidateForSelection;
186
187 // get Step limit proposed by the process
188 G4double steplength = GetContinuousStepLimit(track,previousStepSize,currentMinimumStep, currentSafety);
189
190 // set return value for G4GPILSelection
191 *selection = valueGPILSelection;
192
193 if (verboseLevel>1){
194 G4cout << "G4IVContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength ";
195 G4cout << "[ " << GetProcessName() << "]" <<G4endl;
196 track.GetDynamicParticle()->DumpInfo();
197 G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
198 G4cout << "IntractionLength= " << steplength/CLHEP::cm <<"[cm] " <<G4endl;
199 }
200 return steplength ;
201}
double G4double
Definition: G4Types.hh:64
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
void DumpInfo(G4int mode=0) const
virtual G4double GetContinuousStepLimit(const G4Track &aTrack, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)=0
const G4String & GetName() const
Definition: G4Material.hh:177
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
G4int verboseLevel
Definition: G4VProcess.hh:368
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379

◆ AtRestDoIt()

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

Implements G4VProcess.

Definition at line 104 of file G4IVContinuousDiscreteProcess.hh.

107 {return 0;};

◆ AtRestGetPhysicalInteractionLength()

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

Implements G4VProcess.

Definition at line 98 of file G4IVContinuousDiscreteProcess.hh.

101 { return -1.0; };

◆ GetContinuousStepLimit()

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

◆ GetGPILSelection()

G4GPILSelection G4IVContinuousDiscreteProcess::GetGPILSelection ( ) const
inlineprotected

Definition at line 129 of file G4IVContinuousDiscreteProcess.hh.

129{return valueGPILSelection;};

◆ PostStepDoIt()

G4VParticleChange * G4IVContinuousDiscreteProcess::PostStepDoIt ( const G4Track ,
const G4Step  
)
inlinevirtual

Implements G4VProcess.

Definition at line 156 of file G4IVContinuousDiscreteProcess.hh.

160{
161// clear NumberOfInteractionLengthLeft
163 return pParticleChange;
164}

◆ PostStepGetPhysicalInteractionLength()

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

Implements G4VProcess.

Definition at line 72 of file G4IVContinuousDiscreteProcess.cc.

77{
78 G4double value = DBL_MAX ;
79
80 return value;
81}
#define DBL_MAX
Definition: templates.hh:83

◆ SetGPILSelection()

void G4IVContinuousDiscreteProcess::SetGPILSelection ( G4GPILSelection  selection)
inlineprotected

Definition at line 126 of file G4IVContinuousDiscreteProcess.hh.

127 { valueGPILSelection = selection;};

◆ SubtractNumberOfInteractionLengthLeft()

void G4IVContinuousDiscreteProcess::SubtractNumberOfInteractionLengthLeft ( G4double  previousStepSize)
inlineprotectedvirtual

Definition at line 146 of file G4IVContinuousDiscreteProcess.hh.

150{
151 // dummy routine
152 ;
153}

Member Data Documentation

◆ BIGSTEP

const G4double G4IVContinuousDiscreteProcess::BIGSTEP
protected

Definition at line 139 of file G4IVContinuousDiscreteProcess.hh.

◆ theInverseNlambdaTable

G4PhysicsTable* G4IVContinuousDiscreteProcess::theInverseNlambdaTable
protected

Definition at line 138 of file G4IVContinuousDiscreteProcess.hh.

◆ theNlambdaTable

G4PhysicsTable* G4IVContinuousDiscreteProcess::theNlambdaTable
protected

Definition at line 137 of file G4IVContinuousDiscreteProcess.hh.


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