BLCMDtrace Class Reference

Inheritance diagram for BLCMDtrace:

BLCommand BLManager::SteppingAction BLManager::TrackingAction BLCollectiveComputation

List of all members.


Detailed Description

class BLCMDtrace implements the trace command.

Public Member Functions

 BLCMDtrace ()
 Constructor.
 BLCMDtrace (const BLCMDtrace &r)
G4String commandName ()
int command (BLArgumentVector &argv, BLArgumentMap &namedArgs)
void defineNamedArgs ()
void help (bool detailed)
 help() prints help text.
void newTrace (const G4Track *firstTrack)
void UserSteppingAction (const G4Step *step)
 UserSteppingAction() from BLManager::SteppingAction.
void PreUserTrackingAction (const G4Track *track)
 PreUserTrackingAction() from BLManager::TrackingAction.
void PostUserTrackingAction (const G4Track *track)
 PostUserTrackingAction() from BLManager::TrackingAction.
virtual void beginCollectiveTracking (std::vector< BLTrackData > &v)
 from BLCollectiveComputation
virtual void collectiveStep (std::vector< BLTrackData > &v)
virtual void endCollectiveTracking (std::vector< BLTrackData > &v)

Private Attributes

int nTrace
G4String format
G4int oneNTuple
G4int primaryOnly
G4int traceTune
G4String filename
G4String require
G4String coordinates
BLCoordinateType coordinateType
bool doTrace
TraceNTupletrace
BLEvaluatoreval
bool init
BLManagermanager
BLRunManagerrunManager

Static Private Attributes

static TraceNTupleallTrace = 0
static TraceNTupletuneTrace = 0
static TraceNTuplereferenceTrace = 0

Friends

class TraceNTuple

Constructor & Destructor Documentation

BLCMDtrace::BLCMDtrace (  ) 

Constructor.

References BLCMDTYPE_DATA, BLCOORD_CENTERLINE, coordinates, coordinateType, doTrace, eval, filename, format, init, manager, nTrace, oneNTuple, primaryOnly, BLCommand::registerCommand(), require, runManager, BLCommand::setDescription(), BLCommand::setSynopsis(), trace, and traceTune.

00134                        : BLCommand(), BLManager::SteppingAction(), 
00135             BLManager::TrackingAction(), BLCollectiveComputation(),
00136             format()
00137 {
00138         registerCommand(BLCMDTYPE_DATA);
00139         setSynopsis("Specifies tracing of tracks.");
00140         setDescription("Generates a separate NTuple for each track, with 1\n"
00141                 "row per step, unless oneNTuple is nonzero (in which case "
00142                 "all tracks are put into a single NTuple).\n"
00143                 "So format=ascii generates one file per track with names\n"
00144                 "generated by the pattern in filename (first %d is replaced "
00145                 "by event #, second %d is replaced by trackId); for "
00146                 "oneNTuple, the default filename is AllTracks.txt.\n"
00147                 "\nNote that without a trace command no traces are generated, "
00148                 "so to trace just the tune and reference particles include "
00149                 "a trace command with no arguments.\n\n"
00150                 "In collective tracking mode, oneNTuple must be nonzero, and "
00151                 "the entries will be generated only at collective steps "
00152                 "(usually at a specified deltaT).\n\n"
00153                 "Unlike other NTuple commands, the require expression applies "
00154                 "to entire tracks, not individual entries.\n\n"
00155                 "The standard NTuple fields are:\n"
00156                 "    x,y,z (mm)\n"
00157                 "    Px,Py,Pz (MeV/c)\n"
00158                 "    t (ns)\n"
00159                 "    PDGid (11=e-, 13=mu-, 22=gamma, 211=pi+, 2212=proton, ...)\n"
00160                 "    EventID (may be inexact above 16,777,215)\n"
00161                 "    TrackID\n"
00162                 "    ParentID (0 => primary particle)\n"
00163                 "    Weight (defaults to 1.0)\n"
00164                 "The trace includes the following fields:\n"
00165                 "    Bx, By, Bz (Tesla)\n"
00166                 "    Ex, Ey, Ez (Megavolts/meter)\n\n"
00167                 "The following additional fields are appended for "
00168                 "format=Extended, format=asciiExtended, and "
00169                 "format=rootExtended:\n"
00170                 "    ProperTime (ns)\n"
00171                 "    PathLength (mm)\n"
00172                 "    PolX, PolY, PolZ (polarization)\n"
00173                 "    InitialKE (MeV when track was created)\n\n"
00174                 "Valid Formats (ignore case): ");
00175         nTrace = 0;
00176         format = "";
00177         oneNTuple = 0;
00178         primaryOnly = 0;
00179         traceTune = 1;     //trace  tune track by default (previous behavior)
00180         filename = "";
00181         require = "";
00182         coordinates = "Centerline";
00183         coordinateType = BLCOORD_CENTERLINE;
00184         doTrace = true;
00185         trace = 0;
00186         eval = 0;
00187         init = false;
00188         manager = 0;    // too early to set these -- see command()
00189         runManager = 0;
00190 }

