BLCMDtrackerplane Class Reference

Inheritance diagram for BLCMDtrackerplane:

BLElement BLCommand

List of all members.


Detailed Description

class BLCMDtrackerplane - implements one plane of a tracker.

A tracker plane is a detector consisting of either: A) a large number of parallel wires arranged in a plane B) a counter that gives timing information C) both (A) and (B) For a given track hitting the plane, one and only one wire is hit. Arguments to the command specify the wire spacing, orientation, and offset; the number of wires is implicit from their spacing and orientation and the size of the plane. if wireSpacing=0, no wires are present.

If sigmaTime>0 this plane measures the time of the track. A Gaussian random number with that sigma is added to the true value.

This class works with the BLCMDtracker class to permit the fitting of tracks to the wire hits. modes: true reports the hit wire # to the BLCMDtracker instance (this is G4beamline tracking "true" tracks) fit computes the plane's contribution to Chisq for the track and reports it to the BLCMDtracker instance ignore track is ignored

Public Member Functions

 BLCMDtrackerplane ()
 Default constructor.
virtual ~BLCMDtrackerplane ()
 Destructor.
BLElementclone ()
 clone()
 BLCMDtrackerplane (const BLCMDtrackerplane &r)
 Copy constructor.
G4String commandName ()
 commandName() returns "trackerplane".
int command (BLArgumentVector &argv, BLArgumentMap &namedArgs)
 command() implements the trackerplane command.
void defineNamedArgs ()
 defineNamedArgs() defines the named arguments for the command.
void argChanged ()
 argChanged() does internal computations after some arg changed
virtual void construct (G4RotationMatrix *relativeRotation, G4ThreeVector relativePosition, G4LogicalVolume *parent, G4String parentName, G4RotationMatrix *parentRotation, G4ThreeVector parentPosition)
 construct() will construct the trackerplane. Used for normal placements of a tracker plane object.
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 Attributes

G4String tracker
G4double radius
G4double innerRadius
G4double height
G4double width
G4double length
G4String material
G4String color
G4VSolid * solid
G4double maxStep
G4double theta
G4double wireSpacing
G4double wireOffset
G4String errType
G4double errTheta
G4double errSpacing
G4double errOffset
G4double sigmaTime

Friends

class TrackerPlaneInstance


Constructor & Destructor Documentation

BLCMDtrackerplane::BLCMDtrackerplane (  ) 

Default constructor.

References BLCMDTYPE_DATA, color, errOffset, errSpacing, errTheta, errType, height, innerRadius, length, material, maxStep, radius, BLCommand::registerCommand(), BLCommand::setDescription(), BLCommand::setSynopsis(), sigmaTime, solid, theta, tracker, width, wireOffset, and wireSpacing.

Referenced by clone(), and command().

