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

#include <G4VDiscreteProcess.hh>

Inheritance diagram for G4VDiscreteProcess:
Collaboration diagram for G4VDiscreteProcess:

Public Member Functions

 G4VDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
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
 
- 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 58 of file G4VDiscreteProcess.hh.

Constructor & Destructor Documentation

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

Definition at line 54 of file G4VDiscreteProcess.cc.

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

Here is the call graph for this function:

Here is the caller graph for this function:

G4VDiscreteProcess::G4VDiscreteProcess ( G4VDiscreteProcess right)

Definition at line 66 of file G4VDiscreteProcess.cc.

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

Here is the call graph for this function:

G4VDiscreteProcess::~G4VDiscreteProcess ( )
virtual

Definition at line 62 of file G4VDiscreteProcess.cc.

63 {
64 }

Member Function Documentation

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

Implements G4VProcess.

Reimplemented in G4NuclearStopping.

Definition at line 102 of file G4VDiscreteProcess.hh.

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

Implements G4VProcess.

Reimplemented in G4NuclearStopping.

Definition at line 83 of file G4VDiscreteProcess.hh.

89  { return -1.0; };
virtual G4VParticleChange* G4VDiscreteProcess::AtRestDoIt ( const G4Track ,
const G4Step  
)
inlinevirtual

Implements G4VProcess.

Reimplemented in G4HadronStoppingProcess, G4MuonMinusAtomicCapture, and G4eplusAnnihilation.

Definition at line 97 of file G4VDiscreteProcess.hh.

100  {return 0;};
virtual G4double G4VDiscreteProcess::AtRestGetPhysicalInteractionLength ( const G4Track ,
G4ForceCondition  
)
inlinevirtual

Implements G4VProcess.

Reimplemented in G4eplusAnnihilation, G4HadronStoppingProcess, and G4MuonMinusAtomicCapture.

Definition at line 91 of file G4VDiscreteProcess.hh.

94  { return -1.0; };
G4double G4VDiscreteProcess::PostStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
virtual

Implements G4VProcess.

Reimplemented in G4VEmProcess, G4PolarizedCompton, G4UnknownDecay, G4NeutronKiller, G4HadronStoppingProcess, G4LowECapture, G4MuonMinusAtomicCapture, G4eplusPolarizedAnnihilation, G4ErrorTrackLengthTarget, G4VErrorLimitProcess, G4ErrorStepLengthLimitProcess, and G4ErrorMagFieldLimitProcess.

Definition at line 71 of file G4VDiscreteProcess.cc.

76 {
77  if ( (previousStepSize < 0.0) || (theNumberOfInteractionLengthLeft<=0.0)) {
78  // beggining of tracking (or just after DoIt of this process)
80  } else if ( previousStepSize > 0.0) {
81  // subtract NumberOfInteractionLengthLeft
82  SubtractNumberOfInteractionLengthLeft(previousStepSize);
83  } else {
84  // zero step
85  // DO NOTHING
86  }
87 
88  // condition is set to "Not Forced"
90 
91  // get mean free path
92  currentInteractionLength = GetMeanFreePath(track, previousStepSize, condition);
93 
97  } else {
98  value = DBL_MAX;
99  }
100 #ifdef G4VERBOSE
101  if (verboseLevel>1){
102  G4cout << "G4VDiscreteProcess::PostStepGetPhysicalInteractionLength ";
103  G4cout << "[ " << GetProcessName() << "]" <<G4endl;
104  track.GetDynamicParticle()->DumpInfo();
105  G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
106  G4cout << "InteractionLength= " << value/cm <<"[cm] " <<G4endl;
107  }
108 #endif
109  return value;
110 }
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
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
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
#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: