SND@LHC Software
Loading...
Searching...
No Matches
ShipBFieldMap Class Reference

Class that defines a (3d) magnetic field map (distances in cm, fields in tesla) More...

#include <ShipBFieldMap.h>

Inheritance diagram for ShipBFieldMap:
Collaboration diagram for ShipBFieldMap:

Public Types

typedef std::vector< std::vector< Float_t > > floatArray
 Typedef for a vector containing a vector of floats.
 

Public Member Functions

 ShipBFieldMap (const std::string &label, const std::string &mapFileName, Float_t xOffset=0.0, Float_t yOffset=0.0, Float_t zOffset=0.0, Float_t phi=0.0, Float_t theta=0.0, Float_t psi=0.0, Float_t scale=1.0, Bool_t isSymmetric=kFALSE)
 Constructor.
 
 ShipBFieldMap (const std::string &newName, const ShipBFieldMap &rhs, Float_t newXOffset, Float_t newYOffset, Float_t newZOffset, Float_t newPhi=0.0, Float_t newTheta=0.0, Float_t newPsi=0.0, Float_t newScale=1.0)
 
virtual ~ShipBFieldMap ()
 Destructor.
 
virtual void Field (const Double_t *position, Double_t *B)
 Implementation of evaluating the B field.
 
floatArraygetFieldMap () const
 Retrieve the field map.
 
void SetXOffset (Float_t xValue)
 Set the x global co-ordinate shift.
 
void SetYOffset (Float_t yValue)
 Set the y global co-ordinate shift.
 
void SetZOffset (Float_t zValue)
 Set the z global co-ordinate shift.
 
void SetPhi (Float_t phi)
 Set the first Euler rotation angle phi about the z axis.
 
void SetTheta (Float_t theta)
 Set the second Euler rotation angle theta about the new x axis.
 
void SetPsi (Float_t psi)
 Set the third Euler rotation angle psi about the new z axis.
 
void SetScale (Float_t scale)
 Set the field magnitude scaling factor.
 
void UseSymmetry (Bool_t flag)
 Set the boolean to specify if we have quadrant symmetry.
 
std::string GetMapFileName () const
 Get the name of the map file.
 
Int_t GetNx () const
 Get the number of bins along x.
 
Int_t GetNy () const
 Get the number of bins along y.
 
Int_t GetNz () const
 Get the number of bins along z.
 
Int_t GetNBins () const
 Get the total numer of bins.
 
Float_t GetXMin () const
 Get the minimum value of x for the map.
 
Float_t GetXMax () const
 Get the maximum value of x for the map.
 
Float_t GetdX () const
 Get the bin width along x for the map.
 
Float_t GetXRange () const
 Get the x co-ordinate range for the map.
 
Float_t GetYMin () const
 Get the minimum value of y for the map.
 
Float_t GetYMax () const
 Get the maximum value of y for the map.
 
Float_t GetdY () const
 Get the bin width along y for the map.
 
Float_t GetYRange () const
 Get the y co-ordinate range for the map.
 
Float_t GetZMin () const
 Get the minimum value of z for the map.
 
Float_t GetZMax () const
 Get the maximum value of z for the map.
 
Float_t GetdZ () const
 Get the bin width along z for the map.
 
Float_t GetZRange () const
 Get the z co-ordinate range for the map.
 
Float_t GetXOffset () const
 Get the x offset co-ordinate of the map for global positioning.
 
Float_t GetYOffset () const
 Get the y offset co-ordinate of the map for global positioning.
 
Float_t GetZOffset () const
 Get the z offset co-ordinate of the map for global positioning.
 
Float_t GetPhi () const
 Get the first Euler rotation angle about the z axis for global positioning.
 
Float_t GetTheta () const
 Get the second Euler rotation angle about the new x axis for global positioning.
 
Float_t GetPsi () const
 Get the third Euler rotation angle about the new z axis for global positioning.
 
Float_t GetScale () const
 Get the field magnitude scaling factor.
 
Bool_t HasSymmetry () const
 Get the boolean flag to specify if we have quadrant symmetry.
 
Bool_t IsACopy () const
 Get the boolean flag to specify if we are a "copy".
 
 ClassDef (ShipBFieldMap, 1)
 ClassDef for ROOT.
 

Private Types

enum  CoordAxis { xAxis = 0 , yAxis , zAxis }
 Enumeration to specify the co-ordinate type. More...
 
typedef std::pair< Int_t, Float_t > binPair
 Typedef for an int-double pair.
 

Private Member Functions

 ShipBFieldMap (const ShipBFieldMap &rhs)
 Copy constructor not implemented.
 
ShipBFieldMapoperator= (const ShipBFieldMap &rhs)
 Copy assignment operator not implemented.
 
void initialise ()
 Initialisation.
 
void readMapFile ()
 Read the field map data and store the information internally.
 
void readRootFile ()
 Process the ROOT file containing the field map data.
 
void readTextFile ()
 Process the text file containing the field map data.
 
void setLimits ()
 
Bool_t insideRange (Float_t x, Float_t y, Float_t z)
 Check to see if a point is within the map validity range.
 
binPair getBinInfo (Float_t x, CoordAxis theAxis)
 Get the bin number and fractional distance from the leftmost bin edge.
 
Int_t getMapBin (Int_t iX, Int_t iY, Int_t iZ)
 Find the vector entry of the field map data given the bins iX, iY and iZ.
 
Float_t BInterCalc (CoordAxis theAxis)
 

Private Attributes

floatArrayfieldMap_
 
std::string mapFileName_
 The name of the map file.
 
Bool_t initialised_
 Initialisation boolean.
 
Bool_t isCopy_
 Flag to specify if we are a copy of the field map (just a change of global offsets)
 
Int_t Nx_
 The number of bins along x.
 
Int_t Ny_
 The number of bins along y.
 
Int_t Nz_
 The number of bins along z.
 
Int_t N_
 The total number of bins.
 
Float_t xMin_
 The minimum value of x for the map.
 
Float_t xMax_
 The maximum value of x for the map.
 
Float_t dx_
 The bin width along x.
 
Float_t xRange_
 The co-ordinate range along x.
 
Float_t yMin_
 The minimum value of y for the map.
 
Float_t yMax_
 The maximum value of y for the map.
 
Float_t dy_
 The bin width along y.
 
Float_t yRange_
 The co-ordinate range along y.
 
Float_t zMin_
 The minimum value of z for the map.
 
Float_t zMax_
 The maximum value of z for the map.
 
Float_t dz_
 The bin width along z.
 
Float_t zRange_
 The co-ordinate range along z.
 
Float_t xOffset_
 The x value of the positional offset for the map.
 
Float_t yOffset_
 The y value of the positional offset for the map.
 
