TStepLimiter Class Reference

List of all members.


Detailed Description

class TStepLimiter limits steps based on global time.

The steps in time are not exact, as Geant4 inherently limits steps in space, not time. This class assumes the velocity of the track does not change during the step. This is approximate, but for a 50 MeV/c mu+ in a 20 MegaVolt/meter electric field, it is accurate to 4 ps for steps of 1 ns; for a 200 MeV/c mu+ it is accurate to 0.5 ps for steps of 1 ns. For better timing accuracy, use smaller steps. Obviously the accuracy will be better for tracks with beta near 1. Due to the algorithm used, these timing errors do NOT accumulate.

Public Member Functions

 TStepLimiter ()
virtual ~TStepLimiter ()
G4VProcess * clone ()
virtual G4bool IsApplicable (const G4ParticleDefinition &)
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
virtual G4VParticleChange * AlongStepDoIt (const G4Track &track, const G4Step &stepData)
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
virtual G4VParticleChange * PostStepDoIt (const G4Track &track, const G4Step &aStep)
virtual G4VParticleChange * AtRestDoIt (const G4Track &track, const G4Step &aStep)

Static Public Member Functions

static G4TrackStatus getPrevStatus ()
 getPrevStatus() returns the previous status of the most recently suspended track.

Private Attributes

ParticleChange change
BLRunManagerrunManager

Static Private Attributes

static G4TrackStatus prevStatus = fPostponeToNextEvent

Classes

class  ParticleChange


Constructor & Destructor Documentation

TStepLimiter::TStepLimiter (  )  [inline]

References BLRunManager::getObject(), and runManager.

Referenced by clone().

00034                        : G4VProcess("TStepLimiter",fUserDefined),
00035                 change() { runManager = BLRunManager::getObject(); }

virtual TStepLimiter::~TStepLimiter (  )  [inline, virtual]

00037 { }


Member Function Documentation

G4VProcess* TStepLimiter::clone (  )  [inline]

References TStepLimiter().

00039 { return new TStepLimiter(); }

static G4TrackStatus TStepLimiter::getPrevStatus (  )  [inline, static]

getPrevStatus() returns the previous status of the most recently suspended track.

References prevStatus.

Referenced by BLRunManager::beamOnCollective().

00043 { return prevStatus; }

virtual G4bool TStepLimiter::IsApplicable ( const G4ParticleDefinition &   )  [inline, virtual]

00047                 { return true; }

virtual G4double TStepLimiter::AlongStepGetPhysicalInteractionLength ( const G4Track &  track,
G4double  previousStepSize,
G4double  currentMinimumStep,
G4double &  proposedSafety,
G4GPILSelection *  selection 
) [inline, virtual]

00053                 { return -1.0; }

virtual G4VParticleChange* TStepLimiter::AlongStepDoIt ( const G4Track &  track,
const G4Step &  stepData 
) [inline, virtual]

References change.

00057                 { change.Initialize(track);
00058                   return &change;
00059                 }

virtual G4double TStepLimiter::AtRestGetPhysicalInteractionLength ( const G4Track &  track,
G4ForceCondition *  condition 
) [inline, virtual]

References BLRunManager::getStepTime(), and runManager.

00063                 { *condition = NotForced;
00064                   G4double remaining = runManager->getStepTime() -
00065                                                         track.GetGlobalTime();
00066                   if(remaining <= 0.0) remaining = DBL_MIN;
00067                   return remaining;
00068                 }

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

References BLRunManager::getStepTime(), and runManager.

00073                 { *condition = NotForced;
00074                   G4double remaining = runManager->getStepTime() -
00075                                                         track.GetGlobalTime();
00076                   if(remaining <= 0.0) remaining = DBL_MIN;
00077                   return remaining * track.GetVelocity();
00078                 }

virtual G4VParticleChange* TStepLimiter::PostStepDoIt ( const G4Track &  track,
const G4Step &  aStep 
) [inline, virtual]

References change, and prevStatus.

00082                 { prevStatus = fAlive;
00083                   change.Initialize(track);
00084                   change.ProposeTrackStatus(fSuspend);
00085                   return &change;
00086                 }

virtual G4VParticleChange* TStepLimiter::AtRestDoIt ( const G4Track &  track,
const G4Step &  aStep 
) [inline, virtual]

References change, and prevStatus.

00090                 { prevStatus = fStopButAlive;
00091                   change.Initialize(track);
00092                   change.ProposeTrackStatus(fSuspend);
00093                   return &change;
00094                 }


Member Data Documentation

G4TrackStatus TStepLimiter::prevStatus = fPostponeToNextEvent [static, private]


The documentation for this class was generated from the following file:
g4beamline