G4VRestContinuousProcess Class Reference

#include <G4VRestContinuousProcess.hh>

Inheritance diagram for G4VRestContinuousProcess:

G4VProcess

Public Member Functions

 G4VRestContinuousProcess (const G4String &, G4ProcessType aType=fNotDefined)
 G4VRestContinuousProcess (G4VRestContinuousProcess &)
virtual ~G4VRestContinuousProcess ()
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &, G4double, G4ForceCondition *)
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)

Protected Member Functions

virtual G4double GetContinuousStepLimit (const G4Track &aTrack, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)=0
void SetGPILSelection (G4GPILSelection selection)
G4GPILSelection GetGPILSelection () const
virtual G4double GetMeanLifeTime (const G4Track &aTrack, G4ForceCondition *condition)=0

Detailed Description

Definition at line 56 of file G4VRestContinuousProcess.hh.


Constructor & Destructor Documentation

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

Definition at line 55 of file G4VRestContinuousProcess.cc.

References G4VProcess::enablePostStepDoIt, and G4VRestContinuousProcess().

Referenced by G4VRestContinuousProcess().

00056   : G4VProcess(aName, aType),
00057     valueGPILSelection(CandidateForSelection)
00058 {
00059   enablePostStepDoIt = false;
00060 }

G4VRestContinuousProcess::G4VRestContinuousProcess ( G4VRestContinuousProcess  ) 

Definition at line 66 of file G4VRestContinuousProcess.cc.

References G4VRestContinuousProcess().

00067                   : G4VProcess(right),
00068                     valueGPILSelection(right.valueGPILSelection)
00069 {
00070 }

G4VRestContinuousProcess::~G4VRestContinuousProcess (  )  [virtual]

Definition at line 62 of file G4VRestContinuousProcess.cc.

00063 {
00064 }


Member Function Documentation

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

Implements G4VProcess.

Definition at line 139 of file G4VRestContinuousProcess.cc.

References G4VProcess::pParticleChange.

00143 { 
00144     return pParticleChange;
00145 }

G4double G4VRestContinuousProcess::AlongStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4double  currentMinimumStep,
G4double currentSafety,
G4GPILSelection selection 
) [virtual]

Implements G4VProcess.

Definition at line 72 of file G4VRestContinuousProcess.cc.

References CandidateForSelection, G4DynamicParticle::DumpInfo(), G4cout, G4endl, GetContinuousStepLimit(), G4Track::GetDynamicParticle(), G4Track::GetMaterial(), G4Material::GetName(), G4VProcess::GetProcessName(), and G4VProcess::verboseLevel.

00079 {
00080   // GPILSelection is set to defaule value of CandidateForSelection
00081   valueGPILSelection = CandidateForSelection;
00082 
00083   // get Step limit proposed by the process
00084   G4double steplength = GetContinuousStepLimit(track,previousStepSize,currentMinimumStep, currentSafety);
00085 
00086   // set return value for G4GPILSelection
00087   *selection = valueGPILSelection;
00088 #ifdef G4VERBOSE
00089    if (verboseLevel>1){
00090     G4cout << "G4VRestContinuousProcess::AlongStepGetPhysicalInteractionLength ";
00091     G4cout << "[ " << GetProcessName() << "]" <<G4endl;
00092     track.GetDynamicParticle()->DumpInfo();
00093     G4cout << " in Material  " <<  track.GetMaterial()->GetName() <<G4endl;
00094     G4cout << "IntractionLength= " << steplength/cm <<"[cm] " <<G4endl;
00095   }
00096 #endif
00097    return  steplength ;
00098 }

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

Implements G4VProcess.

Definition at line 128 of file G4VRestContinuousProcess.cc.

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

00132 {
00133 //  clear NumberOfInteractionLengthLeft
00134     ClearNumberOfInteractionLengthLeft();
00135 
00136     return pParticleChange;
00137 }

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

Implements G4VProcess.

Definition at line 100 of file G4VRestContinuousProcess.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.

00104 {
00105   // beggining of tracking 
00106   ResetNumberOfInteractionLengthLeft();
00107 
00108   // condition is set to "Not Forced"
00109   *condition = NotForced;
00110 
00111   // get mean life time
00112   currentInteractionLength = GetMeanLifeTime(track, condition);
00113 
00114 #ifdef G4VERBOSE
00115    if ((currentInteractionLength <0.0) || (verboseLevel>2)){
00116     G4cout << "G4VRestContinuousProcess::AtRestGetPhysicalInteractionLength ";
00117     G4cout << "[ " << GetProcessName() << "]" <<G4endl;
00118     track.GetDynamicParticle()->DumpInfo();
00119     G4cout << " in Material  " << track.GetMaterial()->GetName() <<G4endl;
00120     G4cout << "MeanLifeTime = " << currentInteractionLength/ns << "[ns]" <<G4endl;
00121   }
00122 #endif
00123 
00124   return theNumberOfInteractionLengthLeft * currentInteractionLength;
00125 }

virtual G4double G4VRestContinuousProcess::GetContinuousStepLimit ( const G4Track aTrack,
G4double  previousStepSize,
G4double  currentMinimumStep,
G4double currentSafety 
) [protected, pure virtual]

Referenced by AlongStepGetPhysicalInteractionLength().

G4GPILSelection G4VRestContinuousProcess::GetGPILSelection (  )  const [inline, protected]

Definition at line 125 of file G4VRestContinuousProcess.hh.

00125 {return valueGPILSelection;};

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

Referenced by AtRestGetPhysicalInteractionLength().

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

Implements G4VProcess.

Definition at line 101 of file G4VRestContinuousProcess.hh.

00104                               {return 0;};

virtual G4double G4VRestContinuousProcess::PostStepGetPhysicalInteractionLength ( const G4Track ,
G4double  ,
G4ForceCondition  
) [inline, virtual]

Implements G4VProcess.

Definition at line 94 of file G4VRestContinuousProcess.hh.

00098                              { return -1.0; };

void G4VRestContinuousProcess::SetGPILSelection ( G4GPILSelection  selection  )  [inline, protected]

Definition at line 122 of file G4VRestContinuousProcess.hh.

00123     { valueGPILSelection = selection;};


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