Belle II Software light-2406-ragdoll
geometry
Collaboration diagram for geometry:

Modules

 geometry data objects
 
 geometry modules
 

Namespaces

namespace  Belle2::geometry
 Common code concerning the geometry representation of the detector.
 

Classes

class  BFieldComponent3d
 The BFieldComponent3d class. More...
 
class  BFieldComponentAbs
 The BFieldComponentAbs class. More...
 
class  BFieldComponentBeamline
 The BFieldComponentBeamline class. More...
 
class  BFieldComponentConstant
 The BFieldComponentConstant class. More...
 
class  BFieldComponentKlm1
 The Bfieldcomponentklm1 class. More...
 
class  BFieldComponentQuad
 The BFieldComponentQuad class. More...
 
class  BFieldComponentRadial
 The BFieldComponentRadial class. More...
 
class  BFieldFrameworkInterface
 Simple BFieldComponent to just wrap the existing BFieldMap with the new BFieldManager. More...
 
class  BFieldMap
 This class represents the magnetic field of the Belle II detector. More...
 
class  GeoMagneticField
 The GeoMagneticField class. More...
 
struct  triangle_t
 Triangle structure. More...
 
struct  xy_t
 A simple 2d vector stucture. More...
 
class  TriangularInterpolation
 The TriangularInterpolation class. More...
 
class  BeamlineFieldMapInterpolation
 The BeamlineFieldMapInterpolation class. More...
 
class  GeoComponent
 Describe one component of the Geometry. More...
 
class  GeoConfiguration
 configuration of the geometry More...
 
class  GeoMaterial
 Class to represent a material informaion in the Database. More...
 
class  GeoMaterialComponent
 Component of a material. More...
 
class  GeoMaterialProperty
 Property of a material. More...
 
class  GeoOpticalSurface
 Represent an optical finish of a surface. More...
 
class  MagneticFieldComponent3D
 Describe one component of the Geometry. More...
 
class  MyDBPayloadClass
 Class containing all the parameters needed to create the geometry and suitable to save into a ROOT file to be used from the Database. More...
 
class  MyDBCreator
 Very simple Creator class which actually does not do anything but shows how creators should implement loading the geometry from database. More...
 

Functions

BFieldComponentBeamline ** GetInstancePtr ()
 Static function holding the instance.
 
BFieldComponentRadial::BFieldPoint operator* (const BFieldComponentRadial::BFieldPoint &v, double a)
 multiply a radial bfield point by a real number
 
BFieldComponentRadial::BFieldPoint operator+ (const BFieldComponentRadial::BFieldPoint &u, const BFieldComponentRadial::BFieldPoint &v)
 Add two radial bfield points together.
 
template<class BFIELDCOMP >
BFIELDCOMP & addBFieldComponent ()
 Adds a new BField component to the Belle II magnetic field.
 
ROOT::Math::XYZVector getBField (const ROOT::Math::XYZVector &point) const
 Returns the magnetic field of the Belle II detector at the specified space point.
 
virtual void initialize () override
 Initializes the magnetic field component.
 
virtual ~BFieldComponentBeamline ()
 The BFieldComponentBeamline destructor.
 
bool isInRange (const ROOT::Math::XYZVector &point) const
 Check presence of beamline field at the specific space point in the detector coordinate frame.
 
virtual ROOT::Math::XYZVector calculate (const ROOT::Math::XYZVector &point) const override
 Calculates the magnetic field vector at the specified space point.
 
virtual void terminate () override
 Terminates the magnetic field component.
 
static BFieldComponentBeamlineInstance ()
 BFieldComponentBeamline instance.
 
 BFieldComponentBeamline ()
 The BFieldComponentBeamline constructor.
 
ROOT::Math::XYZVector getField (const ROOT::Math::XYZVector &pos) const override
 return the field assuming we are inside the active region as returned by inside()
 
ROOT::Math::XYZVector interpolate (unsigned int ir, unsigned int iphi, unsigned int iz, double wr, double wphi, double wz) const
 Linear interpolate the magnetic field inside a bin.
 

Detailed Description

Function Documentation

◆ addBFieldComponent()

BFIELDCOMP & addBFieldComponent

Adds a new BField component to the Belle II magnetic field.

The class of the magnetic field component has to inherit from BFieldComponentAbs.

