PillboxField Class Reference

Inheritance diagram for PillboxField:

BLElementField BLManager::SteppingAction

List of all members.


Detailed Description

PillboxField represents one placement of a pillbox.

frequency=0.0 is accepted and yields a pillbox with a constant Ez = maxGradient (useful to verify units are correct).

Public Member Functions

 PillboxField (G4String &_name, BLCoordinateTransform &_global2local, BLCMDpillbox *_pillbox, G4VPhysicalVolume *_timingPV)
 constructor. _zmin/max are global coordinates.
G4String getName ()
 getName() returns the name of this placement.
void addFieldValue (const G4double point[4], G4double field[6]) const
 addFieldValue() adds the field for this pillbox into field[]. point[] is in global coordinates.
void UserSteppingAction (const G4Step *step)
 doStep() handles a step of the tune particle.

Private Types

enum  State { TIMING_UNKNOWN, SETTING_TIMING, TIMING_COMPLETE }

Private Attributes

G4String name
BLCMDpillboxpillbox
G4VPhysicalVolume * timingPV
BLCoordinateTransform global2local
G4RotationMatrix * rotation
G4double timeOffset
G4Track saveTrack
G4int timeCount
G4bool validSaveTrack
State state
BLFieldMapfieldMap

Friends

class BLCMDpillbox


Member Enumeration Documentation

enum PillboxField::State [private]

Enumerator:
TIMING_UNKNOWN 
SETTING_TIMING 
TIMING_COMPLETE 


Constructor & Destructor Documentation

PillboxField::PillboxField ( G4String &  _name,
BLCoordinateTransform _global2local,
BLCMDpillbox _pillbox,
G4VPhysicalVolume *  _timingPV 
)

constructor. _zmin/max are global coordinates.

References BLCMDpillbox::fieldMap, fieldMap, BLCMDpillbox::frequency, BLCoordinateTransform::getGlobal(), BLCoordinateTransform::getRotation(), global2local, BLCMDpillbox::innerLength, BLCMDpillbox::innerRadius, BLCoordinateTransform::isRotated(), isUndetermined(), BLCMDpillbox::maxGradient, name, pillbox, rotation, BLElementField::setGlobalPoint(), state, timeCount, BLCMDpillbox::timeOffset, timeOffset, TIMING_COMPLETE, TIMING_UNKNOWN, timingPV, validSaveTrack, and BLCMDpillbox::wallThick.

00595                                                       :
00596                         BLElementField(), BLManager::SteppingAction()
00597 {
00598         name = _name;
00599         pillbox = _pillbox;
00600         timingPV = _timingPV;
00601         global2local = _global2local;
00602         rotation = 0;
00603         timeOffset = pillbox->timeOffset;
00604         timeCount = 0;
00605         validSaveTrack = false;
00606         state = TIMING_UNKNOWN;
00607         fieldMap = pillbox->fieldMap;
00608 
00609         if(global2local.isRotated()) {
00610                 rotation = new G4RotationMatrix(global2local.getRotation());
00611                 rotation->invert();
00612         }
00613 
00614         if(pillbox->frequency == 0.0 || pillbox->maxGradient == 0.0)
00615                 timeOffset = 0.0;
00616 
00617         // don't tune timeOffset if it is known
00618         if(!isUndetermined(timeOffset))
00619                 state = TIMING_COMPLETE;
00620         else
00621                 timeOffset = 0.0; // don't start tuning at UNDETERMINED
00622 
00623         // set global bounding box
00624         G4double local[4], global[4];
00625         local[3] = 0.0;
00626         G4double dz = pillbox->innerLength/2.0 + pillbox->wallThick;
00627         for(int i=0; i<2; ++i) {
00628                 local[0] = (i==0 ? -1.0 : 1.0) * pillbox->innerRadius;
00629                 for(int j=0; j<2; ++j) {
00630                         local[1] = (j==0 ? -1.0 : 1.0) * pillbox->innerRadius;
00631                         for(int k=0; k<2; ++k) {
00632                                 local[2] = (k==0 ? -1.0 : 1.0) * dz;
00633                                 global2local.getGlobal(local,global);
00634                                 setGlobalPoint(global);
00635                         }
00636                 }
00637         }
00638 }


Member Function Documentation

G4String PillboxField::getName (  )  [inline]

getName() returns the name of this placement.

References name.

Referenced by UserSteppingAction().

00080 { return name; }

void PillboxField::addFieldValue ( const G4double  point[4],
G4double  field[6] 
) const [virtual]

addFieldValue() adds the field for this pillbox into field[]. point[] is in global coordinates.

Implements BLElementField.