BLCMDtrace::BLCMDtrace ( const BLCMDtrace r  ) 

References coordinates, coordinateType, doTrace, eval, filename, format, init, manager, nTrace, oneNTuple, primaryOnly, require, runManager, trace, and traceTune.

00192                                           : BLCommand(r), 
00193             BLManager::SteppingAction(), BLManager::TrackingAction(),
00194             BLCollectiveComputation(), format()
00195 {
00196         nTrace = r.nTrace;
00197         format = r.format;
00198         oneNTuple = r.oneNTuple;
00199         primaryOnly = r.primaryOnly;
00200         traceTune= r.traceTune;
00201         filename = r.filename;
00202         require = r.require;
00203         coordinates = r.coordinates;
00204         coordinateType = r.coordinateType;
00205         doTrace = r.doTrace;
00206         trace = r.trace;
00207         eval = 0;
00208         init = r.init;
00209         manager = r.manager;
00210         runManager = r.runManager;
00211 }


Member Function Documentation

G4String BLCMDtrace::commandName (  )  [inline, virtual]

Implements BLCommand.

00100 { return "trace"; }

int BLCMDtrace::command ( BLArgumentVector argv,
BLArgumentMap namedArgs 
) [virtual]

Implements BLCommand.

References coordinates, coordinateType, eval, filename, format, BLCoordinates::getCoordinateType(), BLRunManager::getObject(), BLManager::getObject(), BLCommand::handleNamedArgs(), init, manager, oneNTuple, BLCommand::print(), BLRunManager::registerCollectiveComputation(), BLManager::registerSteppingAction(), BLManager::registerTrackingAction(), require, and runManager.

00214 {
00215         manager = BLManager::getObject();
00216         runManager = BLRunManager::getObject();
00217 
00218         int retval = handleNamedArgs(namedArgs);
00219 
00220         if(!init) {
00221                 init = true;
00222                 manager->registerSteppingAction(this);
00223                 manager->registerTrackingAction(this);
00224                 runManager->registerCollectiveComputation(this);
00225         }
00226 
00227         if(filename == "")
00228                 filename = (oneNTuple!=0 ? "AllTracks" : "Ev%dTrk%d");
00229 
00230         coordinateType = BLCoordinates::getCoordinateType(coordinates);
00231         // requir applies per track, not per step; handled in newTrace()
00232         if(require != "") eval = new BLEvaluator();
00233 
00234         // handle format
00235         for(G4String::size_type i=0; i<format.size(); ++i)
00236                 format[i] = tolower(format[i]);
00237         G4String::size_type j=format.find("trace");
00238         if(j == format.npos) {
00239                 j = format.find("extended");
00240                 if(j == format.npos) format += "Trace";
00241         }
00242 
00243         print("");
00244 
00245         return retval;
00246 }

void BLCMDtrace::defineNamedArgs (  )  [virtual]

Reimplemented from BLCommand.

References BLCommand::argInt(), BLCommand::argString(), coordinates, filename, format, nTrace, oneNTuple, primaryOnly, require, and traceTune.

00249 {
00250         argInt(nTrace,"nTrace","Number of tracks to trace.");
00251         argString(format,"format","Format of the NTuple (see above for list).");
00252         argInt(oneNTuple,"oneNTuple","Nonzero to put all traces into a single NTuple.");
00253         argInt(primaryOnly,"primaryOnly","Nonzero to trace only primary tracks.");
00254         argInt(traceTune,"traceTune","Nonzero to trace tune tracks (default=1).");
00255         argString(filename,"filename","Filename (Ev%dTrk%d.txt or AllTracks.txt).");
00256         argString(filename,"file","synonym for filename.");
00257         argString(require,"require","Expression which must be nonzero to trace the track (default=1)",false);
00258         argString(coordinates,"coordinates","Coordinates: global, centerline, or reference (default=c).");
00259 }

