BLRunManager Class Reference

#include <BLRunManager.hh>

Inheritance diagram for BLRunManager:

BLCallback

List of all members.


Detailed Description

class BLRunManager -- enhanced version of G4RunManager.

There are three basic modes of BLRunManager: normal, collective, and enhanced.

Normal mode: This class inherents from G4RunManager, and can thus be used just like it, for standard control of the Geant4 simulation (via function BeamOn()).

Collective mode:

In collective mode a vector<G4TrackData*> is generated from the input file's beam commands, and then all tracks are stepped for one time step after which the BLCollectiveComputation::collectiveStep() is called; this is continued until all tracks are killed. Any secondaries generated can either be kept or discarded (see setCollectiveMode() below). The TrackID-s of primary tracks are preserved, and TrackIDs for secondary tracks are assigned sequentially from 1001, ignoring EventID.

Collective mode is selected by calling setCollectiveMode(), and then behaving externally just like normal mode (i.e. just call BeamOn()).

Enhanced mode: This class also implements functions that permit enhanced control of the simulation -- in particular there is a processOneTrack() function that simply processes a single track.

NOTE: the word "track" is ambiguous, and is commonly used as both a noun and a verb. Here it is a noun only (referring to a G4Track), and the verb is "process". Processing a track consists of stepping it through the simulation until it is killed or suspended, while holding all secondaries it generates. The secondaries are themselves tracks, and must be dealt with before another call to processOneTrack() (see below).

Enhanced run control: All of the usual G4beamline initialization must occur, so the following code is best put into a BLCallback function registered to replace the main loop -- that means it gets called after the post-Reference callbacks are called, and that after all replace-main-loop callbacks are called the program exits (thus avoiding the main loop).

BLRunManager *runmgr = BLRunManager::getObject(); runmgr->beginRun(int runid); --- loop over events --- runmgr->beginEvent(int evid); --- loop over tracks --- --- get a track runmgr->processOneTrack(track); --- handle secondaries (see below) runmgr->endEvent(); runmgr->endRun(); // write out NTuples, etc. (if appropriate)

There can be a loop over runs, if necessary.

Note that these enhanced run control routines do not create or handle trajectories, and ignore them if present. Hence no visualization of tracks is possible.

Note that secondaries are tracks, but they have not yet been assigned either a trackid or a parentid. Those should be assigned before the track is deferred, as all deferred tracks should have those id-s.

Handling secondaries -- do one of: runmgr->discardAllSecondaries(); or runmgr->processAllSecondariesAndDeferredTracks(); or runmgr->deferAllSecondaries(); // followed sometime later by: runmgr->processAllDeferredTracksAndTheirSecondaries(); or while((track=runmgr->popOneSecondary()) != 0) { // NOTE: you cannot call processOneTrack() here! // -- that would discard any remaining secondaries. if(--- condition ---) runmgr->deferOneTrack(track); else delete track; } // followed sometime later by: runmgr->processAllDeferredTracksAndTheirSecondaries();

Public Member Functions

 BLRunManager ()
 Constructor.
virtual ~BLRunManager ()
 Destructor.
virtual void BeamOn (int nEvents, const char *macroFile=0, G4int n_select=-1)
 BeamOn() will track beam events, in normal or collective mode. See setCollectiveMode() below.
void DoEventLoop (G4int n_event, const char *macroFile, G4int n_select)
 DoEventLoop() is REPLACED, so abondonCurrentEvent() works.
void abandonCurrentEvent ()
 abandonCurrentEvent() will abandon the current event by performing a longjmp into DoEventLoop(). Not guaranteed to work, but probably does, even if called from BLAlarm::signal(). NOTE: If the current event can be abandoned, this function never returns; if the event cannot be abandoned, it returns immediately.
virtual void callback (int type)
 callback() from BLCallback.
void beamOnCollective (int nEvents)
 beamOnCollective() processes an entire run in collective mode. In collective mode, a vector of BLTrackData is constructed, and all tracks in it are stepped one time step at a time. This is essentially an inversion of the usual loop over steps inside loops over tracks and events. All secondaries are immediately added to the end of trackVector, if keepSecondaries is true. The Tune and Reference particles are always tracked in normal mode.
void registerCollectiveComputation (BLCollectiveComputation *compute)
 registerCollectiveComputation() will register a computation to be called after every time step in collective mode. Multiple computations are called in the order they are registered. Must be called before beamOnCollective() or steps may be missed (usually called in the command() function of a BLCommand, so the order of commands in the input file is preserved).
std::vector< BLTrackData > & getTrackVector ()
 getTrackVector() returns the trackVector. Note: the track vector and its tracks can be modified inside the collective computation (takes effect immediately). As BLRunManager is a singleton, the trackVector never changes. No track is ever removed from the track vector, but they can have their status set to fStopAndKill which effectively does that.
int getCurrentTrackIndex ()
 getCurrentTrackIndex() returns the index in the track vector of the track currently being tracked. Returns -1 if not in collective mode, or if tracking is not in progress.
bool getCollectiveMode ()
 getCollectiveMode() returns true if in collective mode.
