GenericQuadField Class Reference

Inheritance diagram for GenericQuadField:

BLElementField

List of all members.


Detailed Description

GenericQuadField represents one placement of a generic quadrupole magnet.

Public Member Functions

 GenericQuadField (BLCoordinateTransform &_global2local, BLCMDgenericquad *quad)
 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

Constructor & Destructor Documentation

GenericQuadField::GenericQuadField ( BLCoordinateTransform _global2local,
BLCMDgenericquad quad 
)

constructor.

References BLCMDgenericquad::apertureRadius, enge, ENGE_BLOCK, BLCMDgenericquad::fieldLength, BLCMDgenericquad::fieldRadius, fieldRadius, BLCMDgenericquad::fringe, FRINGE_ACCURACY, fringeDepth, BLCMDgenericquad::fringeFactor, fringeMaxZ, BLCoordinateTransform::getGlobal(), BLCommand::getList(), BLCoordinateTransform::getRotation(), BLEngeFunction::getType(), global2local, BLCMDgenericquad::gradient, gradient, halflength, BLCommand::printError(), rotation, BLEngeFunction::set(), and BLElementField::setGlobalPoint().

00378                                 : BLElementField(), rotation(), enge(ENGE_QUAD)
00379 {
00380         fieldRadius = quad->fieldRadius;
00381         halflength = quad->fieldLength/2.0;
00382         gradient = quad->gradient;
00383         global2local = _global2local;
00384         rotation = global2local.getRotation().inverse();
00385         fringeDepth = quad->fringeFactor * quad->apertureRadius * 2.0;
00386 
00387         if(quad->fringe != "") {
00388                 std::vector<G4double> v = BLCommand::getList(quad->fringe,',');
00389                 if(v.size() == 1 && v[0] == 0.0)
00390                         enge.set(0,0,0,0,0,0);
00391                 else if(v.size() == 6)
00392                         enge.set(v[0],v[1],v[2],v[3],v[4],v[5]);
00393                 else
00394                         BLCommand::printError("Invalid fringe value '%s'\n",
00395                                                 quad->fringe.c_str());
00396         }
00397 
00398         if(enge.getType() == ENGE_BLOCK) {
00399                 fringeMaxZ = halflength;
00400         } else {
00401                 for(int i=0; i<1000; ++i) {
00402                         fringeMaxZ = i*fieldRadius/10.0 + halflength;
00403                         if(enge((fringeMaxZ-halflength)/fringeDepth) <
00404                                                         FRINGE_ACCURACY)
00405                                 break;
00406                 }
00407         }
00408 
00409         // set global bounding box
00410         G4double local[4], global[4];
00411         local[3] = 0.0;
00412         for(int i=0; i<2; ++i) {
00413                 local[0] = (i==0 ? -1.0 : 1.0) * fieldRadius;
00414                 for(int j=0; j<2; ++j) {
00415                         local[1] = (j==0 ? -1.0 : 1.0) * fieldRadius;
00416                         for(int k=0; k<2; ++k) {
00417                                 local[2] = (k==0 ? -1.0 : 1.0) * fringeMaxZ;
00418                                 global2local.getGlobal(local,global);
00419                                 setGlobalPoint(global);
00420                         }
00421                 }
00422         }
00423 }


Member Function Documentation

void GenericQuadField::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, enge, fieldRadius, fringeDepth, fringeMaxZ, BLCoordinateTransform::getLocal(), global2local, gradient, halflength, BLCoordinateTransform::isRotated(), BLEngeFunction::prime(), and rotation.

00427 {
00428         G4ThreeVector global(point[0],point[1],point[2]);
00429         G4ThreeVector local;
00430 
00431         global2local.getLocal(local,global);
00432         G4double r = sqrt(local[0]*local[0]+local[1]*local[1]);
00433         if(r > fieldRadius || fabs(local[2]) > fringeMaxZ)
00434                 return;
00435 
00436         /* apply enge() to the scalar potential phi=-G0*x*y*enge(z);
00437            B is minus its gradient. Handle both edges properly. */
00438         G4double G0 = gradient;
00439         double fringeZ = (fabs(local[2])-halflength)/fringeDepth;
00440         G4double f = enge(fringeZ);
00441         G4double fp = enge.prime(fringeZ)/fringeDepth;
00442         G4ThreeVector B(G0*f*local[1],G0*f*local[0],G0*fp*local[0]*local[1]);
00443         if(local[2] < 0.0) B[2] = -B[2];
00444 
00445         if(global2local.isRotated())
00446                 B = rotation * B;
00447 
00448         field[0] += B[0];
00449         field[1] += B[1];
00450         field[2] += B[2];
00451 }


Member Data Documentation

G4double GenericQuadField::fieldRadius [private]

Referenced by addFieldValue(), and GenericQuadField().

G4double GenericQuadField::halflength [private]

Referenced by addFieldValue(), and GenericQuadField().

G4double GenericQuadField::gradient [private]

Referenced by addFieldValue(), and GenericQuadField().

Referenced by addFieldValue(), and GenericQuadField().

G4RotationMatrix GenericQuadField::rotation [private]

Referenced by addFieldValue(), and GenericQuadField().

G4double GenericQuadField::fringeMaxZ [private]

Referenced by addFieldValue(), and GenericQuadField().

G4double GenericQuadField::fringeDepth [private]

Referenced by addFieldValue(), and GenericQuadField().

Referenced by addFieldValue(), and GenericQuadField().


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