Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4VRestDiscreteProcess Class Referenceabstract

#include <G4VRestDiscreteProcess.hh>

Inheritance diagram for G4VRestDiscreteProcess:
Collaboration diagram for G4VRestDiscreteProcess:

Public Member Functions

 G4VRestDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VRestDiscreteProcess (G4VRestDiscreteProcess &)
 
virtual ~G4VRestDiscreteProcess ()
 
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 ()
 
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 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 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
 
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 previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- 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

Definition at line 55 of file G4VRestDiscreteProcess.hh.

Constructor & Destructor Documentation

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

Definition at line 55 of file G4VRestDiscreteProcess.cc.

56  : G4VProcess(aName, aType)
57 {
58  enableAlongStepDoIt = false;
59 }
G4VProcess(const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
Definition: G4VProcess.cc:52
G4bool enableAlongStepDoIt
Definition: G4VProcess.hh:351

Here is the call graph for this function:

Here is the caller graph for this function:

G4VRestDiscreteProcess::G4VRestDiscreteProcess ( G4VRestDiscreteProcess right)

Definition at line 65 of file G4VRestDiscreteProcess.cc.

66  : G4VProcess(right)
67 {
68 }
G4VProcess(const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
Definition: G4VProcess.cc:52

Here is the call graph for this function:

G4VRestDiscreteProcess::~G4VRestDiscreteProcess ( )
virtual

Definition at line 61 of file G4VRestDiscreteProcess.cc.

62 {
63 }

Member Function Documentation

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

Implements G4VProcess.

Definition at line 99 of file G4VRestDiscreteProcess.hh.

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

Implements G4VProcess.

Definition at line 90 of file G4VRestDiscreteProcess.hh.

96  { return -1.0; };
G4VParticleChange * G4VRestDiscreteProcess::AtRestDoIt ( const G4Track ,
const G4Step  
)
virtual

Implements G4VProcess.

Reimplemented in G4Scintillation, G4Decay, and G4DecayWithSpin.

Definition at line 150 of file G4VRestDiscreteProcess.cc.

154 {
155 // clear NumberOfInteractionLengthLeft
157 
158  return pParticleChange;
159 }
void ClearNumberOfInteractionLengthLeft()
Definition: G4VProcess.hh:447
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283

Here is the call graph for this function:

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

Implements G4VProcess.

Reimplemented in G4Decay.

Definition at line 122 of file G4VRestDiscreteProcess.cc.

126 {
127  // beggining of tracking
129 
130  // condition is set to "Not Forced"
131  *condition = NotForced;
132 
133  // get mean life time
135 
136 #ifdef G4VERBOSE
137  if ((currentInteractionLength <0.0) || (verboseLevel>2)){
138  G4cout << "G4VRestDiscreteProcess::AtRestGetPhysicalInteractionLength ";
139  G4cout << "[ " << GetProcessName() << "]" <<G4endl;
140  track.GetDynamicParticle()->DumpInfo();
141  G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
142  G4cout << "MeanLifeTime = " << currentInteractionLength/ns << "[ns]" <<G4endl;
143  }
144 #endif
145 
147 }
G4double condition(const G4ErrorSymMatrix &m)
G4int verboseLevel
Definition: G4VProcess.hh:368
const G4DynamicParticle * GetDynamicParticle() const
const G4String & GetName() const
Definition: G4Material.hh:178
void DumpInfo(G4int mode=0) const
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:293
virtual void ResetNumberOfInteractionLengthLeft()
Definition: G4VProcess.cc:95
G4GLOB_DLL std::ostream G4cout
G4double currentInteractionLength
Definition: G4VProcess.hh:297
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4Material * GetMaterial() const
#define G4endl
Definition: G4ios.hh:61
#define ns
Definition: xmlparse.cc:614
virtual G4double GetMeanLifeTime(const G4Track &aTrack, G4ForceCondition *condition)=0

Here is the call graph for this function:

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

Implemented in G4RadioactiveDecay, G4Decay, and G4Scintillation.

Here is the caller graph for this function:

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

Implemented in G4RadioactiveDecay, G4Decay, and G4Scintillation.

Here is the caller graph for this function:

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

Implements G4VProcess.

Reimplemented in G4Scintillation, G4Decay, and G4DecayWithSpin.

Definition at line 111 of file G4VRestDiscreteProcess.cc.

115 {
116 // reset NumberOfInteractionLengthLeft
118 
119  return pParticleChange;
120 }
void ClearNumberOfInteractionLengthLeft()
Definition: G4VProcess.hh:447
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283

Here is the call graph for this function:

Here is the caller graph for this function:

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

Implements G4VProcess.

Reimplemented in G4Decay.

Definition at line 70 of file G4VRestDiscreteProcess.cc.

75 {
76  if ( (previousStepSize < 0.0) || (theNumberOfInteractionLengthLeft<=0.0)) {
77  // beggining of tracking (or just after DoIt of this process)
79  } else if ( previousStepSize > 0.0) {
80  // subtract NumberOfInteractionLengthLeft
81  SubtractNumberOfInteractionLengthLeft(previousStepSize);
82  } else {
83  // zero step
84  // DO NOTHING
85  }
86 
87  // condition is set to "Not Forced"
89 
90  // get mean free path
91  currentInteractionLength = GetMeanFreePath(track, previousStepSize, condition);
92 
96  } else {
97  value = DBL_MAX;
98  }
99 #ifdef G4VERBOSE
100  if (verboseLevel>1){
101  G4cout << "G4VRestDiscreteProcess::PostStepGetPhysicalInteractionLength ";
102  G4cout << "[ " << GetProcessName() << "]" <<G4endl;
103  track.GetDynamicParticle()->DumpInfo();
104  G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
105  G4cout << "InteractionLength= " << value/cm <<"[cm] " <<G4endl;
106  }
107 #endif
108  return value;
109 }
G4double condition(const G4ErrorSymMatrix &m)
G4int verboseLevel
Definition: G4VProcess.hh:368
const G4DynamicParticle * GetDynamicParticle() const
const G4String & GetName() const
Definition: G4Material.hh:178
void DumpInfo(G4int mode=0) const
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:293
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
virtual void ResetNumberOfInteractionLengthLeft()
Definition: G4VProcess.cc:95
G4GLOB_DLL std::ostream G4cout
const XML_Char int const XML_Char * value
Definition: expat.h:331
static constexpr double cm
Definition: G4SIunits.hh:119
G4double currentInteractionLength
Definition: G4VProcess.hh:297
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4Material * GetMaterial() const
#define G4endl
Definition: G4ios.hh:61
void SubtractNumberOfInteractionLengthLeft(G4double previousStepSize)
Definition: G4VProcess.hh:544
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: