BLCMDpillbox Class Reference

Inheritance diagram for BLCMDpillbox:

BLElement BLCommand

List of all members.


Detailed Description

BLCMDpillbox implements a pillbox.

Each placement of a BLCMDpillbox creates a PillboxField that is linked into BLGlobalField. Each PillboxField sets its timeOffset via the tune particle (re-stepping through its TimingVol if necessary).

If fieldMapFile is non-null, it is read as a BLFieldMap and determines the peak field; both B and E are multiplied by maxGradient of the cavity.

Public Member Functions

 BLCMDpillbox ()
 Default constructor. Defines the command, etc.
virtual ~BLCMDpillbox ()
 Destructor.
 BLCMDpillbox (const BLCMDpillbox &r)
 Copy constructor.
BLElementclone ()
 clone()
G4String commandName ()
 commandName() returns "pillbox".
int command (BLArgumentVector &argv, BLArgumentMap &namedArgs)
 command() implements the pillbox command.
void defineNamedArgs ()
 defineNamedArgs() defines the named arguments of the command.
void argChanged ()
 argChanged() does internal computations after some arg changed
void construct (G4RotationMatrix *relativeRotation, G4ThreeVector relativePosition, G4LogicalVolume *parent, G4String parentName, G4RotationMatrix *parentRotation, G4ThreeVector parentPosition)
 construct() - construct the pillbox. Creates a new PillboxField and adds it to BLGlobalField.
G4double getLength ()
 getLength() returns the overallLength of the pillbox
G4double getWidth ()
 getWidth() returns the outer radius of the pillbox
G4double getHeight ()
 getHeight() returns the outer radius of the pillbox
G4bool isOK ()
 isOK() returns true only if all placed pillboxes of this type have had their timing set, and if the tuning (if used) converged.
void generatePoints (int npoints, std::vector< G4ThreeVector > &v)
 generatePoints() from BLElement
G4bool isOutside (G4ThreeVector &local, G4double tolerance)
 isOutside() from BLElement

Private Attributes

G4double maxGradient
G4String color
G4double frequency
G4double innerLength
G4double innerRadius
G4double pipeThick
G4double wallThick
G4double irisRadius
G4double collarRadialThick
G4double collarThick
G4double win1Thick
G4double win1OuterRadius
G4double win2Thick
G4String winMat
G4double phaseAcc
G4double skinDepth
G4double timingTolerance
G4double maxStep
G4String cavityMaterial
G4double timeOffset
G4double timeIncrement
G4String fieldMapFile
G4int kill
G4double overallLength
G4double outerRadius
G4double omega
G4double rFactor
G4double E2Bfactor
BLFieldMapfieldMap
std::vector< PillboxField * > pillboxField

Friends

class PillboxField


Constructor & Destructor Documentation

BLCMDpillbox::BLCMDpillbox (  ) 

Default constructor. Defines the command, etc.

References BLCMDTYPE_ELEMENT, cavityMaterial, collarRadialThick, collarThick, color, E2Bfactor, fieldMap, fieldMapFile, frequency, innerLength, innerRadius, irisRadius, kill, maxGradient, maxStep, omega, outerRadius, overallLength, phaseAcc, pipeThick, BLCommand::registerCommand(), rFactor, BLCommand::setDescription(), BLCommand::setSynopsis(), skinDepth, timeIncrement, timeOffset, timingTolerance, UNDETERMINED, wallThick, win1OuterRadius, win1Thick, win2Thick, and winMat.

Referenced by clone(), and command().

