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

#include <G4VITProcess.hh>

+ Inheritance diagram for G4VITProcess:

Classes

struct  G4ProcessState
 
class  G4ProcessStateBase
 

Public Member Functions

 G4VITProcess (const G4String &name, G4ProcessType type=fNotDefined)
 
 ~G4VITProcess () override
 
 G4VITProcess (const G4VITProcess &other)
 
G4VITProcessoperator= (const G4VITProcess &other)
 
G4bool operator== (const G4VITProcess &right) const
 
G4bool operator!= (const G4VITProcess &right) const
 
size_t GetProcessID () const
 
G4shared_ptr< G4ProcessState_LockGetProcessState ()
 
void SetProcessState (G4shared_ptr< G4ProcessState_Lock > aProcInfo)
 
void ResetProcessState ()
 
void StartTracking (G4Track *) override
 
void BuildPhysicsTable (const G4ParticleDefinition &) override
 
G4double GetInteractionTimeLeft ()
 
void ResetNumberOfInteractionLengthLeft () override
 
G4bool ProposesTimeStep () const
 
- 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
 
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 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 EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
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 &)
 

Static Public Member Functions

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

Protected Member Functions

void RetrieveProcessInfo ()
 
void CreateInfo ()
 
template<typename T >
T * GetState ()
 
virtual void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
virtual void ClearInteractionTimeLeft ()
 
virtual void ClearNumberOfInteractionLengthLeft ()
 
void SetInstantiateProcessState (G4bool flag)
 
G4bool InstantiateProcessState ()
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Protected Attributes

G4shared_ptr< G4ProcessStatefpState
 
G4bool fProposesTimeStep
 
- 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

G4VITProcess inherits from G4VProcess. A G4VITProcess is able to save its current state for a given track into G4IT. This state may be retrieve latter on to be used by the G4VITProcess. Each G4VITProcess is tagged.

Definition at line 98 of file G4VITProcess.hh.

Constructor & Destructor Documentation

◆ G4VITProcess() [1/2]

G4VITProcess::G4VITProcess ( const G4String & name,
G4ProcessType type = fNotDefined )

Definition at line 37 of file G4VITProcess.cc.

37 :
38 G4VProcess(name, type)
39// fpState (0)//,
40//fProcessID(fNbProcess)
41{
42 fpState.reset();
43 if(fNbProcess == nullptr) fNbProcess = new size_t(0);
44 fProcessID = *fNbProcess;
45 (*fNbProcess)++;
47 currentInteractionLength = nullptr;
48 theInteractionTimeLeft = nullptr;
49 theNumberOfInteractionLengthLeft = nullptr;
50 fProposesTimeStep = false;
51}
void SetInstantiateProcessState(G4bool flag)
G4shared_ptr< G4ProcessState > fpState
G4bool fProposesTimeStep

◆ ~G4VITProcess()

G4VITProcess::~G4VITProcess ( )
overridedefault

◆ G4VITProcess() [2/2]

G4VITProcess::G4VITProcess ( const G4VITProcess & other)

Definition at line 66 of file G4VITProcess.cc.

66 :
67 G4VProcess(other), fProcessID(other.fProcessID)
68{
69 //copy ctor
70 //fpState = 0 ;
71 currentInteractionLength = nullptr;
72 theInteractionTimeLeft = nullptr;
73 theNumberOfInteractionLengthLeft = nullptr;
74 fInstantiateProcessState = other.fInstantiateProcessState;
76}

Member Function Documentation

◆ BuildPhysicsTable()

void G4VITProcess::BuildPhysicsTable ( const G4ParticleDefinition & )
inlineoverridevirtual

Reimplemented from G4VProcess.

Definition at line 140 of file G4VITProcess.hh.

141 {
142 }

◆ ClearInteractionTimeLeft()

void G4VITProcess::ClearInteractionTimeLeft ( )
inlineprotectedvirtual

Definition at line 257 of file G4VITProcess.hh.

258{
259 fpState->theInteractionTimeLeft = -1.0;
260}

Referenced by G4DNAElectronHoleRecombination::AtRestDoIt(), G4DNAMolecularDissociation::AtRestDoIt(), and G4VITRestProcess::AtRestDoIt().

◆ ClearNumberOfInteractionLengthLeft()

void G4VITProcess::ClearNumberOfInteractionLengthLeft ( )
inlineprotectedvirtual

◆ CreateInfo()

void G4VITProcess::CreateInfo ( )
protected

◆ GetInteractionTimeLeft()

G4double G4VITProcess::GetInteractionTimeLeft ( )
inline

Definition at line 272 of file G4VITProcess.hh.

273{
274 if (fpState) return fpState->theInteractionTimeLeft;
275
276 return -1;
277}

Referenced by G4ITStepProcessor::DoDefinePhysicalStepLength().

◆ GetMaxProcessIndex()

const size_t & G4VITProcess::GetMaxProcessIndex ( )
inlinestatic

Definition at line 284 of file G4VITProcess.hh.

285{
286 if (fNbProcess == nullptr) fNbProcess = new size_t(0);
287 return *fNbProcess;
288}

Referenced by G4TrackingInformation::GetProcessState().