void BLCMDtrace::help ( bool  detailed  )  [inline, virtual]

help() prints help text.

Reimplemented from BLCommand.

References BLCommand::description, BLTrackNTuple::getFormatList(), and BLCommand::help().

00107                                  {
00108                 if(description[description.size()-2] == ':')
00109                         description += BLTrackNTuple::getFormatList(); 
00110                 BLCommand::help(detailed);
00111         }

void BLCMDtrace::newTrace ( const G4Track *  firstTrack  ) 

References allTrace, TraceNTuple::close(), coordinateType, eval, BLEvaluator::evaluate(), filename, format, BLManager::getExternalTrackID(), BLEvaluator::isOK(), manager, oneNTuple, BLCommand::printError(), referenceTrace, require, BLEvaluator::setTrackVariables(), snprintf, trace, TraceNTuple, and tuneTrace.

Referenced by PreUserTrackingAction().

00262 {
00263         int evNum = G4RunManager::GetRunManager()->GetCurrentEvent()->GetEventID();
00264         int trkId = manager->getExternalTrackID(firstTrack);
00265 
00266         if(oneNTuple != 0) {
00267                 if(!allTrace) allTrace = new TraceNTuple(format,"AllTracks",
00268                                                 filename, coordinateType);
00269         } else if(evNum == -2) {
00270                 if(!tuneTrace) tuneTrace = new TraceNTuple(format,
00271                         "TuneParticle","TuneParticle",coordinateType);
00272         } else if(evNum == -1) {
00273                 if(!referenceTrace) referenceTrace = new TraceNTuple(format,
00274                         "ReferenceParticle","ReferenceParticle",coordinateType);
00275         } else {
00276                 if(trace) trace->close();
00277         }
00278         trace = 0;
00279 
00280         if(eval) {
00281                 eval->setTrackVariables(firstTrack,coordinateType);
00282                 G4double v = eval->evaluate(require);
00283                 if(!eval->isOK()) {
00284                         BLCommand::printError("trace: invalid expression require='%s'\n",
00285                                 require.c_str());
00286                         delete eval;
00287                         eval = 0;
00288                         v = 1.0;
00289                 }
00290                 if(v == 0.0) return;
00291         }
00292 
00293         if(allTrace != 0) {
00294                 trace = allTrace;
00295         } else if(evNum == -2) {
00296                 trace = tuneTrace;
00297         } else if(evNum == -1) {
00298                 trace = referenceTrace;
00299         } else if(evNum >= 0) {
00300                 char tmp[128];
00301                 snprintf(tmp,128,filename.c_str(),evNum,trkId);
00302                 trace = new TraceNTuple(format,tmp,tmp,coordinateType);
00303         }
00304 }

void BLCMDtrace::UserSteppingAction ( const G4Step *  step  )  [virtual]

UserSteppingAction() from BLManager::SteppingAction.

Implements BLManager::SteppingAction.

References TraceNTuple::appendTrace(), doTrace, BLRunManager::getCollectiveMode(), BLManager::getState(), manager, runManager, SPECIAL, and trace.

00308 {
00309         if(manager->getState() == SPECIAL || runManager->getCollectiveMode())
00310                 return;
00311 
00312         if(!trace || !doTrace) return;
00313 
00314         G4Track *track = step->GetTrack();
00315         trace->appendTrace(track);
00316 }

void BLCMDtrace::PreUserTrackingAction ( const G4Track *  track  )  [virtual]

PreUserTrackingAction() from BLManager::TrackingAction.

Implements BLManager::TrackingAction.

References TraceNTuple::annotate(), BLCOORD_REFERENCE, TraceNTuple::close(), coordinateType, doTrace, BLRunManager::getCollectiveMode(), manager, newTrace(), nTrace, oneNTuple, primaryOnly, runManager, trace, and traceTune.

