BLCMDmultipole Class Reference

Inheritance diagram for BLCMDmultipole:

BLElement BLCommand

List of all members.


Detailed Description

BLCMDmultipole implements a generic multipole magnet with a cylindrical field.

The magnetic field is multipole in x and y, uniform along the z axis, modified by the fringe field computation described below. There's no field in the iron, so setting kill=1 is a good idea.

Public Member Functions

 BLCMDmultipole ()
 Default constructor. Defines the command, args, etc.
virtual ~BLCMDmultipole ()
 Destructor.
 BLCMDmultipole (const BLCMDmultipole &r)
 Copy constructor.
BLElementclone ()
 clone()
G4String commandName ()
 commandName() returns "multipole".
int command (BLArgumentVector &argv, BLArgumentMap &namedArgs)
 command() implements the multipole command.
void defineNamedArgs ()
 defineNamedArgs() defines the named arguments for the command.
void construct (G4RotationMatrix *relativeRotation, G4ThreeVector relativePosition, G4LogicalVolume *parent, G4String parentName, G4RotationMatrix *parentRotation, G4ThreeVector parentPosition)
 construct() - construct the multipole magnet
G4double getLength ()
 getLength() returns the ironLength of the multiple magnet
G4double getWidth ()
 getWidth() returns the outer radius of the multiple magnet
G4double getHeight ()
 getHeight() returns the outer radius of the multiple magnet
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 fieldLength
G4double ironLength
G4double ironRadius
G4double apertureRadius
G4String ironMaterial
G4String fieldMaterial
G4double dipole
G4double quadrupole
G4double sextupole
G4double octopole
G4double decapole
G4double dodecapole
G4String ironColor
G4int kill
G4double maxStep
G4String fringe
G4double fringeFactor
G4int openAperture

Friends

class MultipoleField


Constructor & Destructor Documentation

BLCMDmultipole::BLCMDmultipole (  ) 

Default constructor. Defines the command, args, etc.

References apertureRadius, BLCMDTYPE_ELEMENT, decapole, dipole, dodecapole, fieldLength, fieldMaterial, fringe, fringeFactor, ironColor, ironLength, ironMaterial, ironRadius, kill, maxStep, octopole, openAperture, quadrupole, BLCommand::registerCommand(), BLCommand::setDescription(), BLCommand::setSynopsis(), and sextupole.

Referenced by clone(), and command().

00144                                : BLElement()
00145 {
00146         // register the commandName(), and its synopsis and description.
00147         registerCommand(BLCMDTYPE_ELEMENT);
00148         setSynopsis("construct a generic multipole magnet.");
00149         setDescription("Multipole magnetic fields from dipole through dodecapole "
00150                 "are implemented with a cylindrical field region and "
00151                 "optional surrounding iron (ironLength=0 or ironRadius=0 "
00152                 "omits it). "
00153                 "All fields with positive strengths are oriented so in the "
00154                 "X-Z plane for X>0 (beam left) the field is purely By. "
00155                 "Negative strengths are allowed and reverse the field. "
00156                 "The fringe field computation is not implemented.");
00157 
00158         // provide initial values for fields
00159         fieldLength = 0.0;
00160         ironLength = 0.0;
00161         ironRadius = 0.0;
00162         apertureRadius = 0.0;
00163         ironMaterial = "Fe";
00164         fieldMaterial = "Vacuum";
00165         dipole = 0.0;
00166         quadrupole = 0.0;
00167         sextupole = 0.0;
00168         octopole = 0.0;
00169         decapole = 0.0;
00170         dodecapole = 0.0;
00171         ironColor = "1,1,1";
00172         kill = 0;
00173         maxStep = -1.0;
00174         fringe = "";
00175         fringeFactor = 1.0;
00176         openAperture = 0;
00177 }

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

Destructor.

00074 { }

BLCMDmultipole::BLCMDmultipole ( const BLCMDmultipole r  ) 

Copy constructor.

References apertureRadius, decapole, dipole, dodecapole, fieldLength, fieldMaterial, fringe, fringeFactor, ironColor, ironLength, ironMaterial, ironRadius, kill, maxStep, octopole, openAperture, quadrupole, and sextupole.

00180                                                       : BLElement(r) 
00181 {
00182         // copy fields one at a time (transfers default values from the
00183         // default object to this new object).
00184         fieldLength = r.fieldLength;
00185         ironLength = r.ironLength;
00186         ironRadius = r.ironRadius;
00187         apertureRadius = r.apertureRadius;
00188         ironMaterial = r.ironMaterial;
00189         fieldMaterial = r.fieldMaterial;
00190         dipole = r.dipole;
00191         quadrupole = r.quadrupole;
00192         sextupole = r.sextupole;
00193         octopole = r.octopole;
00194         decapole = r.decapole;
00195         dodecapole = r.dodecapole;
00196         ironColor = r.ironColor;
00197         kill = r.kill;
00198         maxStep = r.maxStep;
00199         fringe = r.fringe;
00200         fringeFactor = r.fringeFactor;
00201         openAperture = r.openAperture;
00202 }