Returns
A reference to the added magnetic field component.

Definition at line 98 of file BFieldMap.h.

99 {
100 BFIELDCOMP* newComponent = new BFIELDCOMP;
101 m_components.push_back(newComponent);
102 return *newComponent;
103 }
std::list< BFieldComponentAbs * > m_components
The components of the magnetic field.
Definition: BFieldMap.h:69

◆ BFieldComponentBeamline()

The BFieldComponentBeamline constructor.

Definition at line 707 of file BFieldComponentBeamline.cc.

708 {
709 auto gInstance = GetInstancePtr();
710 if (*gInstance != nullptr) {
711 B2WARNING("BFieldComponentBeamline: object already instantiated");
712 } else {
713 *gInstance = this;
714 }
715 }
BFieldComponentBeamline ** GetInstancePtr()
Static function holding the instance.

◆ 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 Cartesian coordinates (x,y,z) 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 661 of file BFieldComponentBeamline.cc.

662 {
663 ROOT::Math::XYZVector res;
665 ROOT::Math::XYZVector v = -p; // invert coordinates to match ANSYS one
666 double xc = v.X() * c, zs = v.Z() * s, zc = v.Z() * c, xs = v.X() * s;
667 ROOT::Math::XYZVector hv{xc - zs, v.Y(), zc + xs};
668 ROOT::Math::XYZVector lv{xc + zs, v.Y(), zc - xs};
669 ROOT::Math::XYZVector hb = m_her->interpolateField(hv);
670 ROOT::Math::XYZVector lb = m_ler->interpolateField(lv);
671 ROOT::Math::XYZVector rhb{hb.X()* c + hb.Z()* s, hb.Y(), hb.Z()* c - hb.X()* s};
672 ROOT::Math::XYZVector rlb{lb.X()* c - lb.Z()* s, lb.Y(), lb.Z()* c + lb.X()* s};
673
674 double mhb = std::abs(rhb.X()) + std::abs(rhb.Y()) + std::abs(rhb.Z());
675 double mlb = std::abs(rlb.X()) + std::abs(rlb.Y()) + std::abs(rlb.Z());
676
677 if (mhb < 1e-10) res = rlb;
678 else if (mlb < 1e-10) res = rhb;
679 else {
680 res = 0.5 * (rlb + rhb);
681 }
682 return res;
683 }
BeamlineFieldMapInterpolation * m_ler
Actual magnetic field interpolation object for LER.
BeamlineFieldMapInterpolation * m_her
Actual magnetic field interpolation object for HER.
double m_sinBeamCrossAngle
The sin of the crossing angle of the beams
double m_cosBeamCrossAngle
The cos of the crossing angle of the beams
ROOT::Math::XYZVector interpolateField(const ROOT::Math::XYZVector &v) const
Interpolate the magnetic field vector at the specified space point.

◆ getBField()

ROOT::Math::XYZVector getBField ( const ROOT::Math::XYZVector &  point) const
inlineprivate

Returns the magnetic field of the Belle II detector at the specified space point.

The space point is given in Cartesian coordinates (x,y,z) in [cm].

Parameters
pointThe space point in Cartesian coordinates.
Returns
A three vector of the magnetic field in [T] at the specified space point.

Definition at line 106 of file BFieldMap.h.

107 {
108 ROOT::Math::XYZVector magFieldVec(0.0, 0.0, 0.0);
109 //Check that the point makes sense
110 if (std::isnan(point.X()) || std::isnan(point.Y()) || std::isnan(point.Z())) {
111 B2ERROR("Bfield requested for a position containing NaN, returning field 0");
112 return magFieldVec;
113 }
114 //Loop over all magnetic field components and add their magnetic field vectors
115 for (const BFieldComponentAbs* comp : m_components) {
116 magFieldVec += comp->calculate(point);
117 }
118 return magFieldVec;
119 }

◆ getField()

ROOT::Math::XYZVector getField ( const ROOT::Math::XYZVector &  pos) const
overridevirtual

return the field assuming we are inside the active region as returned by inside()

small helper function to calculate the bin index and fraction inside the index given a relative coordinate and the coordinate index (0=R, 1=Phi, 2=Z).

Example: If the z grid starts at 5 and goes to 15 with 6 points (pitch size of 2), then calling getIndexWeight(10.5, 2) will return tuple(2, 0.75) since the 10.5 lies in bin "2" and is already 75% to the next bin