Float_t zOffset_
 The z value of the positional offset for the map.
 
Float_t phi_
 The first Euler rotation angle about the z axis.
 
Float_t theta_
 The second Euler rotation angle about the new x axis.
 
Float_t psi_
 The third Euler rotation angle about the new z axis.
 
Float_t scale_
 The B field magnitude scaling factor.
 
Bool_t isSymmetric_
 The boolean specifying if we have quadrant symmetry.
 
TGeoMatrix * theTrans_
 The combined translation and rotation transformation.
 
Float_t Tesla_
 Double converting Tesla to kiloGauss (for VMC/FairRoot B field units)
 
Int_t binA_
 Bin A for the trilinear interpolation.
 
Int_t binB_
 Bin B for the trilinear interpolation.
 
Int_t binC_
 Bin C for the trilinear interpolation.
 
Int_t binD_
 Bin D for the trilinear interpolation.
 
Int_t binE_
 Bin E for the trilinear interpolation.
 
Int_t binF_
 Bin F for the trilinear interpolation.
 
Int_t binG_
 Bin G for the trilinear interpolation.
 
Int_t binH_
 Bin H for the trilinear interpolation.
 
Float_t xFrac_
 Fractional bin distance along x.
 
Float_t yFrac_
 Fractional bin distance along y.
 
Float_t zFrac_
 Fractional bin distance along z.
 
Float_t xFrac1_
 Complimentary fractional bin distance along x.
 
Float_t yFrac1_
 Complimentary fractional bin distance along y.
 
Float_t zFrac1_
 Complimentary fractional bin distance along z.
 

Detailed Description

Class that defines a (3d) magnetic field map (distances in cm, fields in tesla)

Class that defines a (3d) magnetic field map.

Author
John Back J.J.B.nosp@m.ack@.nosp@m.warwi.nosp@m.ck.a.nosp@m.c.uk

Definition at line 16 of file ShipBFieldMap.h.

Member Typedef Documentation

◆ binPair

typedef std::pair<Int_t, Float_t> ShipBFieldMap::binPair
private

Typedef for an int-double pair.

Definition at line 324 of file ShipBFieldMap.h.

◆ floatArray

typedef std::vector< std::vector<Float_t> > ShipBFieldMap::floatArray

Typedef for a vector containing a vector of floats.

Definition at line 69 of file ShipBFieldMap.h.

Member Enumeration Documentation

◆ CoordAxis

Enumeration to specify the co-ordinate type.

Enumerator
xAxis 
yAxis 
zAxis 

Definition at line 297 of file ShipBFieldMap.h.

Constructor & Destructor Documentation

◆ ShipBFieldMap() [1/3]

ShipBFieldMap::ShipBFieldMap ( const std::string &  label,
const std::string &  mapFileName,
Float_t  xOffset = 0.0,
Float_t  yOffset = 0.0,
Float_t  zOffset = 0.0,
Float_t  phi = 0.0,
Float_t  theta = 0.0,
Float_t  psi = 0.0,
Float_t  scale = 1.0,
Bool_t  isSymmetric = kFALSE 
)

Constructor.

Parameters
[in]labelA descriptive name/title/label for this field
[in]mapFileNameThe name of the field map file (distances in cm, fields in Tesla)
[in]xOffsetThe x global co-ordinate shift to position the field map (cm)
[in]yOffsetThe y global co-ordinate shift to position the field map (cm)
[in]zOffsetThe z global co-ordinate shift to position the field map (cm)
[in]phiThe first Euler rotation angle about the z axis (degrees)
[in]thetaThe second Euler rotation angle about the new x axis (degrees)
[in]psiThe third Euler rotation angle about the new z axis (degrees)
[in]scaleThe field magnitude scaling factor (default = 1.0)
[in]isSymmetricBoolean to specify if we have quadrant symmetry (default = false)

Definition at line 15 of file ShipBFieldMap.cxx.

24 :
25 TVirtualMagField(label.c_str()),
26 fieldMap_(new floatArray()),
27 mapFileName_(mapFileName),
28 initialised_(kFALSE),
29 isCopy_(kFALSE),
30 Nx_(0), Ny_(0), Nz_(0), N_(0),
31 xMin_(0.0), xMax_(0.0),
32 dx_(0.0), xRange_(0.0),
33 yMin_(0.0), yMax_(0.0),
34 dy_(0.0), yRange_(0.0),
35 zMin_(0.0), zMax_(0.0),
36 dz_(0.0), zRange_(0.0),
37 xOffset_(xOffset),
38 yOffset_(yOffset),
39 zOffset_(zOffset),
40 phi_(phi),
41 theta_(theta),
42 psi_(psi),
43 scale_(scale),
44 isSymmetric_(isSymmetric),
45 theTrans_(0),
46 Tesla_(10.0)
47{
48 this->initialise();
49}
Float_t xMax_
The maximum value of x for the map.
Float_t dy_
The bin width along y.
Float_t psi_
The third Euler rotation angle about the new z axis.
Float_t dz_
The bin width along z.
Float_t zRange_
The co-ordinate range along z.
std::vector< std::vector< Float_t > > floatArray
Typedef for a vector containing a vector of floats.
Int_t Nz_
The number of bins along z.
Float_t Tesla_
Double converting Tesla to kiloGauss (for VMC/FairRoot B field units)
Int_t N_
The total number of bins.
floatArray * fieldMap_
Float_t yMin_
The minimum value of y for the map.
Float_t zOffset_
The z value of the positional offset for the map.
std::string mapFileName_
The name of the map file.
Bool_t initialised_
Initialisation boolean.
Int_t Ny_
The number of bins along y.
Float_t theta_
The second Euler rotation angle about the new x axis.
Bool_t isCopy_
Flag to specify if we are a copy of the field map (just a change of global offsets)
Float_t xMin_
The minimum value of x for the map.
void initialise()
Initialisation.
Int_t Nx_
The number of bins along x.
Float_t xRange_
The co-ordinate range along x.
Float_t yOffset_
The y value of the positional offset for the map.
Float_t scale_
The B field magnitude scaling factor.
TGeoMatrix * theTrans_
The combined translation and rotation transformation.
Float_t yMax_
The maximum value of y for the map.
Float_t yRange_
The co-ordinate range along y.
Float_t zMax_
The maximum value of z for the map.
Float_t phi_
The first Euler rotation angle about the z axis.
Float_t zMin_
The minimum value of z for the map.
Float_t dx_
The bin width along x.
Bool_t isSymmetric_
The boolean specifying if we have quadrant symmetry.
Float_t xOffset_
The x value of the positional offset for the map.

◆ ShipBFieldMap() [2/3]

