BLCMDparticlefilter Class Reference

Inheritance diagram for BLCMDparticlefilter:

BLElement BLCommand

List of all members.


Detailed Description

class BLCMDparticlefilter - implements a ParticleFilter Each placement of this element generates a physical volume that will force the decay of any particle named in the "decay" argument, and will kill any particle named in its "kill" argument, or not named in its "keep" argument if that is non-empty.

Note that the standard Decay process is disabled for all particles named in the "decay" agument, so they will only decay when they enter a ParticleFilter volume.

Note also the G4VProcess of the ParticleFilterPlacement is only applied to particles listed in the decay argument.

There is a bug: the decay happens after the first step in the volume; so maxStep is arranged to ensure two steps occur in the volume. This appears to be in the way Geant4 handles physical volumes while tracking....

Public Member Functions

 BLCMDparticlefilter ()
 Default constructor.
virtual ~BLCMDparticlefilter ()
 Destructor.
BLElementclone ()
 clone()
 BLCMDparticlefilter (const BLCMDparticlefilter &r)
 Copy constructor.
G4String commandName ()
 commandName() returns "particlefilter".
int command (BLArgumentVector &argv, BLArgumentMap &namedArgs)
 command() implements the particlefilter command.
void defineNamedArgs ()
 defineNamedArgs() defines the named arguments for the command.
virtual void construct (G4RotationMatrix *relativeRotation, G4ThreeVector relativePosition, G4LogicalVolume *parent, G4String parentName, G4RotationMatrix *parentRotation, G4ThreeVector parentPosition)
 construct() will construct the particlefilter.
virtual G4double getLength ()
 getLength() returns this element's Length along the Z axis.
virtual G4double getWidth ()
 getWidth() returns this element's Width along the X axis.
virtual G4double getHeight ()
 getHeight() returns this element's height along the Y axis.
virtual G4bool isOK ()
 isOK() returns true.
void generatePoints (int npoints, std::vector< G4ThreeVector > &v)
 generatePoints() from BLElement
G4bool isOutside (G4ThreeVector &local, G4double tolerance)
 isOutside() from BLElement

Private Member Functions

void setup ()

Private Attributes

G4double radius
G4double innerRadius
G4double height
G4double width
G4double length
G4String material
G4String color
G4VSolid * solid
G4double maxStep
G4String decay
G4String kill
G4String keep
G4int nWait
G4int referenceWait
G4String require
G4int steppingVerbose
std::set< G4String > decaySet
std::set< G4String > killSet
std::set< G4String > keepSet

Static Private Attributes

static int instance = 0

Friends

class ParticleFilterPlacement


Constructor & Destructor Documentation

BLCMDparticlefilter::BLCMDparticlefilter (  ) 

Default constructor.

References BLCMDTYPE_ELEMENT, color, decay, height, innerRadius, keep, kill, length, material, maxStep, nWait, radius, referenceWait, BLCommand::registerCommand(), require, BLCommand::setDescription(), BLCommand::setSynopsis(), solid, steppingVerbose, and width.

Referenced by clone(), and command().

00164                                          : BLElement(), decaySet(), killSet(),
00165                                                                 keepSet()
00166 {
00167         registerCommand(BLCMDTYPE_ELEMENT);
00168         setSynopsis("Will kill particles from a list, or force particles to decay.");
00169         setDescription("A particlefilter will force the decay of certain particles\n"
00170                 "when they enter the physical volume of the element.\n"
00171                 "The list of affected particles is in the 'decay' argument,\n"
00172                 "and the normal Decay process is disabled for them.\n"
00173                 "In addition, the 'kill' argument is a list of particles to\n"
00174                 "kill when they enter the element, and the 'keep' argument\n"
00175                 "will kill all particles not named, if it is not empty.\n\n"
00176                 "require is an expression invloving track parameters that will "
00177                 "kill the track if it evaluates to zero (use a comparison "
00178                 "operator). The variables available are:\n"
00179                 "    x,y,z,t,Px,Py,Pz,Ptot,PDGid,EventID,TrackID,ParentID\n"
00180                 "Units are mm, ns, MeV/c.\n\n"
00181                 "If nWait > 1, particles will not be killed until they hit "
00182                 "this element nwait times; this can be used to limit the "
00183                 "number of revolutions around a ring. Decays are unaffected "
00184                 "by nWait. referenceWait does the same for the reference "
00185                 "particle\n"
00186                 "The element can be either a cylinder or a box: set length and "
00187                 "radius, or set length and width and height.");
00188         
00189         // initial field values
00190         radius = -1.0;
00191         innerRadius = 0.0;
00192         height = 0.0;
00193         width = 0.0;
00194         length = 1.0*mm;
00195         material = "";
00196         color = "1,1,1";
00197         solid = 0;
00198         maxStep = -99.0;
00199         decay = "";
00200         kill = "";
00201         keep = "";
00202         nWait = 1;
00203         referenceWait = 1;
00204         require = "";
00205         steppingVerbose = 0;
00206 }

