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().

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

BLCMDparticlecolor::~BLCMDparticlecolor (  )  [inline]

Destructor.

00052 { }


Member Function Documentation

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

commandName() returns "particlecolor"

Implements BLCommand.

00055 { 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().

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

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

defineNamedArgs() defines the named arguments for the command.

Reimplemented from BLCommand.

00061 { }

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().

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

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.

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

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.

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

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

BeginOfEventAction from BLManager::EventAction.

Implements BLManager::EventAction.

00289 {
00290 }

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

EndOfEventAction from BLManager::EventAction.

Implements BLManager::EventAction.

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

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


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