Belle II Software  release-08-01-10
BFieldComponent3d Class Reference

The BFieldComponent3d class. More...

#include <BFieldComponent3d.h>

Inheritance diagram for BFieldComponent3d:
Collaboration diagram for BFieldComponent3d:

Public Member Functions

 BFieldComponent3d ()=default
 The BFieldComponent3d constructor.
 
virtual ~BFieldComponent3d ()=default
 The BFieldComponent3d destructor.
 
virtual void initialize () override
 Initializes the magnetic field component. More...
 
virtual ROOT::Math::XYZVector calculate (const ROOT::Math::XYZVector &point) const override
 Calculates the magnetic field vector at the specified space point. More...
 
virtual void terminate () override
 Terminates the magnetic field component. More...
 
void setMapFilename (const std::string &filename)
 Sets the filename of the magnetic field map. More...
 
void setMapSize (int sizeR, int sizePhi, int sizeZ)
 Sets the size of the magnetic field map. More...
 
void setMapRegionZ (double minZ, double maxZ, double offset)
 Sets the size of the magnetic field map. More...
 
void setMapRegionR (double minR, double maxR)
 Sets the size of the magnetic field map. More...
 
void setGridPitch (double pitchR, double pitchPhi, double pitchZ)
 Sets the grid pitch of the magnetic field map. More...
 
void setExcludeRegionZ (double minZ, double maxZ)
 Sets the size of the magnetic field map to exclude. More...
 
void setExcludeRegionR (double minR, double maxR)
 Sets the size of the magnetic field map to exclude. More...
 
void setErrorRegionR (double minR=-1., double maxR=-1., double errBr=0.0, double errBphi=0.0, double errBz=0.0)
 Sets the size of the magnetic field map to apply error on B. More...
 
void mirrorPhi (bool mirrorPhi=0.)
 
void doInterpolation (bool interpolate=true)
 
void enableCoordinate (const std::string &mapEnable="rphiz")
 Option to reduce 3D to 2D map (in coordinates, NOT Br, Bphi, Bz components) More...
 

Private Member Functions

ROOT::Math::XYZVector interpolate (unsigned int ir, unsigned int iphi, unsigned int iz, double wr, double wphi, double wz) const
 Interpolate the value of B-field between (ir, iphi, iz) and (ir+1, iphi+1, iz+1) using weights (wr, wphi, wz)
 

Private Attributes

std::string m_mapFilename {""}
 The filename of the magnetic field map.
 
std::vector< ROOT::Math::XYZVector > m_bmap
 The memory buffer for the magnetic field map.
 
std::string m_mapEnable {"rphiz"}
 Enable different dimension, "rphiz", "rphi", "phiz" or "rz" >
 
bool m_interpolate {true}
 Flag to switch on/off interpolation >
 
int m_mapSize [3]
 The size of the map in r, phi and z.
 
double m_mapRegionZ [2] {0}
 The min and max boundaries of the map region in z.
 
double m_mapOffset {0}
 Offset required because the accelerator group defines the Belle center as zero.
 
double m_mapRegionR [2] {0}
 The min and max boundaries of the map region in r.
 
double m_gridPitch [3] {0}
 The grid pitch in r,phi,z.
 
double m_igridPitch [3] {0}
 The inverted grid pitch in r,phi,z.
 
double m_exRegionZ [2] {0}
 The min and max boundaries of the excluded region in z.
 
double m_exRegionR [2] {0}
 The min and max boundaries of the excluded region in r.
 
bool m_exRegion {true}
 Flag to indicate whether there is a region to exclude. More...
 
bool m_mirrorPhi {true}
 Flag to indicate whether there is a region to exclude. More...
 
double m_errRegionR [2] {0}
 The min and max boundaries of the region in r to apply error.
 
double m_errB [3] {0}
 The error Br, Bphi, Bz as a scale factor (B_new = m_errB * B_old).
 

Detailed Description

The BFieldComponent3d class.

This class represents a 3d magnetic field map. The magnetic field map is stored as a grid in cylindrical coordinates. It is defined by a minimum radius and a maximum radius, a minimum z and a maximum z value (phi is assumed to be in 360 degrees range, option to mirror about the x-y plane exist), a pitch size in both, r, z and phi and the number of grid points. The ZOffset is used to account for the fact that the acceleration group defines 0 to be in the center of the detector, while the detector group defines the IP to be the center. The filename points to the zip file containing the magnetic field map.

Definition at line 36 of file BFieldComponent3d.h.

Member Function Documentation

◆ calculate()

ROOT::Math::XYZVector calculate ( const ROOT::Math::XYZVector &  point) const
overridevirtual

Calculates the magnetic field vector at the specified space point.

Parameters
pointThe space point in [cm] at which the magnetic field vector should be calculated.
Returns
The magnetic field vector at the given space point in [T]. Returns a zero vector XYZVector(0,0,0) if the space point lies outside the region described by the component.

