FairShip
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 coordinate shift.
 
void SetYOffset (Float_t yValue)
 Set the y global coordinate shift.
 
void SetZOffset (Float_t zValue)
 Set the z global coordinate 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 number 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 coordinate 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 coordinate 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 coordinate range for the map.
 
Float_t GetXOffset () const
 Get the x offset coordinate of the map for global positioning.
 
Float_t GetYOffset () const
 Get the y offset coordinate of the map for global positioning.
 
Float_t GetZOffset () const
 Get the z offset coordinate 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
 
Float_t GetPsi () const
 
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 coordinate type. More...
 
typedef std::pair< Int_t, Float_t > binPair
 Typedef for an int-double pair.
 

Private Member Functions

 ShipBFieldMap (const ShipBFieldMap &rhs)=delete
 
ShipBFieldMapoperator= (const ShipBFieldMap &rhs)=delete
 
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_
 
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 coordinate 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 coordinate 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 coordinate 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
John Back J.J.B.nosp@m.ack@.nosp@m.warwi.nosp@m.ck.a.nosp@m.c.uk

Definition at line 20 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 336 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 83 of file ShipBFieldMap.h.

Member Enumeration Documentation

◆ CoordAxis

Enumeration to specify the coordinate type.

Enumerator
xAxis 
yAxis 
zAxis 

Definition at line 309 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.

for FairLogger, MESSAGE_ORIGIN

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 coordinate shift to position the field map (cm)
[in]yOffsetThe y global coordinate shift to position the field map (cm)
[in]zOffsetThe z global coordinate 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 21 of file ShipBFieldMap.cxx.

