BLCMDlilens Class Reference

Inheritance diagram for BLCMDlilens:

BLElement BLCommand

List of all members.


Detailed Description

BLCMDlilens implements a Lithium lens.

This is only the current-carrying Li cylinder, and its magnetic field. The field exists only between the end planes of the cylinder, out to radial infinity.

Public Member Functions

 BLCMDlilens ()
 Default constructor. Defines the command, args, etc.
virtual ~BLCMDlilens ()
 Destructor.
 BLCMDlilens (const BLCMDlilens &r)
 Copy constructor.
BLElementclone ()
 clone()
G4String commandName ()
 commandName() returns "lilens".
int command (BLArgumentVector &argv, BLArgumentMap &namedArgs)
 command() implements the lilens 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() - construct the Li lens
G4double getLength ()
 getLength() returns the length of the Li lens
G4double getWidth ()
 getWidth() returns the width of the Li lens
G4double getHeight ()
 getHeight() returns the height of the Li lens
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

G4double radius
G4double length
G4double current
G4String material
G4String color
G4double maxStep
G4Tubs * lilens


Constructor & Destructor Documentation

BLCMDlilens::BLCMDlilens (  ) 

Default constructor. Defines the command, args, etc.

References BLCMDTYPE_ELEMENT, color, current, length, lilens, material, maxStep, radius, BLCommand::registerCommand(), BLCommand::setDescription(), and BLCommand::setSynopsis().

Referenced by clone(), and command().

00116                          : BLElement()
00117 {
00118         // register the commandName(), and its synopsis and description.
00119         registerCommand(BLCMDTYPE_ELEMENT);
00120         setSynopsis("construct a simple Lithium lens.");
00121         setDescription(" This element consists of a current-carrying cylinder "
00122                 "and its field. The field exists only between the end planes "
00123                 "of the cylinder, out to adial infinity.");
00124         radius = 5.0;
00125         length = 100.0;
00126         current = 100000.0;
00127         material = "Li";
00128         color = "1,1,1";
00129         maxStep = -1.0;
00130         lilens = 0;
00131 }

virtual BLCMDlilens::~BLCMDlilens (  )  [inline, virtual]

Destructor.

00058 { }

BLCMDlilens::BLCMDlilens ( const BLCMDlilens r  ) 

Copy constructor.

References color, current, length, lilens, material, maxStep, and radius.

00133                                              : BLElement(r)
00134 {
00135         radius = r.radius;
00136         length = r.length;
00137         current = r.current;
00138         material = r.material;
00139         color = r.color;
00140         maxStep = r.maxStep;
00141         lilens = r.lilens;
00142 }


Member Function Documentation

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

clone()

Implements BLElement.

References BLCMDlilens().

00064 { return new BLCMDlilens(*this); }

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

commandName() returns "lilens".

Implements BLCommand.

00067 { return "lilens"; }

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

command() implements the lilens command.

Implements BLCommand.

References BLCMDlilens(), BLParam::getDouble(), BLCommand::getMaterial(), BLCommand::handleNamedArgs(), material, maxStep, Param, BLCommand::print(), BLCommand::printError(), and BLElement::setName().

00145 {
00146         if(argv.size() != 1) {
00147                 printError("lilens: Invalid command, must have name");
00148                 return -1;
00149         }
00150 
00151         if(argv[0] == "default") {
00152                 return defaultLiLens.handleNamedArgs(namedArgs);
00153         }
00154 
00155         BLCMDlilens *t = new BLCMDlilens(defaultLiLens);
00156         t->setName(argv[0]);
00157         int retval = t->handleNamedArgs(namedArgs);
00158 
00159         if(t->maxStep < 0.0) t->maxStep = Param.getDouble("maxStep");
00160 
00161         // check material exists
00162         if(t->material.size() > 0) getMaterial(t->material,false);
00163 
00164         t->print(argv[0]);
00165 
00166         return retval;
00167 }

void BLCMDlilens::defineNamedArgs (  )  [virtual]

defineNamedArgs() defines the named arguments for the command.

Reimplemented from BLCommand.

References BLCommand::argDouble(), BLCommand::argString(), color, current, length, material, maxStep, and radius.