00191                            : BLElement(), pillboxField()
00192 {
00193         // register the commandName(), and its synopsis and description.
00194         registerCommand(BLCMDTYPE_ELEMENT);
00195         setSynopsis("Defines a pillbox RF cavity");
00196         setDescription("A Pillbox RF cavity is the basic RF element used to\n"
00197                 "construct a linac. The phaseAcc parameter sets the\n"
00198                 "phase of the tune particle at the center of the\n"
00199                 "cavity, and the timing offset of the cavity is determined\n"
00200                 "from that the first time that the Tune particle is tracked "
00201                 "through the cavity. Zero degrees is the rising zero-crossing "
00202                 "of the Ez field. If timeOffset is specified, it is used "
00203                 "rather than setting it from the Tune particle.\n\n"
00204                 "The Pipe, walls, and collars are always made of copper. "
00205                 "Pipe, wall, collar, win1, and win2 can be omitted by setting\n"
00206                 "their thickness to 0. Common usage is to set the collar "
00207                 "values so by placing multiple pillboxes sequentially the "
00208                 "collars form a beam pipe between them.\n\n"
00209                 "Note that section 4.4 of the User's "
00210                 "Guide has a dimensioned drawing of a pillbox.");
00211 
00212         // provide initial values for fields
00213         maxGradient = UNDETERMINED; // will be converted to MV/m
00214         color = "1.0,0.0,0.0";
00215         frequency = UNDETERMINED;
00216         innerLength = UNDETERMINED;
00217         innerRadius = UNDETERMINED;
00218         pipeThick = 3.0*mm;
00219         wallThick = 3.0*mm;
00220         irisRadius = 11.0*cm;
00221         collarRadialThick = 5.0*mm;
00222         collarThick = 2.5*mm;
00223         win1Thick = 0.200*mm;
00224         win1OuterRadius = 5.0*cm;
00225         win2Thick = 0.500*mm;
00226         winMat = "Be";
00227         phaseAcc = 40.0*deg;
00228         skinDepth = 0.002*mm;
00229         timingTolerance = 0.001*ns;
00230         maxStep = -1.0;
00231         cavityMaterial = "Vacuum";
00232         timeOffset = UNDETERMINED;
00233         timeIncrement = 0.0;
00234         fieldMapFile = "";
00235         kill = 0;
00236         // non-arguments:
00237         overallLength = UNDETERMINED;
00238         outerRadius = UNDETERMINED;
00239         omega = UNDETERMINED;
00240         rFactor = UNDETERMINED;
00241         E2Bfactor = UNDETERMINED;
00242         fieldMap = 0;
00243 }

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

Destructor.

00138 { }

BLCMDpillbox::BLCMDpillbox ( const BLCMDpillbox r  ) 

Copy constructor.

References cavityMaterial, collarRadialThick, collarThick, color, BLTune::copyTunableArg(), E2Bfactor, fieldMap, fieldMapFile, frequency, innerLength, innerRadius, irisRadius, kill, maxGradient, maxStep, omega, outerRadius, overallLength, phaseAcc, pipeThick, rFactor, skinDepth, timeIncrement, timeOffset, timingTolerance, wallThick, win1OuterRadius, win1Thick, win2Thick, and winMat.

00246                                                 : BLElement(r), pillboxField()
00247 {
00248         // copy fields one at a time (transfers default values from the
00249         // default object to this new object).
00250         BLTune::copyTunableArg(&maxGradient,&r.maxGradient);
00251         color = r.color;
00252         frequency = r.frequency;
00253         innerLength = r.innerLength;
00254         innerRadius = r.innerRadius;
00255         pipeThick = r.pipeThick;
00256         wallThick = r.wallThick;
00257         irisRadius = r.irisRadius;
00258         collarRadialThick = r.collarRadialThick;
00259         collarThick = r.collarThick;
00260         win1Thick = r.win1Thick;
00261         win1OuterRadius = r.win1OuterRadius;
00262         win2Thick = r.win2Thick;
00263         winMat = r.winMat;
00264         phaseAcc = r.phaseAcc;
00265         skinDepth = r.skinDepth;
00266         timingTolerance = r.timingTolerance;
00267         maxStep = r.maxStep;
00268         cavityMaterial = r.cavityMaterial;
00269         timeOffset = r.timeOffset;
00270         timeIncrement = r.timeIncrement;
00271         fieldMapFile = r.fieldMapFile;
00272         kill = r.kill;
00273 
00274         // non-arguments:
00275         overallLength = r.overallLength;
00276         outerRadius = r.outerRadius;
00277         omega = r.omega;
00278         rFactor = r.rFactor;
00279         E2Bfactor = r.E2Bfactor;
00280         fieldMap = r.fieldMap;
00281 }


Member Function Documentation

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

clone()

Implements BLElement.

References BLCMDpillbox().

00144 { return new BLCMDpillbox(*this); }

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

commandName() returns "pillbox".

Implements BLCommand.

00147 { return "pillbox"; }

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

command() implements the pillbox command.

Implements BLCommand.