void setCollectiveMode (bool flag=true)
 setCollectiveMode() arranges for calls to BeamOn() to call beamOnCollective() rather than the usual G4RunManager function. This is for BEAM only; TUNE and REFERENCE particles are always tracked normally via G4RunManager::BeamOn().
G4double getStepTime ()
 getStepTime() will get the ste time. NOTE: inside collectiveStep() this function returns the time of the current step; outside collectiveStep() it returns the next step time (i.e. the time at which tracks should be suspended).
void setStepTime (G4double t)
 setStepTime() will set the time for the next step in time. If called from within collectiveStep(), it may confuse other classes collectiveStep()-s. If called from beginCollectiveTracking() it sets the time of the first step.
G4double getDeltaT ()
 getDeltaT() will get the time interval between time steps
void setDeltaT (G4double dt)
 setDeltaT() will set the time interval between time steps. All tracks will be stepped to the same value of global time. Initiall value is 1.0*ns. deltaT <= 0.0 is an error. Can be called at any time; if called within collectiveStep() it will affect the immediately following time step.
void rejectCollectiveStep ()
 rejectCollectiveStep() will cause the run manager to abandon tracking (if in progress) and to discard all of the current tracks, replacing the vector with a copy saved immediately after the previous call to collectiveStep() returned. Thus the current collective step will be repeated as if no tracking had been performed since the previous return from collectiveStep(). For this to make sense, setDeltaT() should have been called with a smaller value -- this is intended for a collective algorithm to dynamically adapt the value of deltaT. Note: the first step is to bring all beam tracks to a common value of global time, and rejecting this step is ignored.
void beginRun (int runid=0)
 beginRun() begins a run in enhanced run-control mode.
void endRun ()
 endRun() ends a run in enhanced run-control mode.
void beginEvent (int evid=0)
 beginEvent() begins an event in enhanced run-control mode.
void endEvent ()
 endEvent() ends an event in enhanced run-control mode.
bool getNextBeamEventAndTrack (G4Event **pevent, G4Track **ptrack)
 getNextBeamEventAndTrack() will get the next beam event and track. returns true if one is returned, false if none are left. Neithe event actions nor track actions are performed. Both *pevent and *ptrack must be deleted.
void processOneTrack (G4Track *track)
 processOneTrack() tracks a single track until it is suspended or killed. Note that before tracking begins it deletes any held secondaries, so if secondaries are not to be discarded they must be deferred or processed before the next call. You may want to call BLManager::setState(SPECIAL) to prevent this track from being entered into NTuples for virtualdetectors (etc.) encountered during processing. This function MUST NOT BE CALLED if any track is being processed via any means (e.g. the normal BeamOn() is executing).
void discardAllSecondaries ()
 discardAllSecondaries() will discard all secondaries.
int deferAllSecondaries (int secondaryid=10000, int parentid=-1)
 deferAllSecondaries() will defer all secondaries until deferred tracks are processed. Returns # tracks deferred.
void deferOneTrack (G4Track *track)
 deferOneTrack() will defer a single track until deferred tracks are processed.
int processAllSecondariesAndDeferredTracks (int secondaryid=10000, int parentid=-1)
 processAllSecondariesAndDeferredTracks() will defer all secondaries and then process all deferred tracks and their secondaries. returns # tracks processed.
int processAllDeferredTracksAndTheirSecondaries (int trackid=10000)
 processAllDeferredTracksAndTheirSecondaries() will process all deferred tracks, including processing their secondaries. Returns # tracks processed.
G4Track * popOneSecondary ()
 popOneSecondary() will return a pointer to one secondary, removing it from the list of secondaries; returns NULL if no more. Intended to be called in a loop immediately after a call to processOneTrack(). The returned pointer should eventually be deleted. Order is LIFO.
G4Track * popOneDeferredTrack ()
 popOneDeferredTrack() will return a pointer to one deferred track, removing it from the list of deferred tracks; returns NULL if no more. The returned pointer should eventually be deleted. Note that a loop calling this function can defer other tracks or secondaries, which simply extends the loop. Order is LIFO.
void RunInitialization ()
void RunTermination ()
G4EventManager * getEventManager ()
G4TrackingManager * getTrackingManager ()
G4StackManager * getStackManager ()
void setCurrentEvent (G4Event *ev)
G4bool getRunAborted ()

Static Public Member Functions

static BLRunManagergetObject ()
 returns the (singleton) BLRunManager. Note that G4RunManager::GetRunManager() returns a pointr to the base class, not this class.

Private Member Functions

void appendTrack (G4Track *track)

Private Attributes

G4TrackingManager * trackManager
G4StackManager * stackManager
jmp_buf jmpBuf
bool validJmpBuf
std::vector< BLTrackDatatrackVector
std::vector
< BLCollectiveComputation * > 
computeVector
bool collectiveMode
bool keepSecondaries
G4double stepTime
G4double deltaT
G4int currentTrackIndex
bool rejected
int nextSecondaryTrackID

Static Private Attributes

static BLRunManagersingleton = 0


Constructor & Destructor Documentation

BLRunManager::BLRunManager (  ) 

Constructor.