Implements BFieldComponentAbs.

Definition at line 187 of file BFieldComponent3d.cc.

188 {
189  auto getPhiIndexWeight = [this](double y, double x, double & wphi) -> unsigned int {
190  double phi = fast_atan2_minimax<4>(y, x);
191  wphi = phi * m_igridPitch[1];
192  auto iphi = static_cast<unsigned int>(wphi);
193  iphi = min(iphi, static_cast<unsigned int>(m_mapSize[1] - 2));
194  wphi -= iphi;
195  return iphi;
196  };
197 
198  ROOT::Math::XYZVector B(0, 0, 0);
199 
200  // If both '3d' and 'Beamline' components are defined in xml file,
201  // '3d' component returns zero field where 'Beamline' component is defined.
202  // If no 'Beamline' component is defined in xml file, the following function will never be called.
203  if (BFieldComponentBeamline::Instance().isInRange(point)) {
204  return B;
205  }
206 
207  double z = point.Z();
208  // Check if the point lies inside the magnetic field boundaries
209  if (z < m_mapRegionZ[0] || z > m_mapRegionZ[1]) return B;
210 
211  double r2 = point.Perp2();
212  // Check if the point lies inside the magnetic field boundaries
213  if (r2 < m_mapRegionR[0]*m_mapRegionR[0] || r2 >= m_mapRegionR[1]*m_mapRegionR[1]) return B;
214  // Check if the point lies in the exclude region
215  if (m_exRegion && (z >= m_exRegionZ[0]) && (z < m_exRegionZ[1]) &&
216  (r2 >= m_exRegionR[0]*m_exRegionR[0]) && (r2 < m_exRegionR[1]*m_exRegionR[1])) return B;
217 
218  // Calculate the lower index of the point in the Z grid
219  // Note that z coordinate is inverted to match ANSYS frame
220  double wz = (m_mapRegionZ[1] - z) * m_igridPitch[2];
221  unsigned int iz = static_cast<int>(wz);
222  iz = min(iz, static_cast<unsigned int>(m_mapSize[2] - 2));
223  wz -= iz;
224 
225  if (r2 > 0) {
226  double r = sqrt(r2);
227 
228  // Calculate the lower index of the point in the R grid
229  double wr = (r - m_mapRegionR[0]) * m_igridPitch[0];
230  auto ir = static_cast<unsigned int>(wr);
231  ir = min(ir, static_cast<unsigned int>(m_mapSize[0] - 2));
232  wr -= ir;
233 
234  // Calculate the lower index of the point in the Phi grid
235  double ay = std::abs(point.Y());
236  double wphi;
237  unsigned int iphi = getPhiIndexWeight(ay, point.X(), wphi);
238 
239  // Get B-field values from map
240  ROOT::Math::XYZVector b = interpolate(ir, iphi, iz, wr, wphi, wz); // in cylindrical system
241  double norm = 1 / r;
242  double s = ay * norm, c = point.X() * norm;
243  // Flip sign of By if y<0
244  const double sgny = (point.Y() >= 0) - (point.Y() < 0);
245  // in cartesian system
246  B.SetXYZ(-(b.X() * c - b.Y() * s), -sgny * (b.X() * s + b.Y() * c), b.Z());
247  } else {
248  // Get B-field values from map in cartesian system assuming phi=0, so Bx = Br and By = Bphi
249  B = interpolate(0, 0, iz, 0, 0, wz);
250  }
251 
252  return B;
253 }
double m_exRegionR[2]
The min and max boundaries of the excluded region in r.
ROOT::Math::XYZVector interpolate(unsigned int ir, unsigned int iphi, unsigned int iz, double wr, double wphi, double wz) const
Interpolate the value of B-field between (ir, iphi, iz) and (ir+1, iphi+1, iz+1) using weights (wr,...
double m_mapRegionZ[2]
The min and max boundaries of the map region in z.
int m_mapSize[3]
The size of the map in r, phi and z.
double m_igridPitch[3]
The inverted grid pitch in r,phi,z.
bool m_exRegion
Flag to indicate whether there is a region to exclude.
double m_mapRegionR[2]
The min and max boundaries of the map region in r.
double m_exRegionZ[2]
The min and max boundaries of the excluded region in z.
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
static BFieldComponentBeamline & Instance()
BFieldComponentBeamline instance.

◆ doInterpolation()

void doInterpolation ( bool  interpolate = true)
inline
Parameters
interpolateFlag to switch on/off interpolation

Definition at line 137 of file BFieldComponent3d.h.

bool m_interpolate
Flag to switch on/off interpolation >

◆ enableCoordinate()

void enableCoordinate ( const std::string &  mapEnable = "rphiz")
inline

Option to reduce 3D to 2D map (in coordinates, NOT Br, Bphi, Bz components)

Parameters
mapEnableList of dimensions to enable: "rphiz", "rphi", "phiz" or "rz"

Definition at line 143 of file BFieldComponent3d.h.

◆ initialize()

void initialize ( void  )
overridevirtual

Initializes the magnetic field component.

This method opens the magnetic field map file.

Reimplemented from BFieldComponentAbs.

Definition at line 27 of file BFieldComponent3d.cc.

◆ mirrorPhi()

void mirrorPhi ( bool  mirrorPhi = 0.)
inline
Parameters
mirrorPhiFlag to enable mirroring in phi about x-z plane

Definition at line 132 of file BFieldComponent3d.h.

◆ setErrorRegionR()

void setErrorRegionR ( double  minR = -1.,
double  maxR = -1.,
double  errBr = 0.0,
double  errBphi = 0.0,
double  errBz = 0.0 
)
inline

Sets the size of the magnetic field map to apply error on B.

Parameters
minRThe left (min) border of the magnetic field map region in r [m] -> [cm].
maxRThe right (max) border of the magnetic field map region in r [m] -> [cm].
errBrThe size of the error on Br, Br_new = errBr*Br [fraction].
errBphiThe size of the error on Bphi, Bphi_new = errBphi*Bphi [fraction].
errBzThe size of the error on Bz, Bz_new = errBz*Bz [fraction].

Definition at line 126 of file BFieldComponent3d.h.

◆ setExcludeRegionR()

void setExcludeRegionR ( double  minR,
double  maxR 
)
inline

Sets the size of the magnetic field map to exclude.

Parameters
minRThe left (min) border of the magnetic field map region in r [m] -> [cm].
maxRThe right (max) border of the magnetic field map region in r [m] -> [cm].

Definition at line 116 of file BFieldComponent3d.h.

◆ setExcludeRegionZ()

void setExcludeRegionZ ( double  minZ,
double  maxZ 
)
inline

Sets the size of the magnetic field map to exclude.

Parameters
minZThe left (min) border of the magnetic field map region in z [m] -> [cm].
maxZThe right (max) border of the magnetic field map region in z [m] -> [cm].

Definition at line 109 of file BFieldComponent3d.h.

◆ setGridPitch()

void setGridPitch ( double  pitchR,
double  pitchPhi,
double  pitchZ 
)
inline

Sets the grid pitch of the magnetic field map.

Parameters
pitchRThe pitch of the grid in r [m] -> [cm].
pitchPhiThe pitch of the grid in phi [degrees] -> [rad] converted back to degrees in the code.
pitchZThe pitch of the grid in z [m] -> [cm].

Definition at line 102 of file BFieldComponent3d.h.

◆ setMapFilename()

void setMapFilename ( const std::string &  filename)
inline

Sets the filename of the magnetic field map.

Parameters
filenameThe filename of the magnetic field map.

Definition at line 71 of file BFieldComponent3d.h.

◆ setMapRegionR()

void setMapRegionR ( double  minR,
double  maxR 
)
inline

Sets the size of the magnetic field map.

Parameters
minRThe left (min) border of the magnetic field map region in r [m] -> [cm].
maxRThe right (max) border of the magnetic field map region in r [m] -> [cm].

Definition at line 94 of file BFieldComponent3d.h.

◆ setMapRegionZ()

void setMapRegionZ ( double  minZ,
double  maxZ,
double  offset 
)
inline

Sets the size of the magnetic field map.

Parameters
minZThe left (min) border of the magnetic field map region in z [m] -> [cm].
maxZThe right (max) border of the magnetic field map region in z [m] -> [cm].
offsetThe offset in z [m] -> [cm] which is required because the accelerator group defines the Belle center as zero.

Definition at line 87 of file BFieldComponent3d.h.

◆ setMapSize()

void setMapSize ( int  sizeR,
int  sizePhi,
int  sizeZ 
)
inline

Sets the size of the magnetic field map.

Parameters
sizeRThe number of points in the r direction.
sizePhiThe number of points in the phi direction.
sizeZThe number of points in the z direction.

Definition at line 79 of file BFieldComponent3d.h.

◆ terminate()

void terminate ( void  )
overridevirtual

Terminates the magnetic field component.

This method closes the magnetic field map file.

Reimplemented from BFieldComponentAbs.

Definition at line 255 of file BFieldComponent3d.cc.

Member Data Documentation

◆ m_exRegion

bool m_exRegion {true}
private

Flag to indicate whether there is a region to exclude.

>

Definition at line 178 of file BFieldComponent3d.h.

◆ m_mirrorPhi

bool m_mirrorPhi {true}
private

Flag to indicate whether there is a region to exclude.

>

Definition at line 180 of file BFieldComponent3d.h.


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