#include <BLManager.hh>
Note it is uses all the G4 user action classes, implemented as classes under BLManager.
TrackID-s: It is non-trivial to preserve TrackID-s from input files, because Geant4 considers TrackID to be an internal variable, and it assigns them in a non-customizable manner. So BLManager keeps a trackIDMap[] that converts from internal (Geant4) to external (User) TrackID-s. It also uses BLTrackInfo to associate the external IDs with the track -- the external trackID is determined either in BLBeam (when the track is created), or in BLManager::PreUserTrackingAction for secondaries.
Public Member Functions | |
void | delayedConstruction () |
delayedConstruction() performs things which must wait until all static initializers have executed (e.g. in Geant4 routines). | |
~BLManager () | |
Destructor. | |
BLManagerState | getState () |
getState() returns the current state. | |
void | setState (BLManagerState _state) |
setState sets the state | |
int | getSteppingVerbose () |
getSteppingVerbose() returns steppingVerbose. NOTE: use this during tracking instead of Param.getInt("steppingVerbose"); | |
void | setSteppingVerbose (int v) |
setSteppingVerbose() updates steppingVerbose. -- NOTE: many other classes relay on the Parameter, not the valud in this class. | |
int | getEventTimeLimit () |
getEventTimeLimit() returns the CPU time limit for events (seconds). -1 mean infinite. | |
void | setEventTimeLimit (int sec) |
setEventTimeLimit() sets the CPU time limit for events (seconds). -1 mean infinite. | |
G4int | getEventID () const |
getEventID() gets the current eventID; | |
void | setEventID (int evId) |
setEventID() sets the current eventID; | |
bool | skipEvent (int EventID) |
skipEvent() determines if this EventID should be skipped. It checks keepEventList and skipEventList. | |
void | setKeepEvent (int EventID) |
setKeepEvent() puts EventID into the keepEventList. (if keepEventList is empty, all events are kept.) | |
void | setSkipEvent (int EventID) |
setSkipEvent() puts EventID into the skipEventList. | |
void | incrEventsProcessed (int eventID) |
incrEventsProcessed() will increment eventsProcessed. For special uses only (e.g.MPI). | |
bool | showAllExceptions (bool value=true) |
showAllExceptions() sets the flag that prevents thinning out multiple exception printouts. Returns the previous value. | |
void | registerSteppingAction (SteppingAction *sa) |
registerSteppingAction() registers a SteppingAction to be called for each step (regardless of state). | |
void | registerSteppingAction (G4VPhysicalVolume *physicalVol, SteppingAction *sa) |
registerSteppingAction() registers a SteppingAction to be called for each step (regardless of state), for every step involving the physicalVol. The callback will be called if physicalVol is either the pre- or post-step physical volume (once if both). If physicalVol==0 it is called every reference particle step. LIMITATION: only one callback can be registered for a given physicalVol (except 0). | |
void | registerTuneParticleStep (G4VPhysicalVolume *physicalVol, SteppingAction *sa) |
registerTuneParticleStep() registers a SteppingAction to be called for every tune particle step involving the physicalVol. The callback will be called if physicalVol is either the pre- or post-step physical volume (once if both). If physicalVol==0 it is called every tune particle step. LIMITATION: only one callback can be registered for a given physicalVol (except 0). | |
void | registerReferenceParticleStep (G4VPhysicalVolume *physicalVol, SteppingAction *sa) |
registerReferenceParticleStep() registers a SteppingAction to be called for every reference particle step involving the physicalVol. The callback will be called if physicalVol is either the pre- or post-step physical volume (once if both). If physicalVol==0 it is called every reference particle step. LIMITATION: only one callback can be registered for a given physicalVol (except 0). | |
void | registerBeamStep (G4VPhysicalVolume *physicalVol, SteppingAction *sa) |
registerBeamStep() registers a SteppingAction to be called for every beam step involving the physicalVol. (beam particles are everything except reference and tune.) The callback will be called if physicalVol is either the pre- or post-step physical volume (once if both). if physicalVol==0 it is called every beam step. LIMITATION: only one callback can be registered for a given physicalVol (except 0). | |
void | registerZStep (G4double z, ZSteppingAction *sa, G4int when=7) |
registerZStep() will force a step to occur near the given z position, and will call the ZSteppingAction for it, interpolating to the desired z value (Centerline coords). when is a bitwise OR of 1=tune, 2=reference, 4=beam. | |
void | registerStackingAction (StackingAction *sa) |
registerStackingAction() registers a StackingAction to be called by the Geant4 stacking manager. | |
void | setSteppingFormat () |
setSteppingFormat() sets the verbose printing format according to parameter "steppingFormat". | |
G4String | getFormatHelp () |
getFormatHelp() returns a string with help text about valid format items. | |
void | appendVerboseFormat (G4String fmt) |
appendVerboseFormat() appends fmt to the format for printing when param steppingVerbose is nonzero | |
void | steppingVerbosePrint (const G4Step *step, const G4Track *track, int header=0) |
steppingVerbosePrint() will print this step according to the current verboseFormat. Prints header if header != 0). | |
void | setSteppingManager (G4SteppingManager *p) |
setSteppingManager() sets the pointer to the current G4SteppingManager. | |
G4SteppingManager * | getSteppingManager () |
getSteppingManager() returns a pointer to the current G4SteppingManager. | |
void | setTrackingManager (G4TrackingManager *p) |
setTrackingManager() sets the pointer to the current G4TrackingManager. | |
void | initialize () |
initialize() will initialize the BLManager object, and the geant4 kernel, thus constructing the geometry in the world group. Note that registerPhysics() must be called before initialize() (normally done by a "physics" command in the input file). | |
void | trackTuneAndReferenceParticles () |
trackTuneAndReferenceParticles() will generate and track the tune particle and then the reference particle. | |
void | trackBeam () |
trackBeam() will generate the defined beam and track each event. | |
void | displayVisual () |
displayVisual() will display the detector visually. | |
void | displayGeometry (G4VPhysicalVolume *phys=0, int level=0) |
displayGeometry() will display the geant4 geometry. This is a hierarchical ASCII listing of all volumes. if phys==0 then use the worldPhysicalVolume. | |
void | registerPhysics (BLPhysics *_physics) |
registerPhysics() registers the BLPhysics object. It also sets the physics list to the BLRunManager, so following commands can find particle by name. | |
BLPhysics * | getPhysics () |
getPhysics returns the registered BLPhysics object. | |
void | registerBeam (BLBeam *_beam) |
registerBeam() registers a BLBeam object for beam generation. Multiple BLBeam objects can be registered, used in order. | |
void | registerReference (BLBeam *_beam) |
registerReference() registers a BLBeam object for reference particle generation. Multiple BLBeam objects can be registered, used in order. | |
int | nReference () |
nReference() returns the number of reference particles registered. | |
void | registerRunAction (RunAction *a, G4bool beamOnly=true) |
registerRunAction() registers a UserRunAction. If beamOnly is true (the default), the callback is made only if state==BEAM. | |
void | registerEventAction (EventAction *a, G4bool beamOnly=true) |
registerEventAction() registers a UserEventAction. If beamOnly is true (the default), the callback is made only if state==BEAM. | |
void | registerTrackingAction (TrackingAction *a, bool front=false) |
registerTrackingAction() registers a UserTackingAction. By default puts new entry at the back end of the vector; if front is true, puts new entry at the front. | |
void | registerUserCode (BLUserCode *instance) |
registerUserCode() registers a BLUserCode instance. | |
std::vector< BLUserCode * > | getUserCodeInstances (G4String type) |
getUserCodeInstances() returns a vector of all registered instances of BLUserCode with the specified type. | |
void | registerCallback (BLCallback *cb, int type) |
registerCallback() registers a BLCallback. type=0 for pre-Tune particle, type=1 for post-Reference (pre-beam tracking), type=2 for post-beam tracking. type=3 for replacing the main program loop. type=4 for visualization. NOTE: if there are type=3 callbacks, when the last one returns the closes up by summarizing NTuples and callng handleCallbacks(2), and then the program exits. This prevents the main program loop from executing. type=3 callbacks are called just after the type=1 callbacks (i.e. after the post-Reference callbacks). | |
void | handleCallbacks (int type) |
handleCallbacks() calls all applicable registered callbacks. type=0 for pre-reference particle, type=1 for post-center (pre-beam tracking), type=2 for post-beam tracking (just before program exit). type=3 for replacing the main program loop. type=4 for visualization. | |
G4VPhysicalVolume * | getWorldPhysicalVolume () |
getWorldPhysicalVolume() returns a pointer to the world PV. Note that it must already have been consructed, so this function returns NULL before construct() is called (by main). | |
void | setPRNGSeedMethod (PRNGSeedMethod method) |
setPRNGSeedMethod() will set the method used to seed the pseudo random number generator at the start of each event. | |
PRNGSeedMethod | getPRNGSeedMethod () |
getPRNGSeedMethod() will return the method used to seed the pseudo random number generator at the start of each event. | |
void | UserSteppingAction (const G4Step *step) |
UserSteppingAction() from G4UserSteppingAction. | |
G4VPhysicalVolume * | Construct () |
Construct() from G4VUserDetectorConstruction. | |
void | PreUserTrackingAction (const G4Track *track) |
PreUserTrackingAction() from G4UserTrackingAction. | |
void | PostUserTrackingAction (const G4Track *track) |
PostUserTrackingAction() from G4UserTrackingAction. | |
void | BeginOfRunAction (const G4Run *run) |
BeginOfRunAction() from G4UserRunAction. | |
void | EndOfRunAction (const G4Run *run) |
EndOfRunAction() from G4UserRunAction. | |
void | BeginOfEventAction (const G4Event *event) |
BeginOfEventAction() from G4UserEventAction. | |
void | EndOfEventAction (const G4Event *event) |
EndOfEventAction() from G4UserEventAction. | |
void | GeneratePrimaries (G4Event *event) |
GeneratePrimaries() from G4VUserPrimaryGeneratorAction. | |
G4ClassificationOfNewTrack | ClassifyNewTrack (const G4Track *) |
ClassifyNewTrack() from G4UserStackingAction. | |
void | NewStage () |
NewStage() from G4userStackingAction. | |
void | PrepareNewEvent () |
PrepareNewEvent() from G4StackingAction. | |
void | clearTrackIDMap () |
clearTrackIDMap() clears the TrackID map. | |
void | setNextSecondaryTrackID (int next) |
setNextSecondaryTrackID() sets the external TrackID for the next secondary track. Automatically incremented for subsequent secondaries. | |
int | getNextSecondaryTrackID () |
int | getExternalTrackID (const G4Track *track) |
getExternalTrackID() returns the external TrackID for the given track. In collective mode, internal and external trackID-s are the same. | |
int | getExternalParentID (const G4Track *track) |
getExternalParentID() returns the external ParentID for the given track. In collective mode, internal and external trackID-s are the same. | |
void | setExternalTrackID (G4Track *track, int trackID, int parentID) |
setExternalTrackID() will set the external trackID and parentID. if trackID<0 uses nextSecondarTrackID++. | |
int | getPrimaryTrackID () |
getPrimaryTrackID() returns the primaryTrackID set by the beam command | |
int | getPrimaryParentID () |
getPrimaryParentID() returns the primaryParentID set by the beam command | |
void | setPrimaryTrackID (int t, int p) |
setPrimaryTrackID() sets track and parent IDs for a primary track. | |
G4bool | Notify (const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description) |
Notify() from G4VExceptionHandler. | |
void | exceptionSummary () |
exceptionSummary() prints a summary of all exceptions | |
Static Public Member Functions | |
static BLManager * | getObject () |
getObject() will return a pointer to the single BLManager object, creating it if necessary. Does only the immediate constructor, not delayedConstruction(). that means it is Ok to register capabilities, but not much else. | |
static bool | isInitialized () |
isInitialized() returns true if the BLManager has been initialized. | |
static void | gsl_error_handler (const char *reason, const char *file, int lineno, int gsl_errno) |
GSL error handler. | |
Private Member Functions | |
BLManager () | |
private constructor -- immediate construction only. | |
void | insertZStep (std::vector< ZStep > &vector, G4double z, ZSteppingAction *action) |
Private Attributes | |
BLRunManager * | runManager |
BLManagerState | state |
G4int | steppingVerbose |
unsigned int | beamIndex |
G4int | histoUpdate |
time_t | startRun |
time_t | startEvent |
G4int | eventTimeLimit |
BLPhysics * | physics |
G4VPhysicalVolume * | worldPhysicalVolume |
G4int | eventID |
G4int | trackID |
const G4Track * | currentTrack |
G4bool | allExceptions |
G4int | fatalExceptions |
G4int | eventsAborted |
G4int | stuckTracks |
G4int | warnings |
G4int | prevEventID |
G4int | eventsProcessed |
bool | endRun |
PRNGSeedMethod | seedMethod |
std::set< int > | keepEventList |
std::set< int > | skipEventList |
G4SteppingManager * | fpSteppingManager |
G4TrackingManager * | fpTrackingManager |
std::vector< BLBeam * > | beamVector |
std::vector< BLBeam * > | referenceVector |
std::vector< RunAction * > | runActionVector |
std::vector< RunAction * > | beamRunActionVector |
std::vector< EventAction * > | eventActionVector |
std::vector< EventAction * > | beamEventActionVector |
std::vector< TrackingAction * > | trackingActionVector |
std::vector< BLCallback * > | preReferenceCallbackVector |
std::vector< BLCallback * > | postReferenceCallbackVector |
std::vector< BLCallback * > | postTrackingCallbackVector |
std::vector< BLCallback * > | replaceMainLoopCallbackVector |
std::vector< BLCallback * > | visualizationCallbackVector |
std::vector< SteppingAction * > | allStepVector |
std::map< G4VPhysicalVolume *, SteppingAction * > | allStepMap |
std::map< G4VPhysicalVolume *, SteppingAction * > | tpStepMap |
std::map< G4VPhysicalVolume *, SteppingAction * > | rpStepMap |
std::vector< SteppingAction * > | tpStepVector |
std::vector< SteppingAction * > | rpStepVector |
std::map< G4VPhysicalVolume *, SteppingAction * > | beamStepMap |
std::vector< SteppingAction * > | beamStepVector |
std::vector< int > | verboseFormat |
G4double | zTolerance |
std::vector< ZStep > | tuneZStep |
std::vector< ZStep > | referenceZStep |
std::vector< ZStep > | beamZStep |
std::vector< ZStep > * | currentZStep |
unsigned | indexZStep |
G4double | prevZ |
G4int | nStuckSteps |
std::vector< StackingAction * > | stackingActionVector |
int | nextSecondaryTrackID |
std::map< G4int, G4int > | trackIDMap |
int | primaryTrackID |
int | primaryParentID |
std::vector< BLUserCode * > | userCodeVector |
std::map< G4String, int > | exceptionCount |
Static Private Attributes | |
static BLManager * | blManager = 0 |
static bool | initialized = false |
Classes | |
class | EventAction |
class | PrimaryGeneratorAction |
class | RunAction |
class | StackingAction |
class | SteppingAction |
class | TrackingAction |
struct | ZStep |
class | ZSteppingAction |
BLManager::BLManager | ( | ) | [private] |
private constructor -- immediate construction only.
References allExceptions, beamIndex, beamZStep, blManager, CL, currentTrack, currentZStep, endRun, EVENT_NUMBER, eventID, eventsAborted, eventsProcessed, eventTimeLimit, fatalExceptions, BLParam::getDouble(), GLOBAL, gsl_error_handler(), histoUpdate, IDLE, indexZStep, BLAlarm::init(), BLSignal::init(), KE, NEWLINE, nextSecondaryTrackID, NSTEP, nStuckSteps, Param, physics, prevEventID, prevZ, primaryParentID, primaryTrackID, PROCESS, referenceZStep, runManager, seedMethod, startEvent, startRun, state, STEP, steppingVerbose, stuckTracks, trackID, tuneZStep, verboseFormat, VOL, warnings, worldPhysicalVolume, and zTolerance.
Referenced by getObject().
00205 : G4VUserDetectorConstruction(), 00206 G4VExceptionHandler(), keepEventList(), skipEventList(), 00207 beamVector(), referenceVector(), 00208 runActionVector(), beamRunActionVector(), 00209 eventActionVector(), beamEventActionVector(), 00210 trackingActionVector(), preReferenceCallbackVector(), 00211 postReferenceCallbackVector(), postTrackingCallbackVector(), 00212 replaceMainLoopCallbackVector(), visualizationCallbackVector(), 00213 allStepVector(), allStepMap(), tpStepMap(), rpStepMap(), 00214 tpStepVector(), rpStepVector(), 00215 beamStepMap(), beamStepVector(), verboseFormat(), 00216 tuneZStep(), referenceZStep(), beamZStep(), 00217 stackingActionVector(), trackIDMap(), userCodeVector(), 00218 exceptionCount() 00219 { 00220 // Parameter definitions (variables are all unused): 00221 BLSetParam unused_1("histoFile","g4beamline", 00222 "Default (Root) NTuple output filename"); 00223 BLSetParam unused_2("histoUpdate","0", 00224 "Output update interval (events)"); 00225 BLSetParam unused_3("viewer","none", 00226 "Visualization driver selected (default=none)"); 00227 BLSetParam unused_4("eventTimeLimit","30", 00228 "CPU Time Limit (sec)"); 00229 BLSetParam unused_5("steppingVerbose","0", 00230 "Set nonzero to print each step"); 00231 BLSetParam unused_6("steppingFormat", "N GLOBAL CL KE STEP VOL PROCESS", 00232 "Format for printing steps"); 00233 BLSetParam unused_7("zTolerance","2.0", 00234 "Tolerance for Z steps (mm)"); 00235 00236 if(blManager) 00237 G4Exception("BLManager","Object Already Exists",FatalException, 00238 ""); 00239 blManager = this; 00240 00241 state = IDLE; 00242 steppingVerbose = 0; 00243 00244 // ensure ZStep vectors have bookends 00245 tuneZStep.push_back(ZStep(-DBL_MAX,0)); 00246 tuneZStep.push_back(ZStep(DBL_MAX,0)); 00247 referenceZStep.push_back(ZStep(-DBL_MAX,0)); 00248 referenceZStep.push_back(ZStep(DBL_MAX,0)); 00249 beamZStep.push_back(ZStep(-DBL_MAX,0)); 00250 beamZStep.push_back(ZStep(DBL_MAX,0)); 00251 currentZStep = &beamZStep; 00252 indexZStep = 1; 00253 zTolerance = Param.getDouble("zTolerance"); 00254 prevZ = -DBL_MAX; 00255 00256 // default verbose format 00257 verboseFormat.push_back(NSTEP); 00258 verboseFormat.push_back(GLOBAL); 00259 verboseFormat.push_back(CL); 00260 verboseFormat.push_back(KE); 00261 verboseFormat.push_back(STEP); 00262 verboseFormat.push_back(VOL); 00263 verboseFormat.push_back(PROCESS); 00264 verboseFormat.push_back(NEWLINE); 00265 00266 runManager = 0; 00267 state = IDLE; 00268 beamIndex = 0; 00269 histoUpdate = 0; 00270 startRun = 0; 00271 startEvent = 0; 00272 eventTimeLimit = 30; 00273 physics = 0; 00274 worldPhysicalVolume = 0; 00275 eventID = 1; 00276 trackID=-9999; 00277 currentTrack = 0; 00278 allExceptions = false; 00279 fatalExceptions = 0; 00280 eventsAborted=0; 00281 stuckTracks = 0; 00282 warnings=0; 00283 prevEventID = 0; 00284 eventsProcessed = 0; 00285 endRun = false; 00286 seedMethod = EVENT_NUMBER; 00287 nStuckSteps = 0; 00288 primaryTrackID = -1; 00289 primaryParentID = -1; 00290 nextSecondaryTrackID = 1001; 00291 00292 BLSignal::init(); 00293 BLAlarm::init(); 00294 00295 G4StateManager * stateManager = G4StateManager::GetStateManager(); 00296 stateManager->SetExceptionHandler(this); 00297 gsl_set_error_handler(gsl_error_handler); 00298 }
BLManager::~BLManager | ( | ) |
Destructor.
References blManager.
00308 { 00309 // cannot delete runManager, as it deletes the detector, 00310 // steppingaction, eventaction, ..., and that becomes multiple 00311 // deletes of this object. But we are exiting, so we don't really 00312 // need to delete anything. 00313 //delete runManager; 00314 00315 blManager = 0; 00316 }
void BLManager::insertZStep | ( | std::vector< ZStep > & | vector, | |
G4double | z, | |||
ZSteppingAction * | action | |||
) | [private] |
Referenced by registerZStep().
00390 { 00391 // the bookends in vector (see BLManager()) ensure this works 00392 std::vector<ZStep>::iterator i; 00393 for(i=vector.begin(); i<vector.end(); ++i) { 00394 if(z < i->z) { // preserve order for == 00395 ZStep zs(z,action); 00396 vector.insert(i,zs); 00397 break; 00398 } 00399 } 00400 }
BLManager * BLManager::getObject | ( | ) | [static] |
getObject() will return a pointer to the single BLManager object, creating it if necessary. Does only the immediate constructor, not delayedConstruction(). that means it is Ok to register capabilities, but not much else.
References BLManager(), and blManager.
Referenced by BLTrackNTuple::appendTrack(), BLRunManager::BeamOn(), BLRunManager::beamOnCollective(), BLCMDmovie::beamStep(), BLRunManager::beginEvent(), BLRunManager::beginRun(), BLVirtualDetectorNTuple::BLVirtualDetectorNTuple(), BLCMDzntuple::callback(), BLCMDtrackermode::callback(), BLCMDtotalenergy::callback(), BLCMDshowmaterial::callback(), BLCMDreweightprocess::callback(), BLCMDgeometry::callback(), BLCMDtrackcuts::ClassifyNewTrack(), BLCMDzntuple::command(), BLCMDusertrackfilter::command(), BLCMDtune::command(), BLCMDtrackermode::command(), BLCMDtracker::command(), BLCMDtrackcuts::command(), BLCMDtrace::command(), BLCMDtotalenergy::command(), BLCMDtimentuple::command(), BLCMDspacechargelw::command(), BLCMDspacecharge::command(), BLCMDshowmaterial::command(), BLCMDsetdecay::command(), BLCMDreweightprocess::command(), BLCMDreference::command(), BLCMDrandomseed::command(), BLCMDprofile::command(), BLCMDprobefield::command(), BLCMDprintfield::command(), BLCMDprintf::command(), BLCMDphysics::command(), BLCMDparticlesource::command(), BLCMDparticlecolor::command(), BLCMDntuple::command(), BLCMDnewparticlentuple::command(), BLCMDmuminuscapturefix::command(), BLCMDmovie::command(), BLCMDlist::command(), BLCMDhelp::command(), BLCMDgeometry::command(), BLCMDg4ui::command(), BLCMDfieldntuple::command(), BLCMDfieldlines::command(), BLCMDcosmicraybeam::command(), BLCMDbug1021::command(), BLCMDbeamlossntuple::command(), BLCMDbeam::command(), BLCMDvirtualdetector::construct(), BLCMDusertrackfilter::construct(), BLCMDtubs::construct(), BLCMDtrap::construct(), BLCMDtrackerplane::construct(), BLCMDtorus::construct(), BLCMDsphere::construct(), BLCMDsolenoid::construct(), BLCMDrfdevice::construct(), BLCMDpolycone::construct(), BLCMDpillbox::construct(), BLCMDparticlefilter::construct(), BLCMDmultipole::construct(), BLCMDidealsectorbend::construct(), BLCMDhelicalharmonic::construct(), BLCMDhelicaldipole::construct(), BLCMDgenericquad::construct(), BLCMDgenericbend::construct(), BLCMDextrusion::construct(), BLCMDcorner::construct(), BLCMDbox::construct(), BLRunManager::endEvent(), BLCMDparticlecolor::EndOfEventAction(), BLRunManager::endRun(), SetDecayInstance::fatalError(), BLCMDtracker::fitTrack(), BLCMDreference::generateReferenceParticle(), BLRunManager::getNextBeamEventAndTrack(), BLCMDparticlecolor::getVisAttributes(), BLCMDtracker::handlePreviousTracks(), BLCMDparam::help(), BLCoordinates::init(), BLCMDlist::listParticles(), main(), MaterialFilter::MaterialFilter(), BLCMDparticlesource::nextBeamEvent(), BLCMDbeam::nextBeamEvent(), BLCMDtracker::operator()(), ParticleFilterPlacement::ParticleFilterPlacement(), BLCoordinates::ReferenceCoordinates::PostUserTrackingAction(), BLCMDtracker::PostUserTrackingAction(), BLCMDreference::PostUserTrackingAction(), BLCMDparticlecolor::PostUserTrackingAction(), BLCMDbeamlossntuple::PostUserTrackingAction(), BLCoordinates::ReferenceCoordinates::PreUserTrackingAction(), BLCMDtracker::PreUserTrackingAction(), BLCMDreference::PreUserTrackingAction(), BLCMDparticlecolor::PreUserTrackingAction(), BLCMDnewparticlentuple::PreUserTrackingAction(), BLCMDeventcuts::readFile(), BLCoordinates::ReferenceCoordinates::ReferenceCoordinates(), BLUserCode::registerUserCode(), BLBeam::setRandomSeedToGenerate(), BLEvaluator::setTrackVariables(), BLCoordinates::start(), steppingVerbosePrint(), Surface::Surface(), TraceNTuple::TraceNTuple(), BLCoordinates::ReferenceCoordinates::UserSteppingAction(), BLVirtualDetectorNTuple::UserSteppingAction(), BLCMDusertrackfilter::UserSteppingAction(), TrackerPlaneInstance::UserSteppingAction(), TimeNTuple::UserSteppingAction(), RFdeviceField::UserSteppingAction(), PillboxField::UserSteppingAction(), ParticleFilterPlacement::UserSteppingAction(), BLCMDzntuple::Entry::UserZSteppingAction(), BLCMDtune::UserZSteppingAction(), BLCMDtracker::UserZSteppingAction(), BLCMDreference::UserZSteppingAction(), BLCMDprofile::Entry::UserZSteppingAction(), and BLCMDprintf::UserZSteppingAction().
void BLManager::delayedConstruction | ( | ) |
delayedConstruction() performs things which must wait until all static initializers have executed (e.g. in Geant4 routines).
References runManager.
Referenced by main().
00301 { 00302 runManager = new BLRunManager(); 00303 G4StateManager * stateManager = G4StateManager::GetStateManager(); 00304 stateManager->SetExceptionHandler(this); 00305 }
BLManagerState BLManager::getState | ( | ) | [inline] |
getState() returns the current state.
References state.
Referenced by BLCMDspacechargelw::addFieldValue(), BLCMDspacecharge::addFieldValue(), BLCMDspacechargelw::beginTrack(), BLCMDspacechargelw::collectiveStep(), BLCMDparticlecolor::EndOfEventAction(), BLCMDtracker::fitTrack(), BLCMDparticlecolor::getVisAttributes(), BLCMDtracker::PostUserTrackingAction(), BLCMDspacechargelw::PostUserTrackingAction(), BLCMDreference::PostUserTrackingAction(), BLCMDparticlecolor::PostUserTrackingAction(), BLCMDmovie::PostUserTrackingAction(), BLCMDbeamlossntuple::PostUserTrackingAction(), BLCMDtracker::PreUserTrackingAction(), BLCMDspacechargelw::PreUserTrackingAction(), BLCMDspacecharge::PreUserTrackingAction(), BLCMDreference::PreUserTrackingAction(), BLCMDparticlecolor::PreUserTrackingAction(), BLCMDnewparticlentuple::PreUserTrackingAction(), BLCMDmovie::PreUserTrackingAction(), BLVirtualDetectorNTuple::UserSteppingAction(), TrackerPlaneInstance::UserSteppingAction(), BLCMDtrace::UserSteppingAction(), TimeNTuple::UserSteppingAction(), BLCMDspacechargelw::UserSteppingAction(), BLCMDmovie::UserSteppingAction(), BLCMDzntuple::Entry::UserZSteppingAction(), BLCMDprofile::Entry::UserZSteppingAction(), and BLCMDprintf::UserZSteppingAction().
00187 { return state; }
void BLManager::setState | ( | BLManagerState | _state | ) | [inline] |
setState sets the state
References state.
Referenced by BLRunManager::beamOnCollective(), BLRunManager::beginRun(), BLCMDspacechargelw::beginTrack(), BLCMDtrackermode::callback(), BLRunManager::endRun(), BLRunManager::getNextBeamEventAndTrack(), and BLCMDtracker::handlePreviousTracks().
00190 { state = _state; }
int BLManager::getSteppingVerbose | ( | ) | [inline] |
getSteppingVerbose() returns steppingVerbose. NOTE: use this during tracking instead of Param.getInt("steppingVerbose");
References steppingVerbose.
Referenced by BLCMDspacechargelw::beginTrack(), BLVirtualDetectorNTuple::BLVirtualDetectorNTuple(), BLCoordinates::init(), and BLCoordinates::start().
00194 { return steppingVerbose; }
void BLManager::setSteppingVerbose | ( | int | v | ) | [inline] |
setSteppingVerbose() updates steppingVerbose. -- NOTE: many other classes relay on the Parameter, not the valud in this class.
References steppingVerbose.
00198 { steppingVerbose = v; }
int BLManager::getEventTimeLimit | ( | ) | [inline] |
getEventTimeLimit() returns the CPU time limit for events (seconds). -1 mean infinite.
References eventTimeLimit.
00202 { return eventTimeLimit; }
void BLManager::setEventTimeLimit | ( | int | sec | ) | [inline] |
setEventTimeLimit() sets the CPU time limit for events (seconds). -1 mean infinite.
References eventTimeLimit.
Referenced by BLRunManager::beamOnCollective().
00206 { eventTimeLimit = sec; }
G4int BLManager::getEventID | ( | ) | const [inline] |
getEventID() gets the current eventID;
References eventID.
Referenced by BLCMDparticlesource::nextBeamEvent(), BLCMDbeam::nextBeamEvent(), and BLCMDusertrackfilter::UserSteppingAction().
00209 { return eventID; }
void BLManager::setEventID | ( | int | evId | ) | [inline] |
setEventID() sets the current eventID;
References eventID, and prevEventID.
Referenced by BLRunManager::beamOnCollective(), displayVisual(), GeneratePrimaries(), BLCMDtracker::handlePreviousTracks(), trackBeam(), and trackTuneAndReferenceParticles().
00212 { eventID=evId; prevEventID=evId-1; }
bool BLManager::skipEvent | ( | int | EventID | ) | [inline] |
skipEvent() determines if this EventID should be skipped. It checks keepEventList and skipEventList.
References keepEventList, and skipEventList.
Referenced by GeneratePrimaries().
00216 { 00217 if(skipEventList.count(EventID) > 0) return true; 00218 return keepEventList.size() > 0 && keepEventList.count(EventID) == 0; 00219 }
void BLManager::setKeepEvent | ( | int | EventID | ) | [inline] |
setKeepEvent() puts EventID into the keepEventList. (if keepEventList is empty, all events are kept.)
References keepEventList.
Referenced by BLCMDeventcuts::readFile().
00223 { keepEventList.insert(EventID); }
void BLManager::setSkipEvent | ( | int | EventID | ) | [inline] |
setSkipEvent() puts EventID into the skipEventList.
References skipEventList.
Referenced by BLCMDeventcuts::readFile().
00226 { skipEventList.insert(EventID); }
void BLManager::incrEventsProcessed | ( | int | eventID | ) |
incrEventsProcessed() will increment eventsProcessed. For special uses only (e.g.MPI).
References eventsProcessed, BLRunManager::getCollectiveMode(), runManager, startRun, and BLTime::time().
Referenced by EndOfEventAction().
01105 { 01106 // print event number, if appropriate 01107 ++eventsProcessed; 01108 if(!runManager->getCollectiveMode()) { 01109 if(eventsProcessed <= 10 || 01110 (eventsProcessed < 100 && eventsProcessed%10 == 0) || 01111 (eventsProcessed < 1000 && eventsProcessed%100 == 0) || 01112 eventsProcessed%1000 == 0) { 01113 printf("Event %d Completed",eventID); 01114 int t = BLTime::time() - startRun; 01115 if(t <= 0) t = 1; 01116 printf(" %d events realTime=%d sec %.1f ev/sec", 01117 eventsProcessed,t,(double)eventsProcessed/t); 01118 printf("\n"); 01119 fflush(stdout); 01120 } 01121 } 01122 }
bool BLManager::showAllExceptions | ( | bool | value = true |
) | [inline] |
showAllExceptions() sets the flag that prevents thinning out multiple exception printouts. Returns the previous value.
References allExceptions.
Referenced by BLCMDreweightprocess::callback().
00235 { bool tmp=allExceptions; allExceptions = value; return tmp; }
void BLManager::registerSteppingAction | ( | SteppingAction * | sa | ) | [inline] |
registerSteppingAction() registers a SteppingAction to be called for each step (regardless of state).
References allStepVector.
Referenced by BLCMDtrackcuts::command(), BLCMDtrace::command(), BLCMDspacechargelw::command(), BLCMDspacecharge::command(), BLCMDusertrackfilter::construct(), BLCMDparticlefilter::construct(), BLCMDcorner::construct(), and MaterialFilter::MaterialFilter().
00240 { allStepVector.push_back(sa); }
void BLManager::registerSteppingAction | ( | G4VPhysicalVolume * | physicalVol, | |
SteppingAction * | sa | |||
) | [inline] |
registerSteppingAction() registers a SteppingAction to be called for each step (regardless of state), for every step involving the physicalVol. The callback will be called if physicalVol is either the pre- or post-step physical volume (once if both). If physicalVol==0 it is called every reference particle step. LIMITATION: only one callback can be registered for a given physicalVol (except 0).
References allStepMap, and allStepVector.
00252 { if(physicalVol != 0) allStepMap[physicalVol] = sa; 00253 else allStepVector.push_back(sa); }
void BLManager::registerTuneParticleStep | ( | G4VPhysicalVolume * | physicalVol, | |
SteppingAction * | sa | |||
) | [inline] |
registerTuneParticleStep() registers a SteppingAction to be called for every tune particle step involving the physicalVol. The callback will be called if physicalVol is either the pre- or post-step physical volume (once if both). If physicalVol==0 it is called every tune particle step. LIMITATION: only one callback can be registered for a given physicalVol (except 0).
References tpStepMap, and tpStepVector.
Referenced by BLCMDrfdevice::construct(), and BLCMDpillbox::construct().
00264 { if(physicalVol != 0) tpStepMap[physicalVol] = sa; 00265 else tpStepVector.push_back(sa); }
void BLManager::registerReferenceParticleStep | ( | G4VPhysicalVolume * | physicalVol, | |
SteppingAction * | sa | |||
) | [inline] |
registerReferenceParticleStep() registers a SteppingAction to be called for every reference particle step involving the physicalVol. The callback will be called if physicalVol is either the pre- or post-step physical volume (once if both). If physicalVol==0 it is called every reference particle step. LIMITATION: only one callback can be registered for a given physicalVol (except 0).
References rpStepMap, and rpStepVector.
Referenced by BLCMDtimentuple::command(), BLCMDmovie::command(), BLCMDvirtualdetector::construct(), BLCMDtrackerplane::construct(), and BLCoordinates::ReferenceCoordinates::ReferenceCoordinates().
00277 { if(physicalVol != 0) rpStepMap[physicalVol] = sa; 00278 else rpStepVector.push_back(sa); }
void BLManager::registerBeamStep | ( | G4VPhysicalVolume * | physicalVol, | |
SteppingAction * | sa | |||
) | [inline] |
registerBeamStep() registers a SteppingAction to be called for every beam step involving the physicalVol. (beam particles are everything except reference and tune.) The callback will be called if physicalVol is either the pre- or post-step physical volume (once if both). if physicalVol==0 it is called every beam step. LIMITATION: only one callback can be registered for a given physicalVol (except 0).
References beamStepMap, and beamStepVector.
Referenced by BLCMDtotalenergy::command(), BLCMDtimentuple::command(), BLCMDmovie::command(), BLCMDvirtualdetector::construct(), and BLCMDtrackerplane::construct().
00290 { if(physicalVol != 0) beamStepMap[physicalVol] = sa; 00291 else beamStepVector.push_back(sa); }
void BLManager::registerZStep | ( | G4double | z, | |
ZSteppingAction * | sa, | |||
G4int | when = 7 | |||
) |
registerZStep() will force a step to occur near the given z position, and will call the ZSteppingAction for it, interpolating to the desired z value (Centerline coords). when is a bitwise OR of 1=tune, 2=reference, 4=beam.
References beamZStep, insertZStep(), referenceZStep, and tuneZStep.
Referenced by BLCMDzntuple::callback(), BLCMDtune::command(), BLCMDtracker::command(), BLCMDreference::command(), BLCMDprofile::command(), and BLCMDprintf::command().
00403 { 00404 if(when & 1) insertZStep(tuneZStep,z,sa); 00405 if(when & 2) insertZStep(referenceZStep,z,sa); 00406 if(when & 4) insertZStep(beamZStep,z,sa); 00407 }
void BLManager::registerStackingAction | ( | StackingAction * | sa | ) | [inline] |
registerStackingAction() registers a StackingAction to be called by the Geant4 stacking manager.
References stackingActionVector.
Referenced by BLCMDtrackcuts::command().
00302 { stackingActionVector.push_back(sa); }
void BLManager::setSteppingFormat | ( | ) |
setSteppingFormat() sets the verbose printing format according to parameter "steppingFormat".
References appendVerboseFormat(), BLParam::getString(), Param, and BLCommand::parseArgs().
Referenced by initialize().
01233 { 01234 G4String fmt = Param.getString("steppingFormat"); 01235 BLArgumentVector args; 01236 BLArgumentMap namedArgs; 01237 01238 // convert "," to " " 01239 G4String s; 01240 for(unsigned i=0; i<fmt.size(); ++i) 01241 s += (fmt[i]==',' ? ' ' : fmt[i]); 01242 01243 BLCommand::parseArgs(s,args,namedArgs); 01244 appendVerboseFormat(""); 01245 for(unsigned i=0; i<args.size(); ++i) 01246 appendVerboseFormat(args[i]); 01247 appendVerboseFormat("NEWLINE"); 01248 }
G4String BLManager::getFormatHelp | ( | ) |
getFormatHelp() returns a string with help text about valid format items.
01251 { 01252 G4String s = 01253 " EXT toggle extended precision (3 more digits)\n" 01254 " TAG print a '>' (useful to grep output)\n" 01255 " N step number\n" 01256 " NSTEP Synonym of N\n" 01257 " GLOBAL X,Y,Z,T in global coords\n" 01258 " XYZT Synonym of GLOBAL\n" 01259 " CL X,Y,Z,dxdz,dydz in CL coords\n" 01260 " CLX extended precision CL\n" 01261 " KE kinetic energy\n" 01262 " STEP step length\n" 01263 " STEPLEN Synonym of STEP\n" 01264 " VOL volume name\n" 01265 " VOLNAME Synonym of VOL\n" 01266 " PROCESS process name\n" 01267 " B magnetic field\n" 01268 " E electric field\n" 01269 " P 3-momentum\n" 01270 " MAT material name\n" 01271 " ID event ID, track ID, parent ID\n" 01272 " PART particle name\n" 01273 " SEG centerline coord segment number\n" 01274 " WT weight\n" 01275 " NL <newline>\n" 01276 " NEWLINE Synonym of NL\n" 01277 " \\n Synonym of NL\n"; 01278 01279 return s; 01280 }
void BLManager::appendVerboseFormat | ( | G4String | fmt | ) |
appendVerboseFormat() appends fmt to the format for printing when param steppingVerbose is nonzero
References B, CL, CLX, E, EXT, GLOBAL, ID, KE, MAT, NEWLINE, NSTEP, P, PART, PROCESS, SEG, STEP, TAG, verboseFormat, VOL, and WT.
Referenced by setSteppingFormat().
01283 { 01284 G4String f; 01285 for(const char *p=fmt.c_str(); *p; ++p) 01286 f += toupper(*p); 01287 if(f == "") verboseFormat.clear(); 01288 else if(f == "TAG") verboseFormat.push_back(TAG); 01289 else if(f == "EXT") verboseFormat.push_back(EXT); 01290 else if(f == "N") verboseFormat.push_back(NSTEP); 01291 else if(f == "NSTEP") verboseFormat.push_back(NSTEP); 01292 else if(f == "GLOBAL") verboseFormat.push_back(GLOBAL); 01293 else if(f == "XYZT") verboseFormat.push_back(GLOBAL); 01294 else if(f == "CL") verboseFormat.push_back(CL); 01295 else if(f == "CLX") verboseFormat.push_back(CLX); 01296 else if(f == "KE") verboseFormat.push_back(KE); 01297 else if(f == "STEP") verboseFormat.push_back(STEP); 01298 else if(f == "STEPLEN") verboseFormat.push_back(STEP); 01299 else if(f == "VOL") verboseFormat.push_back(VOL); 01300 else if(f == "VOLNAME") verboseFormat.push_back(VOL); 01301 else if(f == "PROCESS") verboseFormat.push_back(PROCESS); 01302 else if(f == "B") verboseFormat.push_back(B); 01303 else if(f == "E") verboseFormat.push_back(E); 01304 else if(f == "P") verboseFormat.push_back(P); 01305 else if(f == "MAT") verboseFormat.push_back(MAT); 01306 else if(f == "ID") verboseFormat.push_back(ID); 01307 else if(f == "PART") verboseFormat.push_back(PART); 01308 else if(f == "SEG") verboseFormat.push_back(SEG); 01309 else if(f == "WT") verboseFormat.push_back(WT); 01310 else if(f == "\n") verboseFormat.push_back(NEWLINE); 01311 else if(f == "NL") verboseFormat.push_back(NEWLINE); 01312 else if(f == "NEWLINE") verboseFormat.push_back(NEWLINE); 01313 else printf("BLManager: invalid verbose format '%s'\n", 01314 fmt.c_str()); 01315 }
void BLManager::steppingVerbosePrint | ( | const G4Step * | step, | |
const G4Track * | track, | |||
int | header = 0 | |||
) |
steppingVerbosePrint() will print this step according to the current verboseFormat. Prints header if header != 0).
References B, BLCOORD_CENTERLINE, CL, CLX, E, EXT, BLCoordinates::getCoords(), getExternalParentID(), getExternalTrackID(), BLGlobalField::GetFieldValue(), getObject(), BLGlobalField::getObject(), BLCoordinates::getRotation(), BLCoordinates::getSegmentCL(), GLOBAL, ID, BLCoordinates::isValid(), KE, MAT, NEWLINE, NSTEP, P, PART, PROCESS, runManager, SEG, STEP, TAG, verboseFormat, VOL, and WT.
Referenced by UserSteppingAction().
01319 { 01320 G4StepPoint *prePoint=0; 01321 G4VPhysicalVolume *preVol=0; 01322 G4StepPoint *postPoint=0; 01323 G4double time=0; 01324 G4ThreeVector position; 01325 G4Material *mat=0; 01326 int eventNum = -999; 01327 int nstep=0; 01328 const char *procName=0; 01329 BLCoordinates *coord = (track!=0 ? (BLCoordinates *)track->GetUserInformation() : 0); 01330 if(coord && !coord->isValid()) coord = 0; 01331 01332 if(!header) { 01333 prePoint = step->GetPreStepPoint(); 01334 preVol = prePoint->GetPhysicalVolume(); 01335 postPoint = step->GetPostStepPoint(); 01336 time = track->GetGlobalTime(); 01337 position = track->GetPosition(); 01338 mat = preVol->GetLogicalVolume()->GetMaterial(); 01339 nstep = track->GetCurrentStepNumber(); 01340 procName = (postPoint->GetProcessDefinedStep() ? 01341 postPoint->GetProcessDefinedStep()->GetProcessName().c_str() : 01342 "UserLimit"); 01343 const G4Event* event = runManager->GetCurrentEvent(); 01344 eventNum = event->GetEventID(); 01345 } 01346 01347 bool extended=false; 01348 for(unsigned i=0; i<verboseFormat.size(); ++i) { 01349 switch(verboseFormat[i]) { 01350 case EXT: 01351 extended = !extended; 01352 break; 01353 case TAG: 01354 if(header) printf(" "); 01355 else printf(">"); 01356 break; 01357 case NSTEP: 01358 if(header) printf(" Step"); 01359 else printf(" %4d",nstep); 01360 break; 01361 case GLOBAL: 01362 if(!extended) { 01363 if(header) printf(" X(mm) Y(mm) Z(mm) T(ns)"); 01364 else printf(" %7.1f %7.1f %7.1f %8.2f", 01365 position[0],position[1],position[2], 01366 time); 01367 } else { 01368 if(header) printf(" X(mm) Y(mm) Z(mm) T(ns) "); 01369 else printf(" %10.4f %10.4f %10.4f %11.5f", 01370 position[0],position[1],position[2], 01371 time); 01372 } 01373 break; 01374 case CL: 01375 if(extended) goto clx; 01376 if(header) { 01377 printf(" CL: X Y Z dxdz dydz"); 01378 } else { 01379 G4ThreeVector dir = track->GetMomentumDirection(); 01380 // transform to centerline coordinates 01381 G4ThreeVector clpos; 01382 coord->getCoords(BLCOORD_CENTERLINE,clpos); 01383 dir = coord->getRotation() * dir; 01384 printf(" %8.1f %8.1f %8.1f %7.4f %7.4f", 01385 clpos[0],clpos[1],clpos[2], 01386 dir[0]/dir[2],dir[1]/dir[2]); 01387 } 01388 break; 01389 case CLX: 01390 clx: if(header) { 01391 printf(" CL: X Y Z dxdz dydz "); 01392 } else { 01393 G4ThreeVector dir = track->GetMomentumDirection(); 01394 // transform to centerline coordinates 01395 G4ThreeVector clpos; 01396 coord->getCoords(BLCOORD_CENTERLINE,clpos); 01397 dir = coord->getRotation() * dir; 01398 printf(" %11.4f %11.4f %11.4f %10.7f %10.7f", 01399 clpos[0],clpos[1],clpos[2], 01400 dir[0]/dir[2],dir[1]/dir[2]); 01401 } 01402 break; 01403 case KE: 01404 if(!extended) { 01405 if(header) printf(" KE(MeV)"); 01406 else printf(" %8.1f",track->GetKineticEnergy()); 01407 } else { 01408 if(header) printf(" KE(MeV) "); 01409 else printf(" %11.4f",track->GetKineticEnergy()); 01410 } 01411 break; 01412 case STEP: 01413 if(!extended) { 01414 if(header) printf(" StepLen"); 01415 else printf(" %8.2f",step->GetStepLength()); 01416 } else { 01417 if(header) printf(" StepLen "); 01418 else printf(" %11.5f",step->GetStepLength()); 01419 } 01420 break; 01421 case VOL: 01422 if(header) printf(" This Volume "); 01423 else printf(" %-16s",preVol->GetName().c_str()); 01424 break; 01425 case PROCESS: 01426 if(header) printf(" Process "); 01427 else printf(" %-16s",procName); 01428 break; 01429 case B: 01430 if(header) { 01431 printf(" Bx,By,Bz (Tesla) "); 01432 } else { 01433 G4double point[4], field[6]; 01434 point[0]=position[0],point[1]=position[1],point[2]=position[2]; 01435 point[3]=time; 01436 BLGlobalField::getObject()->GetFieldValue(point,field); 01437 G4ThreeVector B(field[0],field[1],field[2]); 01438 B = coord->getRotation(BLCOORD_CENTERLINE) * B; 01439 printf(" %8.4f %8.4f %8.4f", 01440 B[0]/tesla,B[1]/tesla,B[2]/tesla); 01441 } 01442 break; 01443 case E: 01444 if(header) { 01445 printf(" Ex,Ey,Ez (MV/meter) "); 01446 } else { 01447 G4double point[4], field[6]; 01448 point[0]=position[0],point[1]=position[1],point[2]=position[2]; 01449 point[3]=time; 01450 BLGlobalField::getObject()->GetFieldValue(point,field); 01451 G4ThreeVector E(field[3],field[4],field[5]); 01452 E = coord->getRotation(BLCOORD_CENTERLINE) * E; 01453 printf(" %8.4f %8.4f %8.4f",E[0]/(megavolt/meter), 01454 E[1]/(megavolt/meter),E[2]/(megavolt/meter)); 01455 } 01456 break; 01457 case P: 01458 if(!extended) { 01459 if(header) { 01460 printf(" Px,Py,Pz (MeV/c) "); 01461 } else { 01462 G4ThreeVector P = track->GetMomentum(); 01463 P = coord->getRotation(BLCOORD_CENTERLINE) * P; 01464 printf(" %8.1f %8.1f %8.1f",P[0]/(MeV), 01465 P[1]/(MeV),P[2]/(MeV)); 01466 } 01467 } else { 01468 if(header) { 01469 printf(" Px,Py,Pz (MeV/c) "); 01470 } else { 01471 G4ThreeVector P = track->GetMomentum(); 01472 P = coord->getRotation(BLCOORD_CENTERLINE) * P; 01473 printf(" %11.4f %11.4f %11.4f",P[0]/(MeV), 01474 P[1]/(MeV),P[2]/(MeV)); 01475 } 01476 } 01477 break; 01478 case MAT: 01479 if(header) { 01480 printf(" Material "); 01481 } else { 01482 printf(" %-9s",mat->GetName().c_str()); 01483 } 01484 break; 01485 case ID: 01486 if(header) { 01487 printf(" Event# Trk Prnt"); 01488 } else { 01489 printf(" %6d %3d %4d",eventNum, 01490 BLManager::getObject()-> 01491 getExternalTrackID(track), 01492 BLManager::getObject()-> 01493 getExternalParentID(track)); 01494 } 01495 break; 01496 case PART: 01497 if(header) { 01498 printf(" Particle"); 01499 } else { 01500 printf(" %8s",track->GetDefinition()->GetParticleName().c_str()); 01501 } 01502 break; 01503 case SEG: 01504 if(header) { 01505 printf(" Seg"); 01506 } else { 01507 printf(" %3d",coord->getSegmentCL()); 01508 } 01509 break; 01510 case WT: 01511 if(header) { 01512 printf(" Wt "); 01513 } else { 01514 printf(" %.3f",track->GetWeight()); 01515 } 01516 break; 01517 case NEWLINE: 01518 printf("\n"); 01519 break; 01520 } 01521 } 01522 01523 fflush(stdout); // help debugging location of crashes 01524 }
void BLManager::setSteppingManager | ( | G4SteppingManager * | p | ) | [inline] |
setSteppingManager() sets the pointer to the current G4SteppingManager.
References fpSteppingManager.
Referenced by BLManager_UserSteppingAction::UserSteppingAction().
00324 { fpSteppingManager = p; }
G4SteppingManager* BLManager::getSteppingManager | ( | ) | [inline] |
getSteppingManager() returns a pointer to the current G4SteppingManager.
References fpSteppingManager.
Referenced by RFdeviceField::UserSteppingAction(), PillboxField::UserSteppingAction(), BLCMDtune::UserZSteppingAction(), and BLCMDreference::UserZSteppingAction().
00328 { return fpSteppingManager; }
void BLManager::setTrackingManager | ( | G4TrackingManager * | p | ) | [inline] |
setTrackingManager() sets the pointer to the current G4TrackingManager.
References fpTrackingManager.
Referenced by BLManager_UserTrackingAction::PreUserTrackingAction().
00332 { fpTrackingManager = p; }
void BLManager::initialize | ( | ) |
initialize() will initialize the BLManager object, and the geant4 kernel, thus constructing the geometry in the world group. Note that registerPhysics() must be called before initialize() (normally done by a "physics" command in the input file).
physics->getPhysicsList() is set in registerPhysicsList()
References beamIndex, beamVector, eventTimeLimit, BLParam::getDouble(), BLParam::getInt(), histoUpdate, initialized, Param, physics, referenceVector, runManager, setSteppingFormat(), steppingVerbose, and zTolerance.
Referenced by BLCMDlist::command(), and main().
00319 { 00320 if(initialized) { 00321 G4Exception("BLManager","Already Initialized",FatalException, 00322 ""); 00323 } 00324 initialized = true; 00325 00326 G4UImanager* UI = G4UImanager::GetUIpointer(); 00327 00328 if(!physics) { 00329 G4Exception("BLManager","No physics registered",FatalException, 00330 ""); 00331 } 00332 00333 setSteppingFormat(); 00334 00335 // set initialization and user classes 00336 runManager->SetUserInitialization((G4VUserDetectorConstruction*)this); 00337 /// physics->getPhysicsList() is set in registerPhysicsList() 00338 runManager->SetUserAction(new BLManager_UserSteppingAction(this)); 00339 runManager->SetUserAction(new BLManager_UserTrackingAction(this)); 00340 runManager->SetUserAction(new BLManager_UserRunAction(this)); 00341 runManager->SetUserAction(new BLManager_UserEventAction(this)); 00342 runManager->SetUserAction(new BLManager_PrimaryGeneratorAction(this)); 00343 runManager->SetUserAction(new BLManager_UserStackingAction(this)); 00344 00345 // apply initial range cuts 00346 UI->ApplyCommand("/range/cutG 2 mm"); 00347 UI->ApplyCommand("/range/cutE 2 mm"); 00348 00349 // initialize the RunManager, construct the apparatus, etc. 00350 runManager->Initialize(); 00351 00352 // set initial verbosities 00353 UI->ApplyCommand("/control/verbose 0"); 00354 UI->ApplyCommand("/run/verbose 0"); 00355 UI->ApplyCommand("/event/verbose 0"); 00356 UI->ApplyCommand("/tracking/verbose 0"); 00357 UI->ApplyCommand("/hits/verbose 0"); 00358 UI->ApplyCommand("/material/verbose 0"); 00359 UI->ApplyCommand("/process/setVerbose 0 all"); 00360 UI->ApplyCommand("/process/verbose 0"); 00361 UI->ApplyCommand("/process/eLoss/verbose 0"); 00362 G4HadronicProcessStore::Instance()->SetVerbose(0); 00363 00364 steppingVerbose = Param.getInt("steppingVerbose"); 00365 00366 // get various parameter values (may have changed since constructor) 00367 histoUpdate = Param.getInt("histoUpdate"); 00368 eventTimeLimit = Param.getInt("eventTimeLimit"); 00369 zTolerance = Param.getDouble("zTolerance"); 00370 00371 // initialize the beam 00372 for(beamIndex=0; beamIndex<beamVector.size(); ++beamIndex) 00373 beamVector[beamIndex]->init(); 00374 for(unsigned i=0; i<referenceVector.size(); ++i) 00375 referenceVector[i]->init(); 00376 beamIndex = 0; 00377 00378 // create the ZStepLimiter 00379 new ZStepLimiter(); 00380 }
static bool BLManager::isInitialized | ( | ) | [inline, static] |
isInitialized() returns true if the BLManager has been initialized.
References initialized.
Referenced by BLCMDlist::command().
00342 { return initialized; }
void BLManager::trackTuneAndReferenceParticles | ( | ) |
trackTuneAndReferenceParticles() will generate and track the tune particle and then the reference particle.
References beamIndex, BLRunManager::BeamOn(), FORCE_OFF, BLRunManager::getCollectiveMode(), IDLE, NORMAL, physics, REFERENCE, referenceVector, runManager, BLRunManager::setCollectiveMode(), BLPhysics::setDoStochastics(), setEventID(), state, and TUNE.
Referenced by main().
00410 { 00411 if(referenceVector.size() == 0) 00412 return; 00413 00414 if(!physics) { 00415 G4Exception("BLManager","No physics registered",FatalException, 00416 ""); 00417 } 00418 00419 // Tune and Reference particles cannot use collective mode 00420 bool collectiveMode = runManager->getCollectiveMode(); 00421 runManager->setCollectiveMode(false); 00422 00423 printf("================= Prepare Tune Particle(s) ===========\n"); 00424 physics->setDoStochastics(FORCE_OFF); 00425 runManager->Initialize(); 00426 00427 printf("================= Begin Tune Particle(s) =============\n"); 00428 state = TUNE; 00429 setEventID(-2); 00430 beamIndex = 0; 00431 runManager->BeamOn(referenceVector.size()); 00432 state = IDLE; 00433 00434 // now track center particle 00435 printf("================== Begin Reference Particle(s) ===============\n"); 00436 state = REFERENCE; 00437 setEventID(-1); 00438 beamIndex = 0; 00439 runManager->BeamOn(referenceVector.size()); 00440 state = IDLE; 00441 beamIndex = 0; 00442 00443 physics->setDoStochastics(NORMAL); 00444 runManager->setCollectiveMode(collectiveMode); 00445 }
void BLManager::trackBeam | ( | ) |
trackBeam() will generate the defined beam and track each event.
References BEAM, beamIndex, BLRunManager::BeamOn(), beamVector, IDLE, NORMAL, physics, runManager, BLPhysics::setDoStochastics(), setEventID(), and state.
Referenced by main().
00448 { 00449 if(beamVector.size() == 0) { 00450 G4Exception("BLManager","No beam registered",FatalException, 00451 ""); 00452 } 00453 if(!physics) { 00454 G4Exception("BLManager","No physics registered",FatalException, 00455 ""); 00456 } 00457 00458 printf("================== Prepare Tracking Beam ==================\n"); 00459 physics->setDoStochastics(NORMAL,2); 00460 runManager->Initialize(); 00461 00462 printf("================== Begin Tracking Beam ===============\n"); 00463 state = BEAM; 00464 setEventID(1); 00465 beamIndex = 0; 00466 runManager->BeamOn(0x7FFFFFFF); 00467 for(beamIndex=0; beamIndex<beamVector.size(); ++beamIndex) 00468 beamVector[beamIndex]->summary(); 00469 beamIndex = 0; 00470 state = IDLE; 00471 }
void BLManager::displayVisual | ( | ) |
displayVisual() will display the detector visually.
References beamIndex, BLParam::getString(), handleCallbacks(), IDLE, BLVisManager::init(), NORMAL, Param, physics, runManager, BLPhysics::setDoStochastics(), setEventID(), state, and VISUAL.
Referenced by main().
00474 { 00475 #ifdef G4VIS_USE 00476 printf("================== Prepare Visualization ==================\n"); 00477 state = VISUAL; 00478 physics->setDoStochastics(NORMAL,2); 00479 runManager->Initialize(); 00480 setEventID(1); 00481 beamIndex = 0; 00482 00483 BLVisManager *visManager = new BLVisManager(Param.getString("viewer")); 00484 visManager->init(); 00485 00486 handleCallbacks(4); 00487 00488 printf("To display a run with 10 events in the viewer, type '/run/beamOn 10<cr>'\n"); 00489 G4UIterminal *terminal = new G4UIterminal(0,false); 00490 terminal->SessionStart(); 00491 00492 delete terminal; 00493 delete visManager; 00494 state = IDLE; 00495 beamIndex = 0; 00496 #else 00497 G4Exception("BLManager","G4VIS_USE not defined",FatalException,""); 00498 #endif 00499 }
void BLManager::displayGeometry | ( | G4VPhysicalVolume * | phys = 0 , |
|
int | level = 0 | |||
) |
displayGeometry() will display the geant4 geometry. This is a hierarchical ASCII listing of all volumes. if phys==0 then use the worldPhysicalVolume.
References BLCommand::dumpRotation(), and worldPhysicalVolume.
Referenced by BLCMDgeometry::callback().
00502 { 00503 if(phys == 0) phys = worldPhysicalVolume; 00504 for(int i=0; i<level; ++i) printf(" "); 00505 G4ThreeVector pos = phys->GetObjectTranslation(); 00506 printf("%s x=%.3f y=%.3f z=%.3f\n",phys->GetName().c_str(), 00507 pos[0],pos[1],pos[2]); 00508 G4RotationMatrix rot = phys->GetObjectRotationValue(); 00509 if(!rot.isIdentity()) { 00510 for(int i=0; i<level; ++i) printf(" "); 00511 BLCommand::dumpRotation(&rot,""); 00512 } 00513 G4LogicalVolume *log = phys->GetLogicalVolume(); 00514 int n = log->GetNoDaughters(); 00515 for(int i=0; i<n; ++i) { 00516 G4VPhysicalVolume *p = log->GetDaughter(i); 00517 if(p) 00518 displayGeometry(p,level+1); 00519 } 00520 }
void BLManager::registerPhysics | ( | BLPhysics * | _physics | ) | [inline] |
registerPhysics() registers the BLPhysics object. It also sets the physics list to the BLRunManager, so following commands can find particle by name.
References BLPhysics::getPhysicsList(), physics, and runManager.
Referenced by BLCMDphysics::command().
00363 { physics = _physics; 00364 runManager->SetUserInitialization(physics->getPhysicsList()); 00365 }
BLPhysics* BLManager::getPhysics | ( | ) | [inline] |
getPhysics returns the registered BLPhysics object.
References physics.
Referenced by BLCMDspacechargelw::beginCollectiveTracking(), BLCMDtrackermode::callback(), BLCMDtracker::fitTrack(), BLCMDreference::PostUserTrackingAction(), and BLCMDreference::PreUserTrackingAction().
00368 { return physics; }
void BLManager::registerBeam | ( | BLBeam * | _beam | ) | [inline] |
registerBeam() registers a BLBeam object for beam generation. Multiple BLBeam objects can be registered, used in order.
References beamVector.
Referenced by BLCMDparticlesource::command(), BLCMDcosmicraybeam::command(), and BLCMDbeam::command().
00372 { beamVector.push_back(_beam); }
void BLManager::registerReference | ( | BLBeam * | _beam | ) | [inline] |
registerReference() registers a BLBeam object for reference particle generation. Multiple BLBeam objects can be registered, used in order.
References referenceVector.
Referenced by BLCMDreference::command().
00378 { referenceVector.push_back(_beam); }
int BLManager::nReference | ( | ) | [inline] |
nReference() returns the number of reference particles registered.
References referenceVector.
Referenced by main().
00381 { return referenceVector.size(); }
void BLManager::registerRunAction | ( | RunAction * | a, | |
G4bool | beamOnly = true | |||
) | [inline] |
registerRunAction() registers a UserRunAction. If beamOnly is true (the default), the callback is made only if state==BEAM.
References beamRunActionVector, and runActionVector.
Referenced by BLCMDtune::command(), BLCMDreference::command(), BLCMDhelicalharmonic::construct(), and BLCMDhelicaldipole::construct().
00387 { if(beamOnly) 00388 beamRunActionVector.push_back(a); 00389 else 00390 runActionVector.push_back(a); 00391 }
void BLManager::registerEventAction | ( | EventAction * | a, | |
G4bool | beamOnly = true | |||
) | [inline] |
registerEventAction() registers a UserEventAction. If beamOnly is true (the default), the callback is made only if state==BEAM.
References beamEventActionVector, and eventActionVector.
Referenced by BLCMDparticlecolor::command(), and BLCMDntuple::command().
00397 { if(beamOnly) 00398 beamEventActionVector.push_back(a); 00399 else 00400 eventActionVector.push_back(a); 00401 }
void BLManager::registerTrackingAction | ( | TrackingAction * | a, | |
bool | front = false | |||
) | [inline] |
registerTrackingAction() registers a UserTackingAction. By default puts new entry at the back end of the vector; if front is true, puts new entry at the front.
References trackingActionVector.
Referenced by BLCMDtracker::command(), BLCMDtrace::command(), BLCMDtotalenergy::command(), BLCMDspacechargelw::command(), BLCMDspacecharge::command(), BLCMDreference::command(), BLCMDparticlecolor::command(), BLCMDntuple::command(), BLCMDnewparticlentuple::command(), BLCMDmovie::command(), BLCMDbeamlossntuple::command(), BLCMDparticlefilter::construct(), and BLCoordinates::ReferenceCoordinates::ReferenceCoordinates().
00407 { if(front) 00408 trackingActionVector.insert(trackingActionVector.begin(),a); 00409 else 00410 trackingActionVector.push_back(a); 00411 }
void BLManager::registerUserCode | ( | BLUserCode * | instance | ) | [inline] |
registerUserCode() registers a BLUserCode instance.
References userCodeVector.
Referenced by BLUserCode::registerUserCode().
00415 { userCodeVector.push_back(instance); }
std::vector<BLUserCode*> BLManager::getUserCodeInstances | ( | G4String | type | ) | [inline] |
getUserCodeInstances() returns a vector of all registered instances of BLUserCode with the specified type.
References userCodeVector.
00420 { std::vector<BLUserCode*> ret; 00421 for(unsigned i=0; i<userCodeVector.size(); ++i) { 00422 if(type == userCodeVector[i]->getType()) 00423 ret.push_back(userCodeVector[i]); 00424 } 00425 return ret; 00426 }
void BLManager::registerCallback | ( | BLCallback * | cb, | |
int | type | |||
) | [inline] |
registerCallback() registers a BLCallback. type=0 for pre-Tune particle, type=1 for post-Reference (pre-beam tracking), type=2 for post-beam tracking. type=3 for replacing the main program loop. type=4 for visualization. NOTE: if there are type=3 callbacks, when the last one returns the closes up by summarizing NTuples and callng handleCallbacks(2), and then the program exits. This prevents the main program loop from executing. type=3 callbacks are called just after the type=1 callbacks (i.e. after the post-Reference callbacks).
References postReferenceCallbackVector, postTrackingCallbackVector, preReferenceCallbackVector, replaceMainLoopCallbackVector, and visualizationCallbackVector.
Referenced by BLCMDzntuple::command(), BLCMDusertrackfilter::command(), BLCMDtrackermode::command(), BLCMDtracker::command(), BLCMDtotalenergy::command(), BLCMDshowmaterial::command(), BLCMDsetdecay::command(), BLCMDreweightprocess::command(), BLCMDprofile::command(), BLCMDprobefield::command(), BLCMDprintfield::command(), BLCMDprintf::command(), BLCMDphysics::command(), BLCMDmuminuscapturefix::command(), BLCMDmovie::command(), BLCMDgeometry::command(), BLCMDg4ui::command(), BLCMDfieldntuple::command(), BLCMDfieldlines::command(), BLCMDbug1021::command(), SetDecayInstance::fatalError(), MaterialFilter::MaterialFilter(), ParticleFilterPlacement::ParticleFilterPlacement(), and Surface::Surface().
00440 { 00441 if(type==0) preReferenceCallbackVector.push_back(cb); 00442 else if(type==1) postReferenceCallbackVector.push_back(cb); 00443 else if(type==2) postTrackingCallbackVector.push_back(cb); 00444 else if(type==3) replaceMainLoopCallbackVector.push_back(cb); 00445 else if(type==4) visualizationCallbackVector.push_back(cb); 00446 }
void BLManager::handleCallbacks | ( | int | type | ) |
handleCallbacks() calls all applicable registered callbacks. type=0 for pre-reference particle, type=1 for post-center (pre-beam tracking), type=2 for post-beam tracking (just before program exit). type=3 for replacing the main program loop. type=4 for visualization.
References blManager, BLRunManager::callback(), BLNTuple::closeAll(), exceptionSummary(), g4bl_exit(), BLRunManager::getObject(), postReferenceCallbackVector, postTrackingCallbackVector, preReferenceCallbackVector, replaceMainLoopCallbackVector, BLNTuple::summary(), and visualizationCallbackVector.
Referenced by BLCMDlist::command(), displayVisual(), main(), and Notify().
00537 { 00538 // startup sequencing prevents RunManager from registering 00539 // so just do it 00540 BLRunManager::getObject()->callback(type); 00541 00542 if(type == 0) { 00543 for(unsigned i=0; i<preReferenceCallbackVector.size(); ++i) 00544 preReferenceCallbackVector[i]->callback(type); 00545 } else if(type == 1) { 00546 for(unsigned i=0; i<postReferenceCallbackVector.size(); ++i) 00547 postReferenceCallbackVector[i]->callback(type); 00548 } else if(type == 2) { 00549 for(unsigned i=0; i<postTrackingCallbackVector.size(); ++i) 00550 postTrackingCallbackVector[i]->callback(type); 00551 exceptionSummary(); 00552 } else if(type == 3) { 00553 for(unsigned i=0; i<replaceMainLoopCallbackVector.size(); ++i) 00554 replaceMainLoopCallbackVector[i]->callback(type); 00555 if(replaceMainLoopCallbackVector.size() > 0) { 00556 // we have replaced the main loop, so closeup and exit 00557 BLNTuple::summary(); 00558 BLNTuple::closeAll(); 00559 // handle post-tracking callbacks 00560 handleCallbacks(2); 00561 blManager = 0; 00562 g4bl_exit(0); 00563 } 00564 } else if(type == 4) { 00565 for(unsigned i=0; i<visualizationCallbackVector.size(); ++i) 00566 visualizationCallbackVector[i]->callback(type); 00567 } 00568 }
G4VPhysicalVolume* BLManager::getWorldPhysicalVolume | ( | ) | [inline] |
getWorldPhysicalVolume() returns a pointer to the world PV. Note that it must already have been consructed, so this function returns NULL before construct() is called (by main).
References worldPhysicalVolume.
Referenced by BLCMDtotalenergy::callback(), BLCMDshowmaterial::callback(), and BLCMDmovie::callback().
00460 { return worldPhysicalVolume; }
void BLManager::setPRNGSeedMethod | ( | PRNGSeedMethod | method | ) | [inline] |
setPRNGSeedMethod() will set the method used to seed the pseudo random number generator at the start of each event.
References seedMethod.
Referenced by BLCMDrandomseed::command(), and BLBeam::setRandomSeedToGenerate().
00465 { seedMethod = method; }
PRNGSeedMethod BLManager::getPRNGSeedMethod | ( | ) | [inline] |
getPRNGSeedMethod() will return the method used to seed the pseudo random number generator at the start of each event.
References seedMethod.
00469 { return seedMethod; }
void BLManager::UserSteppingAction | ( | const G4Step * | step | ) |
UserSteppingAction() from G4UserSteppingAction.
References allStepMap, allStepVector, BEAM, beamStepMap, beamStepVector, currentZStep, eventTimeLimit, BLCoordinates::getCLZ(), indexZStep, BLCoordinates::isValid(), nStuckSteps, prevZ, BLSignal::received(), REFERENCE, rpStepMap, rpStepVector, ZStepLimiter::setMaxStep(), snprintf, startEvent, state, steppingVerbose, steppingVerbosePrint(), stuckTracks, BLTime::time(), tpStepMap, tpStepVector, TUNE, BLCoordinates::update(), BLSignal::value(), and zTolerance.
Referenced by BLManager_UserSteppingAction::UserSteppingAction().
00773 { 00774 if(BLSignal::received()) { 00775 char tmp[64]; 00776 snprintf(tmp,sizeof(tmp),"Signal %d received", 00777 BLSignal::value()); 00778 G4Exception("BLManager",tmp,FatalException,""); 00779 } 00780 00781 G4Track *track = step->GetTrack(); 00782 G4StepPoint *prePoint = step->GetPreStepPoint(); 00783 G4StepPoint *postPoint = step->GetPostStepPoint(); 00784 G4VPhysicalVolume *preVol=0; 00785 G4VPhysicalVolume *postVol=0; 00786 if(prePoint) preVol = prePoint->GetPhysicalVolume(); 00787 if(postPoint) postVol = postPoint->GetPhysicalVolume(); 00788 00789 // need to update BLCoordinates before steppingVerbose print. 00790 // So BLCoordinates is not registered with SteppingAction, it 00791 // is handled specially, right here. 00792 BLCoordinates::update(track); 00793 BLCoordinates *coord = (BLCoordinates *)track->GetUserInformation(); 00794 if(!coord || !coord->isValid()) { 00795 G4Exception("BLManager","Coordinates Got Lost",FatalException, 00796 ""); 00797 } 00798 00799 // steppingVerbose print moved to the top of this routine 00800 if(steppingVerbose > 0) { 00801 static bool first=true; 00802 int nstep = track->GetCurrentStepNumber(); 00803 if(first || (nstep <= 1 && !runManager->getCollectiveMode())) 00804 steppingVerbosePrint(0,0,1); 00805 steppingVerbosePrint(step,track,0); 00806 first = false; 00807 } 00808 00809 G4TrackStatus status = track->GetTrackStatus(); 00810 if(status == fStopAndKill || status == fKillTrackAndSecondaries) 00811 goto quit; 00812 00813 /* check for stuck track */ 00814 if(step->GetStepLength() > 0.0001*mm || 00815 step->GetDeltaTime() > 0.0001*ns) { 00816 nStuckSteps = 0; 00817 } else if(status == fAlive && ++nStuckSteps >= 100) { 00818 if(track->GetVelocity()/c_light < 0.010) { 00819 ++stuckTracks; 00820 track->SetTrackStatus(fStopButAlive); 00821 if(steppingVerbose > 0) 00822 printf("BLManager: Stuck Track -- stopped\n"); 00823 } else { 00824 G4Exception("BLManager Stepping Action", 00825 "Stuck Track -- Killed",JustWarning, 00826 "100 steps in a row, each less than 0.1 micron and 0.1 ps"); 00827 track->SetTrackStatus(fStopAndKill); 00828 } 00829 } 00830 00831 // check event time limit 00832 if(eventTimeLimit > 0 && BLTime::time()-startEvent > eventTimeLimit) { 00833 G4Exception("BLManager","Event Time Limit",EventMustBeAborted,""); 00834 goto quit; 00835 } 00836 00837 // call ZStep actions 00838 // NOTE: indexZstep points to the next entry in the +Z direction. 00839 // Note also the bookends constrain indexZStep to always be valid. 00840 if(currentZStep->size() > 2) { // i.e. not just 2 bookends 00841 G4double thisZ = coord->getCLZ(); 00842 if(thisZ == prevZ) goto noZstep; 00843 G4double minZ = (thisZ<prevZ ? thisZ : prevZ); 00844 G4double maxZ = (thisZ>prevZ ? thisZ : prevZ); 00845 // find first entry that might be spanned by the step 00846 while(indexZStep > 0 && minZ <= (*currentZStep)[indexZStep-1].z) 00847 --indexZStep; 00848 while(minZ > (*currentZStep)[indexZStep].z) 00849 ++indexZStep; 00850 int indexPrev = indexZStep-1; 00851 // loop over all entries actually spanned by this step 00852 for( ; indexZStep<currentZStep->size()-1; ++indexZStep) { 00853 G4double z=(*currentZStep)[indexZStep].z; 00854 if(maxZ <= z) break; 00855 if((*currentZStep)[indexZStep].action == 0) continue; 00856 // save current values from track 00857 G4ThreeVector pos = track->GetPosition(); 00858 G4ThreeVector mom = track->GetMomentum(); 00859 G4double time = track->GetGlobalTime(); 00860 G4double ke = track->GetKineticEnergy(); 00861 // interpolate linearly to z, global coords 00862 G4ThreeVector deltaPos=step->GetDeltaPosition(); 00863 G4ThreeVector deltaMom=(postPoint->GetMomentum() - 00864 prePoint->GetMomentum()); 00865 G4double deltaTime=step->GetDeltaTime(); 00866 G4double deltaE=(postPoint->GetKineticEnergy() - 00867 prePoint->GetKineticEnergy()); 00868 G4double f=(z-thisZ)/(prevZ-thisZ); 00869 track->SetPosition(pos-f*deltaPos); 00870 track->SetMomentumDirection((mom-f*deltaMom).unit()); 00871 track->SetGlobalTime(time-f*deltaTime); 00872 track->SetKineticEnergy(ke-f*deltaE); 00873 BLCoordinates::update(track); 00874 (*currentZStep)[indexZStep].action-> 00875 UserZSteppingAction(track); 00876 // restore current values to track 00877 track->SetPosition(pos); 00878 track->SetMomentumDirection(mom.unit()); 00879 track->SetGlobalTime(time); 00880 track->SetKineticEnergy(ke); 00881 BLCoordinates::update(track); 00882 } 00883 double dz = (thisZ>prevZ ? (*currentZStep)[indexZStep].z-thisZ : 00884 thisZ-(*currentZStep)[indexPrev].z); 00885 dz = (dz>zTolerance*2.0 ? dz : zTolerance*2.0); 00886 ZStepLimiter::setMaxStep(dz); 00887 prevZ = thisZ; 00888 } 00889 noZstep: 00890 00891 // call all-state stepping actions before per-state actions. 00892 { std::vector<BLManager::SteppingAction*>::iterator i; 00893 for(i=allStepVector.begin(); i!=allStepVector.end(); ++i) { 00894 (*i)->UserSteppingAction(step); 00895 if(track->GetTrackStatus() != fAlive) 00896 goto quit; 00897 } 00898 if(preVol && allStepMap.count(preVol) > 0) { 00899 allStepMap[preVol]->UserSteppingAction(step); 00900 if(track->GetTrackStatus() != fAlive) 00901 goto quit; 00902 } 00903 if(postVol && preVol != postVol && allStepMap.count(postVol) > 0) { 00904 allStepMap[postVol]->UserSteppingAction(step); 00905 if(track->GetTrackStatus() != fAlive) 00906 goto quit; 00907 } 00908 } 00909 00910 // call Tune Particle stepping actions 00911 if(state == TUNE) { 00912 std::vector<BLManager::SteppingAction*>::iterator i; 00913 for(i=tpStepVector.begin(); i!=tpStepVector.end(); ++i) { 00914 (*i)->UserSteppingAction(step); 00915 if(track->GetTrackStatus() != fAlive) 00916 goto quit; 00917 } 00918 if(preVol && tpStepMap.count(preVol) > 0) { 00919 tpStepMap[preVol]->UserSteppingAction(step); 00920 if(track->GetTrackStatus() != fAlive) 00921 goto quit; 00922 } 00923 if(postVol && preVol != postVol && 00924 tpStepMap.count(postVol) > 0) { 00925 tpStepMap[postVol]->UserSteppingAction(step); 00926 if(track->GetTrackStatus() != fAlive) 00927 goto quit; 00928 } 00929 } 00930 00931 // call Reference Particle stepping actions 00932 if(state == REFERENCE) { 00933 std::vector<BLManager::SteppingAction*>::iterator i; 00934 for(i=rpStepVector.begin(); i!=rpStepVector.end(); ++i) { 00935 (*i)->UserSteppingAction(step); 00936 if(track->GetTrackStatus() != fAlive) 00937 goto quit; 00938 } 00939 if(preVol && rpStepMap.count(preVol) > 0) { 00940 rpStepMap[preVol]->UserSteppingAction(step); 00941 if(track->GetTrackStatus() != fAlive) 00942 goto quit; 00943 } 00944 if(postVol && preVol != postVol && 00945 rpStepMap.count(postVol) > 0) { 00946 rpStepMap[postVol]->UserSteppingAction(step); 00947 if(track->GetTrackStatus() != fAlive) 00948 goto quit; 00949 } 00950 } 00951 00952 // call beam stepping actions 00953 if(state == BEAM) { 00954 std::vector<BLManager::SteppingAction*>::iterator i; 00955 for(i=beamStepVector.begin(); i!=beamStepVector.end(); ++i) { 00956 (*i)->UserSteppingAction(step); 00957 if(track->GetTrackStatus() != fAlive) 00958 goto quit; 00959 } 00960 if(preVol && beamStepMap.count(preVol) > 0) { 00961 beamStepMap[preVol]->UserSteppingAction(step); 00962 if(track->GetTrackStatus() != fAlive) 00963 goto quit; 00964 } 00965 if(postVol && preVol != postVol && 00966 beamStepMap.count(postVol) > 0) { 00967 beamStepMap[postVol]->UserSteppingAction(step); 00968 if(track->GetTrackStatus() != fAlive) 00969 goto quit; 00970 } 00971 } 00972 00973 quit: ; 00974 }
G4VPhysicalVolume * BLManager::Construct | ( | ) |
Construct() from G4VUserDetectorConstruction.
References BLGroup::constructWorld(), BLGlobalField::getObject(), and worldPhysicalVolume.
00523 { 00524 // ensure the global field is initialized 00525 (void)BLGlobalField::getObject(); 00526 00527 // construct the world 00528 worldPhysicalVolume = BLGroup::constructWorld(); 00529 00530 // make it invisible 00531 worldPhysicalVolume->GetLogicalVolume()->SetVisAttributes(G4VisAttributes::GetInvisible()); 00532 00533 return worldPhysicalVolume; 00534 }
void BLManager::PreUserTrackingAction | ( | const G4Track * | track | ) |
PreUserTrackingAction() from G4UserTrackingAction.
References beamZStep, BLAssert, currentTrack, currentZStep, eventID, fpTrackingManager, BLCoordinates::getCLZ(), BLRunManager::getCollectiveMode(), getExternalParentID(), BLTrackInfo::getExternalTrackID(), indexZStep, BLCoordinates::isValid(), nextSecondaryTrackID, nStuckSteps, prevZ, primaryParentID, primaryTrackID, REFERENCE, referenceZStep, runManager, BLTrackInfo::setExternalParentID(), BLTrackInfo::setExternalTrackID(), BLCoordinates::setGlobal(), ZStepLimiter::setMaxStep(), state, steppingVerbose, trackID, trackIDMap, trackingActionVector, TUNE, tuneZStep, and zTolerance.
Referenced by BLManager_UserTrackingAction::PreUserTrackingAction().
00674 { 00675 static bool first=true; 00676 if(first) { 00677 first=false; 00678 std::vector<BLManager::TrackingAction*>::iterator i; 00679 for(i=trackingActionVector.begin(); 00680 i<trackingActionVector.end(); ++i) { 00681 (*i)->SetTrackingManagerPointer(fpTrackingManager); 00682 } 00683 } 00684 00685 nStuckSteps = 0; 00686 00687 // link a BLTrackInfo into the track (called coord -- historical) 00688 BLTrackInfo *coord = (BLTrackInfo *)track->GetUserInformation(); 00689 if(!coord || !coord->isValid()) { 00690 coord = new BLTrackInfo(); 00691 ((G4Track *)track)->SetUserInformation(coord); 00692 } 00693 coord->setGlobal(track->GetPosition(),track->GetGlobalTime()); 00694 00695 // set the external trackID and parentID 00696 trackID = coord->getExternalTrackID(); 00697 if(trackID <= 0) { 00698 if(primaryTrackID >= 0) { 00699 trackID = primaryTrackID; 00700 coord->setExternalTrackID(trackID); 00701 coord->setExternalParentID(primaryParentID); 00702 primaryTrackID = primaryParentID = -1; 00703 } else { 00704 trackID = nextSecondaryTrackID++; 00705 coord->setExternalTrackID(trackID); 00706 coord->setExternalParentID( 00707 trackIDMap[track->GetParentID()]); 00708 } 00709 } 00710 trackIDMap[track->GetTrackID()] = trackID; 00711 currentTrack = track; 00712 00713 // print header, if appropriate 00714 if(steppingVerbose && !runManager->getCollectiveMode()) { 00715 printf("=========== EventID %d TrackID %d %s ParentID %d", 00716 eventID,trackID, 00717 track->GetDefinition()->GetParticleName().c_str(), 00718 getExternalParentID(track)); 00719 if(getExternalParentID(track) == 0) { 00720 printf(" CreatorProcess=Beam"); 00721 } else { 00722 const G4VProcess *proc = track->GetCreatorProcess(); 00723 if(proc != 0) 00724 printf(" CreatorProcess=%s", 00725 proc->GetProcessName().c_str()); 00726 } 00727 printf(" ===========\n"); 00728 } 00729 00730 // set currentZStep 00731 if(state == TUNE) 00732 currentZStep = &tuneZStep; 00733 else if(state == REFERENCE) 00734 currentZStep = &referenceZStep; 00735 else 00736 currentZStep = &beamZStep; 00737 indexZStep = 1; 00738 prevZ = coord->getCLZ(); 00739 while(prevZ > (*currentZStep)[indexZStep].z) { 00740 ++indexZStep; 00741 } 00742 BLAssert(indexZStep < currentZStep->size()); 00743 G4double dz = (*currentZStep)[indexZStep].z - prevZ; 00744 if(dz < zTolerance) dz = zTolerance; 00745 ZStepLimiter::setMaxStep(dz); 00746 00747 // loop over all registered user tracking actions 00748 std::vector<BLManager::TrackingAction*>::iterator i; 00749 for(i=trackingActionVector.begin(); i<trackingActionVector.end(); ++i) { 00750 (*i)->PreUserTrackingAction(track); 00751 } 00752 }
void BLManager::PostUserTrackingAction | ( | const G4Track * | track | ) |
PostUserTrackingAction() from G4UserTrackingAction.
References currentTrack, BLRunManager::getCollectiveMode(), runManager, trackID, and trackingActionVector.
Referenced by BLManager_UserTrackingAction::PostUserTrackingAction().
00755 { 00756 std::vector<BLManager::TrackingAction*>::iterator i; 00757 for(i=trackingActionVector.begin(); i<trackingActionVector.end(); ++i) { 00758 (*i)->PostUserTrackingAction(track); 00759 } 00760 00761 // delete the BLcoordinates 00762 G4VUserTrackInformation *ti = track->GetUserInformation(); 00763 if(ti && !runManager->getCollectiveMode()) { 00764 delete ti; 00765 ((G4Track *)track)->SetUserInformation(0); 00766 } 00767 00768 trackID = -9999; 00769 currentTrack = 0; 00770 }
void BLManager::BeginOfRunAction | ( | const G4Run * | run | ) |
BeginOfRunAction() from G4UserRunAction.
References BEAM, beamRunActionVector, endRun, eventsProcessed, runActionVector, startEvent, startRun, state, and BLTime::time().
Referenced by BLManager_UserRunAction::BeginOfRunAction().
00571 { 00572 startEvent = startRun = BLTime::time(); 00573 eventsProcessed = 0; 00574 endRun = false; 00575 00576 G4UImanager* UI = G4UImanager::GetUIpointer(); 00577 UI->ApplyCommand("/control/verbose 0"); 00578 UI->ApplyCommand("/run/verbose 0"); 00579 UI->ApplyCommand("/event/verbose 0"); 00580 UI->ApplyCommand("/tracking/verbose 0"); 00581 UI->ApplyCommand("/hits/verbose 0"); 00582 UI->ApplyCommand("/material/verbose 0"); 00583 UI->ApplyCommand("/process/setVerbose 0 all"); 00584 UI->ApplyCommand("/process/verbose 0"); 00585 UI->ApplyCommand("/process/eLoss/verbose 0"); 00586 00587 for(unsigned int i=0; i<runActionVector.size(); ++i) 00588 runActionVector[i]->BeginOfRunAction(run); 00589 00590 if(state == BEAM) { 00591 for(unsigned int i=0; i<beamRunActionVector.size(); ++i) 00592 beamRunActionVector[i]->BeginOfRunAction(run); 00593 } 00594 }
void BLManager::EndOfRunAction | ( | const G4Run * | run | ) |
EndOfRunAction() from G4UserRunAction.
References BEAM, beamRunActionVector, eventsProcessed, REFERENCE, runActionVector, startRun, state, BLTime::time(), and TUNE.
Referenced by BLManager_UserRunAction::EndOfRunAction().
00597 { 00598 if(state == BEAM) { 00599 for(unsigned int i=0; i<beamRunActionVector.size(); ++i) 00600 beamRunActionVector[i]->EndOfRunAction(run); 00601 } 00602 00603 for(unsigned int i=0; i<runActionVector.size(); ++i) 00604 runActionVector[i]->EndOfRunAction(run); 00605 00606 if(state == TUNE || state == REFERENCE) 00607 ++eventsProcessed; 00608 00609 printf("Run complete %d Events %ld seconds\n", 00610 eventsProcessed,(long)(BLTime::time()-startRun)); 00611 }
void BLManager::BeginOfEventAction | ( | const G4Event * | event | ) |
BeginOfEventAction() from G4UserEventAction.
References BEAM, beamEventActionVector, endRun, eventActionVector, eventTimeLimit, BLRunManager::getCollectiveMode(), runManager, BLAlarm::set(), startEvent, state, steppingVerbose, and BLTime::time().
Referenced by BLRunManager::beamOnCollective(), BLRunManager::beginEvent(), and BLManager_UserEventAction::BeginOfEventAction().
00614 { 00615 if(endRun || event->IsAborted()) return; 00616 00617 if(state == BEAM) { 00618 00619 int evId = event->GetEventID(); 00620 if(steppingVerbose && !runManager->getCollectiveMode()) 00621 printf("\n\n=================== Event %d ==================\n", 00622 evId); 00623 00624 // call all registered action-s 00625 std::vector<BLManager::EventAction*>::iterator i; 00626 for(i=beamEventActionVector.begin(); i<beamEventActionVector.end(); ++i) { 00627 (*i)->BeginOfEventAction(event); 00628 } 00629 } 00630 00631 // call all registered action-s 00632 std::vector<BLManager::EventAction*>::iterator i; 00633 for(i=eventActionVector.begin(); i<eventActionVector.end(); ++i) { 00634 (*i)->BeginOfEventAction(event); 00635 } 00636 00637 startEvent = BLTime::time(); 00638 // add 10 seconds so the test in UserSteppingAction() can try first 00639 // (it kills 1 event, BLAlarm kills the entire job). 00640 if(eventTimeLimit > 0) BLAlarm::set(eventTimeLimit+10); 00641 }
void BLManager::EndOfEventAction | ( | const G4Event * | event | ) |
EndOfEventAction() from G4UserEventAction.
References BEAM, beamEventActionVector, BLAlarm::clear(), endRun, eventActionVector, eventsProcessed, BLNTuple::flushAll(), histoUpdate, incrEventsProcessed(), BLSignal::sigusr1(), state, and BLNTuple::update().
Referenced by BLRunManager::endEvent(), and BLManager_UserEventAction::EndOfEventAction().
00644 { 00645 BLAlarm::clear(); 00646 00647 if(endRun || event->IsAborted()) return; 00648 00649 // call all registered action-s 00650 std::vector<BLManager::EventAction*>::iterator i; 00651 for(i=eventActionVector.begin(); i<eventActionVector.end(); ++i) { 00652 (*i)->EndOfEventAction(event); 00653 } 00654 00655 if(state != BEAM) return; 00656 00657 // call all registered action-s 00658 for(i=beamEventActionVector.begin(); i<beamEventActionVector.end(); ++i) { 00659 (*i)->EndOfEventAction(event); 00660 } 00661 BLNTuple::update(); 00662 00663 if(endRun) return; // (might have been set in a user action) 00664 00665 incrEventsProcessed(event->GetEventID()); 00666 00667 // update histogram file, if appropriate 00668 if((histoUpdate > 0 && eventsProcessed%histoUpdate == 0 && 00669 eventsProcessed > 0) || BLSignal::sigusr1()) 00670 BLNTuple::flushAll(); 00671 }
void BLManager::GeneratePrimaries | ( | G4Event * | event | ) |
GeneratePrimaries() from G4VUserPrimaryGeneratorAction.
References BEAM, beamIndex, beamVector, endRun, eventID, IDLE, prevEventID, REFERENCE, referenceVector, runManager, setEventID(), skipEvent(), SPECIAL, state, TUNE, and VISUAL.
Referenced by BLCMDtrackermode::callback(), and BLManager_PrimaryGeneratorAction::GeneratePrimaries().
00977 { 00978 if(beamVector.size() == 0) { 00979 G4Exception("BLManager","No beam registered",FatalException, 00980 ""); 00981 } 00982 00983 switch(state) { 00984 case IDLE: 00985 case SPECIAL: 00986 G4Exception("BLManager","Erroneous call to GeneratePrimaries", 00987 FatalException,""); 00988 case TUNE: 00989 setEventID(-2); 00990 event->SetEventID(-2); 00991 while(beamIndex < referenceVector.size()) { 00992 if(referenceVector[beamIndex++]-> 00993 generateReferenceParticle(event)) 00994 return; 00995 } 00996 goto end_run; 00997 case REFERENCE: 00998 setEventID(-1); 00999 event->SetEventID(-1); 01000 while(beamIndex < referenceVector.size()) { 01001 if(referenceVector[beamIndex++]-> 01002 generateReferenceParticle(event)) 01003 return; 01004 } 01005 goto end_run; 01006 case VISUAL: 01007 case BEAM: 01008 for(;;) { 01009 // use tempEvent so we can skip events when necessary 01010 G4Event *tempEvent = new G4Event(*event); 01011 // default event # (nextBeamEvent() can change it) 01012 eventID = prevEventID + 1; 01013 event->SetEventID(eventID); 01014 // generate the tempEvent, indexing through beamVector[] 01015 if(beamIndex >= beamVector.size()) { 01016 goto end_run; 01017 } 01018 while(!beamVector[beamIndex]->nextBeamEvent(tempEvent)){ 01019 if(++beamIndex >= beamVector.size()) 01020 goto end_run; 01021 } 01022 // update eventID (may have changed) 01023 eventID = tempEvent->GetEventID(); 01024 prevEventID = eventID; // (this increments it, too) 01025 // skip events when necessary 01026 if(!skipEvent(eventID)) { 01027 *event = *tempEvent; 01028 delete tempEvent; 01029 break; 01030 } 01031 delete tempEvent; 01032 } 01033 break; 01034 } 01035 return; 01036 end_run: 01037 // RunManager cannot abort the event from inside 01038 // UserGeneratePrimaries(), so we do a soft abort 01039 // to the RunManager, and abort the event ourself. 01040 // The result is the same as a hard abort. 01041 runManager->AbortRun(true); 01042 event->SetEventAborted(); 01043 beamIndex = 0; 01044 endRun = true; 01045 }
G4ClassificationOfNewTrack BLManager::ClassifyNewTrack | ( | const G4Track * | track | ) |
ClassifyNewTrack() from G4UserStackingAction.
References stackingActionVector.
Referenced by BLManager_UserStackingAction::ClassifyNewTrack().
01048 { 01049 std::vector<BLManager::StackingAction*>::iterator i; 01050 for(i=stackingActionVector.begin(); i<stackingActionVector.end(); ++i) { 01051 G4ClassificationOfNewTrack c = (*i)->ClassifyNewTrack(track); 01052 if(c == fKill) return c; 01053 } 01054 return fUrgent; 01055 }
void BLManager::NewStage | ( | ) |
NewStage() from G4userStackingAction.
References stackingActionVector.
Referenced by BLManager_UserStackingAction::NewStage().
01058 { 01059 std::vector<BLManager::StackingAction*>::iterator i; 01060 for(i=stackingActionVector.begin(); i<stackingActionVector.end(); ++i) { 01061 (*i)->NewStage(); 01062 } 01063 }
void BLManager::PrepareNewEvent | ( | ) |
PrepareNewEvent() from G4StackingAction.
References stackingActionVector.
Referenced by BLManager_UserStackingAction::PrepareNewEvent().
01066 { 01067 std::vector<BLManager::StackingAction*>::iterator i; 01068 for(i=stackingActionVector.begin(); i<stackingActionVector.end(); ++i) { 01069 (*i)->PrepareNewEvent(); 01070 } 01071 }
void BLManager::clearTrackIDMap | ( | ) | [inline] |
clearTrackIDMap() clears the TrackID map.
References trackIDMap.
Referenced by BLCMDreference::generateReferenceParticle(), BLCMDparticlesource::nextBeamEvent(), and BLCMDbeam::nextBeamEvent().
00514 { trackIDMap.clear(); }
void BLManager::setNextSecondaryTrackID | ( | int | next | ) | [inline] |
setNextSecondaryTrackID() sets the external TrackID for the next secondary track. Automatically incremented for subsequent secondaries.
References nextSecondaryTrackID.
Referenced by BLCMDreference::generateReferenceParticle(), BLCMDparticlesource::nextBeamEvent(), and BLCMDbeam::nextBeamEvent().
00519 { nextSecondaryTrackID = next; }
int BLManager::getNextSecondaryTrackID | ( | ) | [inline] |
References nextSecondaryTrackID.
Referenced by BLRunManager::beamOnCollective().
00520 { return nextSecondaryTrackID; }
int BLManager::getExternalTrackID | ( | const G4Track * | track | ) |
getExternalTrackID() returns the external TrackID for the given track. In collective mode, internal and external trackID-s are the same.
References BLRunManager::getCollectiveMode(), BLTrackInfo::getExternalTrackID(), BLCoordinates::isValid(), runManager, and trackIDMap.
Referenced by BLTrackNTuple::appendTrack(), BLCMDmovie::beamStep(), BLCMDtrace::newTrace(), steppingVerbosePrint(), and BLCMDtracker::UserZSteppingAction().
01074 { 01075 if(runManager->getCollectiveMode()) return track->GetTrackID(); 01076 BLTrackInfo *p = (BLTrackInfo *)track->GetUserInformation(); 01077 if(p && p->isValid()) 01078 return p->getExternalTrackID(); 01079 return trackIDMap[track->GetTrackID()]; 01080 }
int BLManager::getExternalParentID | ( | const G4Track * | track | ) |
getExternalParentID() returns the external ParentID for the given track. In collective mode, internal and external trackID-s are the same.
References BLRunManager::getCollectiveMode(), BLTrackInfo::getExternalParentID(), BLCoordinates::isValid(), runManager, and trackIDMap.
Referenced by BLTrackNTuple::appendTrack(), BLCMDmovie::beamStep(), PreUserTrackingAction(), steppingVerbosePrint(), and BLCMDtracker::UserZSteppingAction().
01083 { 01084 if(runManager->getCollectiveMode()) return track->GetParentID(); 01085 BLTrackInfo *p = (BLTrackInfo *)track->GetUserInformation(); 01086 if(p && p->isValid()) 01087 return p->getExternalParentID(); 01088 return trackIDMap[track->GetParentID()]; 01089 }
void BLManager::setExternalTrackID | ( | G4Track * | track, | |
int | trackID, | |||
int | parentID | |||
) |
setExternalTrackID() will set the external trackID and parentID. if trackID<0 uses nextSecondarTrackID++.
References BLCoordinates::isValid(), nextSecondaryTrackID, BLTrackInfo::setExternalParentID(), and BLTrackInfo::setExternalTrackID().
01092 { 01093 if(trackID < 0) trackID = nextSecondaryTrackID++; 01094 01095 BLTrackInfo *ti = (BLTrackInfo *)track->GetUserInformation(); 01096 if(!ti || !ti->isValid()) { 01097 ti = new BLTrackInfo(); 01098 track->SetUserInformation(ti); 01099 } 01100 ti->setExternalTrackID(trackID); 01101 ti->setExternalParentID(parentID); 01102 }
int BLManager::getPrimaryTrackID | ( | ) | [inline] |
getPrimaryTrackID() returns the primaryTrackID set by the beam command
References primaryTrackID.
Referenced by BLRunManager::beamOnCollective(), and BLRunManager::getNextBeamEventAndTrack().
00537 { return primaryTrackID; }
int BLManager::getPrimaryParentID | ( | ) | [inline] |
getPrimaryParentID() returns the primaryParentID set by the beam command
References primaryParentID.
Referenced by BLRunManager::getNextBeamEventAndTrack().
00540 { return primaryParentID; }
void BLManager::setPrimaryTrackID | ( | int | t, | |
int | p | |||
) | [inline] |
setPrimaryTrackID() sets track and parent IDs for a primary track.
References primaryParentID, and primaryTrackID.
Referenced by BLCMDtracker::fitTrack(), BLCMDreference::generateReferenceParticle(), BLCMDparticlesource::nextBeamEvent(), BLCMDbeam::nextBeamEvent(), and BLCMDtracker::operator()().
00544 { primaryTrackID = t; primaryParentID = p; }
G4bool BLManager::Notify | ( | const char * | originOfException, | |
const char * | exceptionCode, | |||
G4ExceptionSeverity | severity, | |||
const char * | description | |||
) |
Notify() from G4VExceptionHandler.
References allExceptions, BLAlarm::clear(), BLNTuple::closeAll(), currentTrack, eventID, eventsAborted, exceptionCount, fatalExceptions, g4bl_exit(), handleCallbacks(), BLMPI::isMPI(), runManager, BLNTuple::summary(), trackID, and warnings.
01127 { 01128 FILE *out=(BLMPI::isMPI() ? stderr : stdout); // originally was stderr 01129 fflush(stdout); 01130 01131 G4bool toBeAborted = false; 01132 const char *p="UNKNOWN"; 01133 switch(severity) { 01134 case FatalException: 01135 p = "Fatal Exception"; 01136 toBeAborted = true; 01137 ++fatalExceptions; 01138 break; 01139 case FatalErrorInArgument: 01140 p = "Fatal Error in Argument"; 01141 toBeAborted = true; 01142 ++fatalExceptions; 01143 break; 01144 case RunMustBeAborted: 01145 p = "Run Must Be Aborted"; 01146 runManager->AbortRun(false); 01147 toBeAborted = false; 01148 ++fatalExceptions; 01149 break; 01150 case EventMustBeAborted: 01151 BLAlarm::clear(); 01152 p = "Event must be Aborted"; 01153 runManager->AbortEvent(); 01154 toBeAborted = false; 01155 ++eventsAborted; 01156 break; 01157 default: 01158 p = "Warning"; 01159 toBeAborted = false; 01160 ++warnings; 01161 break; 01162 } 01163 01164 // thin out the printing of many similar exceptions 01165 int n = ++exceptionCount[exceptionCode]; 01166 int interval = 1; 01167 if(n >= 1000) 01168 interval = 1000; 01169 else if(n >= 100) 01170 interval = 100; 01171 else if(n >= 10) 01172 interval = 10; 01173 if(allExceptions) interval = 1; 01174 if(n%interval != 0 && !toBeAborted) 01175 return false; 01176 01177 fprintf(out,"**************************************************************************\n"); 01178 fprintf(out,"*** G4Exception: %s\n",exceptionCode); 01179 fprintf(out,"*** severity: %s\n",p); 01180 fprintf(out,"*** issued by: %s\n",originOfException); 01181 if(strlen(description) > 0) 01182 fprintf(out,"*** description: %s\n",description); 01183 if(strstr(description,"Missing mandatory data")) { 01184 // brief description already printed, just URL is missing 01185 const char *p = getenv("G4BL_DIR"); 01186 if(!p) p = "Installation_Directory"; 01187 fprintf(out, 01188 "*** The URL is http://geant4.cern.ch\n" 01189 "*** If you unpack it into %s/data\n" 01190 "*** then the g4bl and g4blgui scripts will" 01191 " automatically find it.\n",p); 01192 } else if(G4StateManager::GetStateManager()->GetCurrentState() == 01193 G4State_EventProc && currentTrack != 0) { 01194 fprintf(out,"*** EventID: %d TrackID: %d %s KE=%.4g MeV\n", 01195 eventID,trackID, 01196 currentTrack->GetDefinition()-> 01197 GetParticleName().c_str(), 01198 currentTrack->GetKineticEnergy()); 01199 } 01200 if(interval > 1) 01201 fprintf(out,"*** printing: every %d-th occurrence\n", 01202 interval); 01203 fprintf(out,"**************************************************************************\n"); 01204 01205 fflush(out); // (could be stdout) 01206 01207 if(toBeAborted) { 01208 fprintf(out,"g4beamline: attempting to close up after fatal G4Exception...\n"); 01209 BLNTuple::summary(); 01210 BLNTuple::closeAll(); 01211 handleCallbacks(2); 01212 fprintf(out,"g4beamline: simulation aborted\n"); 01213 g4bl_exit(99); 01214 } 01215 01216 return toBeAborted; 01217 }
void BLManager::exceptionSummary | ( | ) |
exceptionSummary() prints a summary of all exceptions
References eventsAborted, exceptionCount, fatalExceptions, stuckTracks, and warnings.
Referenced by handleCallbacks().
01220 { 01221 printf("\nExceptions: %d Fatal, %d Events Aborted, " 01222 "%d Stuck Tracks (stopped), %d Warnings\n", 01223 fatalExceptions,eventsAborted,stuckTracks,warnings); 01224 std::map<G4String,int>::iterator it; 01225 for(it=exceptionCount.begin(); it!=exceptionCount.end(); ++it) { 01226 G4String except = it->first; 01227 int count = it->second; 01228 printf(" %6d times: %s\n",count,except.c_str()); 01229 } 01230 }
void BLManager::gsl_error_handler | ( | const char * | reason, | |
const char * | file, | |||
int | lineno, | |||
int | gsl_errno | |||
) | [static] |
GSL error handler.
References snprintf.
Referenced by BLManager().
01528 { 01529 // ignore underflows 01530 if(gsl_errno == GSL_EUNDRFLW) return; 01531 01532 char tmp[256]; 01533 snprintf(tmp,sizeof(tmp)-1,"file=%s line=%d gsl_errno=%d",file,lineno, 01534 gsl_errno); 01535 G4Exception("BLManager",reason,EventMustBeAborted,tmp); 01536 }
BLManager * BLManager::blManager = 0 [static, private] |
Referenced by BLManager(), getObject(), handleCallbacks(), and ~BLManager().
bool BLManager::initialized = false [static, private] |
Referenced by initialize(), and isInitialized().
BLRunManager* BLManager::runManager [private] |
Referenced by BeginOfEventAction(), BLManager(), delayedConstruction(), displayVisual(), GeneratePrimaries(), getExternalParentID(), getExternalTrackID(), incrEventsProcessed(), initialize(), Notify(), PostUserTrackingAction(), PreUserTrackingAction(), registerPhysics(), steppingVerbosePrint(), trackBeam(), and trackTuneAndReferenceParticles().
BLManagerState BLManager::state [private] |
G4int BLManager::steppingVerbose [private] |
unsigned int BLManager::beamIndex [private] |
Referenced by BLManager(), displayVisual(), GeneratePrimaries(), initialize(), trackBeam(), and trackTuneAndReferenceParticles().
G4int BLManager::histoUpdate [private] |
Referenced by BLManager(), EndOfEventAction(), and initialize().
time_t BLManager::startRun [private] |
Referenced by BeginOfRunAction(), BLManager(), EndOfRunAction(), and incrEventsProcessed().
time_t BLManager::startEvent [private] |
Referenced by BeginOfEventAction(), BeginOfRunAction(), BLManager(), and UserSteppingAction().
G4int BLManager::eventTimeLimit [private] |
Referenced by BeginOfEventAction(), BLManager(), getEventTimeLimit(), initialize(), setEventTimeLimit(), and UserSteppingAction().
BLPhysics* BLManager::physics [private] |
Referenced by BLManager(), displayVisual(), getPhysics(), initialize(), registerPhysics(), trackBeam(), and trackTuneAndReferenceParticles().
G4VPhysicalVolume* BLManager::worldPhysicalVolume [private] |
Referenced by BLManager(), Construct(), displayGeometry(), and getWorldPhysicalVolume().
G4int BLManager::eventID [private] |
Referenced by BLManager(), GeneratePrimaries(), getEventID(), Notify(), PreUserTrackingAction(), and setEventID().
G4int BLManager::trackID [private] |
Referenced by BLManager(), Notify(), PostUserTrackingAction(), and PreUserTrackingAction().
const G4Track* BLManager::currentTrack [private] |
Referenced by BLManager(), Notify(), PostUserTrackingAction(), and PreUserTrackingAction().
G4bool BLManager::allExceptions [private] |
Referenced by BLManager(), Notify(), and showAllExceptions().
G4int BLManager::fatalExceptions [private] |
Referenced by BLManager(), exceptionSummary(), and Notify().
G4int BLManager::eventsAborted [private] |
Referenced by BLManager(), exceptionSummary(), and Notify().
G4int BLManager::stuckTracks [private] |
Referenced by BLManager(), exceptionSummary(), and UserSteppingAction().
G4int BLManager::warnings [private] |
Referenced by BLManager(), exceptionSummary(), and Notify().
G4int BLManager::prevEventID [private] |
Referenced by BLManager(), GeneratePrimaries(), and setEventID().
G4int BLManager::eventsProcessed [private] |
Referenced by BeginOfRunAction(), BLManager(), EndOfEventAction(), EndOfRunAction(), and incrEventsProcessed().
bool BLManager::endRun [private] |
Referenced by BeginOfEventAction(), BeginOfRunAction(), BLManager(), EndOfEventAction(), and GeneratePrimaries().
PRNGSeedMethod BLManager::seedMethod [private] |
Referenced by BLManager(), getPRNGSeedMethod(), and setPRNGSeedMethod().
std::set<int> BLManager::keepEventList [private] |
Referenced by setKeepEvent(), and skipEvent().
std::set<int> BLManager::skipEventList [private] |
Referenced by setSkipEvent(), and skipEvent().
G4SteppingManager* BLManager::fpSteppingManager [private] |
Referenced by getSteppingManager(), and setSteppingManager().
G4TrackingManager* BLManager::fpTrackingManager [private] |
Referenced by PreUserTrackingAction(), and setTrackingManager().
std::vector<BLBeam*> BLManager::beamVector [private] |
Referenced by GeneratePrimaries(), initialize(), registerBeam(), and trackBeam().
std::vector<BLBeam*> BLManager::referenceVector [private] |
Referenced by GeneratePrimaries(), initialize(), nReference(), registerReference(), and trackTuneAndReferenceParticles().
std::vector<RunAction*> BLManager::runActionVector [private] |
Referenced by BeginOfRunAction(), EndOfRunAction(), and registerRunAction().
std::vector<RunAction*> BLManager::beamRunActionVector [private] |
Referenced by BeginOfRunAction(), EndOfRunAction(), and registerRunAction().
std::vector<EventAction*> BLManager::eventActionVector [private] |
Referenced by BeginOfEventAction(), EndOfEventAction(), and registerEventAction().
std::vector<EventAction*> BLManager::beamEventActionVector [private] |
Referenced by BeginOfEventAction(), EndOfEventAction(), and registerEventAction().
std::vector<TrackingAction*> BLManager::trackingActionVector [private] |
Referenced by PostUserTrackingAction(), PreUserTrackingAction(), and registerTrackingAction().
std::vector<BLCallback*> BLManager::preReferenceCallbackVector [private] |
Referenced by handleCallbacks(), and registerCallback().
std::vector<BLCallback*> BLManager::postReferenceCallbackVector [private] |
Referenced by handleCallbacks(), and registerCallback().
std::vector<BLCallback*> BLManager::postTrackingCallbackVector [private] |
Referenced by handleCallbacks(), and registerCallback().
std::vector<BLCallback*> BLManager::replaceMainLoopCallbackVector [private] |
Referenced by handleCallbacks(), and registerCallback().
std::vector<BLCallback*> BLManager::visualizationCallbackVector [private] |
Referenced by handleCallbacks(), and registerCallback().
std::vector<SteppingAction*> BLManager::allStepVector [private] |
Referenced by registerSteppingAction(), and UserSteppingAction().
std::map<G4VPhysicalVolume*,SteppingAction*> BLManager::allStepMap [private] |
Referenced by registerSteppingAction(), and UserSteppingAction().
std::map<G4VPhysicalVolume*,SteppingAction*> BLManager::tpStepMap [private] |
Referenced by registerTuneParticleStep(), and UserSteppingAction().
std::map<G4VPhysicalVolume*,SteppingAction*> BLManager::rpStepMap [private] |
Referenced by registerReferenceParticleStep(), and UserSteppingAction().
std::vector<SteppingAction*> BLManager::tpStepVector [private] |
Referenced by registerTuneParticleStep(), and UserSteppingAction().
std::vector<SteppingAction*> BLManager::rpStepVector [private] |
Referenced by registerReferenceParticleStep(), and UserSteppingAction().
std::map<G4VPhysicalVolume*,SteppingAction*> BLManager::beamStepMap [private] |
Referenced by registerBeamStep(), and UserSteppingAction().
std::vector<SteppingAction*> BLManager::beamStepVector [private] |
Referenced by registerBeamStep(), and UserSteppingAction().
std::vector<int> BLManager::verboseFormat [private] |
Referenced by appendVerboseFormat(), BLManager(), and steppingVerbosePrint().
G4double BLManager::zTolerance [private] |
Referenced by BLManager(), initialize(), PreUserTrackingAction(), and UserSteppingAction().
std::vector<ZStep> BLManager::tuneZStep [private] |
Referenced by BLManager(), PreUserTrackingAction(), and registerZStep().
std::vector<ZStep> BLManager::referenceZStep [private] |
Referenced by BLManager(), PreUserTrackingAction(), and registerZStep().
std::vector<ZStep> BLManager::beamZStep [private] |
Referenced by BLManager(), PreUserTrackingAction(), and registerZStep().
std::vector<ZStep>* BLManager::currentZStep [private] |
Referenced by BLManager(), PreUserTrackingAction(), and UserSteppingAction().
unsigned BLManager::indexZStep [private] |
Referenced by BLManager(), PreUserTrackingAction(), and UserSteppingAction().
G4double BLManager::prevZ [private] |
Referenced by BLManager(), PreUserTrackingAction(), and UserSteppingAction().
G4int BLManager::nStuckSteps [private] |
Referenced by BLManager(), PreUserTrackingAction(), and UserSteppingAction().
std::vector<StackingAction*> BLManager::stackingActionVector [private] |
Referenced by ClassifyNewTrack(), NewStage(), PrepareNewEvent(), and registerStackingAction().
int BLManager::nextSecondaryTrackID [private] |
Referenced by BLManager(), getNextSecondaryTrackID(), PreUserTrackingAction(), setExternalTrackID(), and setNextSecondaryTrackID().
std::map<G4int,G4int> BLManager::trackIDMap [private] |
Referenced by clearTrackIDMap(), getExternalParentID(), getExternalTrackID(), and PreUserTrackingAction().
int BLManager::primaryTrackID [private] |
Referenced by BLManager(), getPrimaryTrackID(), PreUserTrackingAction(), and setPrimaryTrackID().
int BLManager::primaryParentID [private] |
Referenced by BLManager(), getPrimaryParentID(), PreUserTrackingAction(), and setPrimaryTrackID().
std::vector<BLUserCode*> BLManager::userCodeVector [private] |
Referenced by getUserCodeInstances(), and registerUserCode().
std::map<G4String,int> BLManager::exceptionCount [private] |
Referenced by exceptionSummary(), and Notify().