References BLAssert, collectiveMode, currentTrackIndex, deltaT, keepSecondaries, nextSecondaryTrackID, rejected, singleton, stackManager, stepTime, trackManager, and validJmpBuf.

Referenced by getObject().

00113                            : G4RunManager(), trackVector(), computeVector()
00114 {
00115         BLAssert(singleton==0);
00116         singleton = this;
00117 
00118         trackManager = 0;
00119         stackManager = 0;
00120         validJmpBuf = false;
00121         collectiveMode = false;
00122         keepSecondaries = true;
00123         stepTime = 0.0;
00124         deltaT = 1.0*ns;
00125         currentTrackIndex = -1;
00126         rejected = false;
00127         nextSecondaryTrackID = -9999;
00128 }

BLRunManager::~BLRunManager (  )  [virtual]

Destructor.

References computeVector, and trackVector.

00131 {
00132         trackVector.clear();
00133         computeVector.clear();
00134 }


Member Function Documentation

void BLRunManager::appendTrack ( G4Track *  track  )  [private]

References nextSecondaryTrackID, and trackVector.

Referenced by beamOnCollective().

00383 {
00384         G4int evId = currentEvent->GetEventID();
00385         if(track->GetTrackID() <= 0) {
00386                 track->SetTrackID(nextSecondaryTrackID++);
00387         }
00388 
00389         G4Event *tmp = currentEvent;
00390         currentEvent = new G4Event(*currentEvent);
00391 
00392         trackVector.push_back(BLTrackData(currentEvent,track));
00393 
00394         currentEvent = tmp;
00395 }

BLRunManager * BLRunManager::getObject (  )  [static]

returns the (singleton) BLRunManager. Note that G4RunManager::GetRunManager() returns a pointr to the base class, not this class.

References BLRunManager(), and singleton.

Referenced by BLCMDntuple::BeginOfEventAction(), Bunch::Bunch(), BLCMDtrackermode::callback(), BLCMDtrace::command(), BLCMDspacechargelw::command(), BLCMDspacecharge::command(), BLCMDcollective::command(), BLCMDtracker::fitTrack(), BLManager::handleCallbacks(), BLCMDtracker::handlePreviousTracks(), BLCMDtracker::operator()(), BLCMDtracker::PostUserTrackingAction(), BLCMDtracker::PreUserTrackingAction(), BLAlarm::sighandler(), TraceNTuple::TraceNTuple(), and TStepLimiter::TStepLimiter().

00102 {
00103         if(!singleton) {
00104                 new BLRunManager();
00105                 // startup sequencing prevents registering callback, so it
00106                 // is hard-coded in BLManager::handleCallbacks(). This is due
00107                 // to the BLManager constructor calling this function.
00108         }
00109 
00110         return singleton; 
00111 }

void BLRunManager::BeamOn ( int  nEvents,
const char *  macroFile = 0,
G4int  n_select = -1 
) [virtual]

BeamOn() will track beam events, in normal or collective mode. See setCollectiveMode() below.

References BEAM, beamOnCollective(), collectiveMode, and BLManager::getObject().

Referenced by BLManager::trackBeam(), and BLManager::trackTuneAndReferenceParticles().

00136                                                                             {
00137         if(collectiveMode && BLManager::getObject()->getState() == BEAM)
00138                 beamOnCollective(nEvents);
00139         else
00140                 G4RunManager::BeamOn(nEvents,macroFile,n_select);
00141 }

void BLRunManager::DoEventLoop ( G4int  n_event,
const char *  macroFile,
G4int  n_select 
)

DoEventLoop() is REPLACED, so abondonCurrentEvent() works.

References jmpBuf, and validJmpBuf.

00144 {
00145         if(verboseLevel>0) { timer->Start(); }
00146 
00147         G4String msg;
00148         if(macroFile!=0) { 
00149                 if(n_select<0) n_select = n_event;
00150                 msg = "/control/execute ";
00151                 msg += macroFile;
00152         } else {
00153                 n_select = -1;
00154         }
00155 
00156         int local_n_select = n_select; // argument can be clobbered by longjmp
00157 
00158         // Event loop
00159         volatile G4int i_event = 0;
00160         if(setjmp(jmpBuf) != 0) ++i_event;
00161         validJmpBuf = true;
00162         while(i_event<n_event) {
00163                 currentEvent = GenerateEvent(i_event);
00164                 eventManager->ProcessOneEvent(currentEvent);
00165                 AnalyzeEvent(currentEvent);
00166                 UpdateScoring();
00167                 if(i_event<local_n_select)
00168                         G4UImanager::GetUIpointer()->ApplyCommand(msg);
00169                 StackPreviousEvent(currentEvent);
00170                 currentEvent = 0;
00171                 if(runAborted) break;
00172                 ++i_event;
00173         }
00174         validJmpBuf = false;
00175 
00176         if(verboseLevel>0) {
00177                 timer->Stop();
00178                 G4cout << "Run terminated." << G4endl;
00179                 G4cout << "Run Summary" << G4endl;
00180                 if(runAborted) {
00181                         G4cout << "  Run Aborted after " << i_event + 1 <<
00182                                                 " events processed." << G4endl;
00183                 } else {
00184                         G4cout << "  Number of events processed : " << n_event
00185                                                                 << G4endl; 
00186                 }
00187                 G4cout << "  "  << *timer << G4endl;
00188         }
00189 }

