BLCMDtrackcuts Class Reference

Inheritance diagram for BLCMDtrackcuts:

BLCommand BLManager::SteppingAction BLManager::StackingAction

List of all members.


Detailed Description

class BLCMDtrackcuts implements the trackcuts command


Public Member Functions

 BLCMDtrackcuts ()
 Constructor.
 BLCMDtrackcuts (const BLCMDtrackcuts &r)
G4String commandName ()
int command (BLArgumentVector &argv, BLArgumentMap &namedArgs)
void defineNamedArgs ()
void UserSteppingAction (const G4Step *step)
 UserSteppingAction() from BLManager::SteppingAction.
G4ClassificationOfNewTrack ClassifyNewTrack (const G4Track *track)
 ClassifyNewTrack() from BLManager::StackingAction.

Private Attributes

G4double kineticEnergyCut
G4double kineticEnergyMax
G4int steppingVerbose
G4int killSecondaries
G4String kill
G4String keep
G4double maxTime
G4int keepPrimaries
std::set< G4int > killPDG
std::set< G4int > keepPDG
bool stepping


Constructor & Destructor Documentation

BLCMDtrackcuts::BLCMDtrackcuts (  ) 

Constructor.

References BLCMDTYPE_CUTS, keepPrimaries, killSecondaries, kineticEnergyCut, kineticEnergyMax, maxTime, BLCommand::registerCommand(), BLCommand::setDescription(), BLCommand::setSynopsis(), stepping, and steppingVerbose.

Referenced by command().

00067                                : BLCommand(), BLManager::SteppingAction(), 
00068 BLManager::StackingAction(), kill(), keep()
00069 {
00070         registerCommand(BLCMDTYPE_CUTS);
00071         setSynopsis("Specifies per-track cuts.");
00072         setDescription("Applied to each track before tracking, and at each step.");
00073 
00074         kineticEnergyCut = 0.0;
00075         kineticEnergyMax = DBL_MAX;
00076         steppingVerbose = 0;
00077         killSecondaries = 0;
00078         maxTime = 1.0*millisecond;
00079         keepPrimaries = 0;
00080         stepping = false;
00081 }

BLCMDtrackcuts::BLCMDtrackcuts ( const BLCMDtrackcuts r  ) 


Member Function Documentation

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

Implements BLCommand.

00052 { return "trackcuts"; }

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

Implements BLCommand.

References BLCMDtrackcuts(), BLParam::getInt(), BLManager::getObject(), BLCommand::handleNamedArgs(), keep, keepPDG, kill, killPDG, Param, BLCommand::print(), BLCommand::printError(), BLManager::registerStackingAction(), BLManager::registerSteppingAction(), BLCommand::splitString(), steppingVerbose, and BLCommand::type.

00095 {
00096         G4ParticleTable *pt = G4ParticleTable::GetParticleTable();
00097 
00098         if(namedArgs.size() == 0) {
00099                 int n = pt->entries();
00100                 for(int i=0; i<n; ++i) {
00101                         G4ParticleDefinition *d = pt->GetParticle(i);
00102                         G4String name = d->GetParticleName();
00103                         G4String type = d->GetParticleType();
00104                         G4String subtype = d->GetParticleSubType();
00105                         G4int pdg = d->GetPDGEncoding();
00106                         printf("%15.15s  %s/%s  PDGEncoding=%d\n",
00107                                 name.c_str(),type.c_str(),subtype.c_str(),pdg);
00108                 }
00109                 return 0;
00110         }
00111 
00112         BLCMDtrackcuts *t = new BLCMDtrackcuts(defaultTrackCuts);
00113 
00114         t->steppingVerbose = Param.getInt("steppingVerbose");
00115 
00116         int retval = t->handleNamedArgs(namedArgs);
00117 
00118         // get PDGid-s for particles to kill into killPDG
00119         std::vector<G4String> v=splitString(t->kill,",",true);
00120         for(unsigned i=0; i<v.size(); ++i) {
00121                 if(v[i].size() == 0) continue;
00122                 G4ParticleDefinition *d = pt->FindParticle(v[i]);
00123                 if(!d) {
00124                         printError("Cannot find particle '%s' (physics must come first)",v[i].c_str());
00125                         continue;
00126                 }
00127                 t->killPDG.insert(d->GetPDGEncoding());
00128         }
00129 
00130         // get PDGid-s for particles to keep into keepPDG
00131         v=splitString(t->keep,",",true);
00132         for(unsigned i=0; i<v.size(); ++i) {
00133                 if(v[i].size() == 0) continue;
00134                 G4ParticleDefinition *d = pt->FindParticle(v[i]);
00135                 if(!d) {
00136                         printError("Cannot find particle '%s' (physics must come first)",v[i].c_str());
00137                         continue;
00138                 }
00139                 t->keepPDG.insert(d->GetPDGEncoding());
00140         }
00141 
00142         BLManager::getObject()->registerSteppingAction(t);
00143         BLManager::getObject()->registerStackingAction(t);
00144 
00145         t->print("");
00146 
00147         return retval;
00148 }

void BLCMDtrackcuts::defineNamedArgs (  )  [virtual]

Reimplemented from BLCommand.

