G4VRestDiscreteProcess Class Reference

#include <G4VRestDiscreteProcess.hh>

Inheritance diagram for G4VRestDiscreteProcess:

G4VProcess G4Decay G4RadioactiveDecay G4Scintillation G4DecayWithSpin G4PionDecayMakeSpin

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 &)

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

Definition at line 55 of file G4VRestDiscreteProcess.hh.


Constructor & Destructor Documentation

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

Definition at line 55 of file G4VRestDiscreteProcess.cc.

References G4VProcess::enableAlongStepDoIt, and G4VRestDiscreteProcess().

Referenced by G4VRestDiscreteProcess().

00056                   : G4VProcess(aName, aType)
00057 {
00058   enableAlongStepDoIt  = false;
00059 }

G4VRestDiscreteProcess::G4VRestDiscreteProcess ( G4VRestDiscreteProcess  ) 

Definition at line 65 of file G4VRestDiscreteProcess.cc.

References G4VRestDiscreteProcess().

00066                   : G4VProcess(right)
00067 {
00068 }

G4VRestDiscreteProcess::~G4VRestDiscreteProcess (  )  [virtual]

Definition at line 61 of file G4VRestDiscreteProcess.cc.

00062 {
00063 }


Member Function Documentation

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

Implements G4VProcess.

Definition at line 99 of file G4VRestDiscreteProcess.hh.

00102                               {return 0;};

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

Implements G4VProcess.

Definition at line 90 of file G4VRestDiscreteProcess.hh.

00096                              { return -1.0; };

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

Implements G4VProcess.

Reimplemented in G4Decay, and G4Scintillation.

Definition at line 150 of file G4VRestDiscreteProcess.cc.

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

00154 {
00155 //  clear NumberOfInteractionLengthLeft
00156     ClearNumberOfInteractionLengthLeft();
00157 
00158     return pParticleChange;
00159 }

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

Implements G4VProcess.

Reimplemented in G4Decay.

Definition at line 122 of file G4VRestDiscreteProcess.cc.

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

00126 {
00127   // beggining of tracking 
00128   ResetNumberOfInteractionLengthLeft();
00129 
00130   // condition is set to "Not Forced"
00131   *condition = NotForced;
00132 
00133   // get mean life time
00134   currentInteractionLength = GetMeanLifeTime(track, condition);
00135 
00136 #ifdef G4VERBOSE
00137    if ((currentInteractionLength <0.0) || (verboseLevel>2)){
00138     G4cout << "G4VRestDiscreteProcess::AtRestGetPhysicalInteractionLength ";
00139     G4cout << "[ " << GetProcessName() << "]" <<G4endl;
00140     track.GetDynamicParticle()->DumpInfo();
00141     G4cout << " in Material  " << track.GetMaterial()->GetName() <<G4endl;
00142     G4cout << "MeanLifeTime = " << currentInteractionLength/ns << "[ns]" <<G4endl;
00143   }
00144 #endif
00145 
00146   return theNumberOfInteractionLengthLeft * currentInteractionLength;
00147 }

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

Implemented in G4Decay, G4Scintillation, and G4RadioactiveDecay.

Referenced by PostStepGetPhysicalInteractionLength().

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

Implemented in G4Decay, G4Scintillation, and G4RadioactiveDecay.

Referenced by AtRestGetPhysicalInteractionLength().

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

Implements G4VProcess.

Reimplemented in G4Decay, and G4Scintillation.

Definition at line 111 of file G4VRestDiscreteProcess.cc.

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

Referenced by G4Scintillation::PostStepDoIt().

00115 { 
00116 //  reset NumberOfInteractionLengthLeft
00117     ClearNumberOfInteractionLengthLeft();
00118 
00119     return pParticleChange;
00120 }

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

Implements G4VProcess.

Reimplemented in G4Decay.

Definition at line 70 of file G4VRestDiscreteProcess.cc.

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

00075 {
00076   if ( (previousStepSize < 0.0) || (theNumberOfInteractionLengthLeft<=0.0)) {
00077     // beggining of tracking (or just after DoIt of this process)
00078     ResetNumberOfInteractionLengthLeft();
00079   } else if ( previousStepSize > 0.0) {
00080     // subtract NumberOfInteractionLengthLeft 
00081     SubtractNumberOfInteractionLengthLeft(previousStepSize);
00082   } else {
00083     // zero step
00084     //  DO NOTHING
00085   }
00086 
00087   // condition is set to "Not Forced"
00088   *condition = NotForced;
00089 
00090   // get mean free path
00091   currentInteractionLength = GetMeanFreePath(track, previousStepSize, condition);
00092 
00093   G4double value;
00094   if (currentInteractionLength <DBL_MAX) {
00095     value = theNumberOfInteractionLengthLeft * currentInteractionLength;
00096   } else {
00097     value = DBL_MAX;
00098   }
00099 #ifdef G4VERBOSE
00100    if (verboseLevel>1){
00101     G4cout << "G4VRestDiscreteProcess::PostStepGetPhysicalInteractionLength ";
00102     G4cout << "[ " << GetProcessName() << "]" <<G4endl;
00103     track.GetDynamicParticle()->DumpInfo();
00104     G4cout << " in Material  " <<  track.GetMaterial()->GetName() <<G4endl;
00105     G4cout << "InteractionLength= " << value/cm <<"[cm] " <<G4endl;
00106   }
00107 #endif
00108   return value;
00109 }


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