void BLRunManager::abandonCurrentEvent (  ) 

abandonCurrentEvent() will abandon the current event by performing a longjmp into DoEventLoop(). Not guaranteed to work, but probably does, even if called from BLAlarm::signal(). NOTE: If the current event can be abandoned, this function never returns; if the event cannot be abandoned, it returns immediately.

References jmpBuf, and validJmpBuf.

Referenced by BLAlarm::sighandler().

00192 {
00193         if(validJmpBuf) {
00194                 G4Exception("BLRunManager","Event Time Limit Alarm",
00195                                                         EventMustBeAborted,"");
00196                 G4StateManager::GetStateManager()->
00197                                                 SetNewState(G4State_GeomClosed);
00198                 longjmp(jmpBuf,1);
00199         }
00200 }

void BLRunManager::callback ( int  type  )  [virtual]

callback() from BLCallback.

Reimplemented from BLCallback.

References collectiveMode, and stepTime.

Referenced by BLManager::handleCallbacks().

00203 {
00204         // Only perform actions before Tune particle
00205         if(type != 0) return;
00206 
00207         // in collectiveMode, add a TStepLimiter to every particle
00208         if(collectiveMode) {
00209                 printf("BLRunManager: adding TStepLimiter processes to all particles\n");
00210                 G4ParticleTable::G4PTblDicIterator *theParticleIterator
00211                         = G4ParticleTable::GetParticleTable()->GetIterator();
00212                 theParticleIterator->reset();
00213                 while((*theParticleIterator)()) {
00214                         G4ParticleDefinition *pd = 
00215                                         theParticleIterator->value();
00216                         if(pd->IsShortLived()) continue;
00217                         G4ProcessManager *pmgr = pd->GetProcessManager();
00218                         if(!pmgr) continue;
00219                         TStepLimiter *tsl = new TStepLimiter();
00220                         pmgr->AddProcess(tsl,ordDefault,-1,ordDefault);
00221                 }
00222                 stepTime = DBL_MAX;
00223         }
00224 }

void BLRunManager::beamOnCollective ( int  nEvents  ) 

beamOnCollective() processes an entire run in collective mode. In collective mode, a vector of BLTrackData is constructed, and all tracks in it are stepped one time step at a time. This is essentially an inversion of the usual loop over steps inside loops over tracks and events. All secondaries are immediately added to the end of trackVector, if keepSecondaries is true. The Tune and Reference particles are always tracked in normal mode.

References appendTrack(), BEAM, BLManager::BeginOfEventAction(), BLAssert, computeVector, currentTrackIndex, deltaT, BLTrackData::event, BLManager::getNextSecondaryTrackID(), BLManager::getObject(), TStepLimiter::getPrevStatus(), BLManager::getPrimaryTrackID(), IDLE, keepSecondaries, nextSecondaryTrackID, processOneTrack(), rejected, RunInitialization(), RunTermination(), BLManager::setEventID(), BLManager::setEventTimeLimit(), BLManager::setState(), stackManager, stepTime, BLTrackData::track, trackManager, and trackVector.

Referenced by BeamOn().