References argChanged(), BLCMDpillbox(), cavityMaterial, E, E2Bfactor, fieldMap, fieldMapFile, frequency, g4bl_exit(), BLCommand::getMaterial(), BLCommand::handleNamedArgs(), innerRadius, maxGradient, BLCommand::print(), BLCommand::printError(), BLFieldMap::readFile(), and BLElement::setName().

00284 {
00285         // If the name of the cavity is "testUnits", perform the test and exit.
00286         if(argv[0] == "testUnits") {
00287                 // set the test case values, from Al Moretti's SuperFish example
00288                 maxGradient = 1.0 * megavolt/meter;
00289                 frequency = 0.956188/ns;
00290                 double Hcorrect = 1545.0; // A/m
00291                 double BcorrectTesla = 4.0*pi*1E-7 * Hcorrect; // MKS
00292                 argChanged();
00293                 // the first maximum of J1() has a value of 0.5819
00294                 double Bmax = maxGradient*E2Bfactor*0.5819;
00295                 printf("pillobx testUnits:\n           "         
00296                         "Freq(GHz)  Emax(MV/m)  Bmax(Tesla)  Radius(mm)\n");
00297                 printf("Values:  %10.4f %10.4f %12.4e %11.4f\n",
00298                         frequency/(1e9*hertz), maxGradient/(megavolt/meter),
00299                         Bmax/tesla, innerRadius/mm);
00300                 // these values already are in the correct units:
00301                 printf("Correct: %10.4f %10.4f %12.4e %8.1f\n",
00302                         0.956188,1.0,BcorrectTesla, 120.0);
00303                 printf("pillbox testUnits: exiting\n");
00304                 g4bl_exit(0);
00305         }
00306 
00307         if(argv.size() != 1) {
00308                 printError("pillbox: Invalid command, must have name");
00309                 return -1;
00310         }
00311 
00312         if(argv[0] == "default") {
00313                 return defaultPillbox.handleNamedArgs(namedArgs);
00314         }
00315 
00316         BLCMDpillbox *p = new BLCMDpillbox(defaultPillbox);
00317         p->setName(argv[0]);
00318         int retval = p->handleNamedArgs(namedArgs);
00319 
00320         if(p->fieldMapFile != "") {
00321                 p->fieldMap = new BLFieldMap();
00322                 p->fieldMap->readFile(p->fieldMapFile);
00323         }
00324 
00325         // check material exists
00326         if(p->cavityMaterial.size() > 0) getMaterial(p->cavityMaterial,false);
00327 
00328         p->print(argv[0]);
00329 
00330         return retval;
00331 }

void BLCMDpillbox::defineNamedArgs (  )  [virtual]

defineNamedArgs() defines the named arguments of the command.

Reimplemented from BLCommand.

References BLCommand::argDouble(), BLCommand::argInt(), BLCommand::argString(), BLCommand::argTunable(), cavityMaterial, collarRadialThick, collarThick, color, fieldMapFile, frequency, innerLength, innerRadius, irisRadius, kill, maxGradient, maxStep, phaseAcc, pipeThick, skinDepth, timeIncrement, timeOffset, timingTolerance, wallThick, win1OuterRadius, win1Thick, win2Thick, and winMat.

00334 {
00335         argTunable(maxGradient,"maxGradient","The peak gradient of the cavity (MV/m)",megavolt/meter);
00336         argString(color,"color","The color of the cavity");
00337         argDouble(frequency,"frequency","The frequency of the cavity (GHz)",1/ns,"",false);
00338         argDouble(innerLength,"innerLength","The inside length of the cavity (mm)",mm,"",false);
00339         argDouble(innerRadius,"innerRadius","The inside radius of the cavity (mm)",mm,"",false);
00340         argDouble(pipeThick,"pipeThick","The thickness of the pipe wall (mm)",mm,"",false);
00341         argDouble(wallThick,"wallThick","The thickness of the cavity walls (mm)",mm,"",false);
00342         argDouble(irisRadius,"irisRadius","The radius of the iris (mm)",mm,"",false);
00343         argDouble(collarRadialThick,"collarRadialThick","The radial thickness of the collar (mm)",mm,"",false);
00344         argDouble(collarThick,"collarThick","The thickness of the collar along z(mm)",mm,"",false);
00345         argDouble(win1Thick,"win1Thick","The thickness of the central portion of the\n"
00346                 "windows; zero for no window (mm)",mm,"",false);
00347         argDouble(win1OuterRadius,"win1OuterRadius","The radius of the central portion of\n"
00348                 "the windows (mm)",mm,"",false);
00349         argDouble(win2Thick,"win2Thick","The thickness of the outer portion of the\n"
00350                 "windows; zero for no window (mm)",mm,"",false);
00351         argString(winMat,"winMat","The material of the windows");
00352         argDouble(phaseAcc,"phaseAcc","The reference phase of the cavity (degrees)",deg);
00353         argDouble(skinDepth,"skinDepth","The skin depth (mm)",mm,"",false);
00354         argDouble(timingTolerance,"timingTolerance","Tolerance for timing tuning (ns)",ns);
00355         argDouble(maxStep,"maxStep","The maximum stepsize in the element (mm).",mm);
00356         argString(cavityMaterial,"cavityMaterial","Material of cavity volume (Vacuum).");
00357         argDouble(timeOffset,"timeOffset","Time offset for cavity (default: tuned by tune particle) (ns).",ns);
00358         argDouble(timeIncrement,"timeIncrement","Increment to timeOffset, applied AFTER tuning. (ns).",ns);
00359         argString(fieldMapFile,"fieldMapFile","Filename for BLFieldMap (pillbox if null).",false);
00360         argInt(kill,"kill","Set nonzero to kill tracks that hit the pipe, walls, or collars (0).");
00361 }