01263                                      : BLElement()
01264 {
01265         registerCommand(BLCMDTYPE_DATA);
01266         setSynopsis("Construct a tracker plane.");
01267         setDescription("A trackerplane belongs to a specific tracker, and "
01268                 "represents one measuring element of the tracker. "
01269                 "While the term 'wire' is used, a trackerplane can model "
01270                 "any planar measuring device that measures one dimension "
01271                 "using equally spaced detectors that are either on or off "
01272                 "(hit or not hit). A trackerplane can be "
01273                 "circular (specify radius and possibly innerRadius) or "
01274                 "rectangular (specify height and width).\n\n"
01275                 "The wires are at angle theta from the vertical, so theta=0 "
01276                 "means vertical wires that measure x; theta=90 means "
01277                 "horizontal wires that measure y; theta=180 also "
01278                 "measures x, but increasing x means decreasing wire #.\n\n"
01279                 "Setting wireSpacing=0 means there are no wires, which is "
01280                 "useful for a trackerplane that models a scintillator used "
01281                 "for track timing (set sigmaTime >= 0 to indicate that).\n\n"
01282                 "sigmaTime>0 means this plane can measure the time of the "
01283                 "track with that resolution. A Gaussian random number is added "
01284                 "to the true track's time at this plane when reading the "
01285                 "TrackHit NTuple.\n\n"
01286                 "Survey errors can be modeled using the err-arguments -- "
01287                 "the values they specify are applied during trackfitting "
01288                 "(but not for true tracks). "
01289                 "errType: 'fixed' means error values given are the actual values, "
01290                 "'gaussian' means error values are the sigma of a Gaussian random "
01291                 "number, 'rect' means error values are the half-width of a uniform "
01292                 "random number. The random number is picked before the run "
01293                 "begins; the random-number seed is set from the clock so "
01294                 "every run will have different random errors.\n\n"
01295                 "trackerplane has 3 modes:\n"
01296                 "  true    The trackerplane reports the hit wire and time to\n"
01297                 "          the tracker;\n"
01298                 "  fit     The trackerplane reports the Chisq contribution of\n"
01299                 "          the fit track distance to the true track's hit wire\n"
01300                 "          center, plus survey errors (if any). The track time\n"
01301                 "          also contributes to the chisq.\n"
01302                 "  ignore  Any track is ignored.\n\n"
01303                 "NOTE: the trackerplane-s of a tracker MUST be placed in the "
01304                 "order that particles will hit them; the code does not sort "
01305                 "them. Usually this means that each place command of a "
01306                 "trackerplane must have a larger z value than the previous "
01307                 "place command. They must also come after the trackerZ value "
01308                 "of the tracker command."
01309         );
01310         
01311         // initial field values
01312         tracker = "";
01313         radius = 0.0;
01314         innerRadius = 0.0;
01315         height = 0.0;
01316         width = 0.0;
01317         length = 1.0*mm;
01318         material = "";
01319         color = "1,1,1";
01320         solid = 0;
01321         maxStep = -99.0;
01322         theta = 0.0;
01323         wireSpacing = -1.0;
01324         wireOffset = 0.0;
01325         errType = "fixed";
01326         errTheta = 0.0;
01327         errSpacing = 0.0;
01328         errOffset = 0.0;
01329         sigmaTime = -1.0;
01330 }

BLCMDtrackerplane::~BLCMDtrackerplane (  )  [virtual]

Destructor.

References solid.

01333 {
01334         if(solid) delete solid;
01335 }

BLCMDtrackerplane::BLCMDtrackerplane ( const BLCMDtrackerplane r  ) 

Copy constructor.

References color, errOffset, errSpacing, errTheta, errType, height, innerRadius, length, material, maxStep, radius, sigmaTime, solid, theta, tracker, width, wireOffset, and wireSpacing.

01337                                                                : BLElement(r)
01338 {
01339         tracker = r.tracker;
01340         radius = r.radius;
01341         innerRadius = r.innerRadius;
01342         height = r.height;
01343         width = r.width;
01344         length = r.length;
01345         material = r.material;
01346         color = r.color;
01347         solid = r.solid;
01348         maxStep = r.maxStep;
01349         theta = r.theta;
01350         wireSpacing = r.wireSpacing;
01351         wireOffset = r.wireOffset;
01352         errType = r.errType;
01353         errTheta = r.errTheta;
01354         errSpacing = r.errSpacing;
01355         errOffset = r.errOffset;
01356         sigmaTime = r.sigmaTime;
01357 }


Member Function Documentation

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

clone()

Implements BLElement.

References BLCMDtrackerplane().

00309 { return new BLCMDtrackerplane(*this); }

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

commandName() returns "trackerplane".

Implements BLCommand.

00315 { return "trackerplane"; }

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

command() implements the trackerplane command.

Implements BLCommand.

References BLCMDtrackerplane(), errOffset, errSpacing, errTheta, errType, BLCMDtracker::findTracker(), BLParam::getDouble(), BLCommand::getMaterial(), BLCommand::handleNamedArgs(), material, maxStep, Param, BLCommand::print(), BLCommand::printError(), BLElement::setName(), and tracker.