00227 {
00228         BLManager *manager = BLManager::getObject();
00229         G4TransportationManager *transportManager =
00230                         G4TransportationManager::GetTransportationManager();
00231         G4StateManager* stateManager = G4StateManager::GetStateManager();
00232         numberOfEventToBeProcessed = nEvents;
00233 
00234         if(!ConfirmBeamOnCondition()) {
00235                 G4Exception("BLRunManager","Cannot run beam",FatalException,"");
00236         }
00237         printf("================== In Collective Mode ==================\n");
00238         ConstructScoringWorlds();
00239         RunInitialization();
00240         G4ApplicationState currentState = stateManager->GetCurrentState();
00241         if(currentState != G4State_GeomClosed) {
00242                 G4Exception("BLRunManager","Geometry not closed",FatalException,
00243                                                                         "");
00244         }
00245 
00246         manager->setState(BEAM);
00247         stateManager->SetNewState(G4State_EventProc);
00248         trackManager = 0; // unused in collective mode
00249         stackManager = eventManager->GetStackManager();
00250         runAborted = false;
00251 
00252         // the per-event time limit makes no sense in collective mode
00253         manager->setEventTimeLimit(-1); // infinite
00254 
00255         // event loop to create trackVector
00256         G4PrimaryTransformer *transformer = new G4PrimaryTransformer();
00257         G4Navigator* navigator = transportManager->GetNavigatorForTracking();
00258         stepTime = -DBL_MAX;
00259         int nev;
00260         for(nev=0; nev<nEvents; ++nev) {
00261                 currentEvent = GenerateEvent(nev);
00262                 if(!currentEvent || runAborted) break;
00263                 stackManager->PrepareNewEvent();
00264                 manager->BeginOfEventAction(currentEvent);
00265                 G4TrackVector *tv=transformer->GimmePrimaries(currentEvent,0);
00266                 int evId=currentEvent->GetEventID();
00267                 if(evId <= 0) {
00268                         evId = nev;
00269                         currentEvent->SetEventID(evId);
00270                 }
00271                 nextSecondaryTrackID = manager->getNextSecondaryTrackID();
00272                 for(unsigned j=0; j<tv->size(); ++j) {
00273                         G4Track *track = (*tv)[j];
00274                         BLAssert(j==0); // cannot handle multiple primaries 
00275                         track->SetTrackID(manager->getPrimaryTrackID());
00276                         appendTrack(track);
00277                         G4double t = track->GetGlobalTime();
00278                         if(t > stepTime) stepTime = t;
00279                 }
00280                 if(tv->size() == 0) delete currentEvent;
00281                 tv->clear();
00282                 if(runAborted) break;
00283                 transportManager->SetNavigatorForTracking(navigator);
00284                 manager->setEventID(evId+1);
00285         }
00286 
00287         printf("=========== Collective: %d Events, %ld Tracks ==============\n",
00288                                 nev,(long)trackVector.size());
00289 
00290         for(unsigned i=0; i<computeVector.size(); ++i)
00291                 computeVector[i]->beginCollectiveTracking(trackVector);
00292 
00293         // Step loop
00294         runAborted = false;
00295         int nActive;
00296         bool first=true;
00297         G4double saveStepTime = stepTime;
00298         std::vector<G4Track*> saveVector;
00299         do {
00300                 if(runAborted) goto RunAborted;
00301                 nActive = 0;
00302                 rejected = false;
00303                 // note: secondaries can be appended to trackVector during loop
00304                 for(unsigned i=0; i<trackVector.size(); ++i) {
00305                         if(rejected && !first) break;
00306                         BLTrackData *td = &trackVector[i];
00307                         G4Track *track = td->track;
00308                         G4TrackStatus trackStatus = track->GetTrackStatus();
00309 //printf("Collective Track loop: ev=%d trk=%d status=%d fAlive=%d fSuspend=%d\n",td->event->GetEventID(),track->GetTrackID(),(int)trackStatus, fAlive,fSuspend);
00310                         if(trackStatus != fAlive && 
00311                                                 trackStatus != fStopButAlive)
00312                                 continue;
00313                         ++nActive;
00314                         currentEvent = td->event;
00315                         manager->setEventID(td->event->GetEventID());
00316                         currentTrackIndex = i;
00317                         for(unsigned j=0; j<computeVector.size(); ++j)
00318                                 computeVector[j]->beginTrack(trackVector,i);
00319                         if(runAborted) goto RunAborted;
00320                         processOneTrack(track);
00321                         currentTrackIndex = -1;
00322                         if(td->event->IsAborted())
00323                                 track->SetTrackStatus(fKillTrackAndSecondaries);
00324                         if(runAborted) goto RunAborted;
00325                         trackStatus = track->GetTrackStatus();
00326                         if(trackStatus == fSuspend)
00327                                 track->SetTrackStatus(TStepLimiter::getPrevStatus());
00328                         G4TrackVector *tv = trackManager->GimmeSecondaries();
00329                         if(keepSecondaries && trackStatus != fKillTrackAndSecondaries) {
00330                                 for(unsigned j=0; j<tv->size(); ++j)
00331                                         appendTrack((*tv)[j]);
00332                                 tv->clear();
00333                         } else {
00334                                 for(unsigned j=0; j<tv->size(); ++j)
00335                                         delete (*tv)[j];
00336                                 tv->clear();
00337                         }
00338                 }
00339                 if(runAborted) goto RunAborted;
00340                 // don't call collectiveStep() if rejected during tracking
00341                 if(!rejected || first) {
00342                         for(unsigned i=0; i<computeVector.size(); ++i)
00343                                 computeVector[i]->collectiveStep(trackVector);
00344                 }
00345                 if(rejected && !first) {
00346                         // discard current tracks, and restore saved ones
00347                         unsigned n=saveVector.size();
00348                         // first, discard any added secondaries
00349                         while(trackVector.size() > n) {
00350                                 delete trackVector.back().track;
00351                                 delete trackVector.back().event;
00352                                 trackVector.pop_back();
00353                         }
00354                         BLAssert(trackVector.size() == saveVector.size());
00355                         for(unsigned i=0; i<n; ++i) {
00356                                 delete trackVector[i].track;
00357                                 trackVector[i].track = saveVector[i];
00358                         }
00359                         stepTime = saveStepTime;
00360                 } else {
00361                         // delete saved tracks, and save current ones
00362                         for(unsigned i=0; i<saveVector.size(); ++i)
00363                                 delete saveVector[i];
00364                         saveVector.clear();
00365                         for(unsigned i=0; i<trackVector.size(); ++i)
00366                                 saveVector.push_back(
00367                                         new G4Track(*trackVector[i].track));
00368                         saveStepTime = stepTime;
00369                         stepTime += deltaT;
00370                 }
00371                 first = false;
00372         } while(nActive > 0);
00373 
00374 RunAborted:
00375         for(unsigned i=0; i<computeVector.size(); ++i)
00376                 computeVector[i]->endCollectiveTracking(trackVector);
00377 
00378         RunTermination();
00379         manager->setState(IDLE);
00380 }