void BLCMDpillbox::argChanged (  )  [virtual]

argChanged() does internal computations after some arg changed

Reimplemented from BLCommand.

References collarThick, E2Bfactor, frequency, BLParam::getDouble(), innerLength, innerRadius, maxStep, omega, outerRadius, overallLength, Param, pipeThick, BLCommand::printError(), rFactor, UNDETERMINED, wallThick, win1Thick, and win2Thick.

Referenced by command().

00364 {
00365         // compute non-argument variables
00366         omega = 2.0 * pi * frequency;
00367         if(innerRadius == UNDETERMINED) {
00368                 if(frequency > 0.0) {
00369                         // the first zero of J0(r) is at r=2.405
00370                         innerRadius =  2.405 / (omega/c_light);
00371                 } else {
00372                         printError("pillbox: innerRadius is undefined");
00373                 }
00374         }
00375         double wall=collarThick;
00376         if(wallThick > wall) wall = wallThick;
00377         if(win1Thick > wall) wall = win1Thick;
00378         if(win2Thick > wall) wall = win2Thick;
00379         overallLength = innerLength + 2.0 * wall;
00380         outerRadius = innerRadius + pipeThick;
00381         rFactor = omega/c_light;
00382         // Previous versions had the incorrect value from the BeamTools:
00383         //     E2Bfactor = 1.0e-9 / c_light;   // MeV/m, MKS units
00384         // C.f. Eq. 8.77 in J.D.Jackson, _Classical_Electrodynamics_,
00385         //       first edition (1962). His value is 1/c in vacuum.
00386         E2Bfactor = 1.0/c_light;   // Geant4 units
00387         if(maxStep < 0.0) maxStep = Param.getDouble("maxStep");
00388 }

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

construct() - construct the pillbox. Creates a new PillboxField and adds it to BLGlobalField.

Implements BLElement.

References BLGlobalField::addElementField(), cavityMaterial, collarRadialThick, collarThick, BLCommand::getMaterial(), BLElement::getName(), BLGlobalField::getObject(), BLManager::getObject(), BLCoordinateTransform::getPosition(), PillboxField::global2local, innerLength, innerRadius, irisRadius, BLCoordinateTransform::isRotated(), kill, maxStep, PillboxField::pillbox, pillboxField, PillboxField, pipeThick, BLManager::registerTuneParticleStep(), PillboxField::timingPV, wallThick, win1OuterRadius, win1Thick, win2Thick, and winMat.