References B, E, BLCMDpillbox::E2Bfactor, fieldMap, BLCMDpillbox::frequency, BLFieldMap::getFieldValue(), BLCoordinateTransform::getLocal(), global2local, BLCMDpillbox::innerLength, BLCMDpillbox::innerRadius, BLCMDpillbox::maxGradient, BLCMDpillbox::omega, pillbox, BLCMDpillbox::rFactor, rotation, BLCMDpillbox::skinDepth, timeOffset, and BLCMDpillbox::wallThick.

00641 {
00642         G4double local[4];
00643         global2local.getLocal(local,point);
00644 
00645         if(fieldMap) {
00646                 G4double myField[6];
00647                 fieldMap->getFieldValue(local,myField,pillbox->maxGradient,
00648                                                         pillbox->maxGradient);
00649                 // frequency == 0.0 means a d.c. electric field
00650                 G4double timeE=1.0, timeB=0.0;
00651                 if(pillbox->frequency != 0) {
00652                         double tArg = pillbox->omega * (point[3] - timeOffset);
00653                         timeE = cos(tArg);
00654                         timeB = sin(tArg);
00655                 }
00656                 if(rotation) {
00657                         G4ThreeVector B(myField[0],myField[1],myField[2]);
00658                         G4ThreeVector E(myField[3],myField[4],myField[5]);
00659                         B = *rotation * B;
00660                         E = *rotation * E;
00661                         field[0] += B[0] * timeB;
00662                         field[1] += B[1] * timeB;
00663                         field[2] += B[2] * timeB;
00664                         field[3] += E[0] * timeE;
00665                         field[4] += E[1] * timeE;
00666                         field[5] += E[2] * timeE;
00667                 } else {
00668                         field[0] += myField[0] * timeB;
00669                         field[1] += myField[1] * timeB;
00670                         field[2] += myField[2] * timeB;
00671                         field[3] += myField[3] * timeE;
00672                         field[4] += myField[4] * timeE;
00673                         field[5] += myField[5] * timeE;
00674                 }
00675                 return;
00676         }
00677 
00678         G4double x = local[0];
00679         G4double y = local[1];
00680         G4double z = local[2];
00681         G4double r = sqrt(x*x + y*y);
00682         G4double dz = pillbox->innerLength/2.0 + pillbox->wallThick;
00683         if(z < -dz || z > dz || r > pillbox->innerRadius) return;
00684 
00685         // frequency == 0.0 means a d.c. electric field
00686         if(pillbox->frequency == 0.0) {
00687                 // zmin and zmax include the walls, so re-test z position
00688                 if(z >= -pillbox->innerLength/2.0 && 
00689                    z <=  pillbox->innerLength/2.0) {
00690                         if(rotation) {
00691                                 G4ThreeVector E(0.0, 0.0, pillbox->maxGradient);
00692                                 E = *rotation * E;
00693                                 field[3] += E[0];
00694                                 field[4] += E[1];
00695                                 field[5] += E[2];
00696                         } else {
00697                                 field[5] += pillbox->maxGradient;
00698                         }
00699                 }
00700                 return;
00701         }
00702 
00703         // compute normal RF pillbox field
00704         double arg = pillbox->rFactor*r;
00705         if(fabs(arg) < 1.0e-20) arg = 0.0;
00706         double j0 = gsl_sf_bessel_J0(arg);
00707         double j1 = gsl_sf_bessel_J1(arg);
00708         double tArg = pillbox->omega * (point[3] - timeOffset);
00709         double ez = pillbox->maxGradient * j0 * cos(tArg);
00710         double bphi = -pillbox->maxGradient*pillbox->E2Bfactor*j1*sin(tArg);
00711         // handle skin depth if in wall or window
00712         double f = 1.0;
00713         if(z > pillbox->innerLength/2.0)
00714                 f=exp(-(z-pillbox->innerLength/2.0)/pillbox->skinDepth);
00715         else if(z < -pillbox->innerLength/2.0)
00716                 f=exp(-(-z-pillbox->innerLength/2.0)/pillbox->skinDepth);
00717         //@ double phi = atan2(point[1],point[0]);
00718         //@ TJR 3/29/2011 - this surel should be local, not point
00719         double phi = atan2(local[1],local[0]);
00720         if(rotation) {
00721                 G4ThreeVector B(-bphi*sin(phi)*f, bphi*cos(phi)*f, 0.0);
00722                 G4ThreeVector E(0.0, 0.0, ez*f);
00723                 B = *rotation * B;
00724                 E = *rotation * E;
00725                 field[0] += B[0];
00726                 field[1] += B[1];
00727                 field[2] += B[2];
00728                 field[3] += E[0];
00729                 field[4] += E[1];
00730                 field[5] += E[2];
00731         } else {
00732                 field[0] += -bphi * sin(phi) * f;
00733                 field[1] += bphi * cos(phi) * f;
00734                 field[5] += ez * f;
00735         }
00736 }