void BLRunManager::registerCollectiveComputation ( BLCollectiveComputation compute  )  [inline]

registerCollectiveComputation() will register a computation to be called after every time step in collective mode. Multiple computations are called in the order they are registered. Must be called before beamOnCollective() or steps may be missed (usually called in the command() function of a BLCommand, so the order of commands in the input file is preserved).

References computeVector.

Referenced by BLCMDtrace::command(), BLCMDspacechargelw::command(), BLCMDspacecharge::command(), and BLCMDcollective::command().

00192                 { computeVector.push_back(compute); }

std::vector<BLTrackData>& BLRunManager::getTrackVector (  )  [inline]

getTrackVector() returns the trackVector. Note: the track vector and its tracks can be modified inside the collective computation (takes effect immediately). As BLRunManager is a singleton, the trackVector never changes. No track is ever removed from the track vector, but they can have their status set to fStopAndKill which effectively does that.

References trackVector.

00200 { return trackVector; }

int BLRunManager::getCurrentTrackIndex (  )  [inline]

getCurrentTrackIndex() returns the index in the track vector of the track currently being tracked. Returns -1 if not in collective mode, or if tracking is not in progress.

References currentTrackIndex.

00205 { return currentTrackIndex; }

bool BLRunManager::getCollectiveMode (  )  [inline]

void BLRunManager::setCollectiveMode ( bool  flag = true  )  [inline]

setCollectiveMode() arranges for calls to BeamOn() to call beamOnCollective() rather than the usual G4RunManager function. This is for BEAM only; TUNE and REFERENCE particles are always tracked normally via G4RunManager::BeamOn().

References collectiveMode.

Referenced by BLCMDspacechargelw::command(), BLCMDspacecharge::command(), BLCMDcollective::command(), and BLManager::trackTuneAndReferenceParticles().

00214 { collectiveMode=flag; }

G4double BLRunManager::getStepTime (  )  [inline]

getStepTime() will get the ste time. NOTE: inside collectiveStep() this function returns the time of the current step; outside collectiveStep() it returns the next step time (i.e. the time at which tracks should be suspended).

References stepTime.

Referenced by Bunch::addFieldValue(), TStepLimiter::AtRestGetPhysicalInteractionLength(), BLCMDcollective::collectiveStep(), Bunch::computeField(), Bunch::init(), TStepLimiter::PostStepGetPhysicalInteractionLength(), and Bunch::resizeGrid().

00220 { return stepTime; }

void BLRunManager::setStepTime ( G4double  t  )  [inline]

setStepTime() will set the time for the next step in time. If called from within collectiveStep(), it may confuse other classes collectiveStep()-s. If called from beginCollectiveTracking() it sets the time of the first step.

References stepTime.

00226 { stepTime = t; }

G4double BLRunManager::getDeltaT (  )  [inline]

getDeltaT() will get the time interval between time steps

References deltaT.

Referenced by Bunch::addFieldValue(), Bunch::computeField(), and Bunch::resizeGrid().

00229 { return deltaT; }

void BLRunManager::setDeltaT ( G4double  dt  )  [inline]

setDeltaT() will set the time interval between time steps. All tracks will be stepped to the same value of global time. Initiall value is 1.0*ns. deltaT <= 0.0 is an error. Can be called at any time; if called within collectiveStep() it will affect the immediately following time step.

References BLAssert, and deltaT.

Referenced by BLCMDspacechargelw::beginCollectiveTracking(), BLCMDspacecharge::beginCollectiveTracking(), BLCMDspacechargelw::command(), BLCMDspacecharge::command(), and BLCMDcollective::command().

00236                                     { BLAssert(dt>0.0); deltaT = dt; 
00237         }

void BLRunManager::rejectCollectiveStep (  )  [inline]

rejectCollectiveStep() will cause the run manager to abandon tracking (if in progress) and to discard all of the current tracks, replacing the vector with a copy saved immediately after the previous call to collectiveStep() returned. Thus the current collective step will be repeated as if no tracking had been performed since the previous return from collectiveStep(). For this to make sense, setDeltaT() should have been called with a smaller value -- this is intended for a collective algorithm to dynamically adapt the value of deltaT. Note: the first step is to bring all beam tracks to a common value of global time, and rejecting this step is ignored.

References rejected.

00251 { rejected = true; }

void BLRunManager::beginRun ( int  runid = 0  ) 

beginRun() begins a run in enhanced run-control mode.

References BEAM, BLManager::getObject(), RunInitialization(), and BLManager::setState().

Referenced by BLCMDtrackermode::callback(), and BLCMDtracker::handlePreviousTracks().

00398 {
00399         BLManager::getObject()->setState(BEAM);
00400         ConstructScoringWorlds();
00401         RunInitialization();
00402         SetRunIDCounter(runid);
00403         runAborted = false;
00404 }

void BLRunManager::endRun (  ) 

endRun() ends a run in enhanced run-control mode.

