G4VITRestDiscreteProcess Class Reference

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

#include <G4VITRestDiscreteProcess.hh>

Inheritance diagram for G4VITRestDiscreteProcess:

G4VITProcess G4VProcess

Public Member Functions

 G4VITRestDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 G4VITRestDiscreteProcess (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 &)

Protected Member Functions

virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
virtual G4double GetMeanLifeTime (const G4Track &aTrack, G4ForceCondition *condition)=0

Detailed Description

Identical to G4VRestDiscreteProcess with dependency from G4VITProcess.

Identical to G4VRestDiscreteProcess with dependency from G4VITProcess

Definition at line 52 of file G4VITRestDiscreteProcess.hh.


Constructor & Destructor Documentation

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

Definition at line 45 of file G4VITRestDiscreteProcess.cc.

References G4VProcess::enableAlongStepDoIt, and G4VITRestDiscreteProcess().

Referenced by G4VITRestDiscreteProcess().

00046                   : G4VITProcess(aName, aType)
00047 {
00048   enableAlongStepDoIt  = false;
00049 }

G4VITRestDiscreteProcess::G4VITRestDiscreteProcess ( G4VITRestDiscreteProcess  ) 

Definition at line 55 of file G4VITRestDiscreteProcess.cc.

References G4VITRestDiscreteProcess().

00056                   : G4VITProcess(right)
00057 {
00058 }

G4VITRestDiscreteProcess::~G4VITRestDiscreteProcess (  )  [virtual]

Definition at line 51 of file G4VITRestDiscreteProcess.cc.

00052 {
00053 }


Member Function Documentation

virtual G4VParticleChange* G4VITRestDiscreteProcess::AlongStepDoIt ( const G4Track ,
const G4Step  
) [inline, virtual]

Implements G4VProcess.

Definition at line 96 of file G4VITRestDiscreteProcess.hh.

00099                               {return 0;}

virtual G4double G4VITRestDiscreteProcess::AlongStepGetPhysicalInteractionLength ( const G4Track ,
G4double  ,
G4double  ,
G4double ,
G4GPILSelection  
) [inline, virtual]

Implements G4VProcess.

Definition at line 87 of file G4VITRestDiscreteProcess.hh.

00093                              { return -1.0; }

G4VParticleChange * G4VITRestDiscreteProcess::AtRestDoIt ( const G4Track ,
const G4Step  
) [inline, virtual]

Implements G4VProcess.

Definition at line 207 of file G4VITRestDiscreteProcess.hh.

References G4VITProcess::ClearNumberOfInteractionLengthLeft(), and G4VProcess::pParticleChange.

00211 {
00212 //  clear NumberOfInteractionLengthLeft
00213     ClearNumberOfInteractionLengthLeft();
00214 
00215     return pParticleChange;
00216 }

G4double G4VITRestDiscreteProcess::AtRestGetPhysicalInteractionLength ( const G4Track ,
G4ForceCondition  
) [inline, virtual]

Implements G4VProcess.

Definition at line 179 of file G4VITRestDiscreteProcess.hh.

References G4DynamicParticle::DumpInfo(), G4VITProcess::fpState, G4cout, G4endl, G4Track::GetDynamicParticle(), G4Track::GetMaterial(), GetMeanLifeTime(), G4Material::GetName(), G4VProcess::GetProcessName(), NotForced, ns, G4VITProcess::ResetNumberOfInteractionLengthLeft(), and G4VProcess::verboseLevel.