00396 {
00397         G4String thisname = parentName+getName();
00398 
00399         G4Material *cu = getMaterial("Cu");
00400         G4Material *wm = getMaterial(winMat);
00401         G4Material *cavityMat = getMaterial(cavityMaterial);
00402         G4UserLimits *userLimits = new G4UserLimits(maxStep);
00403 
00404         // find enclosing Tubs
00405         G4double halflen = innerLength/2+wallThick;
00406         halflen = std::max(halflen,innerLength/2+collarThick);
00407         halflen = std::max(halflen,innerLength/2+win1Thick);
00408         halflen = std::max(halflen,innerLength/2+win2Thick);
00409         G4Tubs *tubs = new G4Tubs(getName(),0.0,innerRadius+pipeThick,
00410                                         halflen,0.0,2.0*pi);
00411         G4LogicalVolume *pillbox = new G4LogicalVolume(tubs,cavityMat,
00412                                         parentName+getName(),0,0,userLimits);
00413 
00414         G4ThreeVector loc;
00415         G4VPhysicalVolume *pv;
00416 
00417         // The outer pipe
00418         G4LogicalVolume *pipe=0;
00419         if(pipeThick > 0.001*mm) {
00420                 tubs = new G4Tubs(getName()+"pipe",innerRadius,
00421                                         innerRadius+pipeThick,
00422                                         innerLength/2+wallThick,0.0,2.0*pi);
00423                 pipe = new G4LogicalVolume(tubs,cu,
00424                                         getName()+"pipeLog",0,0,userLimits);
00425                 pv = new G4PVPlacement(0,loc,pipe,parentName+getName()+"Pipe",
00426                                         pillbox,false,0);
00427                 if(kill) BLManager::getObject()->
00428                         registerSteppingAction(pv,new BLKillTrack(thisname));
00429         }
00430 
00431         // the wall (us and ds) - radially: collarOuterRadius to innerRadius
00432         G4LogicalVolume *wall=0;
00433         if(wallThick > 0.001*mm) {
00434                 tubs = new G4Tubs(getName()+"wall",irisRadius+collarRadialThick,
00435                                         innerRadius,wallThick/2,0.0,2.0*pi);
00436                 wall = new G4LogicalVolume(tubs,cu,
00437                                         getName()+"wallLog",0,0,userLimits);
00438                 G4double wallCenterZ = innerLength/2.0 + wallThick/2.0;
00439                 loc.setZ(-wallCenterZ);
00440                 pv = new G4PVPlacement(0,loc,wall,parentName+getName()+"UsWall",
00441                                         pillbox,false,0);
00442                 if(kill) BLManager::getObject()->
00443                         registerSteppingAction(pv,new BLKillTrack(thisname));
00444                 loc.setZ(+wallCenterZ);
00445                 pv = new G4PVPlacement(0,loc,wall,parentName+getName()+"DsWall",
00446                                         pillbox,false,0);
00447                 if(kill) BLManager::getObject()->
00448                         registerSteppingAction(pv,new BLKillTrack(thisname));
00449         }
00450         
00451         // the collar (us and ds) - includes wall thickness (i.e. the inner
00452         // edge of the collar is the same as the inner edge of the wall)
00453         G4LogicalVolume *collar=0;
00454         if(collarThick > 0.001*mm) {
00455                 tubs = new G4Tubs(getName()+"collar",irisRadius,
00456                                         irisRadius+collarRadialThick,
00457                                         collarThick/2,0.0,2.0*pi);
00458                 collar = new G4LogicalVolume(tubs,cu,
00459                                         getName()+"collarLog",0,0,userLimits);
00460                 G4double collarCenterZ = innerLength/2.0 + collarThick/2.0;
00461                 loc.setZ(-collarCenterZ);
00462                 pv = new G4PVPlacement(0,loc,collar,parentName+getName()+"UsCollar",
00463                                         pillbox,false,0);
00464                 if(kill) BLManager::getObject()->
00465                         registerSteppingAction(pv,new BLKillTrack(thisname));
00466                 loc.setZ(+collarCenterZ);
00467                 pv = new G4PVPlacement(0,loc,collar,parentName+getName()+"DsCollar",
00468                                         pillbox,false,0);
00469                 if(kill) BLManager::getObject()->
00470                         registerSteppingAction(pv,new BLKillTrack(thisname));
00471         }
00472         
00473         // win1 (the inner window)
00474         G4LogicalVolume *win1=0;
00475         if(win1Thick > 0.001*mm) {
00476                 tubs = new G4Tubs(getName()+"win1",0.0,win1OuterRadius,
00477                                         win1Thick/2,0.0,2.0*pi);
00478                 win1 = new G4LogicalVolume(tubs,wm,
00479                                         getName()+"win1Log",0,0,userLimits);
00480                 G4double win1CenterZ = innerLength/2.0 + win1Thick/2.0;
00481                 loc.setZ(-win1CenterZ);
00482                 new G4PVPlacement(0,loc,win1,parentName+getName()+"UsWin1",
00483                                         pillbox,false,0);
00484                 loc.setZ(+win1CenterZ);
00485                 new G4PVPlacement(0,loc,win1,parentName+getName()+"DsWin1",
00486                                         pillbox,false,0);
00487         }
00488 
00489         // win2 (the outer window)
00490         G4LogicalVolume *win2=0;
00491         if(win2Thick > 0.001*mm) {
00492                 tubs = new G4Tubs(getName()+"win2",win1OuterRadius,irisRadius,
00493                                         win2Thick/2,0.0,2.0*pi);
00494                 win2 = new G4LogicalVolume(tubs,wm,
00495                                         getName()+"win2Log",0,0,userLimits);
00496                 G4double win2CenterZ = innerLength/2.0 + win2Thick/2.0;
00497                 loc.setZ(-win2CenterZ);
00498                 new G4PVPlacement(0,loc,win2,parentName+getName()+"UsWin2",
00499                                         pillbox,false,0);
00500                 loc.setZ(+win2CenterZ);
00501                 new G4PVPlacement(0,loc,win2,parentName+getName()+"DsWin2",
00502                                         pillbox,false,0);
00503         }
00504 
00505         // half of the interior
00506         tubs = new G4Tubs(getName()+"half",0.0,innerRadius,
00507                                         innerLength/4,0.0,2.0*pi);
00508         G4LogicalVolume *half = new G4LogicalVolume(tubs,cavityMat,getName()+"halfLog",
00509                                         0,0,userLimits);
00510         G4double halfCenterZ = innerLength/4.0;
00511         loc.setZ(+halfCenterZ);
00512         new G4PVPlacement(0,loc,half,parentName+getName()+"DsHalf",
00513                                                         pillbox,false,0);
00514         loc.setZ(-halfCenterZ);
00515         G4VPhysicalVolume *timingPV = new G4PVPlacement(0,loc,half,
00516                                 parentName+getName()+"TimingVol",pillbox,false,0);
00517 
00518 #ifdef G4VIS_USE
00519         G4VisAttributes *c1 = new G4VisAttributes(true,G4Color(0.5,0.0,0.0));
00520         G4VisAttributes *c2 = new G4VisAttributes(true,G4Color(0.7,0.0,0.0));
00521         G4VisAttributes *c3 = new G4VisAttributes(true,G4Color(1.0,0.0,0.0));
00522         const G4VisAttributes *invisible = &G4VisAttributes::GetInvisible();
00523         if(pipe) pipe->SetVisAttributes(c1);
00524         if(wall) wall->SetVisAttributes(c1);
00525         if(collar) collar->SetVisAttributes(c1);
00526         if(win1) win1->SetVisAttributes(c3);
00527         if(win2) win2->SetVisAttributes(c2);
00528         half->SetVisAttributes(invisible);
00529         pillbox->SetVisAttributes(invisible);
00530 #endif
00531 
00532         // geant4 rotation convention is backwards from g4beamline
00533         G4RotationMatrix *g4rot = 0;
00534         if(relativeRotation)
00535                 g4rot = new G4RotationMatrix(relativeRotation->inverse());
00536 
00537         // place the pillbox into the parent
00538         new G4PVPlacement(g4rot,relativePosition,pillbox,
00539                                 parentName+getName(),parent,false,0);
00540 
00541         // get globalRotation and globalPosition
00542         G4RotationMatrix *globalRotation = 0;
00543         if(relativeRotation && parentRotation) {
00544                 globalRotation = 
00545                     new G4RotationMatrix(*parentRotation * *relativeRotation);
00546         } else if(relativeRotation) {
00547                 globalRotation = relativeRotation;
00548         } else if(parentRotation) {
00549                 globalRotation = parentRotation;
00550         }
00551         G4ThreeVector globalPosition(relativePosition + parentPosition);
00552         if(parentRotation)
00553                 globalPosition = *parentRotation * relativePosition +
00554                                 parentPosition;
00555 
00556         BLCoordinateTransform global2local(globalRotation,globalPosition);
00557 
00558         // add this pillbox field to the GlobalField, and to BLManager
00559         G4double dz = innerLength/2.0 + wallThick;
00560         // if rotated, make an overestimate of the field occupancy along z
00561         if(global2local.isRotated())
00562                 dz += innerRadius;
00563         PillboxField *pf = 
00564                 new PillboxField(thisname,global2local,this,timingPV);
00565         BLGlobalField::getObject()->addElementField(pf);
00566         pillboxField.push_back(pf);
00567         BLManager::getObject()->registerTuneParticleStep(timingPV,pf);
00568 
00569         printf("BLCMDpillbox::construct %s parent=%s relZ=%.1f globZ=%.1f\n"
00570                         "\tzmin=%.1f zmax=%.1f\n",
00571                 thisname.c_str(),parentName.c_str(),relativePosition[2],
00572                 global2local.getPosition()[2], global2local.getPosition()[2]-dz,
00573                 global2local.getPosition()[2]+dz);
00574 }

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