ShipBFieldMap::ShipBFieldMap ( const std::string &  newName,
const ShipBFieldMap rhs,
Float_t  newXOffset,
Float_t  newYOffset,
Float_t  newZOffset,
Float_t  newPhi = 0.0,
Float_t  newTheta = 0.0,
Float_t  newPsi = 0.0,
Float_t  newScale = 1.0 
)

Copy constructor with a new global transformation. Use this if you want to reuse the same field map information elsewhere in the geometry

Parameters
[in]rhsThe ShipBFieldMap object to be copied (retaining any symmetry)
[in]newNameThe new description or title of the field
[in]newXOffsetThe new global offset x co-ordinate (cm)
[in]newYOffsetThe new global offset y co-ordinate (cm)
[in]newZOffsetThe new global offset z co-ordinate (cm)
[in]newPhiThe first Euler rotation angle about the z axis (degrees)
[in]newThetaThe second Euler rotation angle about the new x axis (degrees)
[in]newPsiThe third Euler rotation angle about the new z axis (degrees)
[in]newScaleThe field magnitude scaling factor (default = 1.0)
Returns
a copy of the field map object "rhs", keeping the same fieldMap pointer

Definition at line 63 of file ShipBFieldMap.cxx.

65 :
66 TVirtualMagField(newName.c_str()),
69 initialised_(kFALSE),
70 isCopy_(kTRUE),
71 Nx_(rhs.Nx_),
72 Ny_(rhs.Ny_),
73 Nz_(rhs.Nz_),
74 N_(rhs.N_),
75 xMin_(rhs.xMin_),
76 xMax_(rhs.xMax_),
77 dx_(rhs.dx_),
78 xRange_(rhs.xRange_),
79 yMin_(rhs.yMin_),
80 yMax_(rhs.yMax_),
81 dy_(rhs.dy_),
82 yRange_(rhs.yRange_),
83 zMin_(rhs.zMin_),
84 zMax_(rhs.zMax_),
85 dz_(rhs.dz_),
86 zRange_(rhs.zRange_),
87 xOffset_(newXOffset),
88 yOffset_(newYOffset),
89 zOffset_(newZOffset),
90 phi_(newPhi),
91 theta_(newTheta),
92 psi_(newPsi),
93 scale_(newScale),
95 theTrans_(0),
96 Tesla_(10.0)
97{
98 // Copy constructor with new label and different global offset, which uses
99 // the same field map data (pointer) and distance units as the rhs object
100 this->initialise();
101}
std::string GetMapFileName() const
Get the name of the map file.

◆ ~ShipBFieldMap()

ShipBFieldMap::~ShipBFieldMap ( )
virtual

Destructor.

Definition at line 51 of file ShipBFieldMap.cxx.

52{
53 // Delete the internal vector storing the field map values
54 if (fieldMap_ && isCopy_ == kFALSE) {
55 delete fieldMap_; fieldMap_ = 0;
56 }
57
58 if (theTrans_) {delete theTrans_; theTrans_ = 0;}
59
60}

◆ ShipBFieldMap() [3/3]

ShipBFieldMap::ShipBFieldMap ( const ShipBFieldMap rhs)
private

Copy constructor not implemented.

Member Function Documentation

◆ BInterCalc()

Float_t ShipBFieldMap::BInterCalc ( CoordAxis  theAxis)
private

Calculate the magnetic field component using trilinear interpolation. This function uses the various "binX" integers and "uFrac" variables

Parameters
[in]theAxisThe co-ordinate axis (CoordAxis enumeration for x, y or z)
Returns
the magnetic field component for the given axis

Definition at line 486 of file ShipBFieldMap.cxx.

487{
488
489 // Find the magnetic field component along theAxis using trilinear
490 // interpolation based on the current position and neighbouring bins
491 Float_t result(0.0);
492
493 Int_t iAxis(0);
494
495 if (theAxis == CoordAxis::yAxis) {
496 iAxis = 1;
497 } else if (theAxis == CoordAxis::zAxis) {
498 iAxis = 2;
499 }
500
501 if (fieldMap_) {
502
503 // Get the field component values for the neighbouring bins
504 Float_t A = (*fieldMap_)[binA_][iAxis];
505 Float_t B = (*fieldMap_)[binB_][iAxis];
506 Float_t C = (*fieldMap_)[binC_][iAxis];
507 Float_t D = (*fieldMap_)[binD_][iAxis];
508 Float_t E = (*fieldMap_)[binE_][iAxis];
509 Float_t F = (*fieldMap_)[binF_][iAxis];
510 Float_t G = (*fieldMap_)[binG_][iAxis];
511 Float_t H = (*fieldMap_)[binH_][iAxis];
512
513 // Perform linear interpolation along x
514 Float_t F00 = A*xFrac1_ + B*xFrac_;
515 Float_t F10 = C*xFrac1_ + D*xFrac_;
516 Float_t F01 = E*xFrac1_ + F*xFrac_;
517 Float_t F11 = G*xFrac1_ + H*xFrac_;
518
519 // Linear interpolation along y
520 Float_t F0 = F00*yFrac1_ + F10*yFrac_;
521 Float_t F1 = F01*yFrac1_ + F11*yFrac_;
522
523 // Linear interpolation along z
524 result = F0*zFrac1_ + F1*zFrac_;
525
526 }
527
528 return result;
529
530}
Float_t zFrac_
Fractional bin distance along z.
Int_t binD_
Bin D for the trilinear interpolation.
Float_t zFrac1_
Complimentary fractional bin distance along z.
Int_t binE_
Bin E for the trilinear interpolation.
Float_t yFrac1_
Complimentary fractional bin distance along y.
Float_t xFrac_
Fractional bin distance along x.
Float_t yFrac_
Fractional bin distance along y.
Int_t binC_
Bin C for the trilinear interpolation.
Int_t binG_
Bin G for the trilinear interpolation.
Int_t binH_
Bin H for the trilinear interpolation.
Int_t binB_
Bin B for the trilinear interpolation.
Int_t binF_
Bin F for the trilinear interpolation.
Float_t xFrac1_
Complimentary fractional bin distance along x.
Int_t binA_
Bin A for the trilinear interpolation.

◆ ClassDef()

ShipBFieldMap::ClassDef ( ShipBFieldMap  ,
 
)

ClassDef for ROOT.

◆ Field()

void ShipBFieldMap::Field ( const Double_t *  position,
Double_t *  B 
)
virtual

Implementation of evaluating the B field.

Parameters
[in]positionThe x,y,z global co-ordinates of the point (cm)
[out]BThe x,y,z components of the magnetic field (kGauss = 0.1 tesla)

Definition at line 103 of file ShipBFieldMap.cxx.

