BLCMDgenericquad Class Reference

Inheritance diagram for BLCMDgenericquad:

BLElement BLCommand

List of all members.


Detailed Description

BLCMDgenericquad implements an generic quadrupole magnet with a cylindrical field.

The magnetic field is quadrupole 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.

Set apertureRadius to get a circular aperture. Set poleTipRadius, coilRadius, and coilHalfwidth to get a "rounded + sign" aperture, using circles for the pole tips.

Public Member Functions

 BLCMDgenericquad ()
 Default constructor. Defines the command, args, etc.
virtual ~BLCMDgenericquad ()
 Destructor.
 BLCMDgenericquad (const BLCMDgenericquad &r)
 Copy constructor.
BLElementclone ()
 clone()
G4String commandName ()
 commandName() returns "genericquad".
int command (BLArgumentVector &argv, BLArgumentMap &namedArgs)
 command() implements the genericquad 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 generic quadrupole magnet
G4double getLength ()
 getLength() returns the ironLength of the quad
G4double getWidth ()
 getWidth() returns the outer radius of the quad
G4double getHeight ()
 getHeight() returns the outer radius of the quad
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
G4double poleTipRadius
G4double coilRadius
G4double coilHalfwidth
G4double gradient
G4String ironMaterial
G4String fieldMaterial
G4String ironColor
G4int kill
G4double maxStep
G4double fieldRadius
G4String fringe
G4double fringeFactor
G4int openAperture

Friends

class GenericQuadField


Constructor & Destructor Documentation

BLCMDgenericquad::BLCMDgenericquad (  ) 

Default constructor. Defines the command, args, etc.

References apertureRadius, BLCMDTYPE_ELEMENT, coilHalfwidth, coilRadius, fieldLength, fieldMaterial, fieldRadius, fringe, fringeFactor, gradient, ironColor, ironLength, ironMaterial, ironRadius, kill, maxStep, openAperture, poleTipRadius, BLCommand::registerCommand(), BLCommand::setDescription(), and BLCommand::setSynopsis().

Referenced by clone(), and command().

00147                                    : BLElement()
00148 {
00149         // register the commandName(), and its synopsis and description.
00150         registerCommand(BLCMDTYPE_ELEMENT);
00151         setSynopsis("construct a generic quadrupole magnet.");
00152         setDescription("The field region is a tubs with gradient specified.\n"
00153                 "A positive gradient yields a horizontally-focusing\n"
00154                 "quad for positive particles.\n"
00155                 "If apertureRadius>0 the quad has a circular aperture.\n"
00156                 "For a 'rounded +' aperture using circles for the poles,\n"
00157                 "set poleTipRadius, coilRadius, coilHalfwidth.\n"
00158                 "Due to visualization bugs, in the latter case you cannot\n"
00159                 "see through the aperture; it is solid black.\n"
00160                 "A fringe field computation based on the method of COSY\n"
00161                 "INFINITY is included by default, extending the field region.\n"
00162                 "This is first order only, "
00163                 "and the fringe field extends outside of "
00164                 "the magnet aperture only in a cylinder extending the "
00165                 "aperture straight along local z. As the fringe field is first "
00166                 "order only, it is slightly non-Maxwellian.");
00167 
00168         // provide initial values for fields
00169         fieldLength = 0.0;
00170         ironLength = 0.0;
00171         ironRadius = 0.0;
00172         apertureRadius = 0.0;
00173         poleTipRadius = 0.0;
00174         coilRadius = 0.0;
00175         coilHalfwidth = 0.0;
00176         gradient = 0.0;
00177         ironMaterial = "Fe";
00178         fieldMaterial = "Vacuum";
00179         ironColor = "1,1,1";
00180         kill = 0;
00181         maxStep = -1.0;
00182         fieldRadius = 0.0;
00183         fringe = "";
00184         fringeFactor = 1.0;
00185         openAperture = 0;
00186 }

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

Destructor.

00077 { }

BLCMDgenericquad::BLCMDgenericquad ( const BLCMDgenericquad r  ) 

Copy constructor.

