DD4hep  1.30.0
Detector Description Toolkit for High Energy Physics
FieldTypes.cpp
Go to the documentation of this file.
1 //==========================================================================
2 // AIDA Detector description implementation
3 //--------------------------------------------------------------------------
4 // Copyright (C) Organisation europeenne pour la Recherche nucleaire (CERN)
5 // All rights reserved.
6 //
7 // For the licensing terms see $DD4hepINSTALL/LICENSE.
8 // For the list of contributors see $DD4hepINSTALL/doc/CREDITS.
9 //
10 // Author : M.Frank
11 //
12 //==========================================================================
13 
14 #include <DD4hep/FieldTypes.h>
15 #include <DD4hep/detail/Handle.inl>
16 #include <cmath>
17 
18 using namespace dd4hep;
19 
20 #ifndef INFINITY
21 #define INFINITY (numeric_limits<double>::max())
22 #endif
23 
24 // fallthrough only exists from c++17
25 #if defined __has_cpp_attribute
26 #if __has_cpp_attribute(fallthrough)
27 #define ATTR_FALLTHROUGH [[fallthrough]]
28 #else
29 #define ATTR_FALLTHROUGH
30 #endif
31 #else
32 #define ATTR_FALLTHROUGH
33 #endif
34 
39 
41 void ConstantField::fieldComponents(const double* /* pos */, double* field) {
42  field[0] += direction.X();
43  field[1] += direction.Y();
44  field[2] += direction.Z();
45 }
46 
49  : innerField(0), outerField(0), minZ(-INFINITY), maxZ(INFINITY), innerRadius(0), outerRadius(INFINITY)
50 {
52 }
53 
55 void SolenoidField::fieldComponents(const double* pos, double* field) {
56  double z = pos[2] ;
57  // std::cout << " field z=" << z << " maxZ=" << maxZ << " minZ = " << minZ << std::endl ;
58  if( z > minZ && z < maxZ ){
59  double radius = std::sqrt(pos[0] * pos[0] + pos[1] * pos[1]);
60  if (radius < innerRadius)
61  field[2] += innerField;
62  else if (radius < outerRadius)
63  field[2] += outerField;
64  }
65 }
66 
70 }
71 
73 void DipoleField::fieldComponents(const double* pos, double* field) {
74  double z = pos[2], r = std::sqrt(pos[0] * pos[0] + pos[1] * pos[1]);
75  // Check if z coordinate is within dipole z bounds.
76  if (z > zmin && z < zmax && r < rmax) {
77  // Apply all coefficients to this z coordinate.
78  double pp = 1.0;
79  double abs_z = fabs(z);
80  double bx = coefficents[0];
81  for (size_t i = 1; i < coefficents.size(); ++i) {
82  pp *= abs_z;
83  bx += coefficents[i] * pp;
84  }
85  // Flip sign for negative z.
86  if (z < 0)
87  bx = -bx;
88  // Add Bx to the input field.
89  field[0] += bx;
90  }
91 }
92 
93 namespace {
94  constexpr static unsigned char FIELD_INITIALIZED = 1<<0;
95  constexpr static unsigned char FIELD_IDENTITY = 1<<1;
96  constexpr static unsigned char FIELD_ROTATION_ONLY = 1<<2;
97  constexpr static unsigned char FIELD_POSITION_ONLY = 1<<3;
98 }
99 
103 }
104 
106 void MultipoleField::fieldComponents(const double* pos, double* field) {
107  if ( 0 == flag ) {
108  constexpr static double eps = 1e-10;
109  double xx, xy, xz, dx, yx, yy, yz, dy, zx, zy, zz, dz;
110  transform.GetComponents(xx, xy, xz, dx, yx, yy, yz, dy, zx, zy, zz, dz);
111  flag |= FIELD_INITIALIZED;
112  if ( (xx + yy + zz) < (3e0 - eps) ) {
113  flag |= FIELD_ROTATION_ONLY;
114  }
115  else {
116  flag |= FIELD_POSITION_ONLY;
117  if ( (std::abs(dx) + std::abs(dy) + std::abs(dz)) < eps )
118  flag |= FIELD_IDENTITY;
119  }
120  this->inverse = this->transform.Inverse();
121  this->transform.GetRotation(this->rotation);
122  this->transform.GetTranslation(this->translation);
123  }
124  Transform3D::Point p, p0(pos[0],pos[1],pos[2]);
125  if ( flag&FIELD_IDENTITY ) p = std::move(p0);
126  else if ( flag&FIELD_POSITION_ONLY ) p = p0 - this->translation;
127  else p = this->inverse * p0;
128 
129  double x = p.X(), y = p.Y(), z = p.Z();
130  double coord[3] = {x, y, z};
131  // Debug printout:
132  //::printf("Pos: %+15.8e %+15.8e %+15.8e Inverse: %+15.8e %+15.8e %+15.8e\n",
133  // pos[0]/dd4hep::cm,pos[1]/dd4hep::cm,pos[2]/dd4hep::cm, p.X()/dd4hep::cm, p.Y()/dd4hep::cm, p.Z()/dd4hep::cm);
134 
135  if ( 0 == volume.ptr() || volume->Contains(coord) ) {
136  double bx = 0.0;
137  double by = 0.0;
138  double xy = x*y;
139  double x2 = x*x;
140  double y2 = y*y;
141  switch(coefficents.size()) {
142  case 4: // Octupole momentum
143  by += (1./6.) * ( coefficents[3] * (x2*x - 3.0*x*y2) + skews[3]*(y2*y - 3.0*x2*y) );
144  bx += (1./6.) * ( coefficents[3] * (3.0*x2*y - y2*y) + skews[3]*(x2*x - 3.0*x*y2) );
146  case 3: // Sextupole momentum:
147  by += (1./2.) * ( coefficents[2] * (x2 - y2) - skews[2] * 2.0 * xy );
148  bx += (1./2.) * ( coefficents[2] * 2.0 * xy + skews[2] * (x2 - y2) );
150  case 2: // Quadrupole momentum:
151  bx += coefficents[1] * y + skews[1]*x;
152  by += coefficents[1] * x - skews[1]*y;
154  case 1: // Dipole momentum:
155  bx += skews[0];
156  by += coefficents[0];
158  case 0: // Nothing, but still valid
159  break;
160  default: // Error condition
161  throw std::runtime_error("Invalid multipole field definition!");
162  }
163  Transform3D::Point f = this->rotation * Transform3D::Point(bx, by, B_z);
164  field[0] += f.X();
165  field[1] += f.Y();
166  field[2] += f.Z();
167  }
168 }
INFINITY
#define INFINITY
Definition: FieldTypes.cpp:21
dd4hep::DipoleField::zmin
double zmin
Definition: FieldTypes.h:80
dd4hep::DipoleField::rmax
double rmax
Definition: FieldTypes.h:81
dd4hep::MultipoleField
Implementation object of a Multipole magnetic field.
Definition: FieldTypes.h:152
dd4hep::CartesianField::MAGNETIC
@ MAGNETIC
Definition: Fields.h:43
dd4hep::DipoleField
Implementation object of a dipole magnetic field.
Definition: FieldTypes.h:77
dd4hep::ConstantField::fieldComponents
virtual void fieldComponents(const double *, double *field)
Call to access the field components at a given location.
Definition: FieldTypes.cpp:41
dd4hep::SolenoidField::maxZ
double maxZ
Definition: FieldTypes.h:58
dd4hep::ConstantField
Implementation object of a field with constant strength.
Definition: FieldTypes.h:34
dd4hep::DipoleField::zmax
double zmax
Definition: FieldTypes.h:79
dd4hep::SolenoidField::innerRadius
double innerRadius
Definition: FieldTypes.h:59
dd4hep::SolenoidField::fieldComponents
virtual void fieldComponents(const double *pos, double *field)
Call to access the field components at a given location.
Definition: FieldTypes.cpp:55
DD4HEP_INSTANTIATE_HANDLE
DD4HEP_INSTANTIATE_HANDLE(ConstantField)
dd4hep::DipoleField::coefficents
Coefficents coefficents
Definition: FieldTypes.h:82
dd4hep::DipoleField::fieldComponents
virtual void fieldComponents(const double *pos, double *field)
Call to access the field components at a given location.
Definition: FieldTypes.cpp:73
dd4hep::MultipoleField::B_z
double B_z
Constant Z field overlay.
Definition: FieldTypes.h:167
dd4hep::SolenoidField::minZ
double minZ
Definition: FieldTypes.h:57
dd4hep::ConstantField::direction
Direction direction
Field direction.
Definition: FieldTypes.h:37
dd4hep::MultipoleField::transform
Transform3D transform
Position transformation of the field. Only stored here for reference.
Definition: FieldTypes.h:161
dd4hep::DipoleField::DipoleField
DipoleField()
Initializing constructor.
Definition: FieldTypes.cpp:68
dd4hep::SolenoidField::SolenoidField
SolenoidField()
Initializing constructor.
Definition: FieldTypes.cpp:48
dd4hep::MultipoleField::fieldComponents
virtual void fieldComponents(const double *pos, double *field)
Call to access the field components at a given location.
Definition: FieldTypes.cpp:106
dd4hep::MultipoleField::rotation
Rotation3D rotation
The rotation part of the transformation. Need to rotate the field.
Definition: FieldTypes.h:165
dd4hep::CartesianField::TypedObject::field_type
int field_type
Field type.
Definition: Fields.h:56
dd4hep::MultipoleField::MultipoleField
MultipoleField()
Initializing constructor.
Definition: FieldTypes.cpp:101
dd4hep::SolenoidField::outerField
double outerField
Definition: FieldTypes.h:56
dd4hep::MultipoleField::skews
Coefficents skews
Multi-pole skews.
Definition: FieldTypes.h:157
dd4hep::Handle::ptr
T * ptr() const
Access to the held object.
Definition: Handle.h:153
ATTR_FALLTHROUGH
#define ATTR_FALLTHROUGH
Definition: FieldTypes.cpp:32
dd4hep::MultipoleField::inverse
Transform3D inverse
Inverse position transformation of the field.
Definition: FieldTypes.h:163
dd4hep::MultipoleField::volume
Solid volume
Boundary volume (optional)
Definition: FieldTypes.h:159
dd4hep
Namespace for the AIDA detector description toolkit.
Definition: AlignmentsCalib.h:28
dd4hep::MultipoleField::translation
Transform3D::Point translation
Translation of the transformation.
Definition: FieldTypes.h:173
dd4hep::MultipoleField::coefficents
Coefficents coefficents
Multi-pole coefficients.
Definition: FieldTypes.h:155
dd4hep::SolenoidField::outerRadius
double outerRadius
Definition: FieldTypes.h:60
FieldTypes.h
dd4hep::SolenoidField::innerField
double innerField
Definition: FieldTypes.h:55
dd4hep::MultipoleField::flag
unsigned char flag
The access to the field will be optimized. Remember properties.
Definition: FieldTypes.h:171
dd4hep::SolenoidField
Implementation object of a solenoidal magnetic field.
Definition: FieldTypes.h:53