void PillboxField::UserSteppingAction ( const G4Step *  step  )  [virtual]

doStep() handles a step of the tune particle.

Implements BLManager::SteppingAction.

References BLCMDpillbox::frequency, getName(), BLManager::getObject(), BLManager::getSteppingManager(), BLCMDpillbox::innerLength, ITERATION_LIMIT, BLCMDpillbox::maxGradient, BLCMDpillbox::omega, BLCMDpillbox::phaseAcc, pillbox, saveTrack, SETTING_TIMING, state, timeCount, BLCMDpillbox::timeIncrement, timeOffset, TIMING_COMPLETE, timingPV, BLCMDpillbox::timingTolerance, and validSaveTrack.

00740 {
00741     // do nothing if d.c. or no field or timing is complete
00742     if(pillbox->frequency <= 0.0 || pillbox->maxGradient == 0.0) 
00743         state = TIMING_COMPLETE;
00744     if(state == TIMING_COMPLETE) return;
00745 
00746     G4Track *track = step->GetTrack();
00747     G4StepPoint *prePoint = step->GetPreStepPoint();
00748     if(!prePoint) return;
00749     G4VPhysicalVolume *prePV = prePoint->GetPhysicalVolume();
00750     G4StepPoint *postPoint = step->GetPostStepPoint();
00751     if(!postPoint) return;
00752     G4VPhysicalVolume *postPV = postPoint->GetPhysicalVolume();
00753     G4SteppingManager *steppingMgr =
00754                                 BLManager::getObject()->getSteppingManager();
00755 
00756     if(prePV == postPV) return;     // neither entering nor leaving
00757 
00758     if(postPV == timingPV) {        // entering timingPV
00759         // save track for setting timeOffset
00760         saveTrack.CopyTrackInfo(*step->GetTrack());
00761         saveTrack.SetUserInformation(0);
00762         validSaveTrack = true;
00763         state = SETTING_TIMING;
00764         // estimate timingOffset by neglecting acceleration
00765         double v = saveTrack.GetVelocity();
00766         double dist = pillbox->innerLength/2.0;
00767         double arrival = saveTrack.GetGlobalTime() + dist/v;
00768         // pi/2 because E is a cos() and we want phaseAcc=0
00769         // to be the rising zero-crossing
00770         timeOffset = arrival - (pillbox->phaseAcc-pi/2.0)/pillbox->omega;
00771     } else if(prePV == timingPV) {          // leaving timingPV
00772         if(++timeCount > ITERATION_LIMIT)
00773                 G4Exception("pillbox","Iteration Limit",FatalException,"");
00774         // check timing offset
00775         double t = track->GetGlobalTime();
00776         double phase = pillbox->omega * (t-timeOffset) + pi/2.0;
00777         double dt = (phase-pillbox->phaseAcc)/pillbox->omega;
00778         timeOffset = timeOffset + dt;
00779         if(fabs(dt) <= pillbox->timingTolerance) {
00780             // Success!
00781             printf("pillbox %s: Time OK  timeOffset=%.4f ns, incremented to "
00782                 "%.4f ns\n",getName().c_str(),timeOffset,
00783                                         timeOffset+pillbox->timeIncrement);
00784             timeOffset += pillbox->timeIncrement;
00785             state = TIMING_COMPLETE;
00786             return;
00787         }
00788         // restore saved data (i.e. jump back to when the track 
00789         // enteered TimingVol)
00790         if(!validSaveTrack)
00791                 G4Exception("pillbox","Invalid Step",FatalException,"");
00792         steppingMgr->GetfSecondary()->push_back(new G4Track(saveTrack));
00793         track->SetTrackStatus(fStopAndKill);
00794     }
00795 }


Friends And Related Function Documentation

friend class BLCMDpillbox [friend]


Member Data Documentation

G4String PillboxField::name [private]

Referenced by getName(), and PillboxField().

G4VPhysicalVolume* PillboxField::timingPV [private]

G4RotationMatrix* PillboxField::rotation [private]

Referenced by addFieldValue(), and PillboxField().

G4double PillboxField::timeOffset [private]

G4Track PillboxField::saveTrack [private]

Referenced by UserSteppingAction().

G4int PillboxField::timeCount [private]

Referenced by PillboxField(), and UserSteppingAction().

G4bool PillboxField::validSaveTrack [private]

Referenced by PillboxField(), and UserSteppingAction().

Referenced by PillboxField(), and UserSteppingAction().

Referenced by addFieldValue(), and PillboxField().


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