BLRunManager.hh

Go to the documentation of this file.
00001 //      BLRunManager.hh
00002 /*
00003 This source file is part of G4beamline, http://g4beamline.muonsinc.com
00004 Copyright (C) 2003,2004,2005,2006 by Tom Roberts, all rights reserved.
00005 
00006 This program is free software; you can redistribute it and/or
00007 modify it under the terms of the GNU General Public License
00008 as published by the Free Software Foundation; either version 2
00009 of the License, or (at your option) any later version.
00010 
00011 This program is distributed in the hope that it will be useful,
00012 but WITHOUT ANY WARRANTY; without even the implied warranty of
00013 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00014 GNU General Public License for more details.
00015 
00016 http://www.gnu.org/copyleft/gpl.html
00017 */
00018 
00019 #ifndef BLRUNMANAGER_H
00020 #define BLRUNMANAGER_H
00021 
00022 #include <vector>
00023 #include <map>
00024 #include <setjmp.h>
00025 
00026 #include "G4RunManager.hh"
00027 #include "G4EventManager.hh"
00028 #include "G4TrackingManager.hh"
00029 #include "G4StackManager.hh"
00030 #include "G4SteppingManager.hh"
00031 #include "G4Navigator.hh"
00032 #include "G4PropagatorInField.hh"
00033 
00034 #include "BLAssert.hh"
00035 #include "BLCollectiveComputation.hh"
00036 #include "BLCallback.hh"
00037 
00038 
00039 /**     class BLRunManager -- enhanced version of G4RunManager.
00040  *
00041  *      There are three basic modes of BLRunManager: normal, collective,
00042  *      and enhanced.
00043  *
00044  *      Normal mode:
00045  *      This class inherents from G4RunManager, and can thus be used
00046  *      just like it, for standard control of the Geant4 simulation
00047  *      (via function BeamOn()).
00048  *
00049  *      Collective mode:
00050  *
00051  *      In collective mode a vector<G4TrackData*> is generated from the
00052  *      input file's beam commands, and then all tracks are stepped for
00053  *      one time step after which the BLCollectiveComputation::collectiveStep()
00054  *      is called; this is continued until all tracks are killed. Any
00055  *      secondaries generated can either be kept or discarded (see
00056  *      setCollectiveMode() below). The TrackID-s of primary tracks are
00057  *      preserved, and TrackIDs for secondary tracks are assigned sequentially
00058  *      from 1001, ignoring EventID.
00059  *
00060  *      Collective mode is selected by calling setCollectiveMode(), and
00061  *      then behaving externally just like normal mode (i.e. just call
00062  *      BeamOn()).
00063  *
00064  *      Enhanced mode:
00065  *      This class also implements functions that permit enhanced control
00066  *      of the simulation -- in particular there is a processOneTrack()
00067  *      function that simply processes a single track.
00068  *
00069  *      NOTE: the word "track" is ambiguous, and is commonly used as both
00070  *      a noun and a verb. Here it is a noun only (referring to a G4Track),
00071  *      and the verb is "process". Processing a track consists of stepping
00072  *      it through the simulation until it is killed or suspended, while
00073  *      holding all secondaries it generates. The secondaries are themselves
00074  *      tracks, and must be dealt with before another call to processOneTrack()
00075  *      (see below).
00076  *
00077  *      Enhanced run control: All of the usual G4beamline initialization
00078  *      must occur, so the following code is best put into a BLCallback
00079  *      function registered to replace the main loop -- that means it gets
00080  *      called after the post-Reference callbacks are called, and that after
00081  *      all replace-main-loop callbacks are called the program exits (thus
00082  *      avoiding the main loop).
00083  *
00084  *              BLRunManager *runmgr = BLRunManager::getObject();
00085  *              runmgr->beginRun(int runid);
00086  *              --- loop over events ---
00087  *                      runmgr->beginEvent(int evid);
00088  *                      --- loop over tracks ---
00089  *                              --- get a track
00090  *                              runmgr->processOneTrack(track);
00091  *                              --- handle secondaries (see below)
00092  *                      runmgr->endEvent();
00093  *              runmgr->endRun();
00094  *              // write out NTuples, etc. (if appropriate)
00095  *
00096  *      There can be a loop over runs, if necessary.
00097  *
00098  *      Note that these enhanced run control routines do not create or handle
00099  *      trajectories, and ignore them if present. Hence no visualization
00100  *      of tracks is possible.
00101  *
00102  *      Note that secondaries are tracks, but they have not yet been assigned
00103  *      either a trackid or a parentid. Those should be assigned before the
00104  *      track is deferred, as all deferred tracks should have those id-s.
00105  *
00106  *      Handling secondaries -- do one of:
00107  *              runmgr->discardAllSecondaries();
00108  *      or
00109  *              runmgr->processAllSecondariesAndDeferredTracks();
00110  *      or
00111  *              runmgr->deferAllSecondaries();
00112  *              // followed sometime later by:
00113  *              runmgr->processAllDeferredTracksAndTheirSecondaries();
00114  *      or
00115  *              while((track=runmgr->popOneSecondary()) != 0) {
00116  *                      // NOTE: you cannot call processOneTrack() here!
00117  *                      // -- that would discard  any remaining secondaries.
00118  *                      if(--- condition ---)
00119  *                              runmgr->deferOneTrack(track);
00120  *                      else
00121  *                              delete track;
00122  *              }
00123  *              // followed sometime later by:
00124  *              runmgr->processAllDeferredTracksAndTheirSecondaries();
00125  **/
00126 class BLRunManager : public G4RunManager, public BLCallback {
00127         static BLRunManager *singleton;
00128         G4TrackingManager *trackManager;
00129         G4StackManager *stackManager;
00130         jmp_buf jmpBuf;
00131         bool validJmpBuf;
00132         std::vector<BLTrackData> trackVector;
00133         std::vector<BLCollectiveComputation *> computeVector;
00134         bool collectiveMode;
00135         bool keepSecondaries;
00136         G4double stepTime; // collective mode only
00137         G4double deltaT;       // collective mode only
00138         G4int currentTrackIndex; // collective mode only
00139         bool rejected;          // collective mode only
00140         int nextSecondaryTrackID; // collective mode only
00141         void appendTrack(G4Track *track); // collective mode only
00142 public:
00143         /// returns the (singleton) BLRunManager.
00144         /// Note that G4RunManager::GetRunManager() returns a pointr to
00145         /// the base class, not this class.
00146         static BLRunManager *getObject();
00147 
00148         /// Constructor
00149         BLRunManager();
00150 
00151         /// Destructor
00152         virtual ~BLRunManager();
00153 
00154         /// BeamOn() will track beam events, in normal or collective mode.
00155         /// See setCollectiveMode() below.
00156         virtual void BeamOn(int nEvents, const char *macroFile=0,
00157                                                         G4int n_select=-1);
00158 
00159         /// DoEventLoop() is REPLACED, so abondonCurrentEvent() works.
00160         void DoEventLoop(G4int n_event,const char* macroFile,G4int n_select);
00161 
00162         /// abandonCurrentEvent() will abandon the current event by performing
00163         /// a longjmp into DoEventLoop(). Not guaranteed to work, but probably
00164         /// does, even if called from BLAlarm::signal().
00165         /// NOTE: If the current event can be abandoned, this function never
00166         /// returns; if the event cannot be abandoned, it returns immediately.
00167         void abandonCurrentEvent();
00168 
00169         /// callback() from BLCallback.
00170         virtual void callback(int type);
00171 
00172         // Collective run control:
00173         // NONE of these routines are called in normal or enhanced mode.
00174 
00175         /// beamOnCollective() processes an entire run in collective mode.
00176         /// In collective mode, a vector of BLTrackData is constructed, and
00177         /// all tracks in it are stepped one time step at a time. This is 
00178         /// essentially an inversion of the usual loop over steps inside
00179         /// loops over tracks and events.
00180         /// All secondaries are immediately added to the end of trackVector,
00181         /// if keepSecondaries is true.
00182         /// The Tune and Reference particles are always tracked in normal mode.
00183         void beamOnCollective(int nEvents);
00184 
00185         /// registerCollectiveComputation() will register a computation to
00186         /// be called after every time step in collective mode. Multiple
00187         /// computations are called in the order they are registered.
00188         /// Must be called before beamOnCollective() or steps may be missed
00189         /// (usually called in the command() function of a BLCommand, so the
00190         /// order of commands in the input file is preserved).
00191         void registerCollectiveComputation(BLCollectiveComputation *compute)
00192                 { computeVector.push_back(compute); }
00193 
00194         /// getTrackVector() returns the trackVector.
00195         /// Note: the track vector and its tracks can be modified inside the
00196         /// collective computation (takes effect immediately).
00197         /// As BLRunManager is a singleton, the trackVector never changes.
00198         /// No track is ever removed from the track vector, but they can
00199         /// have their status set to fStopAndKill which effectively does that.
00200         std::vector<BLTrackData> &getTrackVector() { return trackVector; }
00201 
00202         /// getCurrentTrackIndex() returns the index in the track vector
00203         /// of the track currently being tracked. Returns -1 if not in
00204         /// collective mode, or if tracking is not in progress.
00205         int getCurrentTrackIndex() { return currentTrackIndex; }
00206 
00207         /// getCollectiveMode() returns true if in collective mode.
00208         bool getCollectiveMode() { return collectiveMode; }
00209 
00210         /// setCollectiveMode() arranges for calls to BeamOn() to call
00211         /// beamOnCollective() rather than the usual G4RunManager function.
00212         /// This is for BEAM only; TUNE and REFERENCE particles are always
00213         /// tracked normally via G4RunManager::BeamOn().
00214         void setCollectiveMode(bool flag=true) { collectiveMode=flag; }
00215 
00216         /// getStepTime() will get the ste time.
00217         /// NOTE: inside collectiveStep() this function returns the time
00218         /// of the current step; outside collectiveStep() it returns the
00219         /// next step time (i.e. the time at which tracks should be suspended).
00220         G4double getStepTime() { return stepTime; }
00221 
00222         /// setStepTime() will set the time for the next step in time.
00223         /// If called from within collectiveStep(), it may confuse other
00224         /// classes collectiveStep()-s. If called from beginCollectiveTracking()
00225         /// it sets the time of the first step.
00226         void setStepTime(G4double t) { stepTime = t; }
00227 
00228         /// getDeltaT() will get the time interval between time steps
00229         G4double getDeltaT() { return deltaT; }
00230 
00231         /// setDeltaT() will set the time interval between time steps.
00232         /// All tracks will be stepped to the same value of global time.
00233         /// Initiall value is 1.0*ns. deltaT <= 0.0 is an error.
00234         /// Can be called at any time; if called within collectiveStep()
00235         /// it will affect the immediately following time step.
00236         void setDeltaT(G4double dt) { BLAssert(dt>0.0); deltaT = dt; 
00237         }
00238 
00239         /// rejectCollectiveStep() will cause the run manager to abandon
00240         /// tracking (if in progress) and to discard all of the current
00241         /// tracks, replacing the vector with a copy saved immediately
00242         /// after the previous call to collectiveStep() returned.
00243         /// Thus the current collective step will be repeated as if
00244         /// no tracking had been performed since the previous return from
00245         /// collectiveStep(). For this to make sense, setDeltaT() 
00246         /// should have been called with a smaller value
00247         /// -- this is intended for a collective algorithm to dynamically
00248         /// adapt the value of deltaT.
00249         /// Note: the first step is to bring all beam tracks to a common
00250         /// value of global time, and rejecting this step is ignored.
00251         void rejectCollectiveStep() { rejected = true; }
00252 
00253 
00254 
00255         // Enhanced run control:
00256         // NONE of these routines are called in normal run-control mode.
00257         // NONE of these routines are called in collective run-control mode.
00258         // In most cases these routines should be called from a callback()
00259         // function registered to replace the main loop (exception: calling
00260         // processOneTrack() as "pre-tracking" for a collective computation).
00261 
00262         /// beginRun() begins a run in enhanced run-control mode.
00263         void beginRun(int runid=0);
00264 
00265         /// endRun() ends a run in enhanced run-control mode.
00266         void endRun();
00267 
00268         /// beginEvent() begins an event in enhanced run-control mode.
00269         void beginEvent(int evid=0);
00270 
00271         /// endEvent() ends an event in enhanced run-control mode.
00272         void endEvent();
00273 
00274         /// getNextBeamEventAndTrack() will get the next beam event and track.
00275         /// returns true if one is returned, false if none are left.
00276         /// Neithe event actions nor track actions are performed.
00277         /// Both *pevent and *ptrack must be deleted.
00278         bool getNextBeamEventAndTrack(G4Event **pevent, G4Track **ptrack);
00279 
00280         /// processOneTrack() tracks a single track until it is suspended or 
00281         /// killed.
00282         /// Note that before tracking begins it deletes any held secondaries,
00283         /// so if secondaries are not to be discarded they must be deferred
00284         /// or processed before the next call.
00285         /// You may want to call BLManager::setState(SPECIAL) to prevent
00286         /// this track from being entered into NTuples for virtualdetectors
00287         /// (etc.) encountered during processing.
00288         /// This function MUST NOT BE CALLED if any track is being processed
00289         /// via any means (e.g. the normal BeamOn() is executing).
00290         void processOneTrack(G4Track *track);
00291 
00292         /// discardAllSecondaries() will discard all secondaries.
00293         void discardAllSecondaries();
00294 
00295         /// deferAllSecondaries() will defer all secondaries until
00296         /// deferred tracks are processed. Returns # tracks deferred.
00297         int deferAllSecondaries(int secondaryid=10000, int parentid=-1);
00298 
00299         /// deferOneTrack() will defer a single track until deferred tracks
00300         /// are processed.
00301         void deferOneTrack(G4Track *track);
00302 
00303         /// processAllSecondariesAndDeferredTracks() will defer all secondaries
00304         /// and then process all deferred tracks and their secondaries.
00305         /// returns # tracks processed.
00306         int processAllSecondariesAndDeferredTracks(int secondaryid=10000,
00307                                                         int parentid=-1) {
00308                 int first = secondaryid;
00309                 secondaryid += deferAllSecondaries(secondaryid,parentid);
00310                 secondaryid +=
00311                        processAllDeferredTracksAndTheirSecondaries(secondaryid);
00312                 return secondaryid-first;
00313         }
00314 
00315         /// processAllDeferredTracksAndTheirSecondaries() will process all
00316         /// deferred tracks, including processing their secondaries.
00317         /// Returns # tracks processed.
00318         int processAllDeferredTracksAndTheirSecondaries(int trackid=10000);
00319 
00320         /// popOneSecondary() will return a pointer to one secondary, removing
00321         /// it from the list of secondaries; returns NULL if no more.
00322         /// Intended to be called in a loop immediately after a call to
00323         /// processOneTrack(). The returned pointer should eventually be
00324         /// deleted. Order is LIFO.
00325         G4Track *popOneSecondary();
00326 
00327         /// popOneDeferredTrack() will return a pointer to one deferred
00328         /// track, removing it from the list of deferred tracks; returns
00329         /// NULL if no more. The returned pointer should eventually be
00330         /// deleted. Note that a loop calling this function can defer other
00331         /// tracks or secondaries, which simply extends the loop.
00332         /// Order is LIFO.
00333         G4Track *popOneDeferredTrack();
00334 
00335         // protected functions made public (for testing and special purposes):
00336         void RunInitialization() { G4RunManager::RunInitialization(); }
00337         void RunTermination() { G4RunManager::RunTermination(); }
00338         G4EventManager *getEventManager() { return eventManager; }
00339         G4TrackingManager *getTrackingManager() { return trackManager; }
00340         G4StackManager *getStackManager() { return stackManager; }
00341         void setCurrentEvent(G4Event *ev) { currentEvent = ev; }
00342         G4bool getRunAborted() { return runAborted; }
00343 };
00344 
00345 #endif // BLRUNMANAGER_H
g4beamline