00319 {
00320         int evNum = runManager->GetCurrentEvent()->GetEventID();
00321 
00322         doTrace = false;
00323 
00324         if(evNum == -2) {
00325                 doTrace = (traceTune && coordinateType != BLCOORD_REFERENCE);
00326                 if(doTrace) newTrace(track);
00327         } else if(evNum < 0) {
00328                 // cannot use reference coordinates for tune or reference
00329                 doTrace = (coordinateType != BLCOORD_REFERENCE);
00330                 if(doTrace) newTrace(track);
00331         } else if(runManager->getCollectiveMode() && nTrace > 0) {
00332                 if(oneNTuple == 0) {
00333                         G4Exception("trace","Collective mode, oneNTuple=0",
00334                                 FatalException,"oneNTuple must be nonzero");
00335                 }
00336                 newTrace(track);
00337         } else if(nTrace > 0) {
00338                 doTrace = ((primaryOnly == 0) || (track->GetParentID() == 0));
00339                 if(doTrace) {
00340                         --nTrace;
00341                         newTrace(track);
00342                 }
00343         } else {
00344                 if(trace) {
00345                         trace->close();
00346                         trace = 0;
00347                 }
00348         }
00349 
00350         if(trace && doTrace && 
00351                         !runManager->getCollectiveMode()) {
00352                 char tmp[64];
00353                 sprintf(tmp,"# Event %d Track %d",evNum,manager->
00354                                 getExternalTrackID(track));
00355                 trace->annotate(tmp);
00356         }
00357 }

void BLCMDtrace::PostUserTrackingAction ( const G4Track *  track  )  [virtual]

PostUserTrackingAction() from BLManager::TrackingAction.

Implements BLManager::TrackingAction.

References TraceNTuple::annotate(), TraceNTuple::close(), doTrace, BLRunManager::getCollectiveMode(), oneNTuple, runManager, and trace.

00360 {
00361         if(!trace) return;
00362 
00363         if(doTrace && !runManager->getCollectiveMode())
00364                 trace->annotate("");
00365 
00366         int evNum = G4RunManager::GetRunManager()->GetCurrentEvent()->
00367                                                                 GetEventID();
00368         if(oneNTuple == 0 && evNum >= 0) {
00369                 trace->close();
00370                 trace = 0;
00371         }
00372 }

void BLCMDtrace::beginCollectiveTracking ( std::vector< BLTrackData > &  v  )  [virtual]

from BLCollectiveComputation

Reimplemented from BLCollectiveComputation.

00375 {
00376         (void)v;
00377 }

void BLCMDtrace::collectiveStep ( std::vector< BLTrackData > &  v  )  [virtual]

Implements BLCollectiveComputation.

References TraceNTuple::appendTrace(), primaryOnly, runManager, BLRunManager::setCurrentEvent(), and trace.

00380 {
00381         if(trace) {
00382                 for(unsigned i=0; i<v.size(); ++i) {
00383                         G4Track *track = v[i].track;
00384                         G4Event *event = v[i].event;
00385                         runManager->setCurrentEvent(event);
00386                         G4TrackStatus trackStatus = track->GetTrackStatus();
00387                         if(trackStatus != fAlive && 
00388                                                 trackStatus != fStopButAlive)
00389                                 continue;
00390                         if(primaryOnly != 0 && track->GetParentID() > 0)
00391                                 continue;
00392                         trace->appendTrace(track);
00393                 }
00394         }
00395 }

void BLCMDtrace::endCollectiveTracking ( std::vector< BLTrackData > &  v  )  [virtual]

Reimplemented from BLCollectiveComputation.

00398 {
00399         (void)v;
00400 }


Friends And Related Function Documentation

friend class TraceNTuple [friend]

Referenced by newTrace().


Member Data Documentation

int BLCMDtrace::nTrace [private]

G4String BLCMDtrace::format [private]

G4int BLCMDtrace::oneNTuple [private]

G4int BLCMDtrace::primaryOnly [private]

G4int BLCMDtrace::traceTune [private]

G4String BLCMDtrace::filename [private]

G4String BLCMDtrace::require [private]

G4String BLCMDtrace::coordinates [private]

bool BLCMDtrace::doTrace [private]

Referenced by BLCMDtrace(), command(), and newTrace().

bool BLCMDtrace::init [private]

Referenced by BLCMDtrace(), and command().

TraceNTuple * BLCMDtrace::allTrace = 0 [static, private]

Referenced by newTrace().

TraceNTuple * BLCMDtrace::tuneTrace = 0 [static, private]

Referenced by newTrace().

TraceNTuple * BLCMDtrace::referenceTrace = 0 [static, private]

Referenced by newTrace().


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