BLCMDparticlefilter::~BLCMDparticlefilter (  )  [virtual]

Destructor.

00209 {
00210 }

BLCMDparticlefilter::BLCMDparticlefilter ( const BLCMDparticlefilter r  ) 

Copy constructor.

References color, decay, height, innerRadius, keep, kill, length, material, maxStep, nWait, radius, referenceWait, require, solid, steppingVerbose, and width.

00212                                                                      : BLElement(r),
00213                 decaySet(r.decaySet), killSet(r.killSet), keepSet(r.keepSet)
00214 {
00215         radius = r.radius;
00216         innerRadius = r.innerRadius;
00217         height = r.height;
00218         width = r.width;
00219         length = r.length;
00220         material = r.material;
00221         color = r.color;
00222         solid = 0;
00223         maxStep = r.maxStep;
00224         decay = r.decay;
00225         kill = r.kill;
00226         keep = r.keep;
00227         nWait = r.nWait;
00228         referenceWait = r.referenceWait;
00229         require = r.require;
00230         steppingVerbose = r.steppingVerbose;
00231 }


Member Function Documentation

void BLCMDparticlefilter::setup (  )  [private]

References decay, decaySet, keep, keepSet, kill, killSet, and BLCommand::splitString().

Referenced by command().

00262 {
00263         // get set of particle names into decaySet
00264         std::vector<G4String> v=splitString(decay,',',true);
00265         for(unsigned i=0; i<v.size(); ++i) {
00266                 if(v[i].size() == 0) continue;
00267                 decaySet.insert(v[i]);
00268         }
00269 /***
00270         unsigned int place, next;
00271         for(place=next=0; next<decay.size(); place=next+1) {
00272                 next = decay.find(",",place);
00273                 G4String p;
00274                 if(next < decay.size())
00275                         p = decay.substr(place,next-place);
00276                 else
00277                         p = decay.substr(place);
00278                 if(p.size() == 0) break;
00279                 decaySet.insert(p);
00280         }
00281 ***/
00282         // get set of particle names into killSet
00283         v=splitString(kill,',',true);
00284         for(unsigned i=0; i<v.size(); ++i) {
00285                 if(v[i].size() == 0) continue;
00286                 killSet.insert(v[i]);
00287         }
00288 /***
00289         for(place=next=0; next<kill.size(); place=next+1) {
00290                 next = kill.find(",",place);
00291                 G4String p;
00292                 if(next < kill.size())
00293                         p = kill.substr(place,next-place);
00294                 else
00295                         p = kill.substr(place);
00296                 if(p.size() == 0) break;
00297                 killSet.insert(p);
00298         }
00299 ***/
00300         // get set of particle names into keepSet
00301         v=splitString(keep,',',true);
00302         for(unsigned i=0; i<v.size(); ++i) {
00303                 if(v[i].size() == 0) continue;
00304                 keepSet.insert(v[i]);
00305         }
00306 /***
00307         for(place=next=0; next<keep.size(); place=next+1) {
00308                 next = keep.find(",",place);
00309                 G4String p;
00310                 if(next < keep.size())
00311                         p = keep.substr(place,next-place);
00312                 else
00313                         p = keep.substr(place);
00314                 if(p.size() == 0) break;
00315                 keepSet.insert(p);
00316         }
00317 ***/
00318         // if keep is nonempty, add decaySet to it (they will decay)
00319         if(keepSet.size() > 0) {
00320                 std::set<G4String>::iterator i;
00321                 for(i=decaySet.begin(); i!=decaySet.end(); ++i)
00322                         keepSet.insert(*i);
00323         }
00324 }

BLElement* BLCMDparticlefilter::clone (  )  [inline, virtual]

clone()

Implements BLElement.

References BLCMDparticlefilter().

00085 { return new BLCMDparticlefilter(*this); }

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

commandName() returns "particlefilter".

Implements BLCommand.

00091 { return "particlefilter"; }

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

command() implements the particlefilter command.

Implements BLCommand.

References BLCMDparticlefilter(), BLParam::getDouble(), BLParam::getInt(), BLCommand::getMaterial(), BLCommand::handleNamedArgs(), material, maxStep, Param, BLCommand::print(), BLCommand::printError(), BLElement::setName(), setup(), and steppingVerbose.

