BLCMDparticlecolor Class Reference

Inheritance diagram for BLCMDparticlecolor:

BLCommand BLManager::TrackingAction BLManager::EventAction

List of all members.


Detailed Description

class BLCMDparticlecolor -- implement the particlecolor command

Public Member Functions

 BLCMDparticlecolor ()
 Constructor.
 ~BLCMDparticlecolor ()
 Destructor.
G4String commandName ()
 commandName() returns "particlecolor"
int command (BLArgumentVector &argv, BLArgumentMap &namedArgs)
 command() executes the command associated with this element.
void defineNamedArgs ()
 defineNamedArgs() defines the named arguments for the command.
const G4VisAttributes * getVisAttributes (const G4ParticleDefinition *pd) const
 getVisAttributes() will return the G4VisAttributes for a particle.
void PreUserTrackingAction (const G4Track *track)
 PreUserTrackingAction() from BLManager::TrackingAction.
void PostUserTrackingAction (const G4Track *track)
 PostUserTrackingAction from BLManager::TrackingAction.
void BeginOfEventAction (const G4Event *event)
 BeginOfEventAction from BLManager::EventAction.
void EndOfEventAction (const G4Event *event)
 EndOfEventAction from BLManager::EventAction.

Private Attributes

const G4VisAttributes * plus
const G4VisAttributes * minus
const G4VisAttributes * neutral
const G4VisAttributes * reference

Static Private Attributes

static std::map< const
G4ParticleDefinition *, const
G4VisAttributes * > 
pd2va

Constructor & Destructor Documentation

BLCMDparticlecolor::BLCMDparticlecolor (  ) 

Constructor.

References BLCMDTYPE_AUX, BLCommand::getVisAttrib(), minus, neutral, plus, reference, BLCommand::registerCommand(), BLCommand::setDescription(), and BLCommand::setSynopsis().

00183                                        : BLCommand(), BLManager::TrackingAction(), BLManager::EventAction()
00184 {
00185         registerCommand(BLCMDTYPE_AUX);
00186         setSynopsis("Set the colors for particle tracks.");
00187         setDescription("Arguments are of the form 'name=1,1,0', where name is\n"
00188                         "the standard name of a particle, and 1,1,0 is the\n"
00189                         "R,G,B value desired for its color ('' for invisible)\n"
00190                         "The special names plus, minus, and neutral will set\n"
00191                         "colors for unnamed particles of each charge.\n"
00192                         "The name reference will apply to the reference\n"
00193                         "track (defaults to invisible).");
00194         
00195         G4String red("1,0,0"), green("0,1,0"), blue("0,0,1");
00196         plus = getVisAttrib(blue);
00197         minus = getVisAttrib(red);
00198         neutral = getVisAttrib(green);
00199         reference = getVisAttrib("Invisible");
00200 }

BLCMDparticlecolor::~BLCMDparticlecolor (  )  [inline]

Destructor.

00053 { }


Member Function Documentation

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

commandName() returns "particlecolor"

Implements BLCommand.

00056 { return "particlecolor"; }

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

command() executes the command associated with this element.

Implements BLCommand.

References BLManager::getObject(), BLCommand::getVisAttrib(), minus, neutral, pd2va, plus, BLCommand::print(), BLCommand::printError(), reference, BLManager::registerEventAction(), and BLManager::registerTrackingAction().

00203 {
00204         int retval = 0;
00205 
00206         if(argv.size() > 0) {
00207                 printError("particlecolor: unnamed arguments are invalid.");
00208                 retval = 1;
00209         }
00210 
00211         BLArgumentMap::iterator i;
00212         for(i=namedArgs.begin(); i!=namedArgs.end(); ++i) {
00213                 G4String name = i->first;
00214                 G4String value = i->second;
00215                 const G4VisAttributes *va = getVisAttrib(value);
00216                 if(name == "plus") {
00217                         plus = va;
00218                         continue;
00219                 } else if(name == "minus") {
00220                         minus = va;
00221                         continue;
00222                 } else if(name == "neutral") {
00223                         neutral = va;
00224                         continue;
00225                 } else if(name == "reference") {
00226                         reference = va;
00227                         continue;
00228                 }
00229                 G4ParticleDefinition *pd = 
00230                         G4ParticleTable::GetParticleTable()->FindParticle(name);
00231                 if(!pd) {
00232                         printError("particlecolor: particle '%s' not found",
00233                                                                 name.c_str());
00234                         retval = 1;
00235                         continue;
00236                 }
00237                 pd2va[pd] = va;
00238         }
00239 
00240         BLManager::getObject()->registerTrackingAction(this);
00241         BLManager::getObject()->registerEventAction(this,false);
00242 
00243         print("",namedArgs);
00244 
00245         return retval;
00246 }

