Public Member Functions | |
MultipoleField (BLCoordinateTransform &_global2local, BLCMDmultipole *mp) | |
constructor. | |
void | addFieldValue (const G4double point[4], G4double field[6]) const |
addFieldValue() adds the field for this solenoid into field[]. point[] is in global coordinates. | |
Private Attributes | |
G4double | fieldRadius |
G4double | halflength |
G4double | gradient |
BLCoordinateTransform | global2local |
G4RotationMatrix | rotation |
G4double | fringeMaxZ |
G4double | fringeDepth |
BLEngeFunction | enge |
BLCMDmultipole * | mp |
MultipoleField::MultipoleField | ( | BLCoordinateTransform & | _global2local, | |
BLCMDmultipole * | mp | |||
) |
constructor.
References BLCMDmultipole::apertureRadius, enge, BLCMDmultipole::fieldLength, fieldRadius, BLCMDmultipole::fringe, fringeDepth, BLCMDmultipole::fringeFactor, fringeMaxZ, BLCoordinateTransform::getGlobal(), BLCoordinateTransform::getRotation(), global2local, halflength, mp, BLCommand::printError(), rotation, BLEngeFunction::set(), and BLElementField::setGlobalPoint().
00324 : BLElementField(), rotation(), enge(ENGE_QUAD) 00325 { 00326 mp = _mp; 00327 fieldRadius = mp->apertureRadius; 00328 halflength = mp->fieldLength/2.0; 00329 global2local = _global2local; 00330 rotation = global2local.getRotation().inverse(); 00331 fringeDepth = mp->fringeFactor * fieldRadius * 2.0; 00332 00333 // @@@ fringe field not supported 00334 if(mp->fringe != "" && mp->fringe != "0") { 00335 BLCommand::printError("Invalid fringe value\n"); 00336 } 00337 enge.set(0,0,0,0,0,0); 00338 fringeMaxZ = halflength; 00339 00340 // set global bounding box 00341 G4double local[4], global[4]; 00342 local[3] = 0.0; 00343 for(int i=0; i<2; ++i) { 00344 local[0] = (i==0 ? -1.0 : 1.0) * fieldRadius; 00345 for(int j=0; j<2; ++j) { 00346 local[1] = (j==0 ? -1.0 : 1.0) * fieldRadius; 00347 for(int k=0; k<2; ++k) { 00348 local[2] = (k==0 ? -1.0 : 1.0) * fringeMaxZ; 00349 global2local.getGlobal(local,global); 00350 setGlobalPoint(global); 00351 } 00352 } 00353 } 00354 }
void MultipoleField::addFieldValue | ( | const G4double | point[4], | |
G4double | field[6] | |||
) | const [virtual] |
addFieldValue() adds the field for this solenoid into field[]. point[] is in global coordinates.
Implements BLElementField.
References B, BLCMDmultipole::decapole, BLCMDmultipole::dipole, BLCMDmultipole::dodecapole, fieldRadius, fringeMaxZ, BLCoordinateTransform::getLocal(), global2local, BLCoordinateTransform::isRotated(), mp, BLCMDmultipole::octopole, BLCMDmultipole::quadrupole, rotation, and BLCMDmultipole::sextupole.
00358 { 00359 G4ThreeVector global(point[0],point[1],point[2]); 00360 G4ThreeVector local; 00361 00362 global2local.getLocal(local,global); 00363 G4double r = sqrt(local[0]*local[0]+local[1]*local[1]); 00364 if(r > fieldRadius || fabs(local[2]) > fringeMaxZ) 00365 return; 00366 00367 G4ThreeVector B(0.0,0.0,0.0); 00368 G4double phi = atan2(local[1],local[0]); 00369 00370 if(mp->dipole != 0.0) { 00371 B[1] += mp->dipole; 00372 } 00373 if(mp->quadrupole != 0.0) { 00374 B[0] += mp->quadrupole * r * sin(phi); 00375 B[1] += mp->quadrupole * r * cos(phi); 00376 } 00377 if(mp->sextupole != 0.0) { 00378 B[0] += mp->sextupole * r*r * sin(2.0*phi); 00379 B[1] += mp->sextupole * r*r * cos(2.0*phi); 00380 } 00381 if(mp->octopole != 0.0) { 00382 B[0] += mp->octopole * r*r*r * sin(3.0*phi); 00383 B[1] += mp->octopole * r*r*r * cos(3.0*phi); 00384 } 00385 if(mp->decapole != 0.0) { 00386 B[0] += mp->decapole * r*r*r*r * sin(4.0*phi); 00387 B[1] += mp->decapole * r*r*r*r * cos(4.0*phi); 00388 } 00389 if(mp->dodecapole != 0.0) { 00390 B[0] += mp->dodecapole * r*r*r*r*r * sin(5.0*phi); 00391 B[1] += mp->dodecapole * r*r*r*r*r * cos(5.0*phi); 00392 } 00393 00394 if(global2local.isRotated()) 00395 B = rotation * B; 00396 00397 field[0] += B[0]; 00398 field[1] += B[1]; 00399 field[2] += B[2]; 00400 }
G4double MultipoleField::fieldRadius [private] |
Referenced by addFieldValue(), and MultipoleField().
G4double MultipoleField::halflength [private] |
Referenced by MultipoleField().
G4double MultipoleField::gradient [private] |
Referenced by addFieldValue(), and MultipoleField().
G4RotationMatrix MultipoleField::rotation [private] |
Referenced by addFieldValue(), and MultipoleField().
G4double MultipoleField::fringeMaxZ [private] |
Referenced by addFieldValue(), and MultipoleField().
G4double MultipoleField::fringeDepth [private] |
Referenced by MultipoleField().
BLEngeFunction MultipoleField::enge [private] |
Referenced by MultipoleField().
BLCMDmultipole* MultipoleField::mp [private] |
Referenced by addFieldValue(), and MultipoleField().