00170 {
00171         argDouble(radius,"radius","The radius of the current-carrying cylinder (5 mm)", mm);
00172         argDouble(length,"length","The length of the cylinder (100 mm).", mm);
00173         argDouble(current,"current","The current in the cylinder (100000 Amp).");
00174         argString(material,"material","The material of the cylinder (Li).");
00175         argString(color,"color","The color of the tube or cylinder (''=invisible)");
00176         argDouble(maxStep,"maxStep","The maximum stepsize in the element (mm).",
00177                         mm);
00178 }

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

construct() - construct the Li lens

Implements BLElement.

References color, current, BLParam::getDouble(), getLength(), BLCommand::getMaterial(), BLElement::getName(), BLCommand::getVisAttrib(), length, lilens, M_PI, material, maxStep, Param, and radius.

00186 {
00187         G4String thisname = parentName+getName();
00188 
00189         if(!lilens)
00190                 lilens = new G4Tubs(thisname+"Tubs", 0.0, radius,
00191                                                 length/2.0, 0.0,2.0*M_PI);
00192         G4Material *mat = getMaterial(material);
00193         G4LogicalVolume *lv = new G4LogicalVolume(lilens,mat,thisname+"LogVol");
00194         lv->SetVisAttributes(getVisAttrib(color));
00195         if(maxStep < 0.0) maxStep = Param.getDouble("maxStep");
00196         lv->SetUserLimits(new G4UserLimits(maxStep));
00197 
00198         // geant4 rotation convention is backwards from g4beamline
00199         G4RotationMatrix *g4rot = 0;
00200         if(relativeRotation)
00201                 g4rot = new G4RotationMatrix(relativeRotation->inverse());
00202 
00203         G4VPhysicalVolume *pv = new G4PVPlacement(g4rot, relativePosition,lv,
00204                                         thisname, parent,false,0);
00205 
00206         // get globalRotation and globalPosition
00207         G4RotationMatrix *globalRotation = 0;
00208         if(relativeRotation && parentRotation) {
00209                 globalRotation = 
00210                     new G4RotationMatrix(*parentRotation * *relativeRotation);
00211         } else if(relativeRotation) {
00212                 globalRotation = relativeRotation;
00213         } else if(parentRotation) {
00214                 globalRotation = parentRotation;
00215         }
00216         G4ThreeVector globalPosition(relativePosition + parentPosition);
00217         if(parentRotation)
00218                 globalPosition = *parentRotation * relativePosition +
00219                                 parentPosition;
00220 
00221         new LiLensField(radius,length,current,globalRotation,globalPosition);
00222 
00223         printf("BLCMDlilens::Construct %s parent=%s relZ=%.1f globZ=%.1f\n"
00224                         "\tzmin=%.1f zmax=%.1f\n",
00225                 thisname.c_str(),parentName.c_str(),relativePosition[2],
00226                 globalPosition[2],
00227                 globalPosition[2]-getLength()/2.0,
00228                 globalPosition[2]+getLength()/2.0);
00229 }

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

getLength() returns the length of the Li lens

Implements BLElement.

References length.

Referenced by construct().

00084 { return length; }

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

getWidth() returns the width of the Li lens

Implements BLElement.

References radius.

00087 { return 2.0*radius; }

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

getHeight() returns the height of the Li lens

Implements BLElement.

References radius.

00090 { return 2.0*radius; }

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

isOK() returns true.

Implements BLElement.

00093 { return true; }

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

generatePoints() from BLElement

Implements BLElement.

References BLElement::generateTubs(), length, M_PI, and radius.

00232 {
00233         generateTubs(npoints, 0.0, radius, 0.0, 2.0*M_PI, length, v);
00234 }

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

isOutside() from BLElement

Implements BLElement.

References length, and radius.

00237 {
00238         G4double r = sqrt(local[0]*local[0]+local[1]*local[1]);
00239         return r > radius-tolerance || fabs(local[2]) > length/2.0-tolerance;
00240 }


Member Data Documentation

G4double BLCMDlilens::radius [private]

G4double BLCMDlilens::length [private]

G4double BLCMDlilens::current [private]

G4String BLCMDlilens::material [private]

G4String BLCMDlilens::color [private]

G4double BLCMDlilens::maxStep [private]

G4Tubs* BLCMDlilens::lilens [private]

Referenced by BLCMDlilens(), and construct().


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