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

#include <G4VRestDiscreteProcess.hh>

+ Inheritance diagram for G4VRestDiscreteProcess:

Public Member Functions

 G4VRestDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VRestDiscreteProcess (G4VRestDiscreteProcess &)
 
virtual ~G4VRestDiscreteProcess ()
 
G4VRestDiscreteProcessoperator= (const G4VRestDiscreteProcess &)=delete
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4VParticleChangeAlongStepDoIt (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
 
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 Member Functions

virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
virtual G4double GetMeanLifeTime (const G4Track &aTrack, G4ForceCondition *condition)=0
 
- 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 G4VRestDiscreteProcess.hh.

Constructor & Destructor Documentation

◆ G4VRestDiscreteProcess() [1/2]

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

Definition at line 51 of file G4VRestDiscreteProcess.cc.

53 : G4VProcess(aName, aType)
54{
55 enableAlongStepDoIt = false;
56}
G4bool enableAlongStepDoIt
Definition: G4VProcess.hh:364

◆ G4VRestDiscreteProcess() [2/2]

G4VRestDiscreteProcess::G4VRestDiscreteProcess ( G4VRestDiscreteProcess right)

Definition at line 64 of file G4VRestDiscreteProcess.cc.

66 : G4VProcess(right)
67{
68}

◆ ~G4VRestDiscreteProcess()

G4VRestDiscreteProcess::~G4VRestDiscreteProcess ( )
virtual

Definition at line 59 of file G4VRestDiscreteProcess.cc.

60{
61}

Member Function Documentation

◆ AlongStepDoIt()

virtual G4VParticleChange * G4VRestDiscreteProcess::AlongStepDoIt ( const G4Track ,
const G4Step  
)
inlinevirtual

Implements G4VProcess.

Definition at line 87 of file G4VRestDiscreteProcess.hh.

90 { return 0; }

◆ AlongStepGetPhysicalInteractionLength()

virtual G4double G4VRestDiscreteProcess::AlongStepGetPhysicalInteractionLength ( const G4Track ,
G4double  ,
G4double  ,
G4double ,
G4GPILSelection  
)
inlinevirtual

Implements G4VProcess.

Definition at line 78 of file G4VRestDiscreteProcess.hh.

84 { return -1.0; }

◆ AtRestDoIt()

G4VParticleChange * G4VRestDiscreteProcess::AtRestDoIt ( const G4Track ,
const G4Step  
)
virtual

Implements G4VProcess.

Reimplemented in G4Decay, G4DecayWithSpin, G4Scintillation, and G4MuonicAtomDecay.

Definition at line 167 of file G4VRestDiscreteProcess.cc.

169{
171
172 return pParticleChange;
173}
void ClearNumberOfInteractionLengthLeft()
Definition: G4VProcess.hh:428
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:325

◆ AtRestGetPhysicalInteractionLength()

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

Implements G4VProcess.

Reimplemented in G4MuonicAtomDecay, and G4Decay.

Definition at line 132 of file G4VRestDiscreteProcess.cc.

135{
136 // beginning of tracking
138
139 // condition is set to "Not Forced"
141
142 // get mean life time
144
147
148#ifdef G4VERBOSE
149 if ((currentInteractionLength <0.0) || (verboseLevel>2))
150 {
153 G4cout << "G4VRestDiscreteProcess::AtRestGetPhysicalInteractionLength() - ";
154 G4cout << "[ " << GetProcessName() << "]" << G4endl;
155 track.GetDynamicParticle()->DumpInfo();
156 G4cout << " in Material " << track.GetMaterial()->GetName() << G4endl;
157 G4cout << "MeanLifeTime = " << t << " [ns]"
158 << G4endl;
159 }
160#endif
161
162 return time;
163}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void DumpInfo(G4int mode=0) const
const G4String & GetName() const
Definition: G4Material.hh:172
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
G4double currentInteractionLength
Definition: G4VProcess.hh:339
virtual void ResetNumberOfInteractionLengthLeft()
Definition: G4VProcess.cc:80
G4int verboseLevel
Definition: G4VProcess.hh:360
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:335
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386
virtual G4double GetMeanLifeTime(const G4Track &aTrack, G4ForceCondition *condition)=0
#define DBL_MAX
Definition: templates.hh:62
#define ns(x)
Definition: xmltok.c:1649

◆ GetMeanFreePath()

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

◆ GetMeanLifeTime()

virtual G4double G4VRestDiscreteProcess::GetMeanLifeTime ( const G4Track aTrack,
G4ForceCondition condition 
)
protectedpure virtual

◆ operator=()

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

◆ PostStepDoIt()

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

Implements G4VProcess.

Reimplemented in G4Decay, G4DecayWithSpin, G4Scintillation, and G4MuonicAtomDecay.

Definition at line 123 of file G4VRestDiscreteProcess.cc.

125{
127
128 return pParticleChange;
129}

Referenced by G4Scintillation::PostStepDoIt().

◆ PostStepGetPhysicalInteractionLength()

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

Implements G4VProcess.

Reimplemented in G4MuonicAtomDecay, and G4Decay.

Definition at line 71 of file G4VRestDiscreteProcess.cc.

75{
76 if ( (previousStepSize < 0.0) || (theNumberOfInteractionLengthLeft<=0.0))
77 {
78 // beginning of tracking (or just after DoIt() of this process)
80 }
81 else if ( previousStepSize > 0.0)
82 {
83 // subtract NumberOfInteractionLengthLeft
85 }
86 else
87 {
88 // zero step
89 // DO NOTHING
90 }
91
92 // condition is set to "Not Forced"
94
95 // get mean free path
96 currentInteractionLength = GetMeanFreePath(track,previousStepSize,condition);
97
98 G4double value;
100 {
102 }
103 else
104 {
105 value = DBL_MAX;
106 }
107#ifdef G4VERBOSE
108 if (verboseLevel>1)
109 {
110 G4double length = (value < DBL_MAX) ? value/cm : value;
111 G4cout << "G4VRestDiscreteProcess::PostStepGetPhysicalInteractionLength() - ";
112 G4cout << "[ " << GetProcessName() << "]" << G4endl;
113 track.GetDynamicParticle()->DumpInfo();
114 G4cout << " in Material " << track.GetMaterial()->GetName() << G4endl;
115 G4cout << "InteractionLength= " << length <<"[cm] " << G4endl;
116 }
117#endif
118 return value;
119}
void SubtractNumberOfInteractionLengthLeft(G4double prevStepSize)
Definition: G4VProcess.hh:528
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0

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