Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4VITRestDiscreteProcess Class Referenceabstract

Identical to G4VRestDiscreteProcess with dependency from G4VITProcess. More...

#include <G4VITRestDiscreteProcess.hh>

Inheritance diagram for G4VITRestDiscreteProcess:
Collaboration diagram for G4VITRestDiscreteProcess:

Public Member Functions

 G4VITRestDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VITRestDiscreteProcess (const G4VITRestDiscreteProcess &)
 
virtual ~G4VITRestDiscreteProcess ()
 
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 G4VITProcess
 G4VITProcess (const G4String &name, G4ProcessType type=fNotDefined)
 
virtual ~G4VITProcess ()
 
 G4VITProcess (const G4VITProcess &other)
 
G4VITProcessoperator= (const G4VITProcess &other)
 
G4int operator== (const G4VITProcess &right) const
 
G4int operator!= (const G4VITProcess &right) const
 
size_t GetProcessID () const
 
G4shared_ptr< G4ProcessState_LockGetProcessState ()
 
void SetProcessState (G4shared_ptr< G4ProcessState_Lock > aProcInfo)
 
void ResetProcessState ()
 
virtual void StartTracking (G4Track *)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
G4double GetInteractionTimeLeft ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4bool ProposesTimeStep () const
 
- 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
 
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 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
 
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
 
 G4VITRestDiscreteProcess ()
 
G4VITRestDiscreteProcessoperator= (const G4VITRestDiscreteProcess &right)
 
- Protected Member Functions inherited from G4VITProcess
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 previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VITProcess
static const size_t & GetMaxProcessIndex ()
 
- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Attributes inherited from G4VITProcess
G4shared_ptr< G4ProcessStatefpState
 
G4bool fProposesTimeStep
 
- 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
 

Detailed Description

Identical to G4VRestDiscreteProcess with dependency from G4VITProcess.

Identical to G4VRestDiscreteProcess with dependency from G4VITProcess

Definition at line 63 of file G4VITRestDiscreteProcess.hh.

Constructor & Destructor Documentation

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

Definition at line 45 of file G4VITRestDiscreteProcess.cc.

46  : G4VITProcess(aName, aType)
47 {
48  enableAlongStepDoIt = false;
49 }
G4VITProcess(const G4String &name, G4ProcessType type=fNotDefined)
Definition: G4VITProcess.cc:35
G4bool enableAlongStepDoIt
Definition: G4VProcess.hh:351
G4VITRestDiscreteProcess::G4VITRestDiscreteProcess ( const G4VITRestDiscreteProcess right)

Definition at line 55 of file G4VITRestDiscreteProcess.cc.

56  : G4VITProcess(right)
57 {
58 }
G4VITProcess(const G4String &name, G4ProcessType type=fNotDefined)
Definition: G4VITProcess.cc:35
G4VITRestDiscreteProcess::~G4VITRestDiscreteProcess ( )
virtual

Definition at line 51 of file G4VITRestDiscreteProcess.cc.

52 {
53 }
G4VITRestDiscreteProcess::G4VITRestDiscreteProcess ( )
protected

Definition at line 38 of file G4VITRestDiscreteProcess.cc.

39  :G4VITProcess("No Name Discrete Process")
40 {
41  G4Exception("G4VITRestDiscreteProcess::G4VITRestDiscreteProcess","Illegal operation",
42  JustWarning,"default constructor is called");
43 }
G4VITProcess(const G4String &name, G4ProcessType type=fNotDefined)
Definition: G4VITProcess.cc:35
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

Here is the call graph for this function:

Member Function Documentation

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

Implements G4VProcess.

Definition at line 99 of file G4VITRestDiscreteProcess.hh.

100  {
101  return 0;
102  }
virtual G4double G4VITRestDiscreteProcess::AlongStepGetPhysicalInteractionLength ( const G4Track ,
G4double  ,
G4double  ,
G4double ,
G4GPILSelection  
)
inlinevirtual

Implements G4VProcess.

Definition at line 89 of file G4VITRestDiscreteProcess.hh.

94  {
95  return -1.0;
96  }
G4VParticleChange * G4VITRestDiscreteProcess::AtRestDoIt ( const G4Track ,
const G4Step  
)
inlinevirtual

Implements G4VProcess.

Reimplemented in G4DNAMolecularDissociation, and G4DNAElectronHoleRecombination.

Definition at line 218 of file G4VITRestDiscreteProcess.hh.

220 {
221 // clear NumberOfInteractionLengthLeft
223 
224  return pParticleChange;
225 }
virtual void ClearNumberOfInteractionLengthLeft()
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283

Here is the call graph for this function:

G4double G4VITRestDiscreteProcess::AtRestGetPhysicalInteractionLength ( const G4Track track,
G4ForceCondition condition 
)
inlinevirtual

Implements G4VProcess.

Reimplemented in G4DNAMolecularDissociation.

Definition at line 190 of file G4VITRestDiscreteProcess.hh.