References BLManager::getObject(), IDLE, RunTermination(), and BLManager::setState().

Referenced by BLCMDtrackermode::callback(), and BLCMDtracker::handlePreviousTracks().

00407 {
00408         RunTermination();
00409         BLManager::getObject()->setState(IDLE);
00410 }

void BLRunManager::beginEvent ( int  evid = 0  ) 

beginEvent() begins an event in enhanced run-control mode.

References BLManager::BeginOfEventAction(), BLManager::getObject(), stackManager, and trackManager.

Referenced by BLCMDtrackermode::callback(), and BLCMDtracker::handlePreviousTracks().

00413 {
00414         if(currentEvent) delete currentEvent;
00415         currentEvent = new G4Event(evid);
00416         trackManager = eventManager->GetTrackingManager();
00417         stackManager = eventManager->GetStackManager();
00418         BLManager::getObject()->BeginOfEventAction(currentEvent);
00419 }

void BLRunManager::endEvent (  ) 

endEvent() ends an event in enhanced run-control mode.

References BLManager::EndOfEventAction(), BLManager::getObject(), stackManager, and trackManager.

Referenced by BLCMDtrackermode::callback(), and BLCMDtracker::handlePreviousTracks().

00422 {
00423         BLManager::getObject()->EndOfEventAction(currentEvent);
00424         if(currentEvent) delete currentEvent;
00425         currentEvent = 0;
00426         trackManager = 0;
00427         stackManager = 0;
00428 }

bool BLRunManager::getNextBeamEventAndTrack ( G4Event **  pevent,
G4Track **  ptrack 
)

getNextBeamEventAndTrack() will get the next beam event and track. returns true if one is returned, false if none are left. Neithe event actions nor track actions are performed. Both *pevent and *ptrack must be deleted.

References BEAM, BLAssert, BLManager::getObject(), BLManager::getPrimaryParentID(), BLManager::getPrimaryTrackID(), BLManager::setState(), and SPECIAL.

00431 {
00432         static bool more=true;
00433         static int nev=0;
00434         static BLManager *manager=0;
00435         static G4PrimaryTransformer *transformer=0;
00436         if(!manager) {
00437                 manager = BLManager::getObject();
00438                 transformer = new G4PrimaryTransformer();
00439         }
00440 
00441         if(!more) return false;
00442 
00443         manager->setState(BEAM);
00444         currentEvent = GenerateEvent(nev);
00445         if(!currentEvent || runAborted) {
00446                 more = false;
00447                 return false;
00448         }
00449         *pevent = currentEvent;
00450         G4TrackVector *tv=transformer->GimmePrimaries(currentEvent,0);
00451         int evId=currentEvent->GetEventID();
00452         if(evId <= 0) {
00453                 evId = nev;
00454                 currentEvent->SetEventID(evId);
00455         }
00456         for(unsigned j=0; j<tv->size(); ++j) {
00457                 G4Track *track = (*tv)[j];
00458                 BLAssert(j==0); // cannot handle multiple primaries 
00459                 track->SetTrackID(manager->getPrimaryTrackID());
00460                 track->SetParentID(manager->getPrimaryParentID());
00461                 *ptrack = track;
00462         }
00463         tv->clear();
00464         if(runAborted) {
00465                 more = false;
00466                 return false;
00467         }
00468         ++nev;
00469         manager->setState(SPECIAL);
00470         return true;
00471 }

void BLRunManager::processOneTrack ( G4Track *  track  ) 

processOneTrack() tracks a single track until it is suspended or killed. Note that before tracking begins it deletes any held secondaries, so if secondaries are not to be discarded they must be deferred or processed before the next call. You may want to call BLManager::setState(SPECIAL) to prevent this track from being entered into NTuples for virtualdetectors (etc.) encountered during processing. This function MUST NOT BE CALLED if any track is being processed via any means (e.g. the normal BeamOn() is executing).

References BLAssert, stackManager, and trackManager.

Referenced by beamOnCollective(), BLCMDspacechargelw::beginTrack(), BLCMDtracker::fitTrack(), BLCMDtracker::operator()(), and processAllDeferredTracksAndTheirSecondaries().

00474 {
00475         if(!trackManager || !stackManager) {
00476                 trackManager = eventManager->GetTrackingManager();
00477                 stackManager = eventManager->GetStackManager();
00478         }
00479         BLAssert(track != 0);
00480         trackManager->ProcessOneTrack(track);
00481 }

void BLRunManager::discardAllSecondaries (  ) 

discardAllSecondaries() will discard all secondaries.

References trackManager.

00484 {
00485         G4TrackVector *secondaries = trackManager->GimmeSecondaries();
00486         for(unsigned j=0; j<secondaries->size(); ++j)
00487                 delete (*secondaries)[j];
00488         secondaries->clear();
00489 }

int BLRunManager::deferAllSecondaries ( int  secondaryid = 10000,
int  parentid = -1 
)

deferAllSecondaries() will defer all secondaries until deferred tracks are processed. Returns # tracks deferred.

References deferOneTrack(), and trackManager.