104{
105
106 // Set the B field components given the global position co-ordinates
107
108 // Convert the global position into a local one for the volume field.
109 // Initialise the local co-ords, which will get overwritten if the
110 // co-ordinate transformation exists. For a global field, any local
111 // volume transformation is ignored
112 Double_t localCoords[3] = {position[0], position[1], position[2]};
113
114 if (theTrans_) {theTrans_->MasterToLocal(position, localCoords);}
115
116 // The local position co-ordinates
117 Float_t x = localCoords[0];
118 Float_t y = localCoords[1];
119 Float_t z = localCoords[2];
120
121 // Now check to see if we have x-y quadrant symmetry (z has no symmetry):
122 // Bx is antisymmetric in x and y, By is symmetric and Bz has no symmetry
123 // so only Bx can change sign. This can happen if either x < 0 or y < 0:
124 // 1. x > 0, y > 0: Bx = Bx
125 // 2. x > 0, y < 0: Bx = -Bx
126 // 3. x < 0, y > 0: Bx = -Bx
127 // 4. x < 0, y < 0: Bx = Bx
128
129 Float_t BxSign(1.0);
130
131 if (isSymmetric_) {
132
133 // The field map co-ordinates only contain x > 0 and y > 0, i.e. we
134 // are using x-y quadrant symmetry. If the local x or y coordinates
135 // are negative we need to change their sign and keep track of the
136 // adjusted sign of Bx which we use as a multiplication factor at the end
137 if (x < 0.0) {
138 x = -x; BxSign *= -1.0;
139 }
140
141 if (y < 0.0) {
142 y = -y; BxSign *= -1.0;
143 }
144
145 }
146
147 // Initialise the B field components to zero
148 B[0] = 0.0;
149 B[1] = 0.0;
150 B[2] = 0.0;
151
152 // First check to see if we are inside the field map range
153 Bool_t inside = this->insideRange(x, y, z);
154 if (inside == kFALSE) {return;}
155
156 // Find the neighbouring bins for the given point
157 binPair xBinInfo = this->getBinInfo(x, ShipBFieldMap::xAxis);
158 binPair yBinInfo = this->getBinInfo(y, ShipBFieldMap::yAxis);
159 binPair zBinInfo = this->getBinInfo(z, ShipBFieldMap::zAxis);
160
161 Int_t iX = xBinInfo.first;
162 Int_t iY = yBinInfo.first;
163 Int_t iZ = zBinInfo.first;
164
165 // Check that the bins are valid
166 if (iX == -1 || iY == -1 || iZ == -1) {return;}
167
168 // Get the various neighbouring bin entries
169 Int_t iX1(iX + 1);
170 Int_t iY1(iY + 1);
171 Int_t iZ1(iZ + 1);
172
173 binA_ = this->getMapBin(iX, iY, iZ);
174 binB_ = this->getMapBin(iX1, iY, iZ);
175 binC_ = this->getMapBin(iX, iY1, iZ);
176 binD_ = this->getMapBin(iX1, iY1, iZ);
177 binE_ = this->getMapBin(iX, iY, iZ1);
178 binF_ = this->getMapBin(iX1, iY, iZ1);
179 binG_ = this->getMapBin(iX, iY1, iZ1);
180 binH_ = this->getMapBin(iX1, iY1, iZ1);
181
182 // Retrieve the fractional bin distances
183 xFrac_ = xBinInfo.second;
184 yFrac_ = yBinInfo.second;
185 zFrac_ = zBinInfo.second;
186
187 // Set the complimentary fractional bin distances
188 xFrac1_ = 1.0 - xFrac_;
189 yFrac1_ = 1.0 - yFrac_;
190 zFrac1_ = 1.0 - zFrac_;
191
192 // Finally get the magnetic field components using trilinear interpolation
193 // and scale with the appropriate multiplication factor (default = 1.0)
194 B[0] = this->BInterCalc(ShipBFieldMap::xAxis)*scale_*BxSign;
197
198}
std::pair< Int_t, Float_t > binPair
Typedef for an int-double pair.
Float_t BInterCalc(CoordAxis theAxis)
binPair getBinInfo(Float_t x, CoordAxis theAxis)
Get the bin number and fractional distance from the leftmost bin edge.
Bool_t insideRange(Float_t x, Float_t y, Float_t z)
Check to see if a point is within the map validity range.
Int_t getMapBin(Int_t iX, Int_t iY, Int_t iZ)
Find the vector entry of the field map data given the bins iX, iY and iZ.

◆ getBinInfo()

ShipBFieldMap::binPair ShipBFieldMap::getBinInfo ( Float_t  x,
ShipBFieldMap::CoordAxis  theAxis 
)
private

Get the bin number and fractional distance from the leftmost bin edge.

Parameters
[in]xThe co-ordinate component of the point (cm)
[in]theAxisThe co-ordinate axis (CoordAxis enumeration for x, y or z)
Returns
the bin number and fractional distance from the leftmost bin edge as a pair

Definition at line 430 of file ShipBFieldMap.cxx.

431{
432
433 Float_t du(0.0), uMin(0.0), Nu(0);
434
435 if (theAxis == ShipBFieldMap::xAxis) {
436 du = dx_; uMin = xMin_; Nu = Nx_;
437 } else if (theAxis == ShipBFieldMap::yAxis) {
438 du = dy_; uMin = yMin_; Nu = Ny_;
439 } else if (theAxis == ShipBFieldMap::zAxis) {
440 du = dz_; uMin = zMin_; Nu = Nz_;
441 }
442
443 Int_t iBin(-1);
444 Float_t fracL(0.0);
445
446 if (du > 1e-10) {
447
448 // Get the number of fractional bin widths the point is from the first volume bin
449 Float_t dist = (u - uMin)/du;
450 // Get the integer equivalent of this distance, which is the bin number
451 iBin = static_cast<Int_t>(dist);
452 // Get the actual fractional distance of the point from the leftmost bin edge
453 fracL = (dist - iBin*1.0);
454
455 }
456
457 // Check that the bin number is valid
458 if (iBin < 0 || iBin >= Nu) {
459 iBin = -1; fracL = 0.0;
460 }
461
462 // Return the bin number and fractional distance of the point from the leftmost bin edge
463 binPair binInfo(iBin, fracL);
464
465 return binInfo;
466
467}

◆ GetdX()

Float_t ShipBFieldMap::GetdX ( ) const
inline

Get the bin width along x for the map.

Returns
the bin width along x (cm)

Definition at line 172 of file ShipBFieldMap.h.

172{return dx_;}

◆ GetdY()

Float_t ShipBFieldMap::GetdY ( ) const
inline

Get the bin width along y for the map.

Returns
the bin width along y (cm)

Definition at line 196 of file ShipBFieldMap.h.

196{return dy_;}

◆ GetdZ()