References BLCommand::argDouble(), BLCommand::argInt(), BLCommand::argString(), keep, keepPrimaries, kill, killSecondaries, kineticEnergyCut, kineticEnergyMax, maxTime, and steppingVerbose.

00151 {
00152         argString(kill,"kill","List of particles to kill (comma separated).");
00153         argString(keep,"keep","List of particles to keep (kill all others).");
00154         argInt(killSecondaries,"killSecondaries","Set nonzero to kill all secondaries.");
00155         argDouble(kineticEnergyCut,"kineticEnergyCut","Minimum K.E. to track (0 MeV).",MeV);
00156         argDouble(kineticEnergyMax,"kineticEnergyMax","Maximum K.E. to track (infinite MeV).",MeV);
00157         argDouble(maxTime,"maxTime","Maximum lab time to track (1000000 ns).");
00158         argInt(keepPrimaries,"keepPrimaries","Set nonzero to keep tracks with ParentID==0 regardless of other tests.");
00159         argInt(steppingVerbose,"steppingVerbose","Set nonzero to print kills\n"
00160                 "(defaults to parameter value).");
00161 }

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

UserSteppingAction() from BLManager::SteppingAction.

Implements BLManager::SteppingAction.

References ClassifyNewTrack(), and stepping.

00165 {
00166         stepping = true;
00167 
00168         G4Track *track = step->GetTrack();
00169         if(ClassifyNewTrack(track) == fKill)
00170                 ((G4Track*)track)->SetTrackStatus(fStopAndKill);
00171 
00172         stepping = false;
00173 }

G4ClassificationOfNewTrack BLCMDtrackcuts::ClassifyNewTrack ( const G4Track *  track  )  [virtual]

ClassifyNewTrack() from BLManager::StackingAction.

Implements BLManager::StackingAction.

References BLManager::getObject(), keepPDG, keepPrimaries, killPDG, killSecondaries, kineticEnergyCut, kineticEnergyMax, maxTime, stepping, steppingVerbose, and TUNE.

Referenced by UserSteppingAction().

00177 {
00178         // handle keepPrimaries independent of other tests
00179         if(keepPrimaries != 0 && track->GetParentID() == 0)
00180                 return fUrgent;
00181 
00182         // Kill particle if Kinetic Energy is below the threshold
00183         G4double kineticEnergy = track->GetKineticEnergy();
00184         if(kineticEnergy < kineticEnergyCut)  {
00185                 if(steppingVerbose != 0)
00186                         printf("Track killed: K.E. %.6f < %.6f MeV\n",
00187                                 kineticEnergy/MeV, kineticEnergyCut/MeV);
00188                 return fKill;
00189         }
00190 
00191         // Kill particle if Kinetic Energy is above the Maximum
00192         if(kineticEnergy > kineticEnergyMax)  {
00193                 if(steppingVerbose != 0)
00194                         printf("Track killed: K.E. %.3f > %.3f MeV\n",
00195                                 kineticEnergy/MeV, kineticEnergyMax/MeV);
00196                 return fKill;
00197         }
00198 
00199         // Kill particle if time is too large
00200         G4double time = track->GetGlobalTime();
00201         if(time > maxTime)  {
00202                 if(steppingVerbose != 0)
00203                         printf("Track killed: time %.3f > %.3f ns\n",
00204                                 time/ns, maxTime/ns);
00205                 return fKill;
00206         }
00207 
00208         // shortcut while stepping, avoid unnecessary tests
00209         if(stepping) return fUrgent;
00210 
00211         // kill secondaries, unless in tune particle mode (the tune particle
00212         // is frequently re-tracked as a "secondary")
00213         if(BLManager::getObject()->getState() != TUNE &&
00214                         killSecondaries != 0 && track->GetParentID() != 0) {
00215                 if(steppingVerbose != 0) printf("Secondary Track killed\n");
00216                 return fKill;
00217         }
00218 
00219         // kill particles in kill
00220         G4ParticleDefinition *pdg = track->GetDefinition();
00221         if(killPDG.count(pdg->GetPDGEncoding()) > 0) {
00222                 if(steppingVerbose != 0) 
00223                         printf("Track killed because it's a %s\n",
00224                                 pdg->GetParticleName().c_str());
00225                 return fKill;
00226         }
00227 
00228         // kill particles not in keep
00229         if(keepPDG.size() > 0 && keepPDG.count(pdg->GetPDGEncoding()) == 0) {
00230                 if(steppingVerbose != 0) 
00231                         printf("Track killed because it's a %s\n",
00232                                 pdg->GetParticleName().c_str());
00233                 return fKill;
00234         }
00235 
00236         return fUrgent;
00237 }


Member Data Documentation

G4String BLCMDtrackcuts::kill [private]

Referenced by command(), and defineNamedArgs().

G4String BLCMDtrackcuts::keep [private]

Referenced by command(), and defineNamedArgs().

G4double BLCMDtrackcuts::maxTime [private]

std::set<G4int> BLCMDtrackcuts::killPDG [private]

Referenced by ClassifyNewTrack(), and command().

std::set<G4int> BLCMDtrackcuts::keepPDG [private]

Referenced by ClassifyNewTrack(), and command().

bool BLCMDtrackcuts::stepping [private]


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