References apertureRadius, coilHalfwidth, coilRadius, fieldLength, fieldMaterial, fieldRadius, fringe, fringeFactor, gradient, ironColor, ironLength, ironMaterial, ironRadius, kill, maxStep, openAperture, and poleTipRadius.

00189                                                             : BLElement(r) 
00190 {
00191         // copy fields one at a time (transfers default values from the
00192         // default object to this new object).
00193         fieldLength = r.fieldLength;
00194         ironLength = r.ironLength;
00195         ironRadius = r.ironRadius;
00196         apertureRadius = r.apertureRadius;
00197         poleTipRadius = r.poleTipRadius;
00198         coilRadius = r.coilRadius;
00199         coilHalfwidth = r.coilHalfwidth;
00200         gradient = r.gradient;
00201         ironMaterial = r.ironMaterial;
00202         fieldMaterial = r.fieldMaterial;
00203         ironColor = r.ironColor;
00204         kill = r.kill;
00205         maxStep = r.maxStep;
00206         fieldRadius = r.fieldRadius;
00207         fringe = r.fringe;
00208         fringeFactor = r.fringeFactor;
00209         openAperture = r.openAperture;
00210 }


Member Function Documentation

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

clone()

Implements BLElement.

References BLCMDgenericquad().

00083 { return new BLCMDgenericquad(*this); }

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

commandName() returns "genericquad".

Implements BLCommand.

00086 { return "genericquad"; }

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

command() implements the genericquad command.

Implements BLCommand.

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

00213 {
00214         if(argv.size() != 1) {
00215                 printError("genericquad: Invalid command, must have name");
00216                 return -1;
00217         }
00218 
00219         if(argv[0] == "default") {
00220                 return defaultGenericQuad.handleNamedArgs(namedArgs);
00221         }
00222 
00223         BLCMDgenericquad *t = new BLCMDgenericquad(defaultGenericQuad);
00224         t->setName(argv[0]);
00225         int retval = t->handleNamedArgs(namedArgs);
00226 
00227         if(t->maxStep < 0.0) t->maxStep = Param.getDouble("maxStep");
00228 
00229         // check material exists
00230         getMaterial(t->fieldMaterial,false);
00231         getMaterial(t->ironMaterial,false);
00232 
00233         t->print(argv[0]);
00234 
00235         return retval;
00236 }

void BLCMDgenericquad::defineNamedArgs (  )  [virtual]

defineNamedArgs() defines the named arguments for the command.

Reimplemented from BLCommand.

References apertureRadius, BLCommand::argDouble(), BLCommand::argInt(), BLCommand::argString(), coilHalfwidth, coilRadius, fieldLength, fieldMaterial, fringe, fringeFactor, gradient, ironColor, ironLength, ironMaterial, ironRadius, kill, maxStep, openAperture, and poleTipRadius.

00239 {
00240         argDouble(fieldLength,"fieldLength","The length of the field region (mm)");
00241         argDouble(ironLength,"ironLength","The length of the iron (mm)");
00242         argDouble(ironRadius,"ironRadius","The outer radius of the iron (mm)");
00243         argDouble(apertureRadius,"apertureRadius","The radius of the aperture (mm)");
00244         argDouble(poleTipRadius,"poleTipRadius","The inner radius of the pole tips (mm).");
00245         argDouble(coilRadius,"coilRadius","The radius of the inside of the coil (mm).");
00246         argDouble(coilHalfwidth,"coilHalfwidth","The halfwidth of the coil (mm).");
00247         argDouble(coilHalfwidth,"coilHalfWidth","Synonym for coilHalfwidth.");
00248         argDouble(gradient,"gradient","The magnetic field gradient, dBy/dx (Tesla/meter)",tesla/meter);
00249         argDouble(maxStep,"maxStep","The maximum stepsize in the element (mm)");
00250         argString(ironMaterial,"ironMaterial","The material of the iron region.");
00251         argString(fieldMaterial,"fieldMaterial","The material of the field region.");
00252         argString(ironColor,"ironColor","The color of the iron region.");
00253         argInt(kill,"kill","Set nonzero to kill tracks hitting the iron.");
00254         argString(fringe,"fringe","Fringe field computation, set to 0 to disable,"
00255                 " or a comma-separated list of 6 Enge function parameters.");
00256         argDouble(fringeFactor,"fringeFactor","Fringe depth factor (1.0).");
00257         argInt(openAperture,"openAperture","Set nonzero to omit the aperture volume.",false);
00258 }

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