◆ GetProcessID()

◆ GetProcessState()

G4shared_ptr< G4ProcessState_Lock > G4VITProcess::GetProcessState ( )
inline

Definition at line 120 of file G4VITProcess.hh.

121 {
123 }
#define UpcastProcessState(destinationType)

◆ GetState()

template<typename T >
T * G4VITProcess::GetState ( )
inlineprotected

Definition at line 211 of file G4VITProcess.hh.

212 {
213 return fpState->GetState<T>();
214 }

Referenced by G4ITTransportation::IsStepLimitedByGeometry().

◆ InstantiateProcessState()

G4bool G4VITProcess::InstantiateProcessState ( )
inlineprotected

Definition at line 233 of file G4VITProcess.hh.

234 {
235 return fInstantiateProcessState;
236 }

Referenced by StartTracking().

◆ operator!=()

G4bool G4VITProcess::operator!= ( const G4VITProcess & right) const

◆ operator=()

G4VITProcess & G4VITProcess::operator= ( const G4VITProcess & other)

Definition at line 78 of file G4VITProcess.cc.

79{
80 if (this == &rhs) return *this; // handle self assignment
81 //assignment operator
82 return *this;
83}

◆ operator==()

G4bool G4VITProcess::operator== ( const G4VITProcess & right) const

◆ ProposesTimeStep()

G4bool G4VITProcess::ProposesTimeStep ( ) const
inline

◆ ResetNumberOfInteractionLengthLeft()

void G4VITProcess::ResetNumberOfInteractionLengthLeft ( )
inlineoverridevirtual

◆ ResetProcessState()

◆ RetrieveProcessInfo()

void G4VITProcess::RetrieveProcessInfo ( )
protected

◆ SetInstantiateProcessState()

◆ SetProcessState()

◆ StartTracking()

void G4VITProcess::StartTracking ( G4Track * track)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 85 of file G4VITProcess.cc.

86{
87 G4TrackingInformation* trackingInfo = GetIT(track)->GetTrackingInfo();
89 {
90 // fpState = new G4ProcessState();
91 fpState = std::make_shared<G4ProcessState>();
92 }
93
94 theNumberOfInteractionLengthLeft = &(fpState->theNumberOfInteractionLengthLeft );
95 theInteractionTimeLeft = &(fpState->theInteractionTimeLeft );
96 currentInteractionLength = &(fpState->currentInteractionLength );
97 trackingInfo->RecordProcessState(fpState,fProcessID);
98 // fpState = 0;
99 fpState.reset();
100}
G4IT * GetIT(const G4Track *track)
Definition G4IT.cc:48
G4TrackingInformation * GetTrackingInfo()
Definition G4IT.hh:143
void RecordProcessState(G4shared_ptr< G4ProcessState_Lock >, size_t index)
G4bool InstantiateProcessState()

Referenced by G4DNAElectronHoleRecombination::StartTracking(), G4DNAPolyNucleotideReactionProcess::StartTracking(), G4DNAScavengerProcess::StartTracking(), G4DNASecondOrderReaction::StartTracking(), and G4ITTransportation::StartTracking().

◆ SubtractNumberOfInteractionLengthLeft()

void G4VITProcess::SubtractNumberOfInteractionLengthLeft ( G4double previousStepSize)
inlineprotectedvirtual

Definition at line 291 of file G4VITProcess.hh.

292{
293 if (fpState->currentInteractionLength > 0.0)
294 {
295 fpState->theNumberOfInteractionLengthLeft -= previousStepSize
296 / fpState->currentInteractionLength;
297 if (fpState->theNumberOfInteractionLengthLeft < 0.)
298 {
299 fpState->theNumberOfInteractionLengthLeft = CLHEP::perMillion;
300 }
301
302 }
303 else
304 {
305#ifdef G4VERBOSE
306 if (verboseLevel > 0)
307 {
308 G4cerr << "G4VITProcess::SubtractNumberOfInteractionLengthLeft()";
309 G4cerr << " [" << theProcessName << "]" << G4endl;
310 G4cerr << " currentInteractionLength = "
311 << fpState->currentInteractionLength << " [mm]";
312 G4cerr << " previousStepSize = " << previousStepSize << " [mm]";
313 G4cerr << G4endl;
314 }
315#endif
316 G4String msg = "Negative currentInteractionLength for ";
317 msg += theProcessName;
318 G4Exception("G4VITProcess::SubtractNumberOfInteractionLengthLeft()",
319 "ProcMan201",EventMustBeAborted,
320 msg);
321 }
322}
@ EventMustBeAborted
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition G4ios.hh:67
G4int verboseLevel
G4String theProcessName

Referenced by G4DNAPolyNucleotideReactionProcess::PostStepGetPhysicalInteractionLength(), G4DNAScavengerProcess::PostStepGetPhysicalInteractionLength(), G4DNASecondOrderReaction::PostStepGetPhysicalInteractionLength(), G4VITDiscreteProcess::PostStepGetPhysicalInteractionLength(), and G4VITRestDiscreteProcess::PostStepGetPhysicalInteractionLength().

Member Data Documentation

◆ fProposesTimeStep

◆ fpState


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