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

#include <G4WrapperProcess.hh>

+ Inheritance diagram for G4WrapperProcess:

Public Member Functions

 G4WrapperProcess (const G4String &aName="Wrapped", G4ProcessType aType=fNotDefined)
 
 G4WrapperProcess (const G4WrapperProcess &right)
 
virtual ~G4WrapperProcess ()
 
G4WrapperProcessoperator= (const G4WrapperProcess &)=delete
 
G4bool operator== (const G4WrapperProcess &right) const
 
G4bool operator!= (const G4WrapperProcess &right) const
 
virtual void RegisterProcess (G4VProcess *)
 
virtual const G4VProcessGetRegisteredProcess () const
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)
 
virtual G4double PostStepGetPhysicalInteractionLength (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 &directory, G4bool ascii=false)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false)
 
virtual void StartTracking (G4Track *)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
- 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 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 Attributes

G4VProcesspRegProcess = nullptr
 
- 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
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Detailed Description

Definition at line 41 of file G4WrapperProcess.hh.

Constructor & Destructor Documentation

◆ G4WrapperProcess() [1/2]

G4WrapperProcess::G4WrapperProcess ( const G4String aName = "Wrapped",
G4ProcessType  aType = fNotDefined 
)

Definition at line 34 of file G4WrapperProcess.cc.

36 : G4VProcess(aName, aType)
37{
38}

◆ G4WrapperProcess() [2/2]

G4WrapperProcess::G4WrapperProcess ( const G4WrapperProcess right)

Definition at line 41 of file G4WrapperProcess.cc.

42 : G4VProcess(*((G4VProcess*)(&right))), pRegProcess(right.pRegProcess)
43{
44}
G4VProcess * pRegProcess

◆ ~G4WrapperProcess()

G4WrapperProcess::~G4WrapperProcess ( )
virtual

Definition at line 47 of file G4WrapperProcess.cc.

48{
49}

Member Function Documentation

◆ AlongStepDoIt()

G4VParticleChange * G4WrapperProcess::AlongStepDoIt ( const G4Track track,
const G4Step stepData 
)
virtual

Implements G4VProcess.

Definition at line 111 of file G4WrapperProcess.cc.

113{
114 return pRegProcess->AlongStepDoIt( track, stepData );
115}
virtual G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData)=0

◆ AlongStepGetPhysicalInteractionLength()

G4double G4WrapperProcess::AlongStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4double  currentMinimumStep,
G4double proposedSafety,
G4GPILSelection selection 
)
virtual

Implements G4VProcess.

Definition at line 58 of file G4WrapperProcess.cc.

64{
65 return pRegProcess->
67 previousStepSize,
68 currentMinimumStep,
69 proposedSafety,
70 selection );
71}
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)

Referenced by AlongStepGetPhysicalInteractionLength().

◆ AtRestDoIt()

G4VParticleChange * G4WrapperProcess::AtRestDoIt ( const G4Track track,
const G4Step stepData 
)
virtual

Implements G4VProcess.

Definition at line 118 of file G4WrapperProcess.cc.

120{
121 return pRegProcess->AtRestDoIt( track, stepData );
122}
virtual G4VParticleChange * AtRestDoIt(const G4Track &track, const G4Step &stepData)=0

◆ AtRestGetPhysicalInteractionLength()

G4double G4WrapperProcess::AtRestGetPhysicalInteractionLength ( const G4Track track,
G4ForceCondition condition 
)
virtual

Implements G4VProcess.

Definition at line 74 of file G4WrapperProcess.cc.

77{
79}
G4double condition(const G4ErrorSymMatrix &m)
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &track, G4ForceCondition *condition)=0

◆ BuildPhysicsTable()

void G4WrapperProcess::BuildPhysicsTable ( const G4ParticleDefinition particle)
virtual

Reimplemented from G4VProcess.

Definition at line 131 of file G4WrapperProcess.cc.

132{
133 return pRegProcess->BuildPhysicsTable(particle);
134}
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
Definition: G4VProcess.hh:187

◆ EndTracking()

void G4WrapperProcess::EndTracking ( )
virtual

Reimplemented from G4VProcess.

Definition at line 167 of file G4WrapperProcess.cc.

168{
170}
virtual void EndTracking()
Definition: G4VProcess.cc:102

◆ GetProcessManager()

const G4ProcessManager * G4WrapperProcess::GetProcessManager ( )
virtual

Reimplemented from G4VProcess.

Definition at line 98 of file G4WrapperProcess.cc.

99{
101}
virtual const G4ProcessManager * GetProcessManager()
Definition: G4VProcess.hh:498

◆ GetRegisteredProcess()

const G4VProcess * G4WrapperProcess::GetRegisteredProcess ( ) const
virtual

Definition at line 181 of file G4WrapperProcess.cc.

182{
183 return pRegProcess;
184}

◆ IsApplicable()

