BLCoil Class Reference

#include <BLCoil.hh>

Inheritance diagram for BLCoil:

BLCMDcoil

List of all members.


Detailed Description

BLCoil creates and implements a coil (part of a solenoid).

A BLCoil consists only of the parameters describing the coil geometry and a field map for unit current density. It is not a BLElement, and cannot be placed geometrically -- see BLSolenoid for that.

Cylindrical symmetry around z is assumed. The field is approximated by computing the field from nSheets infinitesimally-thin current sheets interpolating from innerRadius to OuterRadius, each length long. The solution to Maxwell's Equations for each sheet is used.

The field has been verified to be correct within 0.2% for both a long solenoid and for a pair of Helmholtz coils.

Effect of nSheets on accuracy: consider the coil innerRadius=330.0*mm, outerRadius=505.0*mm, length=167.0*mm. Compute the maximum error for 1000 random points, excluding the region within 2 cm of the coils; use 99 sheets as the "correct" value. In the region 0<=r<2*outerRadius, 0<z<10*outerRadius, for 0.1% accuracy one needs 11 sheets, and for 1% one needs 4 sheets. In the region 0<=r<innerRadius/2, same z, for 0.1% accuracy one needs 5 sheets, and for 1% one needs 2 sheets. Note that inside the coils this approximation of using thin sheets is not very good; but interpolating on a grid and coordinating sheet positions with the grid is probably about as good as one can do -- at least the values will be smooth (which is not true for the bare computation -- there's a singularity at every sheet, and Bz jumps from +infinity to -infinity across the sheet).

Note that Bz does not decrease to 0.1% of its central value until z>10*(innerRadius+outerRadius)/2; it decreases to 1% at ~5 times that average radius.

BLCoil can also read a mapFile defining the field (rather than doing the above computation). For instance, it could be from opera-2d. The mapFile gives tables of Bz and Br on a grid in the z-r plane. mapFile format: lines beginning with # are comments, and are ignored. The first line contains the parameters in name=value form: z0 initial value of z (mm) dz grid interval along z (mm) nz number of grid points along z dr grid interval along r (mm) nr number of grid points along r reflectZ set nonzero to reflect in the z=0 plane (inverting Br) current the current for this field (A/mm^2) All parameters are required. If reflectZ is nonzero then z0 must be 0.0 (only z>=0 is contained in the file). The point with local coordinates (z=0,r=0) is what gets placed.

Then follow the two maps in any order. Each map begins with a line giving its name (Bz or Br) in column 1, followed by nz lines each containing nr floats (Tesla). The first value in each line is for r=0. Values are separated by one or more spaces (the first value on each line may but need not be preceeded by spaces). The main reason to use a mapFile is to accurately include the effects of iron. That implies that the map may not be accurate for values of current other than that for which the map was generated.

See BLCMDfieldmap for a more general EM field read from a map-file.

BLCoil was separated from BLSolenoid so multiple solenoids with the same geometry (coil) could share a single (large) field map for the coil -- the B field is scaled by current/currentOfMap, where currentOfMap is 1.0 for the computed field, and is the current parameter from a mapFile.

This class does not use BLFieldMap becuase it predates it, and because the map must include the parameters so construct() will know whether or not the map corresponds to the current values. You can, of course, build a solenoid out of a BLTubs and a BLFieldMap.

Public Member Functions

 BLCoil ()
 Default constructor.
 BLCoil (G4String _name, G4double _innerRadius, G4double _outerRadius, G4double _length, G4int _nSheets=0, G4double _tolerance=0.0, G4double _maxR=0.0, G4double _maxZ=0.0, G4double _dR=0.0, G4double _dZ=0.0, G4int _nR=0, G4int _nZ=0, G4String _filename="")
 Complex Constructor. All arguments that default to 0 are related to generating the field map; those that are 0 will be determined by tolerance. To construct the map the minimum info is tolerance and maxR; if either of them is 0, no map will be generated or used. NOTE: without a map, field values inside the coils will fluctuate wildly, because points on the sheets are singular; map construction is designed to minimize this. Without a map, nSheets must be given. tolerance is a fraction of the best value at r=0,z=0; values smaller than about 0.001 are unrealistic (and expensive!).
virtual ~BLCoil ()
 Destructor.
 BLCoil (const BLCoil &r)
 Copy constructor.
void printCoil ()
 printCoil() will print a description
void addField (const G4double point[4], G4double field[6], G4double currentAmpPerMm2) const
 addField() adds the field value for the given point[] into field[]. Note that point[] is relative to the center of the coil. Units for current are Amp/mm^2.
void generateFieldMap ()
 generateFieldMap() will generate the field map (for execution speedup). If just maxR and tolerance were given, all other map parameters will be automatically determined so that the estimated error is less than the specified tolerance. The field map is cached in a disk file which is automatically recreated if anything changes (may take several minutes).
void readMapFile ()
 readMapFile() will read the mapFile.
void displayMapValues (G4double r, G4double z)
 displayMapValues() will display the 4 map values around a point.

Static Public Member Functions

static BLCoilfind (G4String name)
 find() returns a pointer to the named BLCoil.

Private Member Functions

void getUnitField (G4double r, G4double z, G4double &Br, G4double &Bz)
 getUnitField() will return Br and Bz for unit current.
void determineMaxR ()
 determineMaxR() will determine maxR based on tolerance. Uses nSheets=99 to determine center field, and walks out far enough so (Bz/Bz0)<0.5*tolerance.
void determineMaxZ ()
 determineMaxZ() will determine maxZ based on tolerance. Uses nSheets=99 to determine on-axis field, and walks out far enough so (Bz/Bz0)<0.5*tolerance.
void determineNsheets ()
 determineNsheets() will determine nSheets, given maxZ and tolerance. Increases nSheets until maxError < tolerance/4 (the other 3/4 of tolerance is for mapping errors). Checks just a handful of points selected to be representative.
G4double estimateMapError ()
 estimateMapError() will estimate the maximum error in mapping, given maxZ, nSheets, dR, dZ, nR, nZ, and an allocated but uninitialized fieldMap. It generates 1000 random points and returns the maximum error; quits early if tolerance is exceeded. Error is (currentMap - best)/best(r=0,z=0), where best is computed with no map and nSheets=99. For each point it generates the 4 surrounding values of fieldMap (this is faster than generating the entire map first). Returns the maximum error of either Bz or Br, normalized to Bz(r=0,z=0,nSheets=99).
G4double computeError (G4double r, G4double z, BLCoil &best, G4bool makeMap)
 computeError() computes the error for a single point. Error is (currentMap - best)/best(r=0,z=0), where best should be an identical BLCoil with no map and nSheets=99. If makeMap is true, it generates the 4 surrounding values of the fieldMap.
void getSheetField (G4double r, G4double z, G4double &Br, G4double &Bz, G4double sheetR, G4double sheetLen, G4double currentAmpPerMm) const
 getSheetField() returns the B field from an infinitesimally-thin sheet of current. Arguments: r and z are the cylindrical coordinates of the point in question relative to center of the coil (axis along z), mm. Br and Bz return the computed B field at that point (geant4 units) sheetR is the radius of the current sheet (mm) sheetLen is the length of the current sheet (mm) currentAmpPerMm is the current of the sheet (Amp per linear mm)

Private Attributes

G4String name
G4double innerRadius
G4double outerRadius
G4double length
G4int nSheets
G4String material
G4double tolerance
G4double maxR
G4double minZ
G4double maxZ
G4double dR
G4double dZ
G4int nR
G4int nZ
G4String filename
G4String mapFile
float * mapBr
float * mapBz
G4double sheetR0
G4double sheetDR
G4double norm
G4bool goodMap
G4bool reflectZ
G4double currentOfMap

Static Private Attributes

static std::map< G4String,
BLCoil * > 
mapCoil

Friends

class BLCMDcoil
class BLCMDsolenoid
class SolenoidField


Constructor & Destructor Documentation

BLCoil::BLCoil (  ) 

Default constructor.

References currentOfMap, dR, dZ, filename, goodMap, innerRadius, length, mapBr, mapBz, mapFile, material, maxR, maxZ, minZ, name, norm, nR, nSheets, nZ, outerRadius, reflectZ, sheetDR, sheetR0, and tolerance.

00034 {
00035         // provide initial values for fields
00036         name = "";
00037         innerRadius = 0.0;
00038         outerRadius = 0.0;
00039         length = 0.0;
00040         nSheets = 0;
00041         material = "Cu";
00042         tolerance = 0.002;
00043         maxR = 0.0;
00044         minZ = 0.0;
00045         maxZ = 0.0;
00046         dR = 0.0;
00047         dZ = 0.0;
00048         nR = 0;
00049         nZ = 0;
00050         filename = "";
00051         mapFile = "";
00052         mapBr = 0;
00053         mapBz = 0;
00054         sheetR0 = 0.0;
00055         sheetDR = 0.0;
00056         norm = 0.0;
00057         goodMap = false;
00058         reflectZ = true;
00059         currentOfMap = 1.0;
00060 }

BLCoil::BLCoil ( G4String  _name,
G4double  _innerRadius,
G4double  _outerRadius,
G4double  _length,
G4int  _nSheets = 0,
G4double  _tolerance = 0.0,
G4double  _maxR = 0.0,
G4double  _maxZ = 0.0,
G4double  _dR = 0.0,
G4double  _dZ = 0.0,
G4int  _nR = 0,
G4int  _nZ = 0,
G4String  _filename = "" 
) [inline]

Complex Constructor. All arguments that default to 0 are related to generating the field map; those that are 0 will be determined by tolerance. To construct the map the minimum info is tolerance and maxR; if either of them is 0, no map will be generated or used. NOTE: without a map, field values inside the coils will fluctuate wildly, because points on the sheets are singular; map construction is designed to minimize this. Without a map, nSheets must be given. tolerance is a fraction of the best value at r=0,z=0; values smaller than about 0.001 are unrealistic (and expensive!).

References dR, dZ, filename, goodMap, innerRadius, length, mapBr, mapBz, maxR, maxZ, name, norm, nR, nSheets, nZ, outerRadius, sheetDR, sheetR0, and tolerance.

00152                                        {
00153                         name=_name; innerRadius=_innerRadius;
00154                         outerRadius=_outerRadius; length=_length;
00155                         nSheets=_nSheets; tolerance=_tolerance; maxR=_maxR;
00156                         maxZ=_maxZ; dR=_dR; dZ=_dZ; nR=_nR; nZ=_nZ; 
00157                         filename = _filename;
00158                         mapBr=0; mapBz=0; sheetR0=0.0; sheetDR=0.0; norm=0.0; 
00159                         goodMap=false; 
00160                         if(nSheets > 0) { 
00161                                 sheetDR=(outerRadius-innerRadius)/nSheets;
00162                                 sheetR0=innerRadius + sheetDR/2.0;
00163                         }
00164                         if(filename == "") filename = name + ".dat";
00165                 }

BLCoil::~BLCoil (  )  [virtual]

Destructor.

References mapBr, and mapBz.

00063 {
00064         if(mapBr) delete[] mapBr;
00065         mapBr = 0;
00066         if(mapBz) delete[] mapBz;
00067         mapBz = 0;
00068 }

BLCoil::BLCoil ( const BLCoil r  ) 

Copy constructor.

References currentOfMap, dR, dZ, filename, goodMap, innerRadius, length, mapBr, mapBz, mapFile, material, maxR, maxZ, minZ, name, norm, nR, nSheets, nZ, outerRadius, reflectZ, sheetDR, sheetR0, and tolerance.

00071 {
00072         name = "";
00073         innerRadius = r.innerRadius;
00074         outerRadius = r.outerRadius;
00075         length = r.length;
00076         nSheets = r.nSheets;
00077         material = r.material;
00078         tolerance = r.tolerance;
00079         maxR = r.maxR;
00080         minZ = r.minZ;
00081         maxZ = r.maxZ;
00082         dR = r.dR;
00083         dZ = r.dZ;
00084         nR = r.nR;
00085         nZ = r.nZ;
00086         filename = r.filename;
00087         mapFile = r.mapFile;
00088         mapBr = 0;
00089         mapBz = 0;
00090         sheetR0 = r.sheetR0;
00091         sheetDR = r.sheetDR;
00092         norm = 0.0;
00093         goodMap = false;
00094         reflectZ = r.reflectZ;
00095         currentOfMap = r.currentOfMap;
00096 }


Member Function Documentation

void BLCoil::printCoil (  ) 

printCoil() will print a description

References dR, dZ, filename, innerRadius, length, mapFile, material, maxR, maxZ, name, nSheets, outerRadius, and tolerance.

Referenced by BLCMDcoil::command().

00099 {
00100         printf("coil    %-7s ",name.c_str());
00101         printf("innerRadius=%.1f ",innerRadius);
00102         printf("outerRadius=%.1f ",outerRadius);
00103         printf("length=%.1f ",length);
00104         printf("material=%s ",material.c_str());
00105         printf("\n\t\t");
00106         if(mapFile == "") {
00107                 printf("tolerance=%.3f ",tolerance);
00108                 printf("nSheets=%d ",nSheets);
00109                 printf("\n\t\t");
00110                 printf("maxR=%.1f ",maxR);
00111                 printf("maxZ=%.1f ",maxZ);
00112                 printf("dR=%.1f ",dR);
00113                 printf("dZ=%.1f ",dZ);
00114                 printf("filename=%s ",filename.c_str());
00115                 printf("\n");
00116         } else {
00117                 printf("mapFile=%s\n",mapFile.c_str());
00118         }
00119 }

BLCoil * BLCoil::find ( G4String  name  )  [static]

find() returns a pointer to the named BLCoil.

References mapCoil.

Referenced by BLCMDsolenoid::command().

00122 {
00123         if(mapCoil.count(name) > 0)
00124                 return mapCoil[name];
00125         return 0;
00126 }

void BLCoil::addField ( const G4double  point[4],
G4double  field[6],
G4double  currentAmpPerMm2 
) const

addField() adds the field value for the given point[] into field[]. Note that point[] is relative to the center of the coil. Units for current are Amp/mm^2.

References currentOfMap, dR, dZ, getSheetField(), goodMap, innerRadius, length, mapBr, mapBz, maxR, maxZ, minZ, name, nR, nSheets, nZ, outerRadius, reflectZ, sheetDR, and sheetR0.

Referenced by SolenoidField::addFieldValue(), and getUnitField().

00130 {
00131         G4double z = point[2];
00132         G4double r = sqrt(point[0]*point[0] + point[1]*point[1]);
00133 // test for nan
00134 if(r != r) {
00135   printf("r is nan!! r=%.3f point=%.3f %.3f %.3f %.3f\n",r,point[0],point[1],point[2],point[3]);
00136 }
00137 
00138         if(goodMap) {
00139                 assert(mapBr != 0 && mapBz != 0 && dR > 0.0 && dZ > 0.0);
00140                 assert(nR > 0 && nZ > 0);
00141                 if(z < minZ+0.002*mm || z > maxZ-0.002*mm || r > maxR-0.002*mm)
00142                         return;
00143                 G4double sign = 1.0;
00144                 if(reflectZ) {
00145                         // map covers only z >= 0
00146                         if(z < 0.0) {
00147                                 sign = -sign;
00148                                 z = -z;
00149                         }
00150                 } else {
00151                         z -= minZ;
00152                 }
00153                 // 2-d linear average of the 4 surrounding points
00154                 int i = (int)(r/dR);
00155                 int j = (int)(z/dZ);
00156 if(i >= 0 && i < nR-1 && j >= 0 && j < nZ-1) {
00157   ;
00158 } else {
00159   printf("\ni=%d j=%d nR=%d nZ=%d r=%.3f z=%.3f dR=%.3f dZ=%.3f\n",i,j,nR,nZ,r,z,dR,dZ);
00160   printf("z=%.6f maxZ=%.6f r=%.6f maxR=%.5f\n",z,maxZ,r,maxR);
00161   fflush(stdout);
00162 }
00163                 assert(i >= 0 && i < nR-1 && j >= 0 && j < nZ-1);
00164                 float fr = 1.0 - (r - i*dR) / dR;
00165                 assert(fr >= 0.0 && fr <= 1.0);
00166                 float fz = 1.0 - (z - j*dZ) / dZ;
00167                 assert(fz >= 0.0 && fz <= 1.0);
00168         
00169                 G4double Bz = mapBz[j*nR+i]*fr*fz +
00170                         mapBz[j*nR+i+1]*(1.0-fr)*fz +
00171                         mapBz[(j+1)*nR+i]*fr*(1.0-fz) +
00172                         mapBz[(j+1)*nR+i+1]*(1.0-fr)*(1.0-fz);
00173                 G4double Br = mapBr[j*nR+i]*fr*fz +
00174                         mapBr[j*nR+i+1]*(1.0-fr)*fz +
00175                         mapBr[(j+1)*nR+i]*fr*(1.0-fz) +
00176                         mapBr[(j+1)*nR+i+1]*(1.0-fr)*(1.0-fz);
00177                 G4double phi = atan2(point[1], point[0]);
00178                 currentAmpPerMm2 /= currentOfMap;
00179                 field[0] += sign * currentAmpPerMm2 * Br * cos(phi);
00180                 field[1] += sign * currentAmpPerMm2 * Br * sin(phi);
00181                 field[2] += currentAmpPerMm2 * Bz;
00182                 return;
00183         }
00184 
00185         // Computation for no map
00186 if(nSheets > 0 && sheetDR > 0.0 && sheetR0 > 0.0); else printf("%s: nSheets=%d sheetR0=%.3f sheetDR=%.3f innerRadius=%.3f outerRadius=%.3f length=%.3f\n",name.c_str(),nSheets,sheetR0,sheetDR,innerRadius,outerRadius,length);
00187         assert(nSheets > 0 && sheetDR > 0.0 && sheetR0 > 0.0);
00188         G4double Br=0.0, Bz=0.0;
00189         G4double rSheet = sheetR0;
00190         G4double currentAmpPerMm = currentAmpPerMm2 * sheetDR;
00191         for(int i=0; i<nSheets; ++i) {
00192                 G4double br, bz;
00193                 // avoid singularity when r==rSheet
00194                 if(fabs(r-rSheet) < sheetDR/8.0) {
00195                     getSheetField(r+sheetDR/4.0,z,br,bz,rSheet,length,
00196                                                         currentAmpPerMm);
00197                     G4double br1, bz1;
00198                     getSheetField(r-sheetDR/4.0,z,br1,bz1,rSheet,length,
00199                                                         currentAmpPerMm);
00200                     br = (br+br1)/2.0;
00201                     bz = (bz+bz1)/2.0;
00202                 } else {
00203                     getSheetField(r,z,br,bz,rSheet,length,currentAmpPerMm);
00204                 }
00205                 Br += br;
00206                 Bz += bz;
00207                 rSheet += sheetDR;
00208         }
00209         G4double phi = atan2(point[1], point[0]);
00210         field[0] += Br * cos(phi);
00211         field[1] += Br * sin(phi);
00212         field[2] += Bz;
00213 }

void BLCoil::generateFieldMap (  ) 

generateFieldMap() will generate the field map (for execution speedup). If just maxR and tolerance were given, all other map parameters will be automatically determined so that the estimated error is less than the specified tolerance. The field map is cached in a disk file which is automatically recreated if anything changes (may take several minutes).

References currentOfMap, determineMaxR(), determineMaxZ(), determineNsheets(), dR, dZ, estimateMapError(), filename, getUnitField(), goodMap, innerRadius, length, mapBr, mapBz, mapFile, maxR, maxZ, minZ, name, norm, nR, nSheets, nZ, outerRadius, readMapFile(), reflectZ, sheetDR, sheetR0, BLTime::time(), and tolerance.

Referenced by BLCMDcoil::command().

00216 {
00217         if(mapFile != "") {
00218                 readMapFile();
00219                 return;
00220         }
00221 
00222         struct FileHeader {
00223                 long magic, nR, nZ, nSheets;
00224                 double innerRadius, outerRadius, length, dR, dZ, error;
00225                 // followed by nR*nZ  float-s of mapBr
00226                 // followed by nR*nZ  float-s of mapBz
00227         };
00228         // Had serious problems with files written on a different OS that
00229         // passed the check but generated gross tracking errors and/or
00230         // exceptions in malloc/free. Probably different compilers padding
00231         // FileHeader differently.
00232 #ifdef __linux
00233         const int MAGIC = 0x544A5232;   // 'TJR2'
00234 #endif
00235 #ifdef WIN32
00236         const int MAGIC = 0x544A5235;   // 'TJR5'
00237 #endif
00238 #ifdef __APPLE__
00239         const int MAGIC = 0x544A5234;   // 'TJR4'
00240 #endif
00241 #ifdef __CYGWIN__
00242         const int MAGIC = 0x544A5235;   // 'TJR5'
00243 #endif
00244 
00245         printf("coilmap %-7s tolerance=%.5f [fraction of "
00246                                 "Bz(r=0,z=0)]\n",name.c_str(),tolerance);
00247         fflush(stdout);
00248         goodMap = false;
00249 
00250         // read FileHeader from filename (if it exists), and check for 
00251         // consistency with this coil
00252 #ifdef WIN32
00253         FILE *in = fopen(filename.c_str(),"rb");
00254 #else
00255         FILE *in = fopen(filename.c_str(),"r");
00256 #endif
00257         FileHeader fh;
00258         if(in) {
00259                 if(fread(&fh,sizeof(fh),1,in) != 1) {
00260                         printf("BLCoil::generateFieldMap Cannot read header\n");
00261                         if(in) fclose(in);
00262                         in = 0;
00263                 }
00264                 if(fh.magic != MAGIC) {
00265                         printf("BLCoil::generateFieldMap wrong magic\n");
00266                         if(in) fclose(in);
00267                         in = 0;
00268                 }
00269                 if(fh.error > tolerance) {
00270                         printf("BLCoil::generateFieldMap error > tolerance\n");
00271                         if(in) fclose(in);
00272                         in = 0;
00273                 }
00274                 if(fabs(fh.innerRadius-innerRadius) > 0.001*mm) {
00275                         printf("BLCoil::generateFieldMap wrong innerRadius\n");
00276                         if(in) fclose(in);
00277                         in = 0;
00278                 }
00279                 if(fabs(fh.outerRadius-outerRadius) > 0.001*mm) {
00280                         printf("BLCoil::generateFieldMap wrong outerRadius\n");
00281                         if(in) fclose(in);
00282                         in = 0;
00283                 }
00284                 if(fabs(fh.length-length) > 0.001*mm) {
00285                         printf("BLCoil::generateFieldMap wrong length\n");
00286                         if(in) fclose(in);
00287                         in = 0;
00288                 }
00289                 if((dR > 0.001*mm && fabs(fh.dR-dR) > 0.001*mm)) {
00290                         printf("BLCoil::generateFieldMap wrong dR\n");
00291                         if(in) fclose(in);
00292                         in = 0;
00293                 }
00294                 if((dZ > 0.001*mm && fabs(fh.dZ-dZ) > 0.001*mm)) {
00295                         printf("BLCoil::generateFieldMap wrong dZ\n");
00296                         if(in) fclose(in);
00297                         in = 0;
00298                 }
00299                 if((fh.nR-1)*fh.dR < maxR || (fh.nZ-1)*fh.dZ < maxZ) {
00300                         printf("BLCoil::generateFieldMap wrong maxR or maxZ\n");
00301                         if(in) fclose(in);
00302                         in = 0;
00303                 }
00304         }
00305         // setup data from file and read the maps from file
00306         if(in) {
00307                 G4double orgMaxZ = maxZ;
00308                 G4double orgMaxR = maxR;
00309                 G4int orgNsheets = nSheets;
00310                 nR = fh.nR;
00311                 nZ = fh.nZ;
00312                 dR = fh.dR;
00313                 dZ = fh.dZ;
00314                 nSheets = fh.nSheets;
00315                 sheetDR=(outerRadius-innerRadius)/nSheets;
00316                 sheetR0=innerRadius + sheetDR/2.0;
00317                 norm = 0.0;
00318                 maxR = (nR-1)*dR;
00319                 maxZ = (nZ-1)*dZ;
00320                 minZ = -maxZ;
00321                 if(mapBr) delete[] mapBr;
00322                 if(mapBz) delete[] mapBz;
00323                 mapBr = new float[nZ*nR];
00324                 mapBz = new float[nZ*nR];
00325                 if(!mapBr || !mapBz) {
00326                         G4Exception("coil","Out of memory",FatalException,"");
00327                 }
00328                 if(fread(mapBr,sizeof(float),nR*nZ,in) != (size_t)nR*nZ ||
00329                    fread(mapBz,sizeof(float),nR*nZ,in) != (size_t)nR*nZ) {
00330                         fclose(in);
00331                         in = 0;
00332                         printf("BLCoil::generateFieldMap cannot read data\n");
00333                         // restore invalid data that came from file
00334                         maxZ = orgMaxZ;
00335                         minZ = -maxZ;
00336                         reflectZ = true;
00337                         maxR = orgMaxR;
00338                         nSheets = orgNsheets;
00339                 }
00340         }
00341         // if file is correct, use it.
00342         if(in) {
00343                 fclose(in);
00344                 in = 0;
00345                 goodMap = true;
00346                 printf("coilmap %-7s read file '%s'  "
00347                     "dR=%.1f dZ=%.1f\n",name.c_str(),filename.c_str(),dR,dZ);
00348                 minZ = -maxZ;
00349                 reflectZ = true;
00350                 currentOfMap = 1.0;
00351                 return;
00352         }
00353 
00354         // file is not correct - generate map
00355 
00356         if(maxR <= 0.0)
00357                 determineMaxR();
00358         if(maxZ <= 0.0)
00359                 determineMaxZ();
00360         if(nSheets <= 0) {
00361                 determineNsheets();
00362         } else {
00363                 sheetDR = (outerRadius-innerRadius) / nSheets;
00364                 sheetR0 = innerRadius + sheetDR/2.0;
00365         }
00366 
00367         // initial estimate of dR,dZ,nR,nZ
00368         if(dR <= 0.0)
00369                 dR = maxR/10;
00370         nR = (int)floor(maxR/dR) + 2;
00371         //@@ for now, assume dZ = dR....
00372         if(dZ <= 0.0)
00373                 dZ = dR;
00374         nZ = (int)floor(maxZ/dZ) + 2;
00375 
00376         // loop until we achieve the required tolerance. Limit to 8 doublings.
00377         const int MAX_ITER=8;
00378         int iter;
00379         G4double err = 1.0;
00380         for(iter=0; iter<MAX_ITER; ++iter) {
00381                 // for simplicity, discard all previous computations
00382                 if(mapBr) delete[] mapBr;
00383                 if(mapBz) delete[] mapBz;
00384                 mapBr = new float[nZ*nR];
00385                 mapBz = new float[nZ*nR];
00386                 if(!mapBr || !mapBz) {
00387                         G4Exception("coil","Out of memory",FatalException,"");
00388                 }
00389                 err = estimateMapError();
00390                 if(err < tolerance)
00391                         break;
00392                 dR /= 2.0;
00393                 dZ /= 2.0;
00394                 nR = (int)floor(maxR/dR) + 2;
00395                 nZ = (int)floor(maxZ/dZ) + 2;
00396         }
00397         if(iter >= MAX_ITER) {
00398                 G4Exception("coilmap","Failed to achieve required accuracy",
00399                                                         FatalException, "");
00400         }
00401         printf("coilmap %-7s nR=%d nZ=%d dR=%.1f dZ=%.1f "
00402                         "error=%.6f\n",name.c_str(),nR,nZ,dR,dZ,err);
00403 
00404         // fill the maps
00405         //@@ for simplicity, ignore all previous computations
00406         goodMap = false;
00407         time_t start = BLTime::time();
00408         for(int i=0; i<nR; ++i) {
00409                 G4double r = i * dR;
00410                 for(int j=0; j<nZ; ++j) {
00411                         G4double z = j * dZ;
00412                         int k = j * nR + i;
00413                         assert( k >= 0 && k < nZ*nR);
00414                         G4double Br,Bz;
00415                         getUnitField(r,z,Br,Bz);
00416                         mapBr[k] = Br;
00417                         mapBz[k] = Bz;
00418                 }
00419 //              printf("coilmap %-7s  %d%% %ld seconds\r",
00420 //                      name.c_str(),(100*(i+1))/nR,(long)(BLTime::time()-start));
00421 //              fflush(stdout);
00422         }
00423         printf("coilmap %-7s  DONE %ld seconds\n",
00424                 name.c_str(),(long)(BLTime::time()-start)); 
00425         fflush(stdout);
00426 
00427         // write the map
00428         if(filename != "") {
00429 #ifdef WIN32
00430                 FILE *out = fopen(filename.c_str(),"wb");
00431 #else
00432                 FILE *out = fopen(filename.c_str(),"w");
00433 #endif
00434                 if(out) {
00435                         FileHeader fh;
00436                         fh.magic = MAGIC;
00437                         fh.nR = nR;
00438                         fh.nZ = nZ;
00439                         fh.nSheets = nSheets;
00440                         fh.innerRadius = innerRadius;
00441                         fh.outerRadius = outerRadius;
00442                         fh.length = length;
00443                         fh.dR = dR;
00444                         fh.dZ = dZ;
00445                         fh.error = err;
00446                         if(fwrite(&fh,sizeof(fh),1,out) != 1 ||
00447                            fwrite(mapBr,sizeof(float),nR*nZ,out) != (size_t)nR*nZ ||
00448                            fwrite(mapBz,sizeof(float),nR*nZ,out) != (size_t)nR*nZ) {
00449                                 printf("coilmap %-7s ERROR "
00450                                         "writing file '%s'\n",
00451                                         name.c_str(),filename.c_str());
00452                         }
00453                         fclose(out);
00454                         printf("coilmap %-7s wrote file '%s'   %lu bytes\n",
00455                                 name.c_str(),filename.c_str(),
00456                                 (long)(sizeof(fh)+sizeof(float)*2*nR*nZ));
00457                 }
00458         }
00459 
00460         goodMap = true;
00461         minZ = -maxZ;
00462         reflectZ = true;
00463         currentOfMap = 1.0;
00464 }

void BLCoil::readMapFile (  ) 

readMapFile() will read the mapFile.

References currentOfMap, dR, dZ, goodMap, mapBr, mapBz, mapFile, maxR, maxZ, minZ, nR, nZ, BLCommand::parseArgs(), BLCommand::printError(), and reflectZ.

Referenced by generateFieldMap().

00467 {
00468         FILE *in = fopen(mapFile.c_str(),"r");
00469         if(!in) {
00470                 BLCommand::printError("Cannot open mapFile '%s'\n",mapFile.c_str());
00471                 return;
00472 error:
00473                 BLCommand::printError("Syntax error in mapFile '%s'\n",mapFile.c_str());
00474                 return;
00475         }
00476 
00477         char line[512];
00478 
00479         // read first line (parameters of the map)
00480         // units are millimeters, Tesla, and Amperes/mm^2.
00481         do {
00482                 if(!fgets(line,sizeof(line),in)) goto error;
00483         } while(line[0] == '#');
00484 
00485         BLArgumentVector argv;
00486         BLArgumentMap argm;
00487         if(BLCommand::parseArgs(line,argv,argm)) goto error;
00488 
00489         // sample parameter line:
00490         // z0=0.0 dz=2.0 nz=1501 dr=5.0 nr=13 reflectZ=1 current=-101.1
00491         G4double z0 = atof(argm["z0"].c_str())*mm;
00492         dZ = atof(argm["dz"].c_str())*mm;
00493         nZ = atoi(argm["nz"].c_str());
00494         dR = atof(argm["dr"].c_str())*mm;
00495         nR = atoi(argm["nr"].c_str());
00496         reflectZ = (atoi(argm["reflectZ"].c_str()) != 0);
00497         currentOfMap = atof(argm["current"].c_str());
00498         maxR = (nR-1)*dR;
00499         maxZ = z0 + (nZ-1)*dZ;
00500         if(reflectZ) {
00501                 minZ = -maxZ;
00502                 if(z0 != 0.0) goto error;
00503         } else {
00504                 minZ = z0;
00505         }
00506         if(argv.size() > 0 || nZ <= 0 || nR <= 0 || maxR <= 0.0 ||
00507                                         maxZ <= 0.0 || currentOfMap == 0.0)
00508                 goto error;
00509         if(mapBr) delete[] mapBr;
00510         if(mapBz) delete[] mapBz;
00511         mapBr = new float[nZ*nR];
00512         mapBz = new float[nZ*nR];
00513         if(!mapBr || !mapBz) {
00514                 G4Exception("coil","Out of memory",FatalException,"");
00515         }
00516 
00517         // read the map -- there are two maps, Bz and Br, each introduced
00518         // by a line with its name in column 1. Each map consists of nZ
00519         // lines of nR floats.
00520         int i=-1;
00521         float *map=0;
00522         while(fgets(line,sizeof(line),in)) {
00523                 if(line[0] == '#') continue;
00524                 if(strncmp(line,"Bz",2) == 0) {
00525                         if(i != -1 && i != nZ) goto error;
00526                         i = 0;
00527                         map = mapBz;
00528                         continue;
00529                 } else if(strncmp(line,"Br",2) == 0) {
00530                         if(i != -1 && i != nZ) goto error;
00531                         i = 0;
00532                         map = mapBr;
00533                         continue;
00534                 }
00535                 if(i < 0 || i >= nZ || !map) goto error;
00536                 char *p = line;
00537                 for(int j=0; j<nR; ++j) {
00538                         char *q=0;
00539                         if(*p == '\0') goto error;
00540                         *map++ = strtod(p,&q)*tesla;
00541                         if(q == 0 || p == q) goto error;
00542                         p = q;
00543                 }
00544                 ++i;
00545         }
00546         if(i != nZ) goto error;
00547         fclose(in);
00548         goodMap = true; 
00549 }

void BLCoil::displayMapValues ( G4double  r,
G4double  z 
)

displayMapValues() will display the 4 map values around a point.

References computeError(), dR, dZ, getUnitField(), goodMap, innerRadius, length, mapBr, mapBz, maxR, name, norm, nR, nZ, and outerRadius.

00594 {
00595         if(!goodMap) {
00596                 printf("BLCoil '%s' has no good map\n",name.c_str());
00597                 return;
00598         }
00599         if(z < 0.0) {
00600                 printf("z must be >= 0\n");
00601                 return;
00602         }
00603 
00604         BLCoil best("temp",innerRadius,outerRadius,length,99,0.0,maxR,maxZ);
00605         G4double Br,Bz;
00606         getUnitField(r,z,Br,Bz);
00607 
00608         int ii = (int)(r/dR);
00609         int jj = (int)(z/dZ);
00610         if(ii+1 >= nR || jj+1 >= nZ) {
00611                 printf("out of range\n");
00612                 return;
00613         }
00614         printf("r=%.3f z=%.3f Bz=%.6f T  Br=%.6f T  Berr=%.6f T\n",
00615                         r,z,Bz/tesla,Bz/tesla,
00616                         computeError(r,z,best,false)/norm/tesla);
00617         G4double bzmin=Bz, bzmax=Bz, bzavg=0.0;
00618         G4double brmin=Br, brmax=Br, bravg=0.0;
00619         for(int i=ii; i<=ii+1; ++i) {
00620             for(int j=jj; j<=jj+1; ++j) {
00621                 if(bzmin > mapBz[j*nR+i]) bzmin = mapBz[j*nR+i];
00622                 if(bzmax < mapBz[j*nR+i]) bzmax = mapBz[j*nR+i];
00623                 if(brmin > mapBr[j*nR+i]) brmin = mapBr[j*nR+i];
00624                 if(brmax < mapBr[j*nR+i]) brmax = mapBr[j*nR+i];
00625                 bzavg += mapBz[j*nR+i]/4.0;
00626                 bravg += mapBr[j*nR+i]/4.0;
00627                 printf("Grid: r=%.3f z=%.3f Bz=%.6f T  Br=%.6f T\n",
00628                         dR*i,dZ*j,mapBz[j*nR+i]/tesla,mapBr[j*nR+i]/tesla);
00629             }
00630         }
00631         printf("BzVariation=%.4f BrVariation=%.4f (fractions of avg)\n",
00632                 (bzmax-bzmin)/bzavg,(brmax-brmin)/bravg);
00633 }

void BLCoil::getUnitField ( G4double  r,
G4double  z,
G4double &  Br,
G4double &  Bz 
) [private]

getUnitField() will return Br and Bz for unit current.

References addField().

Referenced by computeError(), determineMaxR(), determineMaxZ(), displayMapValues(), and generateFieldMap().

00636 {
00637         G4double point[4];
00638         point[0] = point[1] = point[2] = point[3] = 0.0;
00639         G4double f[6];
00640         f[0]=f[1]=f[2]=f[3]=f[4]=f[5]=f[6]=0.0;
00641         point[0] = r;
00642         point[2] = z;
00643         addField(point,f,1.0);
00644         Br = f[0];
00645         Bz = f[2];
00646 }

void BLCoil::determineMaxR (  )  [private]

determineMaxR() will determine maxR based on tolerance. Uses nSheets=99 to determine center field, and walks out far enough so (Bz/Bz0)<0.5*tolerance.

References getUnitField(), innerRadius, length, maxR, name, norm, outerRadius, and tolerance.

Referenced by generateFieldMap().

00649 {
00650         G4double Br,Bz;
00651         BLCoil best("temp",innerRadius,outerRadius,length,99,0.0,
00652                                         10.0*outerRadius,100.0*outerRadius);
00653         if(norm == 0.0) {
00654                 best.getUnitField(0.0,0.0,Br,Bz);
00655                 norm = 1.0/Bz;
00656         }
00657 
00658         // step out in maxR until 0.5*tolerance is achieved
00659         for(maxR=outerRadius+2.0*innerRadius;
00660             maxR<=outerRadius+32.0*innerRadius;
00661             maxR+=innerRadius) {
00662                 best.getUnitField(maxR,0.0,Br,Bz);
00663                 Bz = fabs(Bz);
00664                 if(Bz*norm < 0.5*tolerance) break;
00665         }
00666         if(maxR <= outerRadius+32.0*innerRadius) {
00667                 printf("BLCoil::determineMaxR '%s'  maxR=%.1f err=%.6f\n",
00668                                         name.c_str(),maxR,Bz*norm);
00669         } else {
00670                 G4Exception("coilmap determineMaxR",
00671                                         "Failed to achieve required accuracy",
00672                                         FatalException, "");
00673         }
00674 }

void BLCoil::determineMaxZ (  )  [private]

determineMaxZ() will determine maxZ based on tolerance. Uses nSheets=99 to determine on-axis field, and walks out far enough so (Bz/Bz0)<0.5*tolerance.

References getUnitField(), innerRadius, length, maxR, name, norm, outerRadius, and tolerance.

Referenced by generateFieldMap().

00677 {
00678         G4double Br,Bz;
00679         BLCoil best("temp",innerRadius,outerRadius,length,99,0.0,maxR,
00680                                                         100.0*outerRadius);
00681         if(norm == 0.0) {
00682                 best.getUnitField(0.0,0.0,Br,Bz);
00683                 norm = 1.0/Bz;
00684         }
00685 
00686         // step out in maxZ until 0.5*tolerance is achieved
00687         for(maxZ=length/2.0+4.0*innerRadius; maxZ<=length/2.0+32.0*innerRadius;
00688                                                         maxZ+=innerRadius) {
00689                 best.getUnitField(0.0,maxZ,Br,Bz);
00690                 if(Bz*norm < 0.5*tolerance) break;
00691         }
00692         if(maxZ <= length/2.0+32.0*innerRadius) {
00693                 printf("BLCoil::determineMaxZ '%s'  maxZ=%.1f err=%.6f\n",
00694                                         name.c_str(),maxZ,Bz*norm);
00695         } else {
00696                 G4Exception("coilmap determineMaxZ",
00697                                         "Failed to achieve required accuracy",
00698                                         FatalException, "");
00699         }
00700 }

void BLCoil::determineNsheets (  )  [private]

determineNsheets() will determine nSheets, given maxZ and tolerance. Increases nSheets until maxError < tolerance/4 (the other 3/4 of tolerance is for mapping errors). Checks just a handful of points selected to be representative.

References computeError(), innerRadius, length, maxR, name, nSheets, outerRadius, sheetDR, sheetR0, and tolerance.

Referenced by generateFieldMap().

00703 {
00704         BLCoil best("temp",innerRadius,outerRadius,length,99,0.0,maxR,maxZ);
00705 
00706         G4double exclude = 0.1 * innerRadius;
00707         G4double err = 0.0;
00708         const int MAX_SHEETS=40;
00709         for(nSheets=2; nSheets<=MAX_SHEETS; nSheets+=2) {
00710                 sheetDR = (outerRadius-innerRadius) / nSheets;
00711                 sheetR0 = innerRadius + sheetDR/2.0;
00712                 err = 0.0;
00713                 for(int i=0; i<50; ++i) {
00714                         G4double r = innerRadius - exclude;
00715                         G4double z = G4UniformRand() * length/2.0;
00716                         if(r >= maxR) r = maxR;
00717                         G4double e = computeError(r,z,best,false);
00718                         if(e > err) err = e;
00719                         r = innerRadius + (outerRadius-innerRadius)*G4UniformRand();
00720                         if(r >= maxR) continue;
00721                         z = length/2.0 + exclude;
00722                         e = computeError(r,z,best,false);
00723                         if(e > err) err = e;
00724                         if(err > tolerance/4.0) break;
00725                 }
00726                 if(err <= tolerance/4.0) break;
00727         }
00728 
00729         if(nSheets <= MAX_SHEETS) {
00730                 printf("BLCoil::determineNsheets '%s'  nSheets=%d err=%.6f\n",
00731                                         name.c_str(),nSheets,err);
00732         } else {
00733                 G4Exception("coilmap determineNsheets",
00734                                         "Failed to achieve required accuracy",
00735                                         FatalException, "");
00736         }
00737 }

G4double BLCoil::estimateMapError (  )  [private]

estimateMapError() will estimate the maximum error in mapping, given maxZ, nSheets, dR, dZ, nR, nZ, and an allocated but uninitialized fieldMap. It generates 1000 random points and returns the maximum error; quits early if tolerance is exceeded. Error is (currentMap - best)/best(r=0,z=0), where best is computed with no map and nSheets=99. For each point it generates the 4 surrounding values of fieldMap (this is faster than generating the entire map first). Returns the maximum error of either Bz or Br, normalized to Bz(r=0,z=0,nSheets=99).

References computeError(), innerRadius, length, maxR, outerRadius, and tolerance.

Referenced by generateFieldMap().

00740 {
00741         BLCoil best("temp",innerRadius,outerRadius,length,99,0.0,maxR,maxZ);
00742 
00743         // generate 1000 random positions, "close"
00744         // to the coil in z, but exclude the coil (and nearby).
00745         G4double exclude = 0.1 * innerRadius;
00746         G4double err = 0.0;
00747         G4double genZ = length;
00748         if(length/2.0+exclude+10.0*mm > genZ) genZ = length/2.0+exclude+10.0*mm;
00749         for(int i=0; i<1000; ++i) {
00750                 G4double r, z;
00751                 do {
00752                         r = G4UniformRand() * maxR;
00753                         z = G4UniformRand() * genZ;
00754                 } while(z < length/2+exclude &&
00755                         r > innerRadius-exclude &&
00756                         r < outerRadius+exclude);
00757                 G4double e = computeError(r,z,best,true);
00758                 if(e > err) err = e;
00759                 if(err > tolerance) {
00760                         break;
00761                 }
00762 //              if(i%10 == 0)
00763 //                      printf("BLCoil::estimateMapError nR=%d nZ=%d %d%%\r",
00764 //                              nR,nZ,(100*i)/1000), fflush(stdout);
00765         }
00766 
00767         return err;
00768 }

G4double BLCoil::computeError ( G4double  r,
G4double  z,
BLCoil best,
G4bool  makeMap 
) [private]

computeError() computes the error for a single point. Error is (currentMap - best)/best(r=0,z=0), where best should be an identical BLCoil with no map and nSheets=99. If makeMap is true, it generates the 4 surrounding values of the fieldMap.

References dR, dZ, getUnitField(), goodMap, mapBr, mapBz, maxR, norm, nR, and nZ.

Referenced by determineNsheets(), displayMapValues(), and estimateMapError().

00772 {
00773         G4double Br,Bz;
00774         if(z < -maxZ || z > maxZ || r > maxR) return 0.0;
00775         G4bool originalGoodMap = goodMap;
00776 
00777         if(norm == 0.0) {
00778                 best.getUnitField(0.0,0.0,Br,Bz);
00779                 norm = 1.0/Bz;
00780         }
00781 
00782         if(makeMap) {
00783                 // generate fieldMap for the 4 surrounding grid points
00784                 assert(mapBz != 0 && mapBr != 0 && nR > 0 && nZ > 0);
00785                 assert(z >= 0.0);
00786                 goodMap = false;
00787                 int i = (int)(r/dR);
00788                 int j = (int)(z/dZ);
00789                 assert(i >= 0 && i < nR-1 && j >= 0 && j < nZ-1);
00790         
00791                 getUnitField(i*dR,j*dZ,Br,Bz);
00792                 mapBz[j*nR+i] = Bz;
00793                 mapBr[j*nR+i] = Br;
00794 
00795                 getUnitField((i+1)*dR,j*dZ,Br,Bz);
00796                 mapBz[j*nR+(i+1)] = Bz;
00797                 mapBr[j*nR+(i+1)] = Br;
00798         
00799                 getUnitField(i*dR,(j+1)*dZ,Br,Bz);
00800                 mapBz[(j+1)*nR+i] = Bz;
00801                 mapBr[(j+1)*nR+i] = Br;
00802 
00803                 getUnitField((i+1)*dR,(j+1)*dZ,Br,Bz);
00804                 mapBz[(j+1)*nR+(i+1)] = Bz;
00805                 mapBr[(j+1)*nR+(i+1)] = Br;
00806 
00807                 goodMap = true;
00808         }
00809 
00810         G4double bestBr,bestBz;
00811         best.getUnitField(r,z,bestBr,bestBz);
00812         getUnitField(r,z,Br,Bz);
00813         G4double eBr = fabs(Br-bestBr)*norm;
00814         G4double eBz = fabs(Bz-bestBz)*norm;
00815 
00816         goodMap = originalGoodMap;
00817         return (eBr>eBz ? eBr : eBz);
00818 }

void BLCoil::getSheetField ( G4double  r,
G4double  z,
G4double &  Br,
G4double &  Bz,
G4double  sheetR,
G4double  sheetLen,
G4double  currentAmpPerMm 
) const [private]

getSheetField() returns the B field from an infinitesimally-thin sheet of current. Arguments: r and z are the cylindrical coordinates of the point in question relative to center of the coil (axis along z), mm. Br and Bz return the computed B field at that point (geant4 units) sheetR is the radius of the current sheet (mm) sheetLen is the length of the current sheet (mm) currentAmpPerMm is the current of the sheet (Amp per linear mm)

Referenced by addField().

00866 {
00867   // TJR: within this routine the sheet extends from z=0 to z=sheetLen,
00868   //      externally it is centered at z=0.
00869   z += sheetLen/2.0;
00870 
00871 //printf("BLCoil::getSheetField r=%.6g z=%.6g sheetR=%.6g sheetLen=%.6g\n",r,z,sheetR,sheetLen); fflush(stdout);
00872   G4double k; 
00873 
00874   G4double rposition = r/meter;
00875   G4double zposition = z/meter;
00876 
00877   G4double rad = sheetR/meter;
00878   G4double len = sheetLen/meter;
00879   G4double cur = currentAmpPerMm;
00880 
00881   //  Left-end of sheet
00882 
00883   G4double rho = rposition;
00884   z = -zposition;
00885   
00886   G4double r1 = sqrt( pow(rad + rho,2) + z*z );
00887   k = 2.0 / r1 * sqrt( rad * rho );
00888   G4double c2 = -4.0 * rad*rho / pow(rad+rho,2);
00889 
00890   G4double f1 = z * rad / ((rad + rho) * r1);
00891   G4double f2 = (rad - rho) / (2.0 * rad);
00892 
00893 //printf("gsl_sf_ellint_E k=%.6g\n",k); fflush(stdout);
00894   G4double elle = gsl_sf_ellint_E(pi/2.0,k,GSL_PREC_DOUBLE);
00895 //printf("gsl_sf_ellint_F k=%.6g\n",k); fflush(stdout);
00896   G4double ellf = gsl_sf_ellint_F(pi/2.0,k,GSL_PREC_DOUBLE);
00897 //printf("gsl_sf_ellint_P k=%.6g c2=%.6g\n",k,c2); fflush(stdout);
00898   G4double ellpi = gsl_sf_ellint_P(pi/2.0,k,c2,GSL_PREC_DOUBLE);
00899 
00900   G4double bz1 = f1 * ( ellf + f2 * ( ellpi - ellf ) );
00901 
00902   G4double br1;
00903 
00904   if( rho != 0.0 ) {
00905     f1 = r1 / (4.*rho);
00906     br1 = f1 * ( 2. * ( ellf - elle ) - k * k * ellf );
00907   }
00908   else {
00909     br1 = 0.0;
00910   }
00911 
00912   G4double bz = -bz1;
00913   G4double br = -br1;      
00914 
00915   //  Right-end of sheet
00916       
00917   z = z + len;
00918   r1 = sqrt( pow(rad + rho,2) + z*z );
00919   k = 2.0 / r1 * sqrt( rad * rho );
00920   c2 = -4.0 * rad * rho / pow(rad + rho,2);
00921 
00922   f1 = z * rad / ( (rad + rho) * r1);
00923   f2 = (rad - rho) / (2.0 * rad);
00924 
00925 //printf("gsl_sf_ellint_E k=%.6g\n",k); fflush(stdout);
00926   elle = gsl_sf_ellint_E(pi/2.0,k,GSL_PREC_DOUBLE);
00927 //printf("gsl_sf_ellint_F k=%.6g\n",k); fflush(stdout);
00928   ellf = gsl_sf_ellint_F(pi/2.0,k,GSL_PREC_DOUBLE);
00929 //printf("gsl_sf_ellint_P k=%.6g c2=%.6g\n",k,c2); fflush(stdout);
00930   ellpi = gsl_sf_ellint_P(pi/2.0,k,c2,GSL_PREC_DOUBLE);
00931 
00932   bz1 = f1 * ( ellf + f2 * ( ellpi - ellf ) );
00933 
00934   if ( rho != 0.0 ) {
00935     f1 = r1 / (4.0 * rho);
00936     br1 = f1 * ( 2.0 * ( ellf - elle ) - k * k * ellf );
00937   }      
00938   else {
00939     br1 = 0.0;
00940   }
00941 
00942   bz = bz + bz1;       // note change of sign
00943   br = br + br1;
00944 
00945   G4double permfactor = 4.0e-7;
00946 
00947   Br = permfactor * cur * br;
00948   Bz = permfactor * cur * bz;
00949 
00950  }


Friends And Related Function Documentation

friend class BLCMDcoil [friend]

friend class BLCMDsolenoid [friend]

friend class SolenoidField [friend]


Member Data Documentation

std::map< G4String, BLCoil * > BLCoil::mapCoil [static, private]

Referenced by BLCMDcoil::command(), and find().

G4String BLCoil::name [private]

G4double BLCoil::innerRadius [private]

G4double BLCoil::outerRadius [private]

G4double BLCoil::length [private]

G4int BLCoil::nSheets [private]

G4String BLCoil::material [private]

G4double BLCoil::tolerance [private]

G4double BLCoil::maxR [private]

G4double BLCoil::minZ [private]

G4double BLCoil::maxZ [private]

G4double BLCoil::dR [private]

G4double BLCoil::dZ [private]

G4int BLCoil::nR [private]

G4int BLCoil::nZ [private]

G4String BLCoil::filename [private]

G4String BLCoil::mapFile [private]

float* BLCoil::mapBr [private]

float* BLCoil::mapBz [private]

G4double BLCoil::sheetR0 [private]

G4double BLCoil::sheetDR [private]

G4double BLCoil::norm [private]

G4bool BLCoil::goodMap [private]

G4bool BLCoil::reflectZ [private]

G4double BLCoil::currentOfMap [private]


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