00183 {
00184   // beggining of tracking
00185   ResetNumberOfInteractionLengthLeft();
00186 
00187   // condition is set to "Not Forced"
00188   *condition = NotForced;
00189 
00190   // get mean life time
00191   fpState->currentInteractionLength = GetMeanLifeTime(track, condition);
00192 
00193 #ifdef G4VERBOSE
00194    if ((fpState->currentInteractionLength <0.0) || (verboseLevel>2)){
00195     G4cout << "G4VITRestDiscreteProcess::AtRestGetPhysicalInteractionLength ";
00196     G4cout << "[ " << GetProcessName() << "]" <<G4endl;
00197     track.GetDynamicParticle()->DumpInfo();
00198     G4cout << " in Material  " << track.GetMaterial()->GetName() <<G4endl;
00199     G4cout << "MeanLifeTime = " << fpState->currentInteractionLength/CLHEP::ns << "[ns]" <<G4endl;
00200   }
00201 #endif
00202 
00203   return fpState->theNumberOfInteractionLengthLeft * (fpState->currentInteractionLength);
00204 }

virtual G4double G4VITRestDiscreteProcess::GetMeanFreePath ( const G4Track aTrack,
G4double  previousStepSize,
G4ForceCondition condition 
) [protected, pure virtual]

Referenced by PostStepGetPhysicalInteractionLength().

virtual G4double G4VITRestDiscreteProcess::GetMeanLifeTime ( const G4Track aTrack,
G4ForceCondition condition 
) [protected, pure virtual]

Referenced by AtRestGetPhysicalInteractionLength().

G4VParticleChange * G4VITRestDiscreteProcess::PostStepDoIt ( const G4Track ,
const G4Step  
) [inline, virtual]

Implements G4VProcess.

Definition at line 168 of file G4VITRestDiscreteProcess.hh.

References G4VITProcess::ClearNumberOfInteractionLengthLeft(), and G4VProcess::pParticleChange.

00172 {
00173 //  reset NumberOfInteractionLengthLeft
00174     ClearNumberOfInteractionLengthLeft();
00175 
00176     return pParticleChange;
00177 }

G4double G4VITRestDiscreteProcess::PostStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
) [inline, virtual]

Implements G4VProcess.

Definition at line 127 of file G4VITRestDiscreteProcess.hh.

References DBL_MAX, G4DynamicParticle::DumpInfo(), G4VITProcess::fpState, G4cout, G4endl, G4Track::GetDynamicParticle(), G4Track::GetMaterial(), GetMeanFreePath(), G4Material::GetName(), G4VProcess::GetProcessName(), NotForced, G4VITProcess::ResetNumberOfInteractionLengthLeft(), G4VITProcess::SubtractNumberOfInteractionLengthLeft(), and G4VProcess::verboseLevel.

00132 {
00133   if ( (previousStepSize < 0.0) || (fpState->theNumberOfInteractionLengthLeft<=0.0)) {
00134     // beggining of tracking (or just after DoIt of this process)
00135     ResetNumberOfInteractionLengthLeft();
00136   } else if ( previousStepSize > 0.0) {
00137     // subtract NumberOfInteractionLengthLeft
00138     SubtractNumberOfInteractionLengthLeft(previousStepSize);
00139   } else {
00140     // zero step
00141     //  DO NOTHING
00142   }
00143 
00144   // condition is set to "Not Forced"
00145   *condition = NotForced;
00146 
00147   // get mean free path
00148   fpState->currentInteractionLength = GetMeanFreePath(track, previousStepSize, condition);
00149 
00150   G4double value;
00151   if (fpState->currentInteractionLength <DBL_MAX) {
00152     value = fpState->theNumberOfInteractionLengthLeft * (fpState->currentInteractionLength);
00153   } else {
00154     value = DBL_MAX;
00155   }
00156 #ifdef G4VERBOSE
00157    if (verboseLevel>1){
00158     G4cout << "G4VITRestDiscreteProcess::PostStepGetPhysicalInteractionLength ";
00159     G4cout << "[ " << GetProcessName() << "]" <<G4endl;
00160     track.GetDynamicParticle()->DumpInfo();
00161     G4cout << " in Material  " <<  track.GetMaterial()->GetName() <<G4endl;
00162     G4cout << "InteractionLength= " << value/CLHEP::cm <<"[cm] " <<G4endl;
00163   }
00164 #endif
00165   return value;
00166 }


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:53:50 2013 for Geant4 by  doxygen 1.4.7