G4bool G4WrapperProcess::IsApplicable ( const G4ParticleDefinition particle)
virtual

Reimplemented from G4VProcess.

Definition at line 125 of file G4WrapperProcess.cc.

126{
127 return pRegProcess->IsApplicable(particle);
128}
virtual G4bool IsApplicable(const G4ParticleDefinition &)
Definition: G4VProcess.hh:182

◆ operator!=()

G4bool G4WrapperProcess::operator!= ( const G4WrapperProcess right) const
inline

Definition at line 176 of file G4WrapperProcess.hh.

177{
178 return (this != &right);
179}

◆ operator=()

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

◆ operator==()

G4bool G4WrapperProcess::operator== ( const G4WrapperProcess right) const
inline

Definition at line 170 of file G4WrapperProcess.hh.

171{
172 return (this == &right);
173}

◆ PostStepDoIt()

G4VParticleChange * G4WrapperProcess::PostStepDoIt ( const G4Track track,
const G4Step stepData 
)
virtual

Implements G4VProcess.

Definition at line 104 of file G4WrapperProcess.cc.

106{
107 return pRegProcess->PostStepDoIt( track, stepData );
108}
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData)=0

◆ PostStepGetPhysicalInteractionLength()

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

Implements G4VProcess.

Definition at line 82 of file G4WrapperProcess.cc.

86{
88 previousStepSize,
89 condition );
90}
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0

◆ PreparePhysicsTable()

void G4WrapperProcess::PreparePhysicsTable ( const G4ParticleDefinition particle)
virtual

Reimplemented from G4VProcess.

Definition at line 137 of file G4WrapperProcess.cc.

138{
139 return pRegProcess->PreparePhysicsTable(particle);
140}
virtual void PreparePhysicsTable(const G4ParticleDefinition &)
Definition: G4VProcess.hh:194

◆ RegisterProcess()

void G4WrapperProcess::RegisterProcess ( G4VProcess process)
virtual

Definition at line 173 of file G4WrapperProcess.cc.

174{
175 pRegProcess = process;
176 theProcessName += process->GetProcessName();
177 theProcessType = process->GetProcessType();
178}
G4ProcessType theProcessType
Definition: G4VProcess.hh:350
G4ProcessType GetProcessType() const
Definition: G4VProcess.hh:392
G4String theProcessName
Definition: G4VProcess.hh:345
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386

◆ ResetNumberOfInteractionLengthLeft()

void G4WrapperProcess::ResetNumberOfInteractionLengthLeft ( )
virtual

Reimplemented from G4VProcess.

Definition at line 52 of file G4WrapperProcess.cc.

53{
55}
virtual void ResetNumberOfInteractionLengthLeft()
Definition: G4VProcess.cc:80

◆ RetrievePhysicsTable()

G4bool G4WrapperProcess::RetrievePhysicsTable ( const G4ParticleDefinition particle,
const G4String directory,
G4bool  ascii = false 
)
virtual

Reimplemented from G4VProcess.

Definition at line 152 of file G4WrapperProcess.cc.

156{
157 return pRegProcess->RetrievePhysicsTable(particle, directory, ascii);
158}
virtual G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &, G4bool)
Definition: G4VProcess.hh:211

◆ SetMasterProcess()

void G4WrapperProcess::SetMasterProcess ( G4VProcess masterP)
virtual

Reimplemented from G4VProcess.

Definition at line 187 of file G4WrapperProcess.cc.

188{
189 G4WrapperProcess* master = static_cast<G4WrapperProcess*>(masterP);
190 // Cannot use getter because it returns "const" and we do not want
191 // to cast away constness explicitly (even if this is the same)
193}
virtual void SetMasterProcess(G4VProcess *masterP)
Definition: G4VProcess.cc:218

◆ SetProcessManager()

void G4WrapperProcess::SetProcessManager ( const G4ProcessManager procMan)
virtual

Reimplemented from G4VProcess.

Definition at line 93 of file G4WrapperProcess.cc.

94{
96}
virtual void SetProcessManager(const G4ProcessManager *)
Definition: G4VProcess.hh:492

◆ StartTracking()

void G4WrapperProcess::StartTracking ( G4Track track)
virtual

Reimplemented from G4VProcess.

Definition at line 161 of file G4WrapperProcess.cc.

162{
164}
virtual void StartTracking(G4Track *)
Definition: G4VProcess.cc:87

◆ StorePhysicsTable()

G4bool G4WrapperProcess::StorePhysicsTable ( const G4ParticleDefinition particle,
const G4String directory,
G4bool  ascii = false 
)
virtual

Reimplemented from G4VProcess.

Definition at line 143 of file G4WrapperProcess.cc.

147{
148 return pRegProcess->StorePhysicsTable(particle, directory, ascii);
149}
virtual G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &, G4bool)
Definition: G4VProcess.hh:206

Member Data Documentation

◆ pRegProcess


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