Belle II Software  release-06-01-15
Helix.cc
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * *
5  * See git log for contributors and copyright holders. *
6  * This file is licensed under LGPL-3.0, see LICENSE.md. *
7  **************************************************************************/
8 #include <framework/dataobjects/Helix.h>
9 
10 #include <framework/gearbox/Const.h>
11 
12 #include <boost/math/special_functions/sign.hpp>
13 #include <boost/math/special_functions/sinc.hpp>
14 #include <boost/math/tools/precision.hpp>
15 
16 #include <TMatrixD.h>
17 
18 #include <cassert>
19 
20 using namespace Belle2;
21 using namespace HelixParameterIndex;
22 
23 Helix::Helix():
24  m_d0(0.0),
25  m_phi0(0.0),
26  m_omega(0.0),
27  m_z0(0.0),
28  m_tanLambda(0.0)
29 {
30 }
31 
32 Helix::Helix(const TVector3& position,
33  const TVector3& momentum,
34  const short int charge,
35  const double bZ)
36 {
37  setCartesian(position, momentum, charge, bZ);
38 }
39 
40 Helix::Helix(const double& d0,
41  const double& phi0,
42  const double& omega,
43  const double& z0,
44  const double& tanLambda)
45  : m_d0(d0),
46  m_phi0(phi0),
47  m_omega(omega),
48  m_z0(z0),
49  m_tanLambda(tanLambda)
50 {
51 }
52 
53 double Helix::getPerigeeX() const
54 {
55  return getD0() * getSinPhi0();
56 }
57 
58 double Helix::getPerigeeY() const
59 {
60  return -getD0() * getCosPhi0();
61 }
62 
63 double Helix::getPerigeeZ() const
64 {
65  return getZ0();
66 }
67 
68 TVector3 Helix::getPerigee() const
69 {
70  return TVector3(getPerigeeX(), getPerigeeY(), getPerigeeZ());
71 }
72 
73 double Helix::getMomentumX(const double bZ) const
74 {
75  return getCosPhi0() * getTransverseMomentum(bZ);
76 }
77 
78 double Helix::getMomentumY(const double bZ) const
79 {
80  return getSinPhi0() * getTransverseMomentum(bZ);
81 }
82 
83 double Helix::getMomentumZ(const double bZ) const
84 {
85  return getTanLambda() * getTransverseMomentum(bZ);
86 }
87 
88 TVector3 Helix:: getMomentum(const double bZ) const
89 {
90  return TVector3(getMomentumX(bZ), getMomentumY(bZ), getMomentumZ(bZ));
91 }
92 
93 TVector3 Helix::getDirection() const
94 {
95  return TVector3(getCosPhi0(), getSinPhi0(), getTanLambda());
96 }
97 
98 double Helix::getTransverseMomentum(const double bZ) const
99 {
100  return 1 / std::fabs(getAlpha(bZ) * getOmega());
101 }
102 
103 double Helix::getKappa(const double bZ) const
104 {
105  return getOmega() * getAlpha(bZ);
106 }
107 
108 double Helix::getAlpha(const double bZ)
109 {
110  return 1.0 / (bZ * Const::speedOfLight) * 1E4;
111 }
112 
113 short Helix::getChargeSign() const
114 {
115  return boost::math::sign(getOmega());
116 }
117 
118 void Helix::reverse()
119 {
120  // All except z0 have to be taken to their opposites
121  m_d0 = -m_d0; //d0
122  m_phi0 = reversePhi(m_phi0);
123  m_omega = -m_omega;
124  m_tanLambda = -m_tanLambda;
125 }
126 
127 double Helix::getArcLength2DAtCylindricalR(const double& cylindricalR) const
128 {
129  // Slight trick here
130  // Since the sought point is on the helix we treat it as the perigee
131  // and the origin as the point to extrapolate to.
132  // We know the distance of the origin to the circle, which is just d0
133  // The direct distance from the origin to the imaginary perigee is just the given cylindricalR.
134  const double dr = getD0();
135  const double deltaCylindricalR = cylindricalR;
136  const double absArcLength2D = calcArcLength2DAtDeltaCylindricalRAndDr(deltaCylindricalR, dr);
137  return absArcLength2D;
138 }
139 
140 double Helix::getArcLength2DAtXY(const double& x, const double& y) const
141 {
142  double dr = 0;
143  double arcLength2D = 0;
144  calcArcLength2DAndDrAtXY(x, y, arcLength2D, dr);
145  return arcLength2D;
146 }
147 
148 double Helix::getArcLength2DAtNormalPlane(const double& byX, const double& byY,
149  const double& nX, const double& nY) const
150 {
151  // Construct the tangential vector to the plan in xy space
152  const double tX = nY;
153  const double tY = -nX;
154 
155  // Fetch the helix parameters
156  const double omega = getOmega();
157  const double cosPhi0 = getCosPhi0();
158  const double sinPhi0 = getSinPhi0();
159  const double d0 = getD0();
160 
161  // Prepare a delta vector, which is the vector from the perigee point to the support point of the plane
162  // Split it in component parallel and a component orthogonal to tangent at the perigee.
163  const double deltaParallel = byX * cosPhi0 + byY * sinPhi0;
164  const double deltaOrthogonal = byY * cosPhi0 - byX * sinPhi0 + d0;
165  const double deltaCylindricalR = hypot(deltaOrthogonal, deltaParallel);
166  const double deltaCylindricalRSquared = deltaCylindricalR * deltaCylindricalR;
167 
168  const double UOrthogonal = 1 + omega * deltaOrthogonal; // called nu in the Karimaki paper.
169  const double UParallel = omega * deltaParallel;
170 
171  // Some commonly used terms - compare Karimaki 1990
172  const double A = 2 * deltaOrthogonal + omega * deltaCylindricalRSquared;
173 
174  const double tParallel = tX * cosPhi0 + tY * sinPhi0;
175  const double tOrthogonal = tY * cosPhi0 - tX * sinPhi0;
176  const double tCylindricalR = (tX * tX + tY * tY);
177 
178  const double c = A / 2;
179  const double b = UParallel * tParallel + UOrthogonal * tOrthogonal;
180  const double a = omega / 2 * tCylindricalR;
181 
182  const double discriminant = ((double)b) * b - 4 * a * c;
183  const double root = sqrt(discriminant);
184  const double bigSum = (b > 0) ? -b - root : -b + root;
185 
186  const double bigOffset = bigSum / 2 / a;
187  const double smallOffset = 2 * c / bigSum;
188 
189  const double distance1 = hypot(byX + bigOffset * tX, byY + bigOffset * tY);
190  const double distance2 = hypot(byX + smallOffset * tX, byY + smallOffset * tY);
191 
192  if (distance1 < distance2) {
193  return getArcLength2DAtXY(byX + bigOffset * tX, byY + bigOffset * tY);
194  } else {
195  return getArcLength2DAtXY(byX + smallOffset * tX, byY + smallOffset * tY);
196  }
197 }
198 
199 
200 TVector3 Helix::getPositionAtArcLength2D(const double& arcLength2D) const
201 {
202  /*
203  / \ / \ / \
204  | x | | cos phi0 -sin phi0 | | - sin(chi) / omega |
205  | | = | | * | |
206  | y | | sin phi0 cos phi0 | | -(1 - cos(chi)) / omega - d0 |
207  \ / \ / \ /
208 
209  and
210 
211  z = tanLambda * arclength + z0;
212 
213  where chi = -arcLength2D * omega
214 
215  Here arcLength2D means the arc length of the circle in the xy projection
216  traversed in the forward direction.
217  */
218 
219  // First calculating the position assuming the circle center lies on the y axes (phi0 = 0)
220  // Rotate to the right phi position afterwards
221  // Using the sinus cardinalis yields expressions that are smooth in the limit of omega -> 0
222 
223  // Do calculations in double
224  const double chi = -arcLength2D * getOmega();
225  const double chiHalf = chi / 2.0;
226 
227  using boost::math::sinc_pi;
228 
229  const double x = arcLength2D * sinc_pi(chi);
230  const double y = arcLength2D * sinc_pi(chiHalf) * sin(chiHalf) - getD0();
231 
232  // const double z = s * tan lambda + z0
233  const double z = fma((double)arcLength2D, getTanLambda(), getZ0());
234 
235  // Unrotated position
236  TVector3 position(x, y, z);
237 
238  // Rotate to the right phi0 position
239  position.RotateZ(getPhi0());
240 
241  return position;
242 }
243 
244 TVector3 Helix::getTangentialAtArcLength2D(const double& arcLength2D) const
245 {
246  const double omega = getOmega();
247  const double phi0 = getPhi0();
248  const double tanLambda = getTanLambda();
249 
250  const double chi = - omega * arcLength2D;
251 
252  const double tx = cos(chi + phi0);
253  const double ty = sin(chi + phi0);
254  const double tz = tanLambda;
255 
256  TVector3 tangential(tx, ty, tz);
257  return tangential;
258 }
259 
260 TVector3 Helix::getUnitTangentialAtArcLength2D(const double& arcLength2D) const
261 {
262  TVector3 unitTangential = getTangentialAtArcLength2D(arcLength2D);
263  const double norm = hypot(1, getTanLambda());
264  const double invNorm = 1 / norm;
265  unitTangential *= invNorm;
266  return unitTangential;
267 }
268 
269 TVector3 Helix::getMomentumAtArcLength2D(const double& arcLength2D, const double& bz) const
270 {
271  TVector3 momentum = getTangentialAtArcLength2D(arcLength2D);
272  const double pr = getTransverseMomentum(bz);
273  momentum *= pr;
274 
275  return momentum;
276 }
277 
278 double Helix::passiveMoveBy(const double& byX,
279  const double& byY,
280  const double& byZ)
281 {
282  // First calculate the distance of the new origin to the helix in the xy projection
283  double new_d0 = 0;
284  double arcLength2D = 0;
285  calcArcLength2DAndDrAtXY(byX, byY, arcLength2D, new_d0);
286 
287  // Third the new phi0 and z0 can be calculated from the arc length
288  double chi = -arcLength2D * getOmega();
289  double new_phi0 = m_phi0 + chi;
290  double new_z0 = getZ0() - byZ + getTanLambda() * arcLength2D;
291 
293  m_d0 = new_d0;
294  m_phi0 = new_phi0;
295  m_z0 = new_z0;
296 
297  return arcLength2D;
298 }
299 
300 TMatrixD Helix::calcPassiveMoveByJacobian(const TVector3& by, const double expandBelowChi) const
301 {
302  TMatrixD jacobian(5, 5);
303  calcPassiveMoveByJacobian(by.X(), by.Y(), jacobian, expandBelowChi);
304  return jacobian;
305 }
306 
307 
308 void Helix::calcPassiveMoveByJacobian(const double& byX,
309  const double& byY,
310  TMatrixD& jacobian,
311  const double expandBelowChi) const
312 {
313  // 0. Preparations
314  // Initialise the return value to a unit matrix
315  jacobian.UnitMatrix();
316  assert(jacobian.GetNrows() == jacobian.GetNcols());
317  assert(jacobian.GetNrows() == 5 or jacobian.GetNrows() == 6);
318 
319  // Fetch the helix parameters
320  const double omega = getOmega();
321  const double cosPhi0 = getCosPhi0();
322  const double sinPhi0 = getSinPhi0();
323  const double d0 = getD0();
324  const double tanLambda = getTanLambda();
325 
326  // Prepare a delta vector, which is the vector from the perigee point to the new origin
327  // Split it in component parallel and a component orthogonal to tangent at the perigee.
328  const double deltaParallel = byX * cosPhi0 + byY * sinPhi0;
329  const double deltaOrthogonal = byY * cosPhi0 - byX * sinPhi0 + d0;
330  const double deltaCylindricalR = hypot(deltaOrthogonal, deltaParallel);
331  const double deltaCylindricalRSquared = deltaCylindricalR * deltaCylindricalR;
332 
333  // Some commonly used terms - compare Karimaki 1990
334  const double A = 2 * deltaOrthogonal + omega * deltaCylindricalRSquared;
335  const double USquared = 1 + omega * A;
336  const double U = sqrt(USquared);
337  const double UOrthogonal = 1 + omega * deltaOrthogonal; // called nu in the Karimaki paper.
338  const double UParallel = omega * deltaParallel;
339  // Note U is a vector pointing from the middle of the projected circle scaled by a factor omega.
340  const double u = 1 + omega * d0; // just called u in the Karimaki paper.
341 
342  // ---------------------------------------------------------------------------------------
343  // 1. Set the parts related to the xy coordinates
344  // Fills the upper left 3x3 corner of the jacobain
345 
346  // a. Calculate the row related to omega
347  // The row related to omega is not different from the unit jacobian initialised above.
348  // jacobian(iOmega, iOmega) = 1.0; // From UnitMatrix above.
349 
350  // b. Calculate the row related to d0
351  const double new_d0 = A / (1 + U);
352 
353  jacobian(iD0, iOmega) = (deltaCylindricalR + new_d0) * (deltaCylindricalR - new_d0) / 2 / U;
354  jacobian(iD0, iPhi0) = - u * deltaParallel / U;
355  jacobian(iD0, iD0) = UOrthogonal / U;
356 
357  // c. Calculate the row related to phi0
358  // Also calculate the derivatives of the arc length
359  // which are need for the row related to z.
360  const double dArcLength2D_dD0 = - UParallel / USquared;
361  const double dArcLength2D_dPhi0 = (omega * deltaCylindricalRSquared + deltaOrthogonal - d0 * UOrthogonal) / USquared;
362 
363  jacobian(iPhi0, iD0) = - dArcLength2D_dD0 * omega;
364  jacobian(iPhi0, iPhi0) = u * UOrthogonal / USquared;
365  jacobian(iPhi0, iOmega) = -deltaParallel / USquared;
366 
367  // For jacobian(iPhi0, iOmega) we have to dig deeper
368  // since the normal equations have a divergence for omega -> 0.
369  // This hinders not only the straight line case but extrapolations between points which are close together,
370  // like it happens when the current perigee is very close the new perigee.
371  // To have a smooth transition in this limit we have to carefully
372  // factor out the divergent quotients and approximate them with their taylor series.
373 
374  const double chi = -atan2(UParallel, UOrthogonal);
375  double arcLength2D = 0;
376  double dArcLength2D_dOmega = 0;
377 
378  if (fabs(chi) < std::min(expandBelowChi, M_PI / 2.0)) {
379  // Never expand for the far side of the circle by limiting the expandBelow to maximally half a pi.
380  // Close side of the circle
381  double principleArcLength2D = deltaParallel / UOrthogonal;
382  const double dPrincipleArcLength2D_dOmega = - principleArcLength2D * deltaOrthogonal / UOrthogonal;
383 
384  const double x = principleArcLength2D * omega;
385  const double f = calcATanXDividedByX(x);
386  const double df_dx = calcDerivativeOfATanXDividedByX(x);
387 
388  arcLength2D = principleArcLength2D * f;
389  dArcLength2D_dOmega = (f + x * df_dx) * dPrincipleArcLength2D_dOmega + principleArcLength2D * df_dx * principleArcLength2D;
390  } else {
391  // Far side of the circle
392  // If the far side of the circle is a well definied concept, omega is high enough that
393  // we can divide by it.
394  // Otherwise nothing can rescue us since the far side of the circle is so far away that no reasonable extrapolation can be made.
395  arcLength2D = - chi / omega;
396  dArcLength2D_dOmega = (-arcLength2D - jacobian(iPhi0, iOmega)) / omega;
397  }
398 
399  // ---------------------------------------------------------------------------------------
400  // 2. Set the parts related to the z coordinate
401  // Since tanLambda stays the same there are no entries for tanLambda in the jacobian matrix
402  // For the new z0' = z0 + arcLength2D * tanLambda
403  jacobian(iZ0, iD0) = tanLambda * dArcLength2D_dD0;
404  jacobian(iZ0, iPhi0) = tanLambda * dArcLength2D_dPhi0;
405  jacobian(iZ0, iOmega) = tanLambda * dArcLength2D_dOmega;
406  jacobian(iZ0, iZ0) = 1.0; // From UnitMatrix above.
407  jacobian(iZ0, iTanLambda) = arcLength2D;
408 
409  if (jacobian.GetNrows() == 6) {
410  // Also write the derivates of arcLength2D to the jacobian if the matrix size has space for them.
411  jacobian(iArcLength2D, iD0) = dArcLength2D_dD0;
412  jacobian(iArcLength2D, iPhi0) = dArcLength2D_dPhi0;
413  jacobian(iArcLength2D, iOmega) = dArcLength2D_dOmega;
414  jacobian(iArcLength2D, iZ0) = 0.0;
415  jacobian(iArcLength2D, iTanLambda) = 0;
416  jacobian(iArcLength2D, iArcLength2D) = 1.0;
417  }
418 }
419 
420 double Helix::reversePhi(const double& phi)
421 {
422  return std::remainder(phi + M_PI, 2 * M_PI);
423 }
424 
425 double Helix::calcArcLength2DFromSecantLength(const double& secantLength) const
426 {
427  return secantLength * calcSecantLengthToArcLength2DFactor(secantLength);
428 }
429 
430 
431 double Helix::calcSecantLengthToArcLength2DFactor(const double& secantLength) const
432 {
433  double x = secantLength * getOmega() / 2.0;
434  return calcASinXDividedByX(x);
435 }
436 
437 double Helix::calcASinXDividedByX(const double& x)
438 {
439  // Approximation of asin(x) / x
440  // Inspired by BOOST's sinc
441 
442  BOOST_MATH_STD_USING;
443 
444  auto const taylor_n_bound = boost::math::tools::forth_root_epsilon<double>();
445 
446  if (abs(x) >= taylor_n_bound) {
447  if (fabs(x) == 1) {
448  return M_PI / 2.0;
449 
450  } else {
451  return asin(x) / x;
452 
453  }
454 
455  } else {
456  // approximation by taylor series in x at 0 up to order 0
457  double result = 1.0;
458 
459  auto const taylor_0_bound = boost::math::tools::epsilon<double>();
460  if (abs(x) >= taylor_0_bound) {
461  double x2 = x * x;
462  // approximation by taylor series in x at 0 up to order 2
463  result += x2 / 6.0;
464 
465  auto const taylor_2_bound = boost::math::tools::root_epsilon<double>();
466  if (abs(x) >= taylor_2_bound) {
467  // approximation by taylor series in x at 0 up to order 4
468  result += x2 * x2 * (3.0 / 40.0);
469  }
470  }
471  return result;
472  }
473 
474 }
475 
476 double Helix::calcATanXDividedByX(const double& x)
477 {
478  // Approximation of atan(x) / x
479  // Inspired by BOOST's sinc
480 
481  BOOST_MATH_STD_USING;
482 
483  auto const taylor_n_bound = boost::math::tools::forth_root_epsilon<double>();
484 
485  if (abs(x) >= taylor_n_bound) {
486  return atan(x) / x;
487 
488  } else {
489  // approximation by taylor series in x at 0 up to order 0
490  double result = 1.0;
491 
492  auto const taylor_0_bound = boost::math::tools::epsilon<double>();
493  if (abs(x) >= taylor_0_bound) {
494  double x2 = x * x;
495  // approximation by taylor series in x at 0 up to order 2
496  result -= x2 / 3.0;
497 
498  auto const taylor_2_bound = boost::math::tools::root_epsilon<double>();
499  if (abs(x) >= taylor_2_bound) {
500  // approximation by taylor series in x at 0 up to order 4
501  result += x2 * x2 * (1.0 / 5.0);
502  }
503  }
504  return result;
505  }
506 }
507 
508 double Helix::calcDerivativeOfATanXDividedByX(const double& x)
509 {
510  // Approximation of atan(x) / x
511  // Inspired by BOOST's sinc
512 
513  BOOST_MATH_STD_USING;
514 
515  auto const taylor_n_bound = boost::math::tools::forth_root_epsilon<double>();
516 
517  const double x2 = x * x;
518  if (abs(x) >= taylor_n_bound) {
519  const double chi = atan(x);
520  return ((1 - chi / x) / x - chi) / (1 + x2);
521 
522  } else {
523  // approximation by taylor series in x at 0 up to order 0
524  double result = 1.0;
525 
526  auto const taylor_0_bound = boost::math::tools::epsilon<double>();
527  if (abs(x) >= taylor_0_bound) {
528  // approximation by taylor series in x at 0 up to order 2
529  result -= 2.0 * x / 3.0;
530 
531  auto const taylor_2_bound = boost::math::tools::root_epsilon<double>();
532  if (abs(x) >= taylor_2_bound) {
533  // approximation by taylor series in x at 0 up to order 4
534  result += x2 * x * (4.0 / 5.0);
535  }
536  }
537  return result;
538  }
539 }
540 
541 
542 
543 
544 void Helix::calcArcLength2DAndDrAtXY(const double& x, const double& y, double& arcLength2D, double& dr) const
545 {
546  // Prepare common variables
547  const double omega = getOmega();
548  const double cosPhi0 = getCosPhi0();
549  const double sinPhi0 = getSinPhi0();
550  const double d0 = getD0();
551 
552  const double deltaParallel = x * cosPhi0 + y * sinPhi0;
553  const double deltaOrthogonal = y * cosPhi0 - x * sinPhi0 + d0;
554  const double deltaCylindricalR = hypot(deltaOrthogonal, deltaParallel);
555  const double deltaCylindricalRSquared = deltaCylindricalR * deltaCylindricalR;
556 
557  const double A = 2 * deltaOrthogonal + omega * deltaCylindricalRSquared;
558  const double U = sqrt(1 + omega * A);
559  const double UOrthogonal = 1 + omega * deltaOrthogonal; // called nu in the Karimaki paper.
560  const double UParallel = omega * deltaParallel;
561  // Note U is a vector pointing from the middle of the projected circle scaled by a factor omega.
562 
563  // Calculate dr
564  dr = A / (1 + U);
565 
566  // Calculate the absolute value of the arc length
567  const double chi = -atan2(UParallel, UOrthogonal);
568 
569  if (fabs(chi) < M_PI / 8) { // Rough guess where the critical zone for approximations begins
570  // Close side of the circle
571  double principleArcLength2D = deltaParallel / UOrthogonal;
572  arcLength2D = principleArcLength2D * calcATanXDividedByX(principleArcLength2D * omega);
573  } else {
574  // Far side of the circle
575  // If the far side of the circle is a well definied concept meaning that we have big enough omega.
576  arcLength2D = -chi / omega;
577  }
578 }
579 
580 double Helix::calcArcLength2DAtDeltaCylindricalRAndDr(const double& deltaCylindricalR, const double& dr) const
581 {
582  const double omega = getOmega();
583  double secantLength = sqrt((deltaCylindricalR + dr) * (deltaCylindricalR - dr) / (1 + dr * omega));
584  return calcArcLength2DFromSecantLength(secantLength);
585 }
586 
587 void Helix::setCartesian(const TVector3& position,
588  const TVector3& momentum,
589  const short int charge,
590  const double bZ)
591 {
592  assert(abs(charge) <= 1); // Not prepared for doubly-charged particles.
593  const double alpha = getAlpha(bZ);
594 
595  // We allow for the case that position, momentum are not given
596  // exactly in the perigee. Therefore we have to transform the momentum
597  // with the position as the reference point and then move the coordinate system
598  // to the origin.
599 
600  const double x = position.X();
601  const double y = position.Y();
602  const double z = position.Z();
603 
604  const double px = momentum.X();
605  const double py = momentum.Y();
606  const double pz = momentum.Z();
607 
608  const double ptinv = 1 / hypot(px, py);
609  const double omega = charge * ptinv / alpha;
610  const double tanLambda = ptinv * pz;
611  const double phi0 = atan2(py, px);
612  const double z0 = z;
613  const double d0 = 0;
614 
615  m_omega = omega;
616  m_phi0 = phi0;
617  m_d0 = d0;
618  m_tanLambda = tanLambda;
619  m_z0 = z0;
620 
621  passiveMoveBy(-x, -y, 0);
622 }
623 
624 namespace Belle2 {
630  std::ostream& operator<<(std::ostream& output, const Helix& helix)
631  {
632  return output
633  << "Helix("
634  << "d0=" << helix.getD0() << ", "
635  << "phi0=" << helix.getPhi0() << ", "
636  << "omega=" << helix.getOmega() << ", "
637  << "z0=" << helix.getZ0() << ", "
638  << "tanLambda=" << helix.getTanLambda() << ")";
639  }
641 }
Helix parameter class.
Definition: Helix.h:48
Helix(const HepPoint3D &pivot, const HepVector &a, const HepSymMatrix &Ea)
Constructor with pivot, helix parameter a, and its error matrix.
Definition: Helix.cc:132
double phi0(void) const
Return helix parameter phi0.
Definition: Helix.h:388
const HepVector & a(void) const
Returns helix parameters.
Definition: Helix.h:428
double dr(void) const
Return helix parameter dr.
Definition: Helix.h:380
Hep3Vector momentum(double dPhi=0.) const
returns momentum vector after rotating angle dPhi in phi direction.
Definition: Helix.cc:259
HepPoint3D x(double dPhi=0.) const
returns position after rotating angle dPhi in phi direction.
Definition: Helix.cc:197
double cosPhi0(void) const
Return cos phi0.
Definition: Helix.h:500
double sinPhi0(void) const
Return sin phi0.
Definition: Helix.h:492
static const double speedOfLight
[cm/ns]
Definition: Const.h:575
std::ostream & operator<<(std::ostream &output, const IntervalOfValidity &iov)
Abstract base class for different kinds of events.