It will also cap the index to be inside the valid grid range

Returns a tuple with the bin index as first element and the weight fraction inside the bin as second element

Implements MagneticFieldComponent.

Definition at line 89 of file MagneticFieldComponent3D.cc.

90 {
91 using std::get;
92
106 auto getIndexWeight = [this](double value, int coordinate) -> std::tuple<unsigned int, double> {
107 double weight = value * m_invgridPitch[coordinate];
108 auto index = static_cast<unsigned int>(weight);
109 index = std::min(index, static_cast<unsigned int>(m_mapSize[coordinate] - 2));
110 weight -= index;
111 return std::make_tuple(index, weight);
112 };
113
114 const double z = pos.Z();
115 const double r2 = pos.Perp2();
116
117 // Calculate the lower index of the point in the Z grid
118 // Note that z coordinate is inverted to match ANSYS frame
119 const auto iz = getIndexWeight(m_maxZ - z, 2);
120
121 // directly on z axis: get B-field values from map in cartesian system assuming phi=0, so Bx = Br and By = Bphi
122 if (r2 == 0) return interpolate(0, 0, get<0>(iz), 0, 0, get<1>(iz));
123
124 const double r = std::sqrt(r2);
125 const auto ir = getIndexWeight(r - m_minR, 0);
126
127 // Calculate the lower index of the point in the Phi grid
128 const double ay = std::abs(pos.Y());
129 const auto iphi = getIndexWeight(fast_atan2_minimax<4>(ay, pos.X()), 1);
130
131 // Get B-field values from map in cylindrical coordinates
132 ROOT::Math::XYZVector b = interpolate(get<0>(ir), get<0>(iphi), get<0>(iz), get<1>(ir), get<1>(iphi), get<1>(iz));
133 // and convert it to cartesian
134 const double norm = 1 / r;
135 const double s = ay * norm;
136 const double c = pos.X() * norm;
137 // Flip sign of By if y<0
138 const double sgny = (pos.Y() >= 0) - (pos.Y() < 0);
139 // in cartesian system
140 return ROOT::Math::XYZVector(-(b.X() * c - b.Y() * s), -sgny * (b.X() * s + b.Y() * c), b.Z());
141 }
int m_mapSize[3]
number of bins in r, phi and z
float m_maxZ
maximal Z for which this field is present
float m_invgridPitch[3]
inverted grid pitch in r, phi and z
float m_minR
minimal R=sqrt(x^2+y^2) for which this field is present
ROOT::Math::XYZVector interpolate(unsigned int ir, unsigned int iphi, unsigned int iz, double wr, double wphi, double wz) const
Linear interpolate the magnetic field inside a bin.

◆ GetInstancePtr()

BFieldComponentBeamline ** GetInstancePtr ( )

Static function holding the instance.

Definition at line 690 of file BFieldComponentBeamline.cc.

691 {
692 static BFieldComponentBeamline* gInstance = nullptr;
693 return &gInstance;
694 }

◆ initialize()

void initialize ( )
overridevirtual

Initializes the magnetic field component.

This method opens the magnetic field map file.

Reimplemented from BFieldComponentAbs.

Definition at line 636 of file BFieldComponentBeamline.cc.

637 {
638 if (!m_ler) m_ler = new BeamlineFieldMapInterpolation;
639 if (!m_her) m_her = new BeamlineFieldMapInterpolation;
642 }
std::string m_mapFilename_ler
The filename of the magnetic field map.
std::string m_interFilename_ler
The filename of the map for interpolation.
std::string m_mapFilename_her
Parameter to set Angle of the beam.
std::string m_interFilename_her
The filename of the map for interpolation.
double m_mapRegionR[2]
The min and max boundaries of the map region in r.
void init(const string &fieldmapname, const string &interpolname, double validRadius)
Initializes the magnetic field component.

◆ Instance()

BFieldComponentBeamline & Instance ( )
static

BFieldComponentBeamline instance.

static function

Returns
reference to the BFieldComponentBeamline instance

Definition at line 697 of file BFieldComponentBeamline.cc.