getLength() returns the overallLength of the pillbox

Implements BLElement.

References overallLength.

00168 { return overallLength; }

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

getWidth() returns the outer radius of the pillbox

Implements BLElement.

References outerRadius.

00171 { return outerRadius*2.0; }

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

getHeight() returns the outer radius of the pillbox

Implements BLElement.

References outerRadius.

00174 { return outerRadius*2.0; }

G4bool BLCMDpillbox::isOK (  )  [virtual]

isOK() returns true only if all placed pillboxes of this type have had their timing set, and if the tuning (if used) converged.

Implements BLElement.

References pillboxField, and PillboxField::TIMING_COMPLETE.

00577 {
00578         G4bool retval = true;
00579 
00580         // verify all placed pillboxField-s have set their timeOffset 
00581         std::vector<PillboxField *>::iterator i;
00582         for(i=pillboxField.begin(); i!=pillboxField.end(); ++i) {
00583                 if((*i)->state != PillboxField::TIMING_COMPLETE) {
00584                         printf("*** BLCMDpillbox %s has not determined timing\n",
00585                                 (*i)->getName().c_str());
00586                         retval = false;
00587                 }
00588         }
00589 
00590         return retval;
00591 }

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

generatePoints() from BLElement

Implements BLElement.

References BLElement::generateTubs(), outerRadius, and overallLength.