00234 {
00235         if(argv.size() != 1) {
00236                 printError("particlefilter: Invalid command, must have name");
00237                 return -1;
00238         }
00239 
00240         if(argv[0] == "default") {
00241                 return defaultParticleFilter.handleNamedArgs(namedArgs);
00242         }
00243 
00244         BLCMDparticlefilter *t = new BLCMDparticlefilter(defaultParticleFilter);
00245         t->setName(argv[0]);
00246         t->steppingVerbose = Param.getInt("steppingVerbose");
00247         int retval = t->handleNamedArgs(namedArgs);
00248 
00249         if(t->maxStep < 0.0) t->maxStep = Param.getDouble("maxStep");
00250 
00251         t->setup();
00252 
00253         // check material exists
00254         if(t->material.size() > 0) getMaterial(t->material,false);
00255 
00256         t->print(argv[0]);
00257 
00258         return retval;
00259 }

void BLCMDparticlefilter::defineNamedArgs (  )  [virtual]

defineNamedArgs() defines the named arguments for the command.

Reimplemented from BLCommand.

References BLCommand::argDouble(), BLCommand::argInt(), BLCommand::argString(), color, decay, height, innerRadius, keep, kill, length, material, maxStep, nWait, radius, referenceWait, require, steppingVerbose, and width.

00327 {
00328         argDouble(radius,"radius","The radius of the cylindrical particlefilter (mm).");
00329         argDouble(innerRadius,"innerRadius","The inner radius of the cylindrical particlefilter (0 mm, solid).");
00330         argDouble(height,"height","The height of the rectangular particlefilter (mm).");
00331         argDouble(width,"width","The width of the rectangular particlefilter (mm).");
00332         argDouble(length,"length","The length of the particlefilter (mm).");
00333         argDouble(maxStep,"maxStep","The maximum stepsize in the element (mm).");
00334         argString(material,"material","The material (default=parent's).");
00335         argString(color,"color","The color of the particlefilter (white).");
00336         argString(decay,"decay","A comma-separated list of particle names to decay.",false);
00337         argString(kill,"kill","A comma-separated list of particle names to kill.",false);
00338         argString(keep,"keep","A comma-separated list of particle names to keep.",false);
00339         argInt(nWait,"nWait","Intersection # to do the kill (default = 1)");
00340         argInt(referenceWait,"referenceWait","Intersection # for reference (default = 1)");
00341         argString(require,"require","Expression which will kill the track if zero.",false);
00342         argInt(steppingVerbose,"steppingVerbose","Nonzero to display track kills.");
00343         argString(decay,"decays","Synonym for decay.",false);
00344 
00345         if(radius > 0.0)
00346                 width = height = 2.0 * radius;
00347 }

void BLCMDparticlefilter::construct ( G4RotationMatrix *  relativeRotation,
G4ThreeVector  relativePosition,
G4LogicalVolume *  parent,
G4String  parentName,
G4RotationMatrix *  parentRotation,
G4ThreeVector  parentPosition 
) [virtual]

construct() will construct the particlefilter.

Implements BLElement.

References color, BLParam::getDouble(), BLCommand::getMaterial(), BLElement::getName(), BLManager::getObject(), BLCommand::getVisAttrib(), height, innerRadius, instance, keepSet, killSet, length, material, maxStep, minStep, Param, ParticleFilterPlacement, BLCommand::printError(), radius, BLManager::registerSteppingAction(), BLManager::registerTrackingAction(), require, solid, and width.