698 {
699 auto gInstance = GetInstancePtr();
700 if (*gInstance == nullptr) {
701 // Constructor creates a new instance, inits gInstance.
703 }
704 return **gInstance;
705 }
BFieldComponentBeamline()
The BFieldComponentBeamline constructor.

◆ interpolate()

ROOT::Math::XYZVector interpolate ( unsigned int  ir,
unsigned int  iphi,
unsigned int  iz,
double  wr,
double  wphi,
double  wz 
) const
private

Linear interpolate the magnetic field inside a bin.

Parameters
irnumber of the bin in r
iphinumber of the bin in phi
iznumber of the bin in z
wrr weight: fraction we are away from the lower r corner relative to the pitch size
wphiphi weight: fraction we are away from the lower phi corner relative to the pitch size
wzphi weight: fraction we are away from the lower z corner relative to the pitch size

Definition at line 143 of file MagneticFieldComponent3D.cc.

145 {
146 const unsigned int strideZ = m_mapSize[0] * m_mapSize[1];
147 const unsigned int strideR = m_mapSize[1];
148
149 const double wz0 = 1 - wz1;
150 const double wr0 = 1 - wr1;
151 const double wphi0 = 1 - wphi1;
152 const unsigned int j000 = iz * strideZ + ir * strideR + iphi;
153 const unsigned int j001 = j000 + 1;
154 const unsigned int j010 = j000 + strideR;
155 const unsigned int j011 = j001 + strideR;
156 const unsigned int j100 = j000 + strideZ;
157 const unsigned int j101 = j001 + strideZ;
158 const unsigned int j110 = j010 + strideZ;
159 const unsigned int j111 = j011 + strideZ;
160 const double w00 = wphi0 * wr0;
161 const double w10 = wphi0 * wr1;
162 const double w01 = wphi1 * wr0;
163 const double w11 = wphi1 * wr1;
164 const std::vector<ROOT::Math::XYZVector>& B = m_bmap;
165 return
166 (B[j000] * w00 + B[j001] * w01 + B[j010] * w10 + B[j011] * w11) * wz0 +
167 (B[j100] * w00 + B[j101] * w01 + B[j110] * w10 + B[j111] * w11) * wz1;
168 }
std::vector< ROOT::Math::XYZVector > m_bmap
magnetic field strength

◆ isInRange()

bool isInRange ( const ROOT::Math::XYZVector &  point) const

Check presence of beamline field at the specific space point in the detector coordinate frame.

Parameters
pointThe space point in Cartesian coordinates (x,y,z) in [cm] at which the magnetic field presence is checked
Returns
true is in case of the magnetic field false – otherwise

Definition at line 650 of file BFieldComponentBeamline.cc.

651 {
652 if (!m_ler || !m_her) return false;
654 ROOT::Math::XYZVector v = -p; // invert coordinates to match ANSYS one
655 double xc = v.X() * c, zs = v.Z() * s, zc = v.Z() * c, xs = v.X() * s;
656 ROOT::Math::XYZVector hv{xc - zs, v.Y(), zc + xs};
657 ROOT::Math::XYZVector lv{xc + zs, v.Y(), zc - xs};
658 return m_ler->inRange(lv) || m_her->inRange(hv);
659 }
bool inRange(const ROOT::Math::XYZVector &v) const
Check the space point if the interpolation exists.

◆ operator*()

BFieldComponentRadial::BFieldPoint operator* ( const BFieldComponentRadial::BFieldPoint v,
double  a 
)
inline

multiply a radial bfield point by a real number

Definition at line 72 of file BFieldComponentRadial.cc.

73 {
74 return {v.r * a, v.z * a};
75 }
double r
Magnetic field in r direction.
double z
Magnetic field in z direction.

◆ operator+()

Add two radial bfield points together.

Definition at line 78 of file BFieldComponentRadial.cc.

80 {
81 return {u.r + v.r, u.z + v.z};
82 }

◆ terminate()

void terminate ( )
overridevirtual

Terminates the magnetic field component.

This method closes the magnetic field map file.

Reimplemented from BFieldComponentAbs.

Definition at line 685 of file BFieldComponentBeamline.cc.

686 {
687 }

◆ ~BFieldComponentBeamline()

The BFieldComponentBeamline destructor.

Definition at line 644 of file BFieldComponentBeamline.cc.

645 {
646 if (m_ler) delete m_ler;
647 if (m_her) delete m_her;
648 }