26 : TVirtualMagField(label.c_str()),
27 fieldMap_(new floatArray()),
28 mapFileName_(mapFileName),
29 initialised_(kFALSE),
30 isCopy_(kFALSE),
31 Nx_(0),
32 Ny_(0),
33 Nz_(0),
34 N_(0),
35 xMin_(0.0),
36 xMax_(0.0),
37 dx_(0.0),
38 xRange_(0.0),
39 yMin_(0.0),
40 yMax_(0.0),
41 dy_(0.0),
42 yRange_(0.0),
43 zMin_(0.0),
44 zMax_(0.0),
45 dz_(0.0),
46 zRange_(0.0),
47 xOffset_(xOffset),
48 yOffset_(yOffset),
49 zOffset_(zOffset),
50 phi_(phi),
51 theta_(theta),
52 psi_(psi),
53 scale_(scale),
54 isSymmetric_(isSymmetric),
55 theTrans_(nullptr),
56 Tesla_(10.0) {
57 this->initialise();
58}
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 coordinate range along z.
Int_t Nz_
The number of bins along z.
Float_t Tesla_
Double converting Tesla to kiloGauss (for VMC/FairRoot B field units)
std::vector< std::vector< Float_t > > floatArray
Typedef for a vector containing a vector of floats.
Definition: ShipBFieldMap.h:83
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.
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 coordinate 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 coordinate 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 coordinate (cm)
[in]newYOffsetThe new global offset y coordinate (cm)
[in]newZOffsetThe new global offset z coordinate (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 73 of file ShipBFieldMap.cxx.

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

◆ ~ShipBFieldMap()

ShipBFieldMap::~ShipBFieldMap ( )
virtual

Destructor.

Definition at line 60 of file ShipBFieldMap.cxx.

60 {
61 // Delete the internal vector storing the field map values
62 if (fieldMap_ && isCopy_ == kFALSE) {
63 delete fieldMap_;
64 fieldMap_ = nullptr;
65 }
66
67 if (theTrans_) {
68 delete theTrans_;
69 theTrans_ = nullptr;
70 }
71}

◆ ShipBFieldMap() [3/3]

ShipBFieldMap::ShipBFieldMap ( const ShipBFieldMap rhs)
privatedelete

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 coordinate axis (CoordAxis enumeration for x, y or z)
Returns
the magnetic field component for the given axis

Definition at line 481 of file ShipBFieldMap.cxx.

481 {
482 // Find the magnetic field component along theAxis using trilinear
483 // interpolation based on the current position and neighbouring bins
484 Float_t result(0.0);
485
486 Int_t iAxis(0);
487
488 if (theAxis == CoordAxis::yAxis) {
489 iAxis = 1;
490 } else if (theAxis == CoordAxis::zAxis) {
491 iAxis = 2;
492 }
493
494 if (fieldMap_) {
495 // Get the field component values for the neighbouring bins
496 Float_t A = (*fieldMap_)[binA_][iAxis];
497 Float_t B = (*fieldMap_)[binB_][iAxis];
498 Float_t C = (*fieldMap_)[binC_][iAxis];
499 Float_t D = (*fieldMap_)[binD_][iAxis];
500 Float_t E = (*fieldMap_)[binE_][iAxis];
501 Float_t F = (*fieldMap_)[binF_][iAxis];
502 Float_t G = (*fieldMap_)[binG_][iAxis];
503 Float_t H = (*fieldMap_)[binH_][iAxis];
504
505 // Perform linear interpolation along x
506 Float_t F00 = A * xFrac1_ + B * xFrac_;
507 Float_t F10 = C * xFrac1_ + D * xFrac_;
508 Float_t F01 = E * xFrac1_ + F * xFrac_;
509 Float_t F11 = G * xFrac1_ + H * xFrac_;
510
511 // Linear interpolation along y
512 Float_t F0 = F00 * yFrac1_ + F10 * yFrac_;
513 Float_t F1 = F01 * yFrac1_ + F11 * yFrac_;
514
515 // Linear interpolation along z
516 result = F0 * zFrac1_ + F1 * zFrac_;
517 }
518
519 return result;
520}
Definition: diagrams_a.h:3
Definition: diagrams_b.h:4
Definition: diagrams_c.h:5
Definition: diagrams_d.h:6
Definition: diagrams_e.h:4
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.
def H(p, theta, mDarkPhoton)

◆ 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 coordinates of the point (cm)
[out]BThe x,y,z components of the magnetic field (kGauss = 0.1 tesla)

Definition at line 114 of file ShipBFieldMap.cxx.

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

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

◆ 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 187 of file ShipBFieldMap.h.

187{ 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 211 of file ShipBFieldMap.h.

211{ 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 235 of file ShipBFieldMap.h.

235{ return dz_; }

◆ getFieldMap()

floatArray * ShipBFieldMap::getFieldMap ( ) const
inline

Retrieve the field map.

Returns
the field map

Definition at line 89 of file ShipBFieldMap.h.

89{ 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 467 of file ShipBFieldMap.cxx.

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

◆ 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 145 of file ShipBFieldMap.h.

145{ return mapFileName_; }

◆ GetNBins()

Int_t ShipBFieldMap::GetNBins ( ) const
inline

Get the total number of bins.

Returns
the total number of bins

Definition at line 169 of file ShipBFieldMap.h.

169{ 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 151 of file ShipBFieldMap.h.

151{ 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 157 of file ShipBFieldMap.h.

157{ 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 163 of file ShipBFieldMap.h.

163{ 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 265 of file ShipBFieldMap.h.

265{ 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 280 of file ShipBFieldMap.h.

280{ 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 286 of file ShipBFieldMap.h.

286{ 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 273 of file ShipBFieldMap.h.

273{ return theta_; }

◆ GetXMax()

Float_t ShipBFieldMap::GetXMax ( ) const
inline

Get the maximum value of x for the map.

Returns
the maximum x coordinate (cm)

Definition at line 181 of file ShipBFieldMap.h.

181{ return xMax_; }

◆ GetXMin()

Float_t ShipBFieldMap::GetXMin ( ) const
inline

Get the minimum value of x for the map.

Returns
the minimum x coordinate (cm)

Definition at line 175 of file ShipBFieldMap.h.

175{ return xMin_; }

◆ GetXOffset()

Float_t ShipBFieldMap::GetXOffset ( ) const
inline

Get the x offset coordinate of the map for global positioning.

Returns
the map's x offset coordinate for global positioning (cm)

Definition at line 247 of file ShipBFieldMap.h.

247{ return xOffset_; }

◆ GetXRange()

Float_t ShipBFieldMap::GetXRange ( ) const
inline

Get the x coordinate range for the map.

Returns
the x coordinate range (cm)

Definition at line 193 of file ShipBFieldMap.h.

193{ return xRange_; }

◆ GetYMax()

Float_t ShipBFieldMap::GetYMax ( ) const
inline

Get the maximum value of y for the map.

Returns
the maximum y coordinate (cm)

Definition at line 205 of file ShipBFieldMap.h.

205{ return yMax_; }

◆ GetYMin()

Float_t ShipBFieldMap::GetYMin ( ) const
inline

Get the minimum value of y for the map.

Returns
the minimum y coordinate (cm)

Definition at line 199 of file ShipBFieldMap.h.

199{ return yMin_; }

◆ GetYOffset()

Float_t ShipBFieldMap::GetYOffset ( ) const
inline

Get the y offset coordinate of the map for global positioning.

Returns
the map's y offset coordinate for global positioning (cm)

Definition at line 253 of file ShipBFieldMap.h.

253{ return yOffset_; }

◆ GetYRange()

Float_t ShipBFieldMap::GetYRange ( ) const
inline

Get the y coordinate range for the map.

Returns
the y coordinate range (cm)

Definition at line 217 of file ShipBFieldMap.h.

217{ return yRange_; }

◆ GetZMax()

Float_t ShipBFieldMap::GetZMax ( ) const
inline

Get the maximum value of z for the map.

Returns
the maximum z coordinate (cm)

Definition at line 229 of file ShipBFieldMap.h.

229{ return zMax_; }

◆ GetZMin()

Float_t ShipBFieldMap::GetZMin ( ) const
inline

Get the minimum value of z for the map.

Returns
the minimum z coordinate (cm)

Definition at line 223 of file ShipBFieldMap.h.

223{ return zMin_; }

◆ GetZOffset()

Float_t ShipBFieldMap::GetZOffset ( ) const
inline

Get the z offset coordinate of the map for global positioning.

Returns
the map's z offset coordinate for global positioning (cm)

Definition at line 259 of file ShipBFieldMap.h.

259{ return zOffset_; }

◆ GetZRange()

Float_t ShipBFieldMap::GetZRange ( ) const
inline

Get the z coordinate range for the map.

Returns
the z coordinate range (cm)

Definition at line 241 of file ShipBFieldMap.h.

241{ 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 292 of file ShipBFieldMap.h.

292{ return isSymmetric_; }

◆ initialise()

void ShipBFieldMap::initialise ( )
private

Initialisation.

Definition at line 216 of file ShipBFieldMap.cxx.

216 {
217 if (initialised_ == kFALSE) {
218 if (isCopy_ == kFALSE) {
219 this->readMapFile();
220 }
221
222 // Set the global coordinate translation and rotation info
223 if (fabs(phi_) > 1e-6 || fabs(theta_) > 1e-6 || fabs(psi_) > 1e-6) {
224 // We have non-zero rotation angles. Create a combined translation and
225 // rotation
226 TGeoTranslation tr("offsets", xOffset_, yOffset_, zOffset_);
227 TGeoRotation rot("angles", phi_, theta_, psi_);
228 theTrans_ = new TGeoCombiTrans(tr, rot);
229
230 } else {
231 // We only need a translation
232 theTrans_ = new TGeoTranslation("offsets", xOffset_, yOffset_, zOffset_);
233 }
234
235 initialised_ = kTRUE;
236 }
237}
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 coordinate of the point (cm)
[in]yThe y coordinate of the point (cm)
[in]zThe z coordinate of the point (cm)
Returns
true/false if the point is inside the field map range

Definition at line 375 of file ShipBFieldMap.cxx.

375 {
376 Bool_t inside(kFALSE);
377
378 if (x >= xMin_ && x <= xMax_ && y >= yMin_ && y <= yMax_ && z >= zMin_ &&
379 z <= zMax_) {
380 inside = kTRUE;
381 }
382
383 return inside;
384}

◆ 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 298 of file ShipBFieldMap.h.

298{ return isCopy_; }

◆ operator=()

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

◆ readMapFile()

void ShipBFieldMap::readMapFile ( )
private

Read the field map data and store the information internally.

Definition at line 239 of file ShipBFieldMap.cxx.

239 {
240 LOG(info) << "ShipBFieldMap::readMapFile() creating field " << this->GetName()
241 << " using file " << mapFileName_;
242
243 // Check to see if we have a ROOT file
244 if (mapFileName_.find(".root") != std::string::npos) {
245 this->readRootFile();
246
247 } else {
248 this->readTextFile();
249 }
250}
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 252 of file ShipBFieldMap.cxx.

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

◆ readTextFile()

void ShipBFieldMap::readTextFile ( )
private

Process the text file containing the field map data.

Definition at line 329 of file ShipBFieldMap.cxx.

329 {
330 std::ifstream getData(mapFileName_.c_str());
331
332 if (!getData.is_open()) {
333 LOG(fatal) << "Error: Cannot open magnetic field map file: "
334 << mapFileName_;
335 }
336
337 std::string tmpString("");
338
339 getData >> tmpString >> xMin_ >> xMax_ >> dx_ >> yMin_ >> yMax_ >> dy_ >>
340 zMin_ >> zMax_ >> dz_;
341
342 this->setLimits();
343
344 // Check to see if this object is a "copy"
345 if (isCopy_ == kFALSE) {
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 getData >> Bx >> By >> Bz;
360
361 // B field values are in Tesla. This means these values are multiplied
362 // by a factor of ten since both FairRoot and the VMC interface use kGauss
363 Bx *= Tesla_;
364 By *= Tesla_;
365 Bz *= Tesla_;
366
367 // Store the B field 3-vector
368 fieldMap_->push_back({Bx, By, Bz});
369 }
370 }
371
372 getData.close();
373}

◆ setLimits()

void ShipBFieldMap::setLimits ( )
private

Definition at line 386 of file ShipBFieldMap.cxx.

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

◆ 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 113 of file ShipBFieldMap.h.

113{ phi_ = phi; }

◆ 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 127 of file ShipBFieldMap.h.

127{ 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 133 of file ShipBFieldMap.h.

133{ scale_ = scale; }

◆ 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 120 of file ShipBFieldMap.h.

120{ theta_ = theta; }

◆ SetXOffset()

void ShipBFieldMap::SetXOffset ( Float_t  xValue)
inline

Set the x global coordinate shift.

Parameters
[in]xValueThe value of the x global coordinate shift (cm)

Definition at line 95 of file ShipBFieldMap.h.

95{ xOffset_ = xValue; }

◆ SetYOffset()

void ShipBFieldMap::SetYOffset ( Float_t  yValue)
inline

Set the y global coordinate shift.

Parameters
[in]yValueThe value of the y global coordinate shift (cm)

Definition at line 101 of file ShipBFieldMap.h.

101{ yOffset_ = yValue; }

◆ SetZOffset()

void ShipBFieldMap::SetZOffset ( Float_t  zValue)
inline

Set the z global coordinate shift.

Parameters
[in]zValueThe value of the z global coordinate shift (cm)

Definition at line 107 of file ShipBFieldMap.h.

107{ 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 139 of file ShipBFieldMap.h.

139{ isSymmetric_ = flag; }

Member Data Documentation

◆ binA_

Int_t ShipBFieldMap::binA_
private

Bin A for the trilinear interpolation.

Definition at line 459 of file ShipBFieldMap.h.

◆ binB_

Int_t ShipBFieldMap::binB_
private

Bin B for the trilinear interpolation.

Definition at line 462 of file ShipBFieldMap.h.

◆ binC_

Int_t ShipBFieldMap::binC_
private

Bin C for the trilinear interpolation.

Definition at line 465 of file ShipBFieldMap.h.

◆ binD_

Int_t ShipBFieldMap::binD_
private

Bin D for the trilinear interpolation.

Definition at line 468 of file ShipBFieldMap.h.

◆ binE_

Int_t ShipBFieldMap::binE_
private

Bin E for the trilinear interpolation.

Definition at line 471 of file ShipBFieldMap.h.

◆ binF_

Int_t ShipBFieldMap::binF_
private

Bin F for the trilinear interpolation.

Definition at line 474 of file ShipBFieldMap.h.

◆ binG_

Int_t ShipBFieldMap::binG_
private

Bin G for the trilinear interpolation.

Definition at line 477 of file ShipBFieldMap.h.

◆ binH_

Int_t ShipBFieldMap::binH_
private

Bin H for the trilinear interpolation.

Definition at line 480 of file ShipBFieldMap.h.

◆ dx_

Float_t ShipBFieldMap::dx_
private

The bin width along x.

Definition at line 399 of file ShipBFieldMap.h.

◆ dy_

Float_t ShipBFieldMap::dy_
private

The bin width along y.

Definition at line 411 of file ShipBFieldMap.h.

◆ dz_

Float_t ShipBFieldMap::dz_
private

The bin width along z.

Definition at line 423 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 368 of file ShipBFieldMap.h.

◆ initialised_

Bool_t ShipBFieldMap::initialised_
private

Initialisation boolean.

Definition at line 374 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 378 of file ShipBFieldMap.h.

◆ isSymmetric_

Bool_t ShipBFieldMap::isSymmetric_
private

The boolean specifying if we have quadrant symmetry.

Definition at line 450 of file ShipBFieldMap.h.

◆ mapFileName_

std::string ShipBFieldMap::mapFileName_
private

The name of the map file.

Definition at line 371 of file ShipBFieldMap.h.

◆ N_

Int_t ShipBFieldMap::N_
private

The total number of bins.

Definition at line 390 of file ShipBFieldMap.h.

◆ Nx_

Int_t ShipBFieldMap::Nx_
private

The number of bins along x.

Definition at line 381 of file ShipBFieldMap.h.

◆ Ny_

Int_t ShipBFieldMap::Ny_
private

The number of bins along y.

Definition at line 384 of file ShipBFieldMap.h.

◆ Nz_

Int_t ShipBFieldMap::Nz_
private

The number of bins along z.

Definition at line 387 of file ShipBFieldMap.h.

◆ phi_

Float_t ShipBFieldMap::phi_
private

The first Euler rotation angle about the z axis.

Definition at line 438 of file ShipBFieldMap.h.

◆ psi_

Float_t ShipBFieldMap::psi_
private

The third Euler rotation angle about the new z axis.

Definition at line 444 of file ShipBFieldMap.h.

◆ scale_

Float_t ShipBFieldMap::scale_
private

The B field magnitude scaling factor.

Definition at line 447 of file ShipBFieldMap.h.

◆ Tesla_

Float_t ShipBFieldMap::Tesla_
private

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

Definition at line 456 of file ShipBFieldMap.h.

◆ theta_

Float_t ShipBFieldMap::theta_
private

The second Euler rotation angle about the new x axis.

Definition at line 441 of file ShipBFieldMap.h.

◆ theTrans_

TGeoMatrix* ShipBFieldMap::theTrans_
private

The combined translation and rotation transformation.

Definition at line 453 of file ShipBFieldMap.h.

◆ xFrac1_

Float_t ShipBFieldMap::xFrac1_
private

Complimentary fractional bin distance along x.

Definition at line 492 of file ShipBFieldMap.h.

◆ xFrac_

Float_t ShipBFieldMap::xFrac_
private

Fractional bin distance along x.

Definition at line 483 of file ShipBFieldMap.h.

◆ xMax_

Float_t ShipBFieldMap::xMax_
private

The maximum value of x for the map.

Definition at line 396 of file ShipBFieldMap.h.

◆ xMin_

Float_t ShipBFieldMap::xMin_
private

The minimum value of x for the map.

Definition at line 393 of file ShipBFieldMap.h.

◆ xOffset_

Float_t ShipBFieldMap::xOffset_
private

The x value of the positional offset for the map.

Definition at line 429 of file ShipBFieldMap.h.

◆ xRange_

Float_t ShipBFieldMap::xRange_
private

The coordinate range along x.

Definition at line 402 of file ShipBFieldMap.h.

◆ yFrac1_

Float_t ShipBFieldMap::yFrac1_
private

Complimentary fractional bin distance along y.

Definition at line 495 of file ShipBFieldMap.h.

◆ yFrac_

Float_t ShipBFieldMap::yFrac_
private

Fractional bin distance along y.

Definition at line 486 of file ShipBFieldMap.h.

◆ yMax_

Float_t ShipBFieldMap::yMax_
private

The maximum value of y for the map.

Definition at line 408 of file ShipBFieldMap.h.

◆ yMin_

Float_t ShipBFieldMap::yMin_
private

The minimum value of y for the map.

Definition at line 405 of file ShipBFieldMap.h.

◆ yOffset_

Float_t ShipBFieldMap::yOffset_
private

The y value of the positional offset for the map.

Definition at line 432 of file ShipBFieldMap.h.

◆ yRange_

Float_t ShipBFieldMap::yRange_
private

The coordinate range along y.

Definition at line 414 of file ShipBFieldMap.h.

◆ zFrac1_

Float_t ShipBFieldMap::zFrac1_
private

Complimentary fractional bin distance along z.

Definition at line 498 of file ShipBFieldMap.h.

◆ zFrac_

Float_t ShipBFieldMap::zFrac_
private

Fractional bin distance along z.

Definition at line 489 of file ShipBFieldMap.h.

◆ zMax_

Float_t ShipBFieldMap::zMax_
private

The maximum value of z for the map.

Definition at line 420 of file ShipBFieldMap.h.

◆ zMin_

Float_t ShipBFieldMap::zMin_
private

The minimum value of z for the map.

Definition at line 417 of file ShipBFieldMap.h.

◆ zOffset_

Float_t ShipBFieldMap::zOffset_
private

The z value of the positional offset for the map.

Definition at line 435 of file ShipBFieldMap.h.

◆ zRange_

Float_t ShipBFieldMap::zRange_
private

The coordinate range along z.

Definition at line 426 of file ShipBFieldMap.h.


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