Float_t ShipBFieldMap::GetdZ ( ) const
inline

Get the bin width along z for the map.

Returns
the bin width along z (cm)

Definition at line 220 of file ShipBFieldMap.h.

220{return dz_;}

◆ getFieldMap()

floatArray * ShipBFieldMap::getFieldMap ( ) const
inline

Retrieve the field map.

Returns
the field map

Definition at line 75 of file ShipBFieldMap.h.

75{return fieldMap_;}

◆ getMapBin()

Int_t ShipBFieldMap::getMapBin ( Int_t  iX,
Int_t  iY,
Int_t  iZ 
)
private

Find the vector entry of the field map data given the bins iX, iY and iZ.

Parameters
[in]iXThe bin along the x axis
[in]iYThe bin along the y axis
[in]iZThe bin along the z axis
Returns
the index entry for the field map data vector

Definition at line 469 of file ShipBFieldMap.cxx.

470{
471
472 // Get the index of the map entry corresponding to the x,y,z bins.
473 // Remember that the map is ordered in ascending z, y, then x
474
475 Int_t index = (iX*Ny_ + iY)*Nz_ + iZ;
476 if (index < 0) {
477 index = 0;
478 } else if (index >= N_) {
479 index = N_-1;
480 }
481
482 return index;
483
484}

◆ GetMapFileName()

std::string ShipBFieldMap::GetMapFileName ( ) const
inline

Get the name of the map file.

Returns
the name of the map file as an STL string

Definition at line 130 of file ShipBFieldMap.h.

130{return mapFileName_;}

◆ GetNBins()

Int_t ShipBFieldMap::GetNBins ( ) const
inline

Get the total numer of bins.

Returns
the total number of bins

Definition at line 154 of file ShipBFieldMap.h.

154{return N_;}

◆ GetNx()

Int_t ShipBFieldMap::GetNx ( ) const
inline

Get the number of bins along x.

Returns
the number of bins along x

Definition at line 136 of file ShipBFieldMap.h.

136{return Nx_;}

◆ GetNy()

Int_t ShipBFieldMap::GetNy ( ) const
inline

Get the number of bins along y.

Returns
the number of bins along y

Definition at line 142 of file ShipBFieldMap.h.

142{return Ny_;}

◆ GetNz()

Int_t ShipBFieldMap::GetNz ( ) const
inline

Get the number of bins along z.

Returns
the number of bins along z

Definition at line 148 of file ShipBFieldMap.h.

148{return Nz_;}

◆ GetPhi()

Float_t ShipBFieldMap::GetPhi ( ) const
inline

Get the first Euler rotation angle about the z axis for global positioning.

Returns
the map's first Euler rotation angle about the z axis (degrees)

Definition at line 250 of file ShipBFieldMap.h.

250{return phi_;}

◆ GetPsi()

Float_t ShipBFieldMap::GetPsi ( ) const
inline

Get the third Euler rotation angle about the new z axis for global positioning.

Returns
the map's third Euler rotation angle about the new z axis (degrees)

Definition at line 262 of file ShipBFieldMap.h.

262{return psi_;}

◆ GetScale()

Float_t ShipBFieldMap::GetScale ( ) const
inline

Get the field magnitude scaling factor.

Returns
the scaling factor for the field magnitude

Definition at line 268 of file ShipBFieldMap.h.

268{return scale_;}

◆ GetTheta()

Float_t ShipBFieldMap::GetTheta ( ) const
inline

Get the second Euler rotation angle about the new x axis for global positioning.

Returns
the map's second Euler rotation angle about the new x axis (degrees)

Definition at line 256 of file ShipBFieldMap.h.

256{return theta_;}

◆ GetXMax()

Float_t ShipBFieldMap::GetXMax ( ) const
inline

Get the maximum value of x for the map.

Returns
the maximum x co-ordinate (cm)

Definition at line 166 of file ShipBFieldMap.h.

166{return xMax_;}

◆ GetXMin()

Float_t ShipBFieldMap::GetXMin ( ) const
inline

Get the minimum value of x for the map.

Returns
the minimum x co-ordinate (cm)

Definition at line 160 of file ShipBFieldMap.h.

160{return xMin_;}

◆ GetXOffset()

Float_t ShipBFieldMap::GetXOffset ( ) const
inline

Get the x offset co-ordinate of the map for global positioning.

Returns
the map's x offset co-ordinate for global positioning (cm)

Definition at line 232 of file ShipBFieldMap.h.

232{return xOffset_;}

◆ GetXRange()

Float_t ShipBFieldMap::GetXRange ( ) const
inline

Get the x co-ordinate range for the map.

Returns
the x co-ordinate range (cm)

Definition at line 178 of file ShipBFieldMap.h.

178{return xRange_;}

◆ GetYMax()

Float_t ShipBFieldMap::GetYMax ( ) const
inline

Get the maximum value of y for the map.

Returns
the maximum y co-ordinate (cm)

Definition at line 190 of file ShipBFieldMap.h.

190{return yMax_;}

◆ GetYMin()

Float_t ShipBFieldMap::GetYMin ( ) const
inline

Get the minimum value of y for the map.

Returns
the minimum y co-ordinate (cm)

Definition at line 184 of file ShipBFieldMap.h.

184{return yMin_;}

◆ GetYOffset()

Float_t ShipBFieldMap::GetYOffset ( ) const
inline

Get the y offset co-ordinate of the map for global positioning.

Returns
the map's y offset co-ordinate for global positioning (cm)

Definition at line 238 of file ShipBFieldMap.h.

238{return yOffset_;}

◆ GetYRange()

Float_t ShipBFieldMap::GetYRange ( ) const
inline

Get the y co-ordinate range for the map.

Returns
the y co-ordinate range (cm)

Definition at line 202 of file ShipBFieldMap.h.

202{return yRange_;}

◆ GetZMax()

Float_t ShipBFieldMap::GetZMax ( ) const
inline

Get the maximum value of z for the map.

Returns
the maximum z co-ordinate (cm)

Definition at line 214 of file ShipBFieldMap.h.

214{return zMax_;}

◆ GetZMin()

Float_t ShipBFieldMap::GetZMin ( ) const
inline

Get the minimum value of z for the map.

Returns
the minimum z co-ordinate (cm)

Definition at line 208 of file ShipBFieldMap.h.

208{return zMin_;}

◆ GetZOffset()

Float_t ShipBFieldMap::GetZOffset ( ) const
inline

Get the z offset co-ordinate of the map for global positioning.

Returns
the map's z offset co-ordinate for global positioning (cm)

Definition at line 244 of file ShipBFieldMap.h.

244{return zOffset_;}

◆ GetZRange()

Float_t ShipBFieldMap::GetZRange ( ) const
inline

Get the z co-ordinate range for the map.

