ZStepLimiter Class Reference

List of all members.


Detailed Description

class ZStepLimiter limits the ZStep

Public Member Functions

 ZStepLimiter ()
G4VProcess * clone ()
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
G4VParticleChange * PostStepDoIt (const G4Track &track, const G4Step &stepData)
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 G4VParticleChange * AtRestDoIt (const G4Track &track, const G4Step &stepData)

Static Public Member Functions

static void setMaxStep (G4double ms)

Private Attributes

ParticleChange change

Static Private Attributes

static G4double maxStep = DBL_MAX

Classes

class  ParticleChange

Constructor & Destructor Documentation

ZStepLimiter::ZStepLimiter (  )  [inline]

Referenced by clone().

00077                        : G4VProcess("ZStepLimiter",fUserDefined), change() {
00078             G4ParticleTable::G4PTblDicIterator *theParticleIterator =
00079                         G4ParticleTable::GetParticleTable()->GetIterator();
00080             theParticleIterator->reset();
00081             while((*theParticleIterator)()) {
00082                 G4ParticleDefinition *pd = theParticleIterator->value();
00083                 G4ProcessManager *pmgr = pd->GetProcessManager();
00084                 if(!pmgr) {
00085                         printf("ZStepLimiter: particle '%s' has no ProcessManager!\n",
00086                                 pd->GetParticleName().c_str());
00087                         continue;
00088                 }
00089                 pmgr->AddProcess(this,-1,-1,ordDefault);
00090             }
00091         }


Member Function Documentation

G4VProcess* ZStepLimiter::clone (  )  [inline]

References ZStepLimiter().

00093 { return new ZStepLimiter(*this); }

static void ZStepLimiter::setMaxStep ( G4double  ms  )  [inline, static]

References maxStep.

Referenced by BLManager::PreUserTrackingAction(), and BLManager::UserSteppingAction().

00095 { maxStep = ms; }

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

References maxStep.

00100                 { *condition = NotForced;
00101                   double limit = maxStep-0.010*mm;
00102                   if(limit < 0.0) limit = DBL_MIN;
00103                   return limit;
00104                 }

G4VParticleChange* ZStepLimiter::PostStepDoIt ( const G4Track &  track,
const G4Step &  stepData 
) [inline]

References change.

00107                 { change.Initialize(track);
00108                   return &change;
00109                 }

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

00114                 { return -1.0; }

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

References change.

00117                 { change.Initialize(track);
00118                   return &change;
00119                 }

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

00122                 { return -1.0; }

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

References change.

00125                 { change.Initialize(track);
00126                   return &change;
00127                 }


Member Data Documentation

G4double ZStepLimiter::maxStep = DBL_MAX [static, private]


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