construct() - construct the generic quadrupole magnet

Implements BLElement.

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

00267 {
00268         G4String thisname = parentName+getName();
00269 
00270         G4LogicalVolume *lviron=0;
00271         G4LogicalVolume *lvfield=0;
00272 
00273         G4Material *ironMat = getMaterial(ironMaterial);
00274         G4Material *fieldMat = getMaterial(fieldMaterial);
00275         if(apertureRadius > 0.0) {
00276                 G4VSolid *ironSolid = new G4Tubs(thisname, apertureRadius,
00277                                 ironRadius,ironLength/2.0,0.0,360.0*deg);
00278                 lviron = new G4LogicalVolume(ironSolid,ironMat,thisname);
00279                 lviron->SetVisAttributes(getVisAttrib(ironColor));
00280                 lviron->SetUserLimits(new G4UserLimits(kill ? 0.1 :maxStep));
00281                 if(openAperture == 0) {
00282                         G4VSolid *fieldSolid = new G4Tubs(thisname, 0,
00283                                 apertureRadius,ironLength/2.0,0.0,360.0*deg);
00284                         lvfield = new G4LogicalVolume(fieldSolid,fieldMat,
00285                                                                 thisname);
00286                         lvfield->SetVisAttributes(&G4VisAttributes::GetInvisible());
00287                         lvfield->SetUserLimits(new G4UserLimits(maxStep));
00288                 }
00289                 fieldRadius = apertureRadius;
00290         } else {
00291                 G4VSolid *box=0;
00292                 G4VSolid *ironSolid = new G4Tubs(thisname, 0.0,
00293                                 ironRadius,ironLength/2.0,0.0,360.0*deg);
00294                 box = new G4Box(thisname,coilRadius,
00295                                         coilRadius, ironLength/2.0+0.02*mm);
00296                 // compute circle approximating the hyperbolic pole tip --
00297                 // center is at (x,x) on the diagonal; (x0,y0) is the pole tip,
00298                 // (x1,y1) is the intersection of pole face with coil.
00299                 G4double x0=poleTipRadius/sqrt(2.0);
00300                 G4double y0=x0;
00301                 G4double x1=coilRadius;
00302                 G4double y1=coilHalfwidth;
00303                 G4double x=0.5*(x1*x1+y1*y1-x0*x0-y0*y0)/(x1+y1-x0-y0);
00304                 G4double r=sqrt(2.0)*x-poleTipRadius;
00305                 G4Tubs *t = new G4Tubs("xxx",0.0,r,ironLength/2.0+0.1,
00306                                                 0.0,360.0*deg);
00307                 // Subtract 4 copies of t from box, to get the "rounded +"
00308                 G4ThreeVector offset(x,x,0.0);
00309                 G4RotationMatrix norot;
00310                 box = new G4SubtractionSolid(thisname,box,
00311                                 t,&norot,offset);
00312                 offset[0] = -x;
00313                 box = new G4SubtractionSolid(thisname,box,
00314                                 t,&norot,offset);
00315                 offset[1] = -x;
00316                 box = new G4SubtractionSolid(thisname,box,
00317                                 t,&norot,offset);
00318                 offset[0] = x;
00319                 box = new G4SubtractionSolid(thisname,box,
00320                                 t,&norot,offset);
00321                 lviron = new G4LogicalVolume(ironSolid,ironMat,thisname);
00322                 lviron->SetVisAttributes(getVisAttrib(ironColor));
00323                 lviron->SetUserLimits(new G4UserLimits(kill ? 0.1 : maxStep));
00324                 // place the "rounded +" as daughter in the iron Tubs; black.
00325                 G4LogicalVolume *lvbox = new G4LogicalVolume(box,fieldMat,
00326                                         thisname);
00327                 lvbox->SetVisAttributes(getVisAttrib("0,0,0"));
00328                 lvbox->SetUserLimits(new G4UserLimits(maxStep));
00329                 G4ThreeVector no_offset(0.0,0.0,0.0);
00330                 new G4PVPlacement(0,no_offset,lvbox,thisname,lviron,false,0);
00331                 fieldRadius = 
00332                         sqrt(coilRadius*coilRadius+coilHalfwidth*coilHalfwidth);
00333                 lvfield = 0;
00334                 // set apertureRadius for fringe field computation
00335                 apertureRadius = poleTipRadius;
00336         }
00337 
00338         // geant4 rotation convention is backwards from g4beamline
00339         G4RotationMatrix *g4rot = 0;
00340         if(relativeRotation)
00341                 g4rot = new G4RotationMatrix(relativeRotation->inverse());
00342 
00343         G4PVPlacement *pv = new G4PVPlacement(g4rot,relativePosition,lviron,
00344                                                 thisname,parent,false,0);
00345         if(lvfield) new G4PVPlacement(g4rot,relativePosition,lvfield,
00346                                                 thisname,parent,false,0);
00347 
00348         if(kill)
00349                 BLManager::getObject()->
00350                         registerSteppingAction(pv,new BLKillTrack(thisname));
00351 
00352         // get globalRotation and globalPosition
00353         G4RotationMatrix *globalRotation = 0;
00354         if(relativeRotation && parentRotation) {
00355                 globalRotation = 
00356                     new G4RotationMatrix(*parentRotation * *relativeRotation);
00357         } else if(relativeRotation) {
00358                 globalRotation = relativeRotation;
00359         } else if(parentRotation) {
00360                 globalRotation = parentRotation;
00361         }
00362         G4ThreeVector globalPosition(relativePosition + parentPosition);
00363         if(parentRotation)
00364                 globalPosition = *parentRotation * relativePosition +
00365                                 parentPosition;
00366 
00367         BLCoordinateTransform global2local(globalRotation,globalPosition);
00368 
00369         GenericQuadField *p = new GenericQuadField(global2local,this);
00370         BLGlobalField::getObject()->addElementField(p);
00371 
00372         printf("BLCMDgenericquad::Construct %s parent=%s relZ=%.1f globZ=%.1f\n",
00373                 thisname.c_str(),parentName.c_str(),relativePosition[2],
00374                 globalPosition[2]);
00375 }

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