void BLCMDparticlecolor::defineNamedArgs (  )  [inline, virtual]

defineNamedArgs() defines the named arguments for the command.

Reimplemented from BLCommand.

00062 { }

const G4VisAttributes * BLCMDparticlecolor::getVisAttributes ( const G4ParticleDefinition *  pd  )  const

getVisAttributes() will return the G4VisAttributes for a particle.

References BLManager::getObject(), BLManager::getState(), minus, neutral, pd2va, plus, reference, REFERENCE, and VISUAL.

Referenced by BLTrajectory::BLTrajectory().

00249 {
00250         BLManagerState state = BLManager::getObject()->getState();
00251         if(state == REFERENCE)
00252                 return reference;
00253         if(state != VISUAL)
00254                 return 0;
00255 
00256         if(pd2va.count(pd) > 0)
00257                 return pd2va[pd];
00258         G4double charge = pd->GetPDGCharge();
00259         if(charge > 0.0)
00260                 return plus;
00261         else if(charge < 0.0)
00262                 return minus;
00263         return neutral;
00264 }

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

PreUserTrackingAction() from BLManager::TrackingAction.

Implements BLManager::TrackingAction.

References BLManager::TrackingAction::fpTrackingManager, BLManager::getObject(), BLManager::getState(), BLParam::getString(), Param, and REFERENCE.

00267 {
00268         if(Param.getString("viewer") == "none")
00269                 return;
00270         if(fpTrackingManager)
00271                 fpTrackingManager->SetTrajectory(new BLTrajectory(track));
00272         BLManagerState state = BLManager::getObject()->getState();
00273         if(state == REFERENCE)
00274                 fpTrackingManager->SetStoreTrajectory(true);
00275 }

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

PostUserTrackingAction from BLManager::TrackingAction.

Implements BLManager::TrackingAction.

References CenterTrajectories, BLManager::TrackingAction::fpTrackingManager, BLManager::getObject(), BLManager::getState(), and REFERENCE.

00278 {
00279         if(!fpTrackingManager) return;
00280         // get Center particle trajectories
00281         BLManagerState state = BLManager::getObject()->getState();
00282         if(state == REFERENCE) {
00283                 G4VTrajectory *t = fpTrackingManager->GimmeTrajectory();
00284                 // this is a deep copy, so let the tracking manager delete t
00285                 if(t) CenterTrajectories.push_back(new BLSaveTrajectory(t));
00286         }
00287 }

void BLCMDparticlecolor::BeginOfEventAction ( const G4Event *  event  )  [virtual]

BeginOfEventAction from BLManager::EventAction.

Implements BLManager::EventAction.

00290 {
00291 }

void BLCMDparticlecolor::EndOfEventAction ( const G4Event *  event  )  [virtual]

EndOfEventAction from BLManager::EventAction.

Implements BLManager::EventAction.

References CenterTrajectories, BLManager::getObject(), BLManager::getState(), and VISUAL.

00294 {
00295         BLManagerState state = BLManager::getObject()->getState();
00296         if(state == VISUAL && CenterTrajectories.size() > 0) {
00297                 // add CenterTrajectories to this event
00298                 G4TrajectoryContainer *tc = event->GetTrajectoryContainer();
00299                 for(unsigned i=0; i<CenterTrajectories.size(); ++i) {
00300                         if(tc)
00301                                 tc->push_back(CenterTrajectories[i]);
00302                 }
00303                 // clear CenterTrajectories so only one copy is drawn, and
00304                 // because the tracking manager is going to delete them.
00305                 CenterTrajectories.clear();
00306         }
00307 }


Member Data Documentation

std::map< const G4ParticleDefinition *, const G4VisAttributes * > BLCMDparticlecolor::pd2va [static, private]

Referenced by command(), and getVisAttributes().

const G4VisAttributes* BLCMDparticlecolor::plus [private]

const G4VisAttributes* BLCMDparticlecolor::minus [private]

const G4VisAttributes* BLCMDparticlecolor::neutral [private]

const G4VisAttributes* BLCMDparticlecolor::reference [private]


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