Returns
the z co-ordinate range (cm)

Definition at line 226 of file ShipBFieldMap.h.

226{return zRange_;}

◆ HasSymmetry()

Bool_t ShipBFieldMap::HasSymmetry ( ) const
inline

Get the boolean flag to specify if we have quadrant symmetry.

Returns
the boolean specifying if we have quadrant symmetry

Definition at line 274 of file ShipBFieldMap.h.

274{return isSymmetric_;}

◆ initialise()

void ShipBFieldMap::initialise ( )
private

Initialisation.

Definition at line 200 of file ShipBFieldMap.cxx.

201{
202
203 if (initialised_ == kFALSE) {
204
205 if (isCopy_ == kFALSE) {this->readMapFile();}
206
207 // Set the global co-ordinate translation and rotation info
208 if (fabs(phi_) > 1e-6 || fabs(theta_) > 1e-6 || fabs(psi_) > 1e-6) {
209
210 // We have non-zero rotation angles. Create a combined translation and rotation
211 TGeoTranslation tr("offsets", xOffset_, yOffset_, zOffset_);
212 TGeoRotation rot("angles", phi_, theta_, psi_);
213 theTrans_ = new TGeoCombiTrans(tr, rot);
214
215 } else {
216
217 // We only need a translation
218 theTrans_ = new TGeoTranslation("offsets", xOffset_, yOffset_, zOffset_);
219
220 }
221
222 initialised_ = kTRUE;
223
224 }
225
226}
void readMapFile()
Read the field map data and store the information internally.

◆ insideRange()

Bool_t ShipBFieldMap::insideRange ( Float_t  x,
Float_t  y,
Float_t  z 
)
private

Check to see if a point is within the map validity range.

Parameters
[in]xThe x co-ordinate of the point (cm)
[in]yThe y co-ordinate of the point (cm)
[in]zThe z co-ordinate of the point (cm)
Returns
true/false if the point is inside the field map range

Definition at line 381 of file ShipBFieldMap.cxx.

382{
383
384 Bool_t inside(kFALSE);
385
386 if (x >= xMin_ && x <= xMax_ && y >= yMin_ &&
387 y <= yMax_ && z >= zMin_ && z <= zMax_) {inside = kTRUE;}
388
389 return inside;
390
391}

◆ IsACopy()

Bool_t ShipBFieldMap::IsACopy ( ) const
inline

Get the boolean flag to specify if we are a "copy".

Returns
the boolean copy flag status

Definition at line 280 of file ShipBFieldMap.h.

280{return isCopy_;}

◆ operator=()

ShipBFieldMap & ShipBFieldMap::operator= ( const ShipBFieldMap rhs)
private

Copy assignment operator not implemented.

◆ readMapFile()

void ShipBFieldMap::readMapFile ( )
private

Read the field map data and store the information internally.

Definition at line 228 of file ShipBFieldMap.cxx.

229{
230
231 std::cout<<"ShipBFieldMap::readMapFile() creating field "<<this->GetName()
232 <<" using file "<<mapFileName_<<std::endl;
233
234 // Check to see if we have a ROOT file
235 if (mapFileName_.find(".root") != std::string::npos) {
236
237 this->readRootFile();
238
239 } else {
240
241 this->readTextFile();
242
243 }
244
245}
void readTextFile()
Process the text file containing the field map data.
void readRootFile()
Process the ROOT file containing the field map data.

◆ readRootFile()

void ShipBFieldMap::readRootFile ( )
private

Process the ROOT file containing the field map data.

Definition at line 247 of file ShipBFieldMap.cxx.

247 {
248
249 TFile* theFile = TFile::Open(mapFileName_.c_str());
250
251 if (!theFile) {
252 std::cout<<"ShipBFieldMap: could not find the file "<<mapFileName_<<std::endl;
253 return;
254 }
255
256 // Coordinate ranges
257 TTree* rTree = dynamic_cast<TTree*>(theFile->Get("Range"));
258 if (!rTree) {
259 std::cout<<"ShipBFieldMap: could not find Range tree in "<<mapFileName_<<std::endl;
260 return;
261 }
262
263 rTree->SetBranchAddress("xMin", &xMin_);
264 rTree->SetBranchAddress("xMax", &xMax_);
265 rTree->SetBranchAddress("dx", &dx_);
266 rTree->SetBranchAddress("yMin", &yMin_);
267 rTree->SetBranchAddress("yMax", &yMax_);
268 rTree->SetBranchAddress("dy", &dy_);
269 rTree->SetBranchAddress("zMin", &zMin_);
270 rTree->SetBranchAddress("zMax", &zMax_);
271 rTree->SetBranchAddress("dz", &dz_);
272
273 // Fill the ranges
274 rTree->GetEntry(0);
275
276 this->setLimits();
277
278 // Make sure we don't have a copy
279 if (isCopy_ == kFALSE) {
280
281 // The data is expected to contain Bx,By,Bz data values
282 // in ascending z,y,x co-ordinate order
283
284 TTree* dTree = dynamic_cast<TTree*>(theFile->Get("Data"));
285 if (!dTree) {
286 std::cout<<"ShipBFieldMap: could not find Data tree in "<<mapFileName_<<std::endl;
287 return;
288 }
289
290 Float_t Bx, By, Bz;
291 // Only enable the field components
292 dTree->SetBranchStatus("*", 0);
293 dTree->SetBranchStatus("Bx", 1);
294 dTree->SetBranchStatus("By", 1);
295 dTree->SetBranchStatus("Bz", 1);
296
297 dTree->SetBranchAddress("Bx", &Bx);
298 dTree->SetBranchAddress("By", &By);
299 dTree->SetBranchAddress("Bz", &Bz);
300
301 Int_t nEntries = dTree->GetEntries();
302 if (nEntries != N_) {
303 std::cout<<"Expected "<<N_<<" field map entries but found "<<nEntries<<std::endl;
304 nEntries = 0;
305 }
306
307 fieldMap_->reserve(nEntries);
308
309 for (Int_t i = 0; i < nEntries; i++) {
310
311 dTree->GetEntry(i);
312
313 // B field values are in Tesla. This means these values are multiplied
314 // by a factor of ten since both FairRoot and the VMC interface use kGauss
315 Bx *= Tesla_;
316 By *= Tesla_;
317 Bz *= Tesla_;
318
319 // Store the B field 3-vector
320 std::vector<Float_t> BVector(3);
321 BVector[0] = Bx; BVector[1] = By; BVector[2] = Bz;
322 fieldMap_->push_back(BVector);
323
324 }
325
326 }
327
328 theFile->Close();
329
330}
int i
Definition ShipAna.py:86

◆ readTextFile()

void ShipBFieldMap::readTextFile ( )
private