getLength() returns the ironLength of the quad

Implements BLElement.

References ironLength.

00103 { return ironLength; }

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

getWidth() returns the outer radius of the quad

Implements BLElement.

References ironRadius.

00106 { return ironRadius*2.0; }

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

getHeight() returns the outer radius of the quad

Implements BLElement.

References ironRadius.

00109 { return ironRadius*2.0; }

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

isOK() returns true.

Implements BLElement.

00112 { return true; }

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

generatePoints() from BLElement

Implements BLElement.

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

00454 {
00455         generateTubs(npoints, apertureRadius, ironRadius, 0.0, 360.0*deg,
00456                         ironLength, v);
00457 }

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

isOutside() from BLElement

Implements BLElement.

References ironLength, and ironRadius.

00460 {
00461         G4double r = sqrt(local[0]*local[0]+local[1]*local[1]);
00462         return r < apertureRadius+tolerance || r > ironRadius-tolerance ||
00463                 fabs(local[2]) > ironLength/2.0-tolerance;
00464 }


Friends And Related Function Documentation

friend class GenericQuadField [friend]

Referenced by construct().


Member Data Documentation

G4double BLCMDgenericquad::fieldLength [private]

G4double BLCMDgenericquad::ironLength [private]

G4double BLCMDgenericquad::ironRadius [private]

G4double BLCMDgenericquad::poleTipRadius [private]

G4double BLCMDgenericquad::coilRadius [private]

G4double BLCMDgenericquad::coilHalfwidth [private]

G4double BLCMDgenericquad::gradient [private]

G4String BLCMDgenericquad::ironMaterial [private]

G4String BLCMDgenericquad::fieldMaterial [private]

G4String BLCMDgenericquad::ironColor [private]

G4int BLCMDgenericquad::kill [private]

G4double BLCMDgenericquad::maxStep [private]

G4double BLCMDgenericquad::fieldRadius [private]

G4String BLCMDgenericquad::fringe [private]

G4double BLCMDgenericquad::fringeFactor [private]


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