Member Function Documentation

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

clone()

Implements BLElement.

References BLCMDmultipole().

00080 { return new BLCMDmultipole(*this); }

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

commandName() returns "multipole".

Implements BLCommand.

00083 { return "multipole"; }

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

command() implements the multipole command.

Implements BLCommand.

References BLCMDmultipole(), fieldMaterial, BLParam::getDouble(), BLCommand::getMaterial(), BLCommand::handleNamedArgs(), ironMaterial, maxStep, Param, BLCommand::print(), BLCommand::printError(), and BLElement::setName().

00205 {
00206         if(argv.size() != 1) {
00207                 printError("multipole: Invalid command, must have name");
00208                 return -1;
00209         }
00210 
00211         if(argv[0] == "default") {
00212                 return defaultMultipole.handleNamedArgs(namedArgs);
00213         }
00214 
00215         BLCMDmultipole *t = new BLCMDmultipole(defaultMultipole);
00216         t->setName(argv[0]);
00217         int retval = t->handleNamedArgs(namedArgs);
00218 
00219         if(t->maxStep < 0.0) t->maxStep = Param.getDouble("maxStep");
00220 
00221         // check material exists
00222         getMaterial(t->fieldMaterial,false);
00223         getMaterial(t->ironMaterial,false);
00224 
00225         t->print(argv[0]);
00226 
00227         return retval;
00228 }

void BLCMDmultipole::defineNamedArgs (  )  [virtual]

defineNamedArgs() defines the named arguments for the command.

Reimplemented from BLCommand.

References apertureRadius, BLCommand::argDouble(), BLCommand::argInt(), BLCommand::argString(), decapole, dipole, dodecapole, fieldLength, fieldMaterial, fringe, fringeFactor, ironColor, ironLength, ironMaterial, ironRadius, kill, maxStep, octopole, openAperture, quadrupole, and sextupole.

00231 {
00232         argDouble(fieldLength,"fieldLength","The length of the field region (mm)",mm);
00233         argDouble(ironLength,"ironLength","The length of the iron (mm)",mm);
00234         argDouble(ironRadius,"ironRadius","The outer radius of the iron (mm)",mm);
00235         argDouble(apertureRadius,"apertureRadius","The radius of the aperture (mm)",mm);
00236         argString(ironMaterial,"ironMaterial","The material of the iron region.");
00237         argString(fieldMaterial,"fieldMaterial","The material of the field region.");
00238         argDouble(dipole,"dipole","Strength of dipole (Tesla)",tesla);
00239         argDouble(quadrupole,"quadrupole","Strength of quadrupole (T/m)",tesla/meter);
00240         argDouble(sextupole,"sextupole","Strength of sextupole (T/m^2)",tesla/meter/meter);
00241         argDouble(octopole,"octopole","Strength of octopole (T/m^3)",tesla/meter/meter/meter);
00242         argDouble(decapole,"decapole","Strength of decapole (T/m^4)",tesla/meter/meter/meter/meter);
00243         argDouble(dodecapole,"dodecapole","Strength of dodecapole (T/m^5)",tesla/meter/meter/meter/meter/meter);
00244         argString(ironColor,"ironColor","The color of the iron region.");
00245         argInt(kill,"kill","Set nonzero to kill tracks hitting the iron.");
00246         argDouble(maxStep,"maxStep","The maximum stepsize in the element (mm)");
00247         argString(fringe,"fringe","Fringe field computation, set to 0 to disable");
00248         argDouble(fringeFactor,"fringeFactor","Fringe depth factor (1.0).");
00249         argInt(openAperture,"openAperture","Set nonzero to omit the aperture volume.",false);
00250 }

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

construct() - construct the multipole magnet

Implements BLElement.

References BLGlobalField::addElementField(), apertureRadius, fieldMaterial, BLCommand::getMaterial(), BLElement::getName(), BLGlobalField::getObject(), BLManager::getObject(), BLCommand::getVisAttrib(), ironColor, ironLength, ironMaterial, ironRadius, kill, maxStep, MultipoleField, and openAperture.