01360 {
01361         if(argv.size() != 1) {
01362                 printError("trackerplane: Invalid command, must have name");
01363                 return -1;
01364         }
01365 
01366         if(argv[0] == "default") {
01367                 return defaultTrackerPlane.handleNamedArgs(namedArgs);
01368         }
01369 
01370         BLCMDtrackerplane *t = new BLCMDtrackerplane(defaultTrackerPlane);
01371         t->setName(argv[0]);
01372         int retval = t->handleNamedArgs(namedArgs);
01373 
01374         if(t->maxStep < 0.0) t->maxStep = Param.getDouble("maxStep");
01375 
01376         if(BLCMDtracker::findTracker(t->tracker) == 0)
01377                 printError("trackerplane: invalid tracker");
01378 
01379         // handle errType
01380         if(t->errType.at(0) == 'g' || t->errType.at(0) == 'G') {
01381                 t->errTheta = t->errTheta*CLHEP::RandGauss::shoot();
01382                 t->errSpacing = t->errSpacing*CLHEP::RandGauss::shoot();
01383                 t->errOffset = t->errOffset*CLHEP::RandGauss::shoot();
01384         } else if(t->errType.at(0) == 'r' || t->errType.at(0) == 'R') {
01385                 t->errTheta = t->errTheta - t->errTheta*2.0*G4UniformRand();
01386                 t->errSpacing = t->errSpacing - t->errSpacing*2.0*G4UniformRand();
01387                 t->errOffset = t->errOffset - t->errOffset*2.0*G4UniformRand();
01388         }
01389 
01390         // check material exists
01391         if(t->material.size() > 0) getMaterial(t->material,false);
01392 
01393         t->print(argv[0]);
01394 
01395         return retval;
01396 }

void BLCMDtrackerplane::defineNamedArgs (  )  [virtual]

defineNamedArgs() defines the named arguments for the command.

Reimplemented from BLCommand.

References BLCommand::argDouble(), BLCommand::argString(), color, errOffset, errSpacing, errTheta, errType, height, innerRadius, length, material, maxStep, radius, sigmaTime, theta, tracker, width, wireOffset, and wireSpacing.

01399 {
01400         argString(tracker,"tracker","The tracker to which this plane belongs; REQUIRED.");
01401         argDouble(radius,"radius","The radius of the circular tracker plane (mm).");
01402         argDouble(innerRadius,"innerRadius","The inner radius of the circular tracker plane (0 mm).");
01403         argDouble(height,"height","The height of the rectangular tracker plane (mm).");
01404         argDouble(width,"width","The width of the rectangular tracker plane (mm).");
01405         argDouble(length,"length","The length of the tracker plane (mm).");
01406         argDouble(theta,"theta","Wire angle in X-Y plane (deg). 0=>x, 90=>y...",deg);
01407         argDouble(wireSpacing,"wireSpacing","Wire spacing (mm).");
01408         argDouble(wireOffset,"wireOffset","Wire # 0 offset (mm).");
01409         argString(errType,"errType","Error type: fixed, gaussian, rect.");
01410         argDouble(errTheta,"errTheta","Error in wire angle (deg).",deg);
01411         argDouble(errSpacing,"errSpacing","Error in wire spacing (mm).");
01412         argDouble(errOffset,"errOffset","Error in wire offset (mm).");
01413         argDouble(sigmaTime,"sigmaTime","Sigma for timing plane (ns), default=-1.");
01414         argDouble(maxStep,"maxStep","The maximum stepsize in the element (mm).");
01415         argString(material,"material","The material of the tracker plane.");
01416         argString(color,"color","The color of the tracker plane (''=invisible).");
01417 }

void BLCMDtrackerplane::argChanged (  )  [virtual]

argChanged() does internal computations after some arg changed

Reimplemented from BLCommand.

References height, radius, and width.

01420 {
01421         if(radius > 0.0) width = height = radius*2.0;
01422 }

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

construct() will construct the trackerplane. Used for normal placements of a tracker plane object.

Implements BLElement.

References color, BLParam::getDouble(), BLCommand::getMaterial(), BLElement::getName(), BLManager::getObject(), BLCommand::getVisAttrib(), height, innerRadius, length, material, maxStep, Param, BLCommand::printError(), radius, BLManager::registerBeamStep(), BLManager::registerReferenceParticleStep(), solid, TrackerPlaneInstance, and width.