192 {
193  // beggining of tracking
195 
196  // condition is set to "Not Forced"
197  *condition = NotForced;
198 
199  // get mean life time
200  fpState->currentInteractionLength = GetMeanLifeTime(track, condition);
201 
202 #ifdef G4VERBOSE
203  if((fpState->currentInteractionLength < 0.0) || (verboseLevel > 2))
204  {
205  G4cout << "G4VITRestDiscreteProcess::AtRestGetPhysicalInteractionLength ";
206  G4cout << "[ " << GetProcessName() << "]" << G4endl;
207  track.GetDynamicParticle()->DumpInfo();
208  G4cout << " in Material " << track.GetMaterial()->GetName() << G4endl;
209  G4cout << "MeanLifeTime = " << fpState->currentInteractionLength / CLHEP::ns
210  << "[ns]" << G4endl;
211  }
212 #endif
213 
214  return fpState->theNumberOfInteractionLengthLeft
215  * (fpState->currentInteractionLength);
216 }
G4double condition(const G4ErrorSymMatrix &m)
G4int verboseLevel
Definition: G4VProcess.hh:368
const G4DynamicParticle * GetDynamicParticle() const
static constexpr double ns
virtual void ResetNumberOfInteractionLengthLeft()
const G4String & GetName() const
Definition: G4Material.hh:178
void DumpInfo(G4int mode=0) const
G4GLOB_DLL std::ostream G4cout
G4shared_ptr< G4ProcessState > fpState
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4Material * GetMaterial() const
virtual G4double GetMeanLifeTime(const G4Track &aTrack, G4ForceCondition *condition)=0
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

Here is the caller graph for this function:

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

Implemented in G4DNAMolecularDissociation, and G4DNAElectronHoleRecombination.

Here is the caller graph for this function:

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

Implemented in G4DNAMolecularDissociation, and G4DNAElectronHoleRecombination.

Here is the caller graph for this function:

G4VITRestDiscreteProcess& G4VITRestDiscreteProcess::operator= ( const G4VITRestDiscreteProcess right)
protected
G4VParticleChange * G4VITRestDiscreteProcess::PostStepDoIt ( const G4Track ,
const G4Step  
)
inlinevirtual

Implements G4VProcess.

Reimplemented in G4DNAMolecularDissociation, and G4DNAElectronHoleRecombination.

Definition at line 181 of file G4VITRestDiscreteProcess.hh.

183 {
184 // reset NumberOfInteractionLengthLeft
186 
187  return pParticleChange;
188 }
virtual void ClearNumberOfInteractionLengthLeft()
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283

Here is the call graph for this function:

G4double G4VITRestDiscreteProcess::PostStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
inlinevirtual

Implements G4VProcess.

Reimplemented in G4DNAMolecularDissociation.

Definition at line 130 of file G4VITRestDiscreteProcess.hh.

133 {
134  if((previousStepSize < 0.0) || (fpState->theNumberOfInteractionLengthLeft
135  <= 0.0))
136  {
137  // beggining of tracking (or just after DoIt of this process)
139  }
140  else if(previousStepSize > 0.0)
141  {
142  // subtract NumberOfInteractionLengthLeft
143  SubtractNumberOfInteractionLengthLeft(previousStepSize);
144  }
145  else
146  {
147  // zero step
148  // DO NOTHING
149  }
150 
151  // condition is set to "Not Forced"
152  *condition = NotForced;
153 
154  // get mean free path
155  fpState->currentInteractionLength = GetMeanFreePath(track,
156  previousStepSize,
157  condition);
158 
159  G4double value;
160  if(fpState->currentInteractionLength < DBL_MAX)
161  {
162  value = fpState->theNumberOfInteractionLengthLeft * (fpState->currentInteractionLength);
163  }
164  else
165  {
166  value = DBL_MAX;
167  }
168 #ifdef G4VERBOSE
169  if(verboseLevel > 1)
170  {
171  G4cout << "G4VITRestDiscreteProcess::PostStepGetPhysicalInteractionLength ";
172  G4cout << "[ " << GetProcessName() << "]" << G4endl;
173  track.GetDynamicParticle()->DumpInfo();
174  G4cout << " in Material " << track.GetMaterial()->GetName() << G4endl;
175  G4cout << "InteractionLength= " << value / CLHEP::cm << "[cm] " << G4endl;
176  }
177 #endif
178  return value;
179 }
G4double condition(const G4ErrorSymMatrix &m)
static constexpr double cm
Definition: SystemOfUnits.h:99
G4int verboseLevel
Definition: G4VProcess.hh:368
const G4DynamicParticle * GetDynamicParticle() const
virtual void ResetNumberOfInteractionLengthLeft()
const G4String & GetName() const
Definition: G4Material.hh:178
void DumpInfo(G4int mode=0) const
virtual void SubtractNumberOfInteractionLengthLeft(G4double previousStepSize)
G4GLOB_DLL std::ostream G4cout
const XML_Char int const XML_Char * value
Definition: expat.h:331
G4shared_ptr< G4ProcessState > fpState
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4Material * GetMaterial() const
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83

Here is the call graph for this function:


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