Referenced by processAllDeferredTracksAndTheirSecondaries(), and processAllSecondariesAndDeferredTracks().

00492 {
00493         int n=0;
00494         G4TrackVector *secondaries = trackManager->GimmeSecondaries();
00495         for(unsigned j=0; j<secondaries->size(); ++j) {
00496                 G4Track *track = (*secondaries)[j];
00497                 track->SetTrackID(secondaryid++);
00498                 track->SetParentID(parentid); 
00499                 deferOneTrack(track);
00500                 ++n;
00501         }
00502         secondaries->clear();
00503 
00504         return n;
00505 }

void BLRunManager::deferOneTrack ( G4Track *  track  ) 

deferOneTrack() will defer a single track until deferred tracks are processed.

References stackManager.

Referenced by deferAllSecondaries().

00508 {
00509         stackManager->PushOneTrack(track);
00510 }

int BLRunManager::processAllSecondariesAndDeferredTracks ( int  secondaryid = 10000,
int  parentid = -1 
) [inline]

processAllSecondariesAndDeferredTracks() will defer all secondaries and then process all deferred tracks and their secondaries. returns # tracks processed.

References deferAllSecondaries(), and processAllDeferredTracksAndTheirSecondaries().

00307                                                                          {
00308                 int first = secondaryid;
00309                 secondaryid += deferAllSecondaries(secondaryid,parentid);
00310                 secondaryid +=
00311                        processAllDeferredTracksAndTheirSecondaries(secondaryid);
00312                 return secondaryid-first;
00313         }

int BLRunManager::processAllDeferredTracksAndTheirSecondaries ( int  trackid = 10000  ) 

processAllDeferredTracksAndTheirSecondaries() will process all deferred tracks, including processing their secondaries. Returns # tracks processed.

References deferAllSecondaries(), processOneTrack(), and stackManager.

Referenced by processAllSecondariesAndDeferredTracks().

00513 {
00514         int first=secondaryid;
00515         for(;;) {
00516                 G4VTrajectory *trajectory=0;
00517                 G4Track *track = stackManager->PopNextTrack(&trajectory);
00518                 if(!track) break;
00519                 processOneTrack(track);
00520                 int parentid = track->GetTrackID();
00521                 delete track;
00522                 secondaryid += deferAllSecondaries(secondaryid,parentid);
00523         }
00524 
00525         return secondaryid-first;
00526 }

G4Track * BLRunManager::popOneSecondary (  ) 

popOneSecondary() will return a pointer to one secondary, removing it from the list of secondaries; returns NULL if no more. Intended to be called in a loop immediately after a call to processOneTrack(). The returned pointer should eventually be deleted. Order is LIFO.

References trackManager.

00529 {
00530         G4TrackVector *secondaries = trackManager->GimmeSecondaries();
00531         int i = secondaries->size();
00532         if(i == 0) return 0;
00533         G4Track *track = secondaries->back();
00534         secondaries->pop_back();
00535         return track;
00536 }

G4Track * BLRunManager::popOneDeferredTrack (  ) 

popOneDeferredTrack() will return a pointer to one deferred track, removing it from the list of deferred tracks; returns NULL if no more. The returned pointer should eventually be deleted. Note that a loop calling this function can defer other tracks or secondaries, which simply extends the loop. Order is LIFO.

References stackManager.

00539 {
00540         G4VTrajectory *trajectory=0;
00541         return stackManager->PopNextTrack(&trajectory);
00542 }

void BLRunManager::RunInitialization (  )  [inline]

Referenced by beamOnCollective(), and beginRun().

00336 { G4RunManager::RunInitialization(); }

void BLRunManager::RunTermination (  )  [inline]

Referenced by beamOnCollective(), and endRun().

00337 { G4RunManager::RunTermination(); }

G4EventManager* BLRunManager::getEventManager (  )  [inline]

Referenced by BLCMDtrackermode::callback().

00338 { return eventManager; }

G4TrackingManager* BLRunManager::getTrackingManager (  )  [inline]

References trackManager.

00339 { return trackManager; }

G4StackManager* BLRunManager::getStackManager (  )  [inline]

References stackManager.

00340 { return stackManager; }

void BLRunManager::setCurrentEvent ( G4Event *  ev  )  [inline]

Referenced by BLCMDtrace::collectiveStep().

00341 { currentEvent = ev; }

G4bool BLRunManager::getRunAborted (  )  [inline]

Referenced by BLCMDtrackermode::callback().

00342 { return runAborted; }


Member Data Documentation

BLRunManager * BLRunManager::singleton = 0 [static, private]

Referenced by BLRunManager(), and getObject().

G4TrackingManager* BLRunManager::trackManager [private]

G4StackManager* BLRunManager::stackManager [private]

jmp_buf BLRunManager::jmpBuf [private]

Referenced by abandonCurrentEvent(), and DoEventLoop().

bool BLRunManager::validJmpBuf [private]

std::vector<BLTrackData> BLRunManager::trackVector [private]

Referenced by beamOnCollective(), and BLRunManager().

G4double BLRunManager::stepTime [private]

G4double BLRunManager::deltaT [private]

bool BLRunManager::rejected [private]


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