00355 {
00356         if(maxStep < 0.0) maxStep = Param.getDouble("maxStep");
00357 
00358         // ensure at least two steps will be taken inside the physical volume
00359         G4double minStep = (width<height ? width : height);
00360         minStep = (minStep<length ? minStep : length);
00361         minStep = minStep/2.0;
00362         maxStep = (minStep<maxStep ? minStep : maxStep);
00363 
00364         G4Material *mat;
00365         if(material != "")
00366                 mat = getMaterial(material);
00367         else
00368                 mat = parent->GetMaterial();
00369 
00370         G4String thisname = parentName+getName();
00371 
00372         if(!solid) {
00373                 if(radius > 0.0) {
00374                         solid = new G4Tubs(thisname+"Tubs", innerRadius, radius,
00375                                         length/2.0, 0.0, 2.0*pi);
00376                 } else if(height > 0.0 && width > 0.0) {
00377                         solid = new G4Box(thisname+"Box",width/2.0,
00378                                         height/2.0,length/2.0);
00379                 } else {
00380                         printError("particlefilter::construct %s INVALID - no "
00381                                 "radius or height&width",thisname.c_str());
00382                         return;
00383                 }
00384         }
00385         G4LogicalVolume *lv = new G4LogicalVolume(solid,mat, thisname+"LogVol");
00386         lv->SetVisAttributes(getVisAttrib(color));
00387         if(maxStep < 0.0) maxStep = Param.getDouble("maxStep");
00388         lv->SetUserLimits(new G4UserLimits(maxStep));
00389 
00390         // geant4 rotation convention is backwards from g4beamline
00391         G4RotationMatrix *g4rot = 0;
00392         if(relativeRotation)
00393                 g4rot = new G4RotationMatrix(relativeRotation->inverse());
00394 
00395         G4VPhysicalVolume *vol = new G4PVPlacement(g4rot, relativePosition,
00396                                         lv,thisname,parent,false,0);
00397 
00398         char processname[64];
00399         sprintf(processname,"BLCMDparticlefilter%d",++instance);
00400         ParticleFilterPlacement *place = 
00401                         new ParticleFilterPlacement(this,vol,processname);
00402 
00403         if(killSet.size() > 0 || keepSet.size() > 0 || require != "") {
00404             BLManager::getObject()->registerTrackingAction(place);
00405             BLManager::getObject()->registerSteppingAction(vol,place);
00406         }
00407 
00408         printf("BLCMDparticlefilter::Construct %s parent=%s\n",
00409                 thisname.c_str(),parentName.c_str());
00410 }

virtual G4double BLCMDparticlefilter::getLength (  )  [inline, virtual]

getLength() returns this element's Length along the Z axis.

Implements BLElement.

References length.

00108 { return length; }

virtual G4double BLCMDparticlefilter::getWidth (  )  [inline, virtual]

getWidth() returns this element's Width along the X axis.

Implements BLElement.

References width.

00111 { return width; }

virtual G4double BLCMDparticlefilter::getHeight (  )  [inline, virtual]

getHeight() returns this element's height along the Y axis.

Implements BLElement.

References height.

00114 { return height; }

virtual G4bool BLCMDparticlefilter::isOK (  )  [inline, virtual]

isOK() returns true.

Implements BLElement.

00117 { return true; }

void BLCMDparticlefilter::generatePoints ( int  npoints,
std::vector< G4ThreeVector > &  v 
) [virtual]

generatePoints() from BLElement

Implements BLElement.

References BLElement::generateBox(), BLElement::generateTubs(), height, innerRadius, length, radius, and width.

00413 {
00414         if(radius > 0.0)
00415                 generateTubs(npoints, innerRadius, radius, 0.0, 360.0*deg, length, v);
00416         else
00417                 generateBox(npoints,width,height,length,v);
00418 }

G4bool BLCMDparticlefilter::isOutside ( G4ThreeVector &  local,
G4double  tolerance 
) [virtual]

isOutside() from BLElement

Implements BLElement.

References height, innerRadius, length, radius, and width.

00421 {
00422         if(radius > 0.0) {
00423                 G4double r = sqrt(local[0]*local[0]+local[1]*local[1]);
00424                 return r > radius-tolerance || r < innerRadius+tolerance ||
00425                         fabs(local[2]) > length/2.0-tolerance;
00426         } else {
00427                 return fabs(local[0]) > width/2.0-tolerance ||
00428                         fabs(local[1]) > height/2.0-tolerance ||
00429                         fabs(local[2]) > length/2.0-tolerance;
00430         }
00431 }


Friends And Related Function Documentation

friend class ParticleFilterPlacement [friend]

Referenced by construct().


Member Data Documentation

G4double BLCMDparticlefilter::radius [private]

G4double BLCMDparticlefilter::height [private]

G4double BLCMDparticlefilter::width [private]

G4double BLCMDparticlefilter::length [private]

G4String BLCMDparticlefilter::material [private]

G4String BLCMDparticlefilter::color [private]

G4VSolid* BLCMDparticlefilter::solid [private]

Referenced by BLCMDparticlefilter(), and construct().

G4double BLCMDparticlefilter::maxStep [private]

G4String BLCMDparticlefilter::decay [private]

G4String BLCMDparticlefilter::kill [private]

G4String BLCMDparticlefilter::keep [private]

G4int BLCMDparticlefilter::nWait [private]

G4String BLCMDparticlefilter::require [private]

std::set<G4String> BLCMDparticlefilter::decaySet [private]

std::set<G4String> BLCMDparticlefilter::killSet [private]

std::set<G4String> BLCMDparticlefilter::keepSet [private]

int BLCMDparticlefilter::instance = 0 [static, private]

Referenced by construct().


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