Process the text file containing the field map data.

Definition at line 332 of file ShipBFieldMap.cxx.

332 {
333
334 std::ifstream getData(mapFileName_.c_str());
335
336 std::string tmpString("");
337
338 getData >> tmpString >> xMin_ >> xMax_ >> dx_
339 >> yMin_ >> yMax_ >> dy_ >> zMin_ >> zMax_ >> dz_;
340
341 this->setLimits();
342
343 // Check to see if this object is a "copy"
344 if (isCopy_ == kFALSE) {
345
346 // Read the field map and store the magnetic field components
347
348 // Second line can be ignored since they are
349 // just labels for the data columns for readability
350 getData >> tmpString >> tmpString >> tmpString;
351
352 // The remaining lines contain Bx,By,Bz data values
353 // in ascending z,y,x co-ord order
354 fieldMap_->reserve(N_);
355
356 Float_t Bx(0.0), By(0.0), Bz(0.0);
357
358 for (Int_t i = 0; i < N_; i++) {
359
360 getData >> Bx >> By >> Bz;
361
362 // B field values are in Tesla. This means these values are multiplied
363 // by a factor of ten since both FairRoot and the VMC interface use kGauss
364 Bx *= Tesla_;
365 By *= Tesla_;
366 Bz *= Tesla_;
367
368 // Store the B field 3-vector
369 std::vector<Float_t> BVector(3);
370 BVector[0] = Bx; BVector[1] = By; BVector[2] = Bz;
371 fieldMap_->push_back(BVector);
372
373 }
374
375 }
376
377 getData.close();
378
379}

◆ setLimits()

void ShipBFieldMap::setLimits ( )
private

Definition at line 393 of file ShipBFieldMap.cxx.

393 {
394
395 // Since the default SHIP distance unit is cm, we do not need to convert
396 // these map limits, i.e. cm = 1 already
397
398 xRange_ = xMax_ - xMin_;
399 yRange_ = yMax_ - yMin_;
400 zRange_ = zMax_ - zMin_;
401
402 // Find the number of bins using the limits and bin sizes. The number of bins
403 // includes both the minimum and maximum values. To ensure correct rounding
404 // up to the nearest integer we need to add 1.5 not 1.0.
405 if (dx_ > 0.0) {
406 Nx_ = static_cast<Int_t>(((xMax_ - xMin_)/dx_) + 1.5);
407 }
408 if (dy_ > 0.0) {
409 Ny_ = static_cast<Int_t>(((yMax_ - yMin_)/dy_) + 1.5);
410 }
411 if (dz_ > 0.0) {
412 Nz_ = static_cast<Int_t>(((zMax_ - zMin_)/dz_) + 1.5);
413 }
414
415 N_ = Nx_*Ny_*Nz_;
416
417 std::cout<<"x limits: "<<xMin_<<", "<<xMax_<<", dx = "<<dx_<<std::endl;
418 std::cout<<"y limits: "<<yMin_<<", "<<yMax_<<", dy = "<<dy_<<std::endl;
419 std::cout<<"z limits: "<<zMin_<<", "<<zMax_<<", dz = "<<dz_<<std::endl;
420
421 std::cout<<"Offsets: x = "<<xOffset_<<", y = "<<yOffset_<<", z = "<<zOffset_<<std::endl;
422 std::cout<<"Angles : phi = "<<phi_<<", theta = "<<theta_<<", psi = "<<psi_<<std::endl;
423
424 std::cout<<"Total number of bins = "<<N_
425 <<"; Nx = "<<Nx_<<", Ny = "<<Ny_<<", Nz = "<<Nz_<<std::endl;
426
427}

◆ SetPhi()

void ShipBFieldMap::SetPhi ( Float_t  phi)
inline

Set the first Euler rotation angle phi about the z axis.

Parameters
[in]phiThe first Euler rotation angle about the z axis (degrees)

Definition at line 99 of file ShipBFieldMap.h.

◆ SetPsi()

void ShipBFieldMap::SetPsi ( Float_t  psi)
inline

Set the third Euler rotation angle psi about the new z axis.

Parameters
[in]psiThe third Euler rotation angle about the new z axis (degrees)

Definition at line 111 of file ShipBFieldMap.h.

111{psi_ = psi;}

◆ SetScale()

void ShipBFieldMap::SetScale ( Float_t  scale)
inline

Set the field magnitude scaling factor.

Parameters
[in]scaleThe scaling factor for the field magnitude

Definition at line 117 of file ShipBFieldMap.h.

◆ SetTheta()

void ShipBFieldMap::SetTheta ( Float_t  theta)
inline

Set the second Euler rotation angle theta about the new x axis.

Parameters
[in]thetaThe second Euler rotation angle about the new x axis (degrees)

Definition at line 105 of file ShipBFieldMap.h.

◆ SetXOffset()

void ShipBFieldMap::SetXOffset ( Float_t  xValue)
inline

Set the x global co-ordinate shift.

Parameters
[in]xValueThe value of the x global co-ordinate shift (cm)

Definition at line 81 of file ShipBFieldMap.h.

81{xOffset_ = xValue;}

◆ SetYOffset()

void ShipBFieldMap::SetYOffset ( Float_t  yValue)
inline

Set the y global co-ordinate shift.

Parameters
[in]yValueThe value of the y global co-ordinate shift (cm)

Definition at line 87 of file ShipBFieldMap.h.

87{yOffset_ = yValue;}

◆ SetZOffset()

void ShipBFieldMap::SetZOffset ( Float_t  zValue)
inline

Set the z global co-ordinate shift.

Parameters
[in]zValueThe value of the z global co-ordinate shift (cm)

Definition at line 93 of file ShipBFieldMap.h.

93{zOffset_ = zValue;}

◆ UseSymmetry()

void ShipBFieldMap::UseSymmetry ( Bool_t  flag)
inline

Set the boolean to specify if we have quadrant symmetry.

Parameters
[in]isSymmetricBoolean to specify if we have quadrant symmetry

Definition at line 123 of file ShipBFieldMap.h.

123{isSymmetric_ = flag;}

Member Data Documentation

◆ binA_

Int_t ShipBFieldMap::binA_
private

Bin A for the trilinear interpolation.

Definition at line 443 of file ShipBFieldMap.h.

◆ binB_

Int_t ShipBFieldMap::binB_
private

Bin B for the trilinear interpolation.

Definition at line 446 of file ShipBFieldMap.h.

◆ binC_

Int_t ShipBFieldMap::binC_
private

Bin C for the trilinear interpolation.

Definition at line 449 of file ShipBFieldMap.h.

◆ binD_

Int_t ShipBFieldMap::binD_
private

Bin D for the trilinear interpolation.