00798 {
00799         generateTubs(npoints, 0.0, outerRadius, 0.0, 360.0*deg,
00800                         overallLength, v);
00801 }

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

isOutside() from BLElement

Implements BLElement.

References outerRadius, and overallLength.

00804 {
00805         G4double r = sqrt(local[0]*local[0]+local[1]*local[1]);
00806         return r > outerRadius-tolerance ||
00807                 fabs(local[2]) > overallLength/2.0-tolerance;
00808 }


Friends And Related Function Documentation

friend class PillboxField [friend]

Referenced by construct().


Member Data Documentation

G4double BLCMDpillbox::maxGradient [private]

G4String BLCMDpillbox::color [private]

Referenced by BLCMDpillbox(), and defineNamedArgs().

G4double BLCMDpillbox::frequency [private]

G4double BLCMDpillbox::innerLength [private]

G4double BLCMDpillbox::innerRadius [private]

G4double BLCMDpillbox::pipeThick [private]

G4double BLCMDpillbox::wallThick [private]

G4double BLCMDpillbox::irisRadius [private]

G4double BLCMDpillbox::collarRadialThick [private]

G4double BLCMDpillbox::collarThick [private]

G4double BLCMDpillbox::win1Thick [private]

G4double BLCMDpillbox::win1OuterRadius [private]

G4double BLCMDpillbox::win2Thick [private]

G4String BLCMDpillbox::winMat [private]

G4double BLCMDpillbox::phaseAcc [private]

G4double BLCMDpillbox::skinDepth [private]

G4double BLCMDpillbox::timingTolerance [private]

G4double BLCMDpillbox::maxStep [private]

G4String BLCMDpillbox::cavityMaterial [private]

G4double BLCMDpillbox::timeOffset [private]

G4double BLCMDpillbox::timeIncrement [private]

G4String BLCMDpillbox::fieldMapFile [private]

G4int BLCMDpillbox::kill [private]

G4double BLCMDpillbox::overallLength [private]

G4double BLCMDpillbox::outerRadius [private]

G4double BLCMDpillbox::omega [private]

G4double BLCMDpillbox::rFactor [private]

G4double BLCMDpillbox::E2Bfactor [private]

std::vector<PillboxField *> BLCMDpillbox::pillboxField [private]

Referenced by construct(), and isOK().


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