Public Member Functions | |
BLCMDhelicalharmonic () | |
Default constructor. Defines the command, args, etc. | |
virtual | ~BLCMDhelicalharmonic () |
Destructor. | |
BLCMDhelicalharmonic (const BLCMDhelicalharmonic &r) | |
Copy constructor. | |
BLElement * | clone () |
clone() | |
G4String | commandName () |
commandName() returns "helicalharmonic". | |
int | command (BLArgumentVector &argv, BLArgumentMap &namedArgs) |
command() implements the helicalharmonic 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 helicalharmonic magnet | |
G4double | getLength () |
getLength() returns the fieldLength of the hh | |
G4double | getWidth () |
getWidth() returns the outer radius of the hh | |
G4double | getHeight () |
getHeight() returns the outer radius of the hh | |
G4bool | isOK () |
isOK() returns true. | |
bool | isOutside (G4ThreeVector &local, G4double tolerance) |
isOutside() from BLElement. | |
void | generatePoints (int npoints, std::vector< G4ThreeVector > &v) |
generatePoints() from BLElement. | |
void | BeginOfRunAction (const G4Run *run) |
BeginOfRunAction() from BLManager::RunAction. | |
void | EndOfRunAction (const G4Run *run) |
EndOfRunAction() from BLManager::RunAction. | |
Private Attributes | |
G4double | radius |
G4double | length |
int | n |
G4double | b |
G4double | lambda |
G4double | phi0 |
BLCoordinateTransform | global2local |
G4Polymarker | markers |
Friends | |
class | HelicalHarmonicField |
BLCMDhelicalharmonic::BLCMDhelicalharmonic | ( | ) |
Default constructor. Defines the command, args, etc.
References b, BLCMDTYPE_ELEMENT, lambda, length, n, phi0, radius, BLCommand::registerCommand(), BLCommand::setDescription(), and BLCommand::setSynopsis().
Referenced by clone(), and command().
00144 : BLElement(), BLManager::RunAction() 00145 { 00146 // register the commandName(), and its synopsis and description. 00147 registerCommand(BLCMDTYPE_ELEMENT); 00148 setSynopsis("construct a helicalharmonic magnet."); 00149 setDescription( 00150 "Creates a cylindrical region containing the field of \n" 00151 "a magnetic helical harmonic of given order [n]. \n" 00152 "The field is defined by the value of the (n-1) order \n" 00153 "derivative [b] of the vertical field component (when \n" 00154 "the initial phase is 0) with respect to the horizontal \n" 00155 "coordinate at the center of the helix: \n" 00156 " b=d^(n-1)B_phi/dr^(n-1) @ [r=0 & phi-k*z+phi0=0], \n" 00157 "where k=2*pi/lambda is the helix's wave number, \n" 00158 "[lambda] is the length of the helix's period, and \n" 00159 "phi0 is the initial phase. \n" 00160 "The field components in the cylindrical frame are given by: \n" 00161 " B_phi=(2/(n*k))^(n-1)*b*(I[n-1](n*k*r)-I[n+1](n*k*r))* \n" 00162 " cos(n*(phi-k*z+phi0)), \n" 00163 " B_r =(2/(n*k))^(n-1)*b*(I[n-1](n*k*r)+I[n+1](n*k*r))* \n" 00164 " sin(n*(phi-k*z+phi0)), \n" 00165 " B_z =-2*(2/(n*k))^(n-1)*b*I[n](n*k*r)*cos(n*(phi-k*z+phi0)), \n" 00166 "where I[n](x) is the modified Bessel function of the first kind \n" 00167 "of order [n]. \n\n" 00168 "Note that this Element generates magnetic field only,\n" 00169 "and only within the cylinder defined by length and radius.\n" 00170 "So it has no solid associated with it, and is invisible.\n"); 00171 00172 // provide initial values for fields 00173 radius = 0.0; 00174 length = 0.0; 00175 n=1; 00176 b = 0.0; 00177 lambda = 1.0; 00178 phi0 = 0.0; 00179 }
BLCMDhelicalharmonic::BLCMDhelicalharmonic | ( | const BLCMDhelicalharmonic & | r | ) |
Copy constructor.
References b, lambda, length, n, phi0, and radius.
00182 : BLElement(r), 00183 BLManager::RunAction(r) 00184 { 00185 // copy fields one at a time (transfers default values from the 00186 // default object to this new object). 00187 radius = r.radius; 00188 length = r.length; 00189 n = r.n; 00190 b = r.b; 00191 lambda = r.lambda; 00192 phi0 = r.phi0; 00193 }
BLElement* BLCMDhelicalharmonic::clone | ( | ) | [inline, virtual] |
Implements BLElement.
References BLCMDhelicalharmonic().
00072 { return new BLCMDhelicalharmonic(*this); }
G4String BLCMDhelicalharmonic::commandName | ( | ) | [inline, virtual] |
int BLCMDhelicalharmonic::command | ( | BLArgumentVector & | argv, | |
BLArgumentMap & | namedArgs | |||
) | [virtual] |
command() implements the helicalharmonic command.
Implements BLCommand.
References BLCMDhelicalharmonic(), BLCommand::handleNamedArgs(), BLCommand::print(), BLCommand::printError(), and BLElement::setName().
00196 { 00197 const char *CoDE="BLCMDhelicalharmonic::command"; 00198 00199 if(argv.size() != 1) { 00200 printError("helicalharmonic: Invalid command, must have name"); 00201 return -1; 00202 } 00203 00204 if(argv[0] == "default") { 00205 return defaultHelicalHarmonic.handleNamedArgs(namedArgs); 00206 } 00207 00208 BLCMDhelicalharmonic *t = new BLCMDhelicalharmonic(defaultHelicalHarmonic); 00209 t->setName(argv[0]); 00210 t->handleNamedArgs(namedArgs); // call it twice to ensure that n is 00211 // set for units of b. 00212 int retval = t->handleNamedArgs(namedArgs); 00213 00214 t->print(argv[0]); 00215 00216 return retval; 00217 }
void BLCMDhelicalharmonic::defineNamedArgs | ( | ) | [virtual] |
defineNamedArgs() defines the named arguments for the command.
Reimplemented from BLCommand.
References BLCommand::argDouble(), BLCommand::argInt(), b, lambda, length, n, phi0, and radius.
00220 { 00221 argDouble(radius,"radius","The radius of the field region (mm)",mm); 00222 argDouble(length,"length","The length of the field region (mm)",mm); 00223 argInt(n,"n","Order of helical harmonic (i.e. n=1 for dipole)"); 00224 // see comment in command() about setting n 00225 argDouble(b,"b","(n-1)-order derivative of the field at the center " 00226 "(T/m^(n-1))",tesla/pow(meter,(n-1))); 00227 argDouble(lambda,"lambda","Helix period along the Z axis (mm).",mm); 00228 argDouble(phi0,"phi0","The phase of the XY field at the entrance (rad)."); 00229 }
void BLCMDhelicalharmonic::construct | ( | G4RotationMatrix * | relativeRotation, | |
G4ThreeVector | relativePosition, | |||
G4LogicalVolume * | parent, | |||
G4String | parentName, | |||
G4RotationMatrix * | parentRotation, | |||
G4ThreeVector | parentPosition | |||
) | [virtual] |
construct() - construct the helicalharmonic magnet
Implements BLElement.
References BLGlobalField::addElementField(), getLength(), BLElement::getName(), BLManager::getObject(), BLGlobalField::getObject(), global2local, HelicalHarmonicField, and BLManager::registerRunAction().
00238 { 00239 G4String thisname = parentName+getName(); 00240 00241 // get globalRotation and globalPosition 00242 G4RotationMatrix *globalRotation = 0; 00243 if(relativeRotation && parentRotation) { 00244 globalRotation = 00245 new G4RotationMatrix(*parentRotation * *relativeRotation); 00246 } else if(relativeRotation) { 00247 globalRotation = relativeRotation; 00248 } else if(parentRotation) { 00249 globalRotation = parentRotation; 00250 } 00251 G4ThreeVector globalPosition(relativePosition + parentPosition); 00252 if(parentRotation) 00253 globalPosition = *parentRotation * relativePosition + 00254 parentPosition; 00255 00256 global2local = BLCoordinateTransform(globalRotation,globalPosition); 00257 00258 G4double zmin = globalPosition[2]-getLength()/2.0; 00259 G4double zmax = globalPosition[2]+getLength()/2.0; 00260 00261 HelicalHarmonicField *p = new HelicalHarmonicField(global2local,this); 00262 BLGlobalField::getObject()->addElementField(p); 00263 00264 printf("BLCMDhelicalharmonic::Construct %s parent=%s relZ=%.1f globZ=%.1f\n" 00265 "\tzmin=%.1f zmax=%.1f\n", 00266 thisname.c_str(),parentName.c_str(),relativePosition[2], 00267 globalPosition[2], zmin,zmax); 00268 00269 BLManager::getObject()->registerRunAction(this,false); 00270 }
G4double BLCMDhelicalharmonic::getLength | ( | ) | [inline, virtual] |
getLength() returns the fieldLength of the hh
Implements BLElement.
References length.
Referenced by construct().
00092 { return length; }
G4double BLCMDhelicalharmonic::getWidth | ( | ) | [inline, virtual] |
getWidth() returns the outer radius of the hh
Implements BLElement.
References radius.
00095 { return radius*2.0; }
G4double BLCMDhelicalharmonic::getHeight | ( | ) | [inline, virtual] |
getHeight() returns the outer radius of the hh
Implements BLElement.
References radius.
00098 { return radius*2.0; }
G4bool BLCMDhelicalharmonic::isOK | ( | ) | [inline, virtual] |
bool BLCMDhelicalharmonic::isOutside | ( | G4ThreeVector & | local, | |
G4double | tolerance | |||
) | [inline, virtual] |
void BLCMDhelicalharmonic::generatePoints | ( | int | npoints, | |
std::vector< G4ThreeVector > & | v | |||
) | [inline, virtual] |
void BLCMDhelicalharmonic::BeginOfRunAction | ( | const G4Run * | run | ) | [virtual] |
BeginOfRunAction() from BLManager::RunAction.
Implements BLManager::RunAction.
References markers.
00273 { 00274 markers.clear(); 00275 }
void BLCMDhelicalharmonic::EndOfRunAction | ( | const G4Run * | run | ) | [virtual] |
EndOfRunAction() from BLManager::RunAction.
Implements BLManager::RunAction.
References BLCoordinateTransform::getGlobal(), global2local, lambda, length, MARKER_SIZE, markers, n, phi0, and radius.
00278 { 00279 G4VVisManager* pVVisManager = G4VVisManager::GetConcreteInstance(); 00280 if (!pVVisManager) return; 00281 00282 #ifdef STUB // omit markers 00283 double dz = lambda/10.0; 00284 int n = (int)(length/dz) + 1; 00285 for(int i=0; i<n; ++i) { 00286 G4double local[4], global[4]; 00287 G4double phi = phi0 + i * dz * 2.0*pi/lambda; 00288 local[0] = radius/2.0*cos(phi); 00289 local[1] = radius/2.0*sin(phi); 00290 local[2] = -length/2.0 + i * dz; 00291 local[3] = 0.0; 00292 global2local.getGlobal(local,global); 00293 G4ThreeVector point(global[0],global[1],global[2]); 00294 markers.push_back(point); 00295 } 00296 markers.SetMarkerType(G4Polymarker::circles); 00297 markers.SetScreenSize(MARKER_SIZE); 00298 markers.SetFillStyle(G4VMarker::filled); 00299 G4VisAttributes va(G4Colour(1.,1.,1.)); // white 00300 markers.SetVisAttributes(&va); 00301 00302 pVVisManager->Draw(markers); 00303 #endif // STUB 00304 }
friend class HelicalHarmonicField [friend] |
Referenced by construct().
G4double BLCMDhelicalharmonic::radius [private] |
G4double BLCMDhelicalharmonic::length [private] |
int BLCMDhelicalharmonic::n [private] |
Referenced by BLCMDhelicalharmonic(), defineNamedArgs(), EndOfRunAction(), and HelicalHarmonicField::HelicalHarmonicField().
G4double BLCMDhelicalharmonic::b [private] |
Referenced by BLCMDhelicalharmonic(), defineNamedArgs(), and HelicalHarmonicField::HelicalHarmonicField().
G4double BLCMDhelicalharmonic::lambda [private] |
Referenced by BLCMDhelicalharmonic(), defineNamedArgs(), EndOfRunAction(), and HelicalHarmonicField::HelicalHarmonicField().
G4double BLCMDhelicalharmonic::phi0 [private] |
Referenced by BLCMDhelicalharmonic(), defineNamedArgs(), EndOfRunAction(), and HelicalHarmonicField::HelicalHarmonicField().
Referenced by construct(), and EndOfRunAction().
G4Polymarker BLCMDhelicalharmonic::markers [private] |
Referenced by BeginOfRunAction(), and EndOfRunAction().