00259 {
00260         G4String thisname = parentName+getName();
00261 
00262         G4LogicalVolume *lviron=0;
00263         G4LogicalVolume *lvfield=0;
00264 
00265         if(ironLength > 0.0 && ironRadius > 0.0) {
00266                 G4Material *ironMat = getMaterial(ironMaterial);
00267                 G4Material *fieldMat = getMaterial(fieldMaterial);
00268                 G4VSolid *ironSolid = new G4Tubs(thisname, apertureRadius,
00269                                 ironRadius,ironLength/2.0,0.0,360.0*deg);
00270                 G4VSolid *fieldSolid = new G4Tubs(thisname, 0,
00271                                 apertureRadius,ironLength/2.0,0.0,360.0*deg);
00272                 lviron = new G4LogicalVolume(ironSolid,ironMat,thisname);
00273                 lviron->SetVisAttributes(getVisAttrib(ironColor));
00274                 lviron->SetUserLimits(new G4UserLimits(kill ? 0.1 :maxStep));
00275                 if(openAperture == 0) {
00276                     lvfield = new G4LogicalVolume(fieldSolid,fieldMat,thisname);
00277                     lvfield->SetVisAttributes(&G4VisAttributes::GetInvisible());
00278                     lvfield->SetUserLimits(new G4UserLimits(maxStep));
00279                 }
00280         }
00281 
00282         // geant4 rotation convention is backwards from g4beamline
00283         G4RotationMatrix *g4rot = 0;
00284         if(relativeRotation)
00285                 g4rot = new G4RotationMatrix(relativeRotation->inverse());
00286 
00287         G4PVPlacement *pv=0;
00288         if(lviron) pv = new G4PVPlacement(g4rot,relativePosition,lviron,
00289                                                 thisname,parent,false,0);
00290         if(lvfield) new G4PVPlacement(g4rot,relativePosition,lvfield,
00291                                                 thisname,parent,false,0);
00292 
00293         if(kill && pv)
00294                 BLManager::getObject()->
00295                         registerSteppingAction(pv,new BLKillTrack(thisname));
00296 
00297         // get globalRotation and globalPosition
00298         G4RotationMatrix *globalRotation = 0;
00299         if(relativeRotation && parentRotation) {
00300                 globalRotation = 
00301                     new G4RotationMatrix(*parentRotation * *relativeRotation);
00302         } else if(relativeRotation) {
00303                 globalRotation = relativeRotation;
00304         } else if(parentRotation) {
00305                 globalRotation = parentRotation;
00306         }
00307         G4ThreeVector globalPosition(relativePosition + parentPosition);
00308         if(parentRotation)
00309                 globalPosition = *parentRotation * relativePosition +
00310                                 parentPosition;
00311 
00312         BLCoordinateTransform global2local(globalRotation,globalPosition);
00313 
00314         MultipoleField *p = new MultipoleField(global2local,this);
00315         BLGlobalField::getObject()->addElementField(p);
00316 
00317         printf("BLCMDmultipole::Construct %s parent=%s relZ=%.1f globZ=%.1f\n",
00318                 thisname.c_str(),parentName.c_str(),relativePosition[2],
00319                 globalPosition[2]);
00320 }

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

getLength() returns the ironLength of the multiple magnet

Implements BLElement.

References ironLength.

00100 { return ironLength; }

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

getWidth() returns the outer radius of the multiple magnet

Implements BLElement.

References ironRadius.

00103 { return ironRadius*2.0; }

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

getHeight() returns the outer radius of the multiple magnet

Implements BLElement.

References ironRadius.

00106 { return ironRadius*2.0; }

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

isOK() returns true.

Implements BLElement.

00109 { return true; }

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

generatePoints() from BLElement

Implements BLElement.

References apertureRadius, BLElement::generateTubs(), ironLength, and ironRadius.

00403 {
00404         generateTubs(npoints, apertureRadius, ironRadius, 0.0, 360.0*deg,
00405                         ironLength, v);
00406 }

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

isOutside() from BLElement

Implements BLElement.

References ironLength, and ironRadius.

00409 {
00410         G4double r = sqrt(local[0]*local[0]+local[1]*local[1]);
00411         return r < apertureRadius+tolerance || r > ironRadius-tolerance ||
00412                 fabs(local[2]) > ironLength/2.0-tolerance;
00413 }


Friends And Related Function Documentation

friend class MultipoleField [friend]

Referenced by construct().


Member Data Documentation

G4double BLCMDmultipole::fieldLength [private]

G4double BLCMDmultipole::ironLength [private]

G4double BLCMDmultipole::ironRadius [private]

G4double BLCMDmultipole::apertureRadius [private]

G4String BLCMDmultipole::ironMaterial [private]

G4String BLCMDmultipole::fieldMaterial [private]

G4double BLCMDmultipole::dipole [private]

G4double BLCMDmultipole::quadrupole [private]

G4double BLCMDmultipole::sextupole [private]

G4double BLCMDmultipole::octopole [private]

G4double BLCMDmultipole::decapole [private]

G4double BLCMDmultipole::dodecapole [private]

G4String BLCMDmultipole::ironColor [private]

G4int BLCMDmultipole::kill [private]

G4double BLCMDmultipole::maxStep [private]

G4String BLCMDmultipole::fringe [private]

G4double BLCMDmultipole::fringeFactor [private]


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