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

#include <G4VContinuousProcess.hh>

+ Inheritance diagram for G4VContinuousProcess:

Public Member Functions

 G4VContinuousProcess (const G4String &aName, G4ProcessType aType=fNotDefined)
 
 G4VContinuousProcess (G4VContinuousProcess &)
 
virtual ~G4VContinuousProcess ()
 
G4VContinuousProcessoperator= (const G4VContinuousProcess &)=delete
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &, G4double, G4ForceCondition *)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
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
 
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 G4VContinuousProcess.hh.

Constructor & Destructor Documentation

◆ G4VContinuousProcess() [1/2]

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

Definition at line 50 of file G4VContinuousProcess.cc.

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

Referenced by G4VContinuousProcess().

◆ G4VContinuousProcess() [2/2]

G4VContinuousProcess::G4VContinuousProcess ( G4VContinuousProcess & right)

Definition at line 64 of file G4VContinuousProcess.cc.

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

◆ ~G4VContinuousProcess()

G4VContinuousProcess::~G4VContinuousProcess ( )
virtual

Definition at line 59 of file G4VContinuousProcess.cc.

60{
61}

Member Function Documentation

◆ AlongStepDoIt()

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

Implements G4VProcess.

Reimplemented in G4AdjointAlongStepWeightCorrection, G4ContinuousGainOfEnergy, and G4ErrorEnergyLoss.

Definition at line 102 of file G4VContinuousProcess.cc.

104{
105 return pParticleChange;
106}
G4VParticleChange * pParticleChange

◆ AlongStepGetPhysicalInteractionLength()

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

Implements G4VProcess.

Definition at line 71 of file G4VContinuousProcess.cc.

77{
78 // GPILSelection is set to defaule value of CandidateForSelection
79 valueGPILSelection = CandidateForSelection;
80
81 // get Step limit proposed by the process
82 G4double steplength = GetContinuousStepLimit(track,previousStepSize,currentMinimumStep,currentSafety);
83
84 // set return value for G4GPILSelection
85 *selection = valueGPILSelection;
86
87#ifdef G4VERBOSE
88 if (verboseLevel>1)
89 {
90 G4cout << "G4VContinuousProcess::AlongStepGetPhysicalInteractionLength() - "
91 << "[ " << GetProcessName() << "]" << G4endl;
93 G4cout << " in Material " << track.GetMaterial()->GetName() << G4endl;
94 G4cout << "IntractionLength= " << steplength/cm <<"[cm] " << G4endl;
95 }
96#endif
97
98 return steplength ;
99}
@ 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 * G4VContinuousProcess::AtRestDoIt ( const G4Track & ,
const G4Step &  )
inlinevirtual

Implements G4VProcess.

Definition at line 85 of file G4VContinuousProcess.hh.

88 { return 0; }

◆ AtRestGetPhysicalInteractionLength()

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

Implements G4VProcess.

Definition at line 78 of file G4VContinuousProcess.hh.

81 { return -1.0; }

◆ GetContinuousStepLimit()

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

◆ GetGPILSelection()

G4GPILSelection G4VContinuousProcess::GetGPILSelection ( ) const
inlineprotected

Definition at line 107 of file G4VContinuousProcess.hh.

108 { return valueGPILSelection; }

◆ operator=()

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

◆ PostStepDoIt()

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

Implements G4VProcess.

Definition at line 90 of file G4VContinuousProcess.hh.

93 { return 0; }

◆ PostStepGetPhysicalInteractionLength()

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

Implements G4VProcess.

Definition at line 72 of file G4VContinuousProcess.hh.

76 { return -1.0; }

◆ SetGPILSelection()

void G4VContinuousProcess::SetGPILSelection ( G4GPILSelection selection)
inlineprotected

Definition at line 104 of file G4VContinuousProcess.hh.

105 { valueGPILSelection = selection; }

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