Definition at line 452 of file ShipBFieldMap.h.

◆ binE_

Int_t ShipBFieldMap::binE_
private

Bin E for the trilinear interpolation.

Definition at line 455 of file ShipBFieldMap.h.

◆ binF_

Int_t ShipBFieldMap::binF_
private

Bin F for the trilinear interpolation.

Definition at line 458 of file ShipBFieldMap.h.

◆ binG_

Int_t ShipBFieldMap::binG_
private

Bin G for the trilinear interpolation.

Definition at line 461 of file ShipBFieldMap.h.

◆ binH_

Int_t ShipBFieldMap::binH_
private

Bin H for the trilinear interpolation.

Definition at line 464 of file ShipBFieldMap.h.

◆ dx_

Float_t ShipBFieldMap::dx_
private

The bin width along x.

Definition at line 383 of file ShipBFieldMap.h.

◆ dy_

Float_t ShipBFieldMap::dy_
private

The bin width along y.

Definition at line 395 of file ShipBFieldMap.h.

◆ dz_

Float_t ShipBFieldMap::dz_
private

The bin width along z.

Definition at line 407 of file ShipBFieldMap.h.

◆ fieldMap_

floatArray* ShipBFieldMap::fieldMap_
private

Store the field map information as a vector of 3 floats. Map data ordering is given by first incrementing z, then y, then x

Definition at line 353 of file ShipBFieldMap.h.

◆ initialised_

Bool_t ShipBFieldMap::initialised_
private

Initialisation boolean.

Definition at line 359 of file ShipBFieldMap.h.

◆ isCopy_

Bool_t ShipBFieldMap::isCopy_
private

Flag to specify if we are a copy of the field map (just a change of global offsets)

Definition at line 362 of file ShipBFieldMap.h.

◆ isSymmetric_

Bool_t ShipBFieldMap::isSymmetric_
private

The boolean specifying if we have quadrant symmetry.

Definition at line 434 of file ShipBFieldMap.h.

◆ mapFileName_

std::string ShipBFieldMap::mapFileName_
private

The name of the map file.

Definition at line 356 of file ShipBFieldMap.h.

◆ N_

Int_t ShipBFieldMap::N_
private

The total number of bins.

Definition at line 374 of file ShipBFieldMap.h.

◆ Nx_

Int_t ShipBFieldMap::Nx_
private

The number of bins along x.

Definition at line 365 of file ShipBFieldMap.h.

◆ Ny_

Int_t ShipBFieldMap::Ny_
private

The number of bins along y.

Definition at line 368 of file ShipBFieldMap.h.

◆ Nz_

Int_t ShipBFieldMap::Nz_
private

The number of bins along z.

Definition at line 371 of file ShipBFieldMap.h.

◆ phi_

Float_t ShipBFieldMap::phi_
private

The first Euler rotation angle about the z axis.

Definition at line 422 of file ShipBFieldMap.h.

◆ psi_

Float_t ShipBFieldMap::psi_
private

The third Euler rotation angle about the new z axis.

Definition at line 428 of file ShipBFieldMap.h.

◆ scale_

Float_t ShipBFieldMap::scale_
private

The B field magnitude scaling factor.

Definition at line 431 of file ShipBFieldMap.h.

◆ Tesla_

Float_t ShipBFieldMap::Tesla_
private

Double converting Tesla to kiloGauss (for VMC/FairRoot B field units)

Definition at line 440 of file ShipBFieldMap.h.

◆ theta_

Float_t ShipBFieldMap::theta_
private

The second Euler rotation angle about the new x axis.

Definition at line 425 of file ShipBFieldMap.h.

◆ theTrans_

TGeoMatrix* ShipBFieldMap::theTrans_
private

The combined translation and rotation transformation.

Definition at line 437 of file ShipBFieldMap.h.

◆ xFrac1_

Float_t ShipBFieldMap::xFrac1_
private

Complimentary fractional bin distance along x.

Definition at line 476 of file ShipBFieldMap.h.

◆ xFrac_

Float_t ShipBFieldMap::xFrac_
private

Fractional bin distance along x.

Definition at line 467 of file ShipBFieldMap.h.

◆ xMax_

Float_t ShipBFieldMap::xMax_
private

The maximum value of x for the map.

Definition at line 380 of file ShipBFieldMap.h.

◆ xMin_

Float_t ShipBFieldMap::xMin_
private

The minimum value of x for the map.

Definition at line 377 of file ShipBFieldMap.h.

◆ xOffset_

Float_t ShipBFieldMap::xOffset_
private

The x value of the positional offset for the map.

Definition at line 413 of file ShipBFieldMap.h.

◆ xRange_

Float_t ShipBFieldMap::xRange_
private

The co-ordinate range along x.

Definition at line 386 of file ShipBFieldMap.h.

◆ yFrac1_

Float_t ShipBFieldMap::yFrac1_
private

Complimentary fractional bin distance along y.

Definition at line 479 of file ShipBFieldMap.h.

◆ yFrac_

Float_t ShipBFieldMap::yFrac_
private

Fractional bin distance along y.

Definition at line 470 of file ShipBFieldMap.h.

◆ yMax_

Float_t ShipBFieldMap::yMax_
private

The maximum value of y for the map.

Definition at line 392 of file ShipBFieldMap.h.

◆ yMin_

Float_t ShipBFieldMap::yMin_
private

The minimum value of y for the map.

Definition at line 389 of file ShipBFieldMap.h.

◆ yOffset_

Float_t ShipBFieldMap::yOffset_
private

The y value of the positional offset for the map.

Definition at line 416 of file ShipBFieldMap.h.

◆ yRange_

Float_t ShipBFieldMap::yRange_
private

The co-ordinate range along y.

Definition at line 398 of file ShipBFieldMap.h.

◆ zFrac1_

Float_t ShipBFieldMap::zFrac1_
private

Complimentary fractional bin distance along z.

Definition at line 482 of file ShipBFieldMap.h.

◆ zFrac_

Float_t ShipBFieldMap::zFrac_
private

Fractional bin distance along z.

Definition at line 473 of file ShipBFieldMap.h.

◆ zMax_

Float_t ShipBFieldMap::zMax_
private

The maximum value of z for the map.

Definition at line 404 of file ShipBFieldMap.h.

◆ zMin_

Float_t ShipBFieldMap::zMin_
private

The minimum value of z for the map.

Definition at line 401 of file ShipBFieldMap.h.

◆ zOffset_

Float_t ShipBFieldMap::zOffset_
private

The z value of the positional offset for the map.

Definition at line 419 of file ShipBFieldMap.h.

◆ zRange_

Float_t ShipBFieldMap::zRange_
private

The co-ordinate range along z.

Definition at line 410 of file ShipBFieldMap.h.


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