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 |
BLCMDpillbox * | pillbox |
G4VPhysicalVolume * | timingPV |
BLCoordinateTransform | global2local |
G4RotationMatrix * | rotation |
G4double | timeOffset |
G4Track | saveTrack |
G4int | timeCount |
G4bool | validSaveTrack |
State | state |
BLFieldMap * | fieldMap |
Friends | |
class | BLCMDpillbox |
enum PillboxField::State [private] |
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.
00594 : 00595 BLElementField(), BLManager::SteppingAction() 00596 { 00597 name = _name; 00598 pillbox = _pillbox; 00599 timingPV = _timingPV; 00600 global2local = _global2local; 00601 rotation = 0; 00602 timeOffset = pillbox->timeOffset; 00603 timeCount = 0; 00604 validSaveTrack = false; 00605 state = TIMING_UNKNOWN; 00606 fieldMap = pillbox->fieldMap; 00607 00608 if(global2local.isRotated()) { 00609 rotation = new G4RotationMatrix(global2local.getRotation()); 00610 rotation->invert(); 00611 } 00612 00613 if(pillbox->frequency == 0.0 || pillbox->maxGradient == 0.0) 00614 timeOffset = 0.0; 00615 00616 // don't tune timeOffset if it is known 00617 if(!isUndetermined(timeOffset)) 00618 state = TIMING_COMPLETE; 00619 00620 // set global bounding box 00621 G4double local[4], global[4]; 00622 local[3] = 0.0; 00623 G4double dz = pillbox->innerLength/2.0 + pillbox->wallThick; 00624 for(int i=0; i<2; ++i) { 00625 local[0] = (i==0 ? -1.0 : 1.0) * pillbox->innerRadius; 00626 for(int j=0; j<2; ++j) { 00627 local[1] = (j==0 ? -1.0 : 1.0) * pillbox->innerRadius; 00628 for(int k=0; k<2; ++k) { 00629 local[2] = (k==0 ? -1.0 : 1.0) * dz; 00630 global2local.getGlobal(local,global); 00631 setGlobalPoint(global); 00632 } 00633 } 00634 } 00635 }
G4String PillboxField::getName | ( | ) | [inline] |
getName() returns the name of this placement.
References name.
Referenced by UserSteppingAction().
00078 { 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.
00638 { 00639 G4double local[4]; 00640 global2local.getLocal(local,point); 00641 00642 if(fieldMap) { 00643 G4double myField[6]; 00644 fieldMap->getFieldValue(local,myField,pillbox->maxGradient, 00645 pillbox->maxGradient); 00646 // frequency == 0.0 means a d.c. electric field 00647 G4double timeE=1.0, timeB=0.0; 00648 if(pillbox->frequency != 0) { 00649 double tArg = pillbox->omega * (point[3] - timeOffset); 00650 timeE = cos(tArg); 00651 timeB = sin(tArg); 00652 } 00653 if(rotation) { 00654 G4ThreeVector B(myField[0],myField[1],myField[2]); 00655 G4ThreeVector E(myField[3],myField[4],myField[5]); 00656 B = *rotation * B; 00657 E = *rotation * E; 00658 field[0] += B[0] * timeB; 00659 field[1] += B[1] * timeB; 00660 field[2] += B[2] * timeB; 00661 field[3] += E[0] * timeE; 00662 field[4] += E[1] * timeE; 00663 field[5] += E[2] * timeE; 00664 } else { 00665 field[0] += myField[0] * timeB; 00666 field[1] += myField[1] * timeB; 00667 field[2] += myField[2] * timeB; 00668 field[3] += myField[3] * timeE; 00669 field[4] += myField[4] * timeE; 00670 field[5] += myField[5] * timeE; 00671 } 00672 return; 00673 } 00674 00675 G4double x = local[0]; 00676 G4double y = local[1]; 00677 G4double z = local[2]; 00678 G4double r = sqrt(x*x + y*y); 00679 G4double dz = pillbox->innerLength/2.0 + pillbox->wallThick; 00680 if(z < -dz || z > dz || r > pillbox->innerRadius) return; 00681 00682 // frequency == 0.0 means a d.c. electric field 00683 if(pillbox->frequency == 0.0) { 00684 // zmin and zmax include the walls, so re-test z position 00685 if(z >= -pillbox->innerLength/2.0 && 00686 z <= pillbox->innerLength/2.0) 00687 field[5] += pillbox->maxGradient; 00688 return; 00689 } 00690 00691 // compute normal RF pillbox field 00692 double arg = pillbox->rFactor*r; 00693 if(fabs(arg) < 1.0e-20) arg = 0.0; 00694 double j0 = gsl_sf_bessel_J0(arg); 00695 double j1 = gsl_sf_bessel_J1(arg); 00696 double tArg = pillbox->omega * (point[3] - timeOffset); 00697 double ez = pillbox->maxGradient * j0 * cos(tArg); 00698 double bphi = -pillbox->maxGradient*pillbox->E2Bfactor*j1*sin(tArg); 00699 // handle skin depth if in wall or window 00700 double f = 1.0; 00701 if(z > pillbox->innerLength/2.0) 00702 f=exp(-(z-pillbox->innerLength/2.0)/pillbox->skinDepth); 00703 else if(z < -pillbox->innerLength/2.0) 00704 f=exp(-(-z-pillbox->innerLength/2.0)/pillbox->skinDepth); 00705 double phi = atan2(point[1],point[0]); 00706 if(rotation) { 00707 G4ThreeVector B(-bphi*sin(phi)*f, bphi*cos(phi)*f, 0.0); 00708 G4ThreeVector E(0.0, 0.0, ez*f); 00709 B = *rotation * B; 00710 E = *rotation * E; 00711 field[0] += B[0]; 00712 field[1] += B[1]; 00713 field[2] += B[2]; 00714 field[3] += E[0]; 00715 field[4] += E[1]; 00716 field[5] += E[2]; 00717 } else { 00718 field[0] += -bphi * sin(phi) * f; 00719 field[1] += bphi * cos(phi) * f; 00720 field[5] += ez * f; 00721 } 00722 }
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.
00726 { 00727 // do nothing if d.c. or no field or timing is complete 00728 if(pillbox->frequency <= 0.0 || pillbox->maxGradient == 0.0) 00729 state = TIMING_COMPLETE; 00730 if(state == TIMING_COMPLETE) return; 00731 00732 G4Track *track = step->GetTrack(); 00733 G4StepPoint *prePoint = step->GetPreStepPoint(); 00734 if(!prePoint) return; 00735 G4VPhysicalVolume *prePV = prePoint->GetPhysicalVolume(); 00736 G4StepPoint *postPoint = step->GetPostStepPoint(); 00737 if(!postPoint) return; 00738 G4VPhysicalVolume *postPV = postPoint->GetPhysicalVolume(); 00739 G4SteppingManager *steppingMgr = 00740 BLManager::getObject()->getSteppingManager(); 00741 00742 if(prePV == postPV) return; // neither entering nor leaving 00743 00744 if(postPV == timingPV) { // entering timingPV 00745 // save track for setting timeOffset 00746 saveTrack.CopyTrackInfo(*step->GetTrack()); 00747 saveTrack.SetUserInformation(0); 00748 validSaveTrack = true; 00749 state = SETTING_TIMING; 00750 // estimate timingOffset by neglecting acceleration 00751 double v = saveTrack.GetVelocity(); 00752 double dist = pillbox->innerLength/2.0; 00753 double arrival = saveTrack.GetGlobalTime() + dist/v; 00754 // pi/2 because E is a cos() and we want phaseAcc=0 00755 // to be the rising zero-crossing 00756 timeOffset = arrival - (pillbox->phaseAcc-pi/2.0)/pillbox->omega; 00757 } else if(prePV == timingPV) { // leaving timingPV 00758 if(++timeCount > ITERATION_LIMIT) 00759 G4Exception("pillbox","Iteration Limit",FatalException,""); 00760 // check timing offset 00761 double t = track->GetGlobalTime(); 00762 double phase = pillbox->omega * (t-timeOffset) + pi/2.0; 00763 double dt = (phase-pillbox->phaseAcc)/pillbox->omega; 00764 timeOffset = timeOffset + dt; 00765 if(fabs(dt) <= pillbox->timingTolerance) { 00766 // Success! 00767 printf("pillbox %s: Time OK timeOffset=%.4f ns, incremented to " 00768 "%.4f ns\n",getName().c_str(),timeOffset, 00769 timeOffset+pillbox->timeIncrement); 00770 timeOffset += pillbox->timeIncrement; 00771 state = TIMING_COMPLETE; 00772 return; 00773 } 00774 // restore saved data (i.e. jump back to when the track 00775 // enteered TimingVol) 00776 if(!validSaveTrack) 00777 G4Exception("pillbox","Invalid Step",FatalException,""); 00778 steppingMgr->GetSecondary()->push_back(new G4Track(saveTrack)); 00779 track->SetTrackStatus(fStopAndKill); 00780 } 00781 }
friend class BLCMDpillbox [friend] |
G4String PillboxField::name [private] |
Referenced by getName(), and PillboxField().
BLCMDpillbox* PillboxField::pillbox [private] |
Referenced by addFieldValue(), BLCMDpillbox::construct(), PillboxField(), and UserSteppingAction().
G4VPhysicalVolume* PillboxField::timingPV [private] |
Referenced by BLCMDpillbox::construct(), PillboxField(), and UserSteppingAction().
Referenced by addFieldValue(), BLCMDpillbox::construct(), and PillboxField().
G4RotationMatrix* PillboxField::rotation [private] |
Referenced by addFieldValue(), and PillboxField().
G4double PillboxField::timeOffset [private] |
Referenced by addFieldValue(), PillboxField(), and UserSteppingAction().
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().
State PillboxField::state [private] |
Referenced by PillboxField(), and UserSteppingAction().
BLFieldMap* PillboxField::fieldMap [private] |
Referenced by addFieldValue(), and PillboxField().