01430 {
01431         G4Material *mat;
01432         if(material != "")
01433                 mat = getMaterial(material);
01434         else
01435                 mat = parent->GetMaterial();
01436 
01437         G4String thisname = parentName+getName();
01438 
01439         if(!solid) {
01440                 if(radius > 0.0) {
01441                         solid = new G4Tubs(thisname+"Tubs", innerRadius, radius,
01442                                         length/2.0, 0.0, 2.0*pi);
01443                 } else if(height > 0.0 && width > 0.0) {
01444                         solid = new G4Box(thisname+"Box",width/2.0,
01445                                         height/2.0,length/2.0);
01446                 } else {
01447                         printError("trackerplane::construct %s INVALID - no "
01448                                 "radius or height&width",thisname.c_str());
01449                         return;
01450                 }
01451         }
01452         G4LogicalVolume *lv = new G4LogicalVolume(solid,mat, thisname+"LogVol");
01453         lv->SetVisAttributes(getVisAttrib(color));
01454         if(maxStep < 0.0) maxStep = Param.getDouble("maxStep");
01455         lv->SetUserLimits(new G4UserLimits(maxStep));
01456 
01457         // geant4 rotation convention is backwards from g4beamline
01458         G4RotationMatrix *g4rot = 0;
01459         if(relativeRotation)
01460                 g4rot = new G4RotationMatrix(relativeRotation->inverse());
01461 
01462         G4VPhysicalVolume *pv = new G4PVPlacement(g4rot, relativePosition,
01463                                         lv,thisname,parent,false,0);
01464 
01465         // get globalRotation and globalPosition
01466         G4RotationMatrix *globalRotation = 0;
01467         if(relativeRotation && parentRotation) {
01468                 globalRotation = 
01469                     new G4RotationMatrix(*parentRotation * *relativeRotation);
01470         } else if(relativeRotation) {
01471                 globalRotation = relativeRotation;
01472         } else if(parentRotation) {
01473                 globalRotation = parentRotation;
01474         }
01475         G4ThreeVector globalPosition(relativePosition + parentPosition);
01476         if(globalRotation)
01477                 globalPosition = *globalRotation * relativePosition +
01478                                 parentPosition;
01479 
01480         TrackerPlaneInstance *tpi = new TrackerPlaneInstance(thisname,pv,this);
01481 
01482         BLManager::getObject()->registerBeamStep(pv,tpi);
01483         // I'm not sure if this needs reference step any more...
01484         BLManager::getObject()->registerReferenceParticleStep(pv,tpi);
01485 
01486         printf("BLCMDtrackerplane::Construct %s parent=%s relZ=%.1f globZ=%.1f\n",
01487                 thisname.c_str(),parentName.c_str(),relativePosition[2],
01488                 globalPosition[2]);
01489 }

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

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

Implements BLElement.

References length.

00336 { return length; }

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

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

Implements BLElement.

References width.

00339 { return width; }

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

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

Implements BLElement.

References height.

00342 { return height; }

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

isOK() returns true.

Implements BLElement.

00345 { return true; }

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

generatePoints() from BLElement

Implements BLElement.

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

01492 {
01493         if(radius > 0.0)
01494                 generateTubs(npoints, innerRadius, radius, 0.0, 360.0*deg,
01495                                                                 length, v);
01496         else
01497                 generateBox(npoints,width,height,length,v);
01498 }

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

isOutside() from BLElement

Implements BLElement.

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

01501 {
01502         if(radius > 0.0) {
01503                 G4double r = sqrt(local[0]*local[0]+local[1]*local[1]);
01504                 return r > radius-tolerance || r < innerRadius+tolerance ||
01505                         fabs(local[2]) > length/2.0-tolerance;
01506         } else {
01507                 return fabs(local[0]) > width/2.0-tolerance ||
01508                         fabs(local[1]) > height/2.0-tolerance ||
01509                         fabs(local[2]) > length/2.0-tolerance;
01510         }
01511 }


Friends And Related Function Documentation

friend class TrackerPlaneInstance [friend]

Referenced by construct().


Member Data Documentation

G4String BLCMDtrackerplane::tracker [private]

G4double BLCMDtrackerplane::radius [private]

G4double BLCMDtrackerplane::innerRadius [private]

G4double BLCMDtrackerplane::height [private]

G4double BLCMDtrackerplane::width [private]

G4double BLCMDtrackerplane::length [private]

G4String BLCMDtrackerplane::material [private]

G4String BLCMDtrackerplane::color [private]

G4VSolid* BLCMDtrackerplane::solid [private]

G4double BLCMDtrackerplane::maxStep [private]

G4double BLCMDtrackerplane::theta [private]

G4double BLCMDtrackerplane::wireSpacing [private]

G4double BLCMDtrackerplane::wireOffset [private]

G4String BLCMDtrackerplane::errType [private]

G4double BLCMDtrackerplane::errTheta [private]

G4double BLCMDtrackerplane::errSpacing [private]

G4double BLCMDtrackerplane::errOffset [private]

G4double BLCMDtrackerplane::sigmaTime [private]


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