8 #include <framework/dataobjects/Helix.h>
10 #include <framework/gearbox/Const.h>
12 #include <boost/math/special_functions/sign.hpp>
13 #include <boost/math/special_functions/sinc.hpp>
14 #include <boost/math/tools/precision.hpp>
21 using namespace HelixParameterIndex;
32 Helix::Helix(
const TVector3& position,
33 const TVector3& momentum,
34 const short int charge,
37 setCartesian(position,
momentum, charge, bZ);
44 const double& tanLambda)
49 m_tanLambda(tanLambda)
53 double Helix::getPerigeeX()
const
55 return getD0() * getSinPhi0();
58 double Helix::getPerigeeY()
const
60 return -getD0() * getCosPhi0();
63 double Helix::getPerigeeZ()
const
68 TVector3 Helix::getPerigee()
const
70 return TVector3(getPerigeeX(), getPerigeeY(), getPerigeeZ());
73 double Helix::getMomentumX(
const double bZ)
const
75 return getCosPhi0() * getTransverseMomentum(bZ);
78 double Helix::getMomentumY(
const double bZ)
const
80 return getSinPhi0() * getTransverseMomentum(bZ);
83 double Helix::getMomentumZ(
const double bZ)
const
85 return getTanLambda() * getTransverseMomentum(bZ);
88 TVector3 Helix:: getMomentum(
const double bZ)
const
90 return TVector3(getMomentumX(bZ), getMomentumY(bZ), getMomentumZ(bZ));
93 TVector3 Helix::getDirection()
const
95 return TVector3(getCosPhi0(), getSinPhi0(), getTanLambda());
98 double Helix::getTransverseMomentum(
const double bZ)
const
100 return 1 / std::fabs(getAlpha(bZ) * getOmega());
103 double Helix::getKappa(
const double bZ)
const
105 return getOmega() * getAlpha(bZ);
108 double Helix::getAlpha(
const double bZ)
113 short Helix::getChargeSign()
const
115 return boost::math::sign(getOmega());
118 void Helix::reverse()
122 m_phi0 = reversePhi(m_phi0);
124 m_tanLambda = -m_tanLambda;
127 double Helix::getArcLength2DAtCylindricalR(
const double& cylindricalR)
const
134 const double dr = getD0();
135 const double deltaCylindricalR = cylindricalR;
136 const double absArcLength2D = calcArcLength2DAtDeltaCylindricalRAndDr(deltaCylindricalR,
dr);
137 return absArcLength2D;
140 double Helix::getArcLength2DAtXY(
const double& x,
const double& y)
const
143 double arcLength2D = 0;
144 calcArcLength2DAndDrAtXY(
x, y, arcLength2D,
dr);
148 double Helix::getArcLength2DAtNormalPlane(
const double& byX,
const double& byY,
149 const double& nX,
const double& nY)
const
152 const double tX = nY;
153 const double tY = -nX;
156 const double omega = getOmega();
157 const double cosPhi0 = getCosPhi0();
158 const double sinPhi0 = getSinPhi0();
159 const double d0 = getD0();
164 const double deltaOrthogonal = byY *
cosPhi0 - byX *
sinPhi0 + d0;
165 const double deltaCylindricalR = hypot(deltaOrthogonal, deltaParallel);
166 const double deltaCylindricalRSquared = deltaCylindricalR * deltaCylindricalR;
168 const double UOrthogonal = 1 + omega * deltaOrthogonal;
169 const double UParallel = omega * deltaParallel;
172 const double A = 2 * deltaOrthogonal + omega * deltaCylindricalRSquared;
176 const double tCylindricalR = (tX * tX + tY * tY);
178 const double c = A / 2;
179 const double b = UParallel * tParallel + UOrthogonal * tOrthogonal;
180 const double a = omega / 2 * tCylindricalR;
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;
186 const double bigOffset = bigSum / 2 /
a;
187 const double smallOffset = 2 * c / bigSum;
189 const double distance1 = hypot(byX + bigOffset * tX, byY + bigOffset * tY);
190 const double distance2 = hypot(byX + smallOffset * tX, byY + smallOffset * tY);
192 if (distance1 < distance2) {
193 return getArcLength2DAtXY(byX + bigOffset * tX, byY + bigOffset * tY);
195 return getArcLength2DAtXY(byX + smallOffset * tX, byY + smallOffset * tY);
200 TVector3 Helix::getPositionAtArcLength2D(
const double& arcLength2D)
const
224 const double chi = -arcLength2D * getOmega();
225 const double chiHalf = chi / 2.0;
227 using boost::math::sinc_pi;
229 const double x = arcLength2D * sinc_pi(chi);
230 const double y = arcLength2D * sinc_pi(chiHalf) * sin(chiHalf) - getD0();
233 const double z = fma((
double)arcLength2D, getTanLambda(), getZ0());
236 TVector3 position(
x, y, z);
239 position.RotateZ(getPhi0());
244 TVector3 Helix::getTangentialAtArcLength2D(
const double& arcLength2D)
const
246 const double omega = getOmega();
247 const double phi0 = getPhi0();
248 const double tanLambda = getTanLambda();
250 const double chi = - omega * arcLength2D;
252 const double tx = cos(chi +
phi0);
253 const double ty = sin(chi +
phi0);
254 const double tz = tanLambda;
256 TVector3 tangential(tx, ty, tz);
260 TVector3 Helix::getUnitTangentialAtArcLength2D(
const double& arcLength2D)
const
262 TVector3 unitTangential = getTangentialAtArcLength2D(arcLength2D);
263 const double norm = hypot(1, getTanLambda());
264 const double invNorm = 1 / norm;
265 unitTangential *= invNorm;
266 return unitTangential;
269 TVector3 Helix::getMomentumAtArcLength2D(
const double& arcLength2D,
const double& bz)
const
271 TVector3
momentum = getTangentialAtArcLength2D(arcLength2D);
272 const double pr = getTransverseMomentum(bz);
278 double Helix::passiveMoveBy(
const double& byX,
284 double arcLength2D = 0;
285 calcArcLength2DAndDrAtXY(byX, byY, arcLength2D, new_d0);
288 double chi = -arcLength2D * getOmega();
289 double new_phi0 = m_phi0 + chi;
290 double new_z0 = getZ0() - byZ + getTanLambda() * arcLength2D;
300 TMatrixD Helix::calcPassiveMoveByJacobian(
const TVector3& by,
const double expandBelowChi)
const
302 TMatrixD jacobian(5, 5);
303 calcPassiveMoveByJacobian(by.X(), by.Y(), jacobian, expandBelowChi);
308 void Helix::calcPassiveMoveByJacobian(
const double& byX,
311 const double expandBelowChi)
const
315 jacobian.UnitMatrix();
316 assert(jacobian.GetNrows() == jacobian.GetNcols());
317 assert(jacobian.GetNrows() == 5 or jacobian.GetNrows() == 6);
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();
329 const double deltaOrthogonal = byY *
cosPhi0 - byX *
sinPhi0 + d0;
330 const double deltaCylindricalR = hypot(deltaOrthogonal, deltaParallel);
331 const double deltaCylindricalRSquared = deltaCylindricalR * deltaCylindricalR;
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;
338 const double UParallel = omega * deltaParallel;
340 const double u = 1 + omega * d0;
351 const double new_d0 = A / (1 + U);
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;
360 const double dArcLength2D_dD0 = - UParallel / USquared;
361 const double dArcLength2D_dPhi0 = (omega * deltaCylindricalRSquared + deltaOrthogonal - d0 * UOrthogonal) / USquared;
363 jacobian(iPhi0, iD0) = - dArcLength2D_dD0 * omega;
364 jacobian(iPhi0, iPhi0) = u * UOrthogonal / USquared;
365 jacobian(iPhi0, iOmega) = -deltaParallel / USquared;
374 const double chi = -atan2(UParallel, UOrthogonal);
375 double arcLength2D = 0;
376 double dArcLength2D_dOmega = 0;
378 if (fabs(chi) < std::min(expandBelowChi, M_PI / 2.0)) {
381 double principleArcLength2D = deltaParallel / UOrthogonal;
382 const double dPrincipleArcLength2D_dOmega = - principleArcLength2D * deltaOrthogonal / UOrthogonal;
384 const double x = principleArcLength2D * omega;
385 const double f = calcATanXDividedByX(
x);
386 const double df_dx = calcDerivativeOfATanXDividedByX(
x);
388 arcLength2D = principleArcLength2D * f;
389 dArcLength2D_dOmega = (f +
x * df_dx) * dPrincipleArcLength2D_dOmega + principleArcLength2D * df_dx * principleArcLength2D;
395 arcLength2D = - chi / omega;
396 dArcLength2D_dOmega = (-arcLength2D - jacobian(iPhi0, iOmega)) / omega;
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;
407 jacobian(iZ0, iTanLambda) = arcLength2D;
409 if (jacobian.GetNrows() == 6) {
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;
420 double Helix::reversePhi(
const double& phi)
422 return std::remainder(phi + M_PI, 2 * M_PI);
425 double Helix::calcArcLength2DFromSecantLength(
const double& secantLength)
const
427 return secantLength * calcSecantLengthToArcLength2DFactor(secantLength);
431 double Helix::calcSecantLengthToArcLength2DFactor(
const double& secantLength)
const
433 double x = secantLength * getOmega() / 2.0;
434 return calcASinXDividedByX(
x);
437 double Helix::calcASinXDividedByX(
const double& x)
442 BOOST_MATH_STD_USING;
444 auto const taylor_n_bound = boost::math::tools::forth_root_epsilon<double>();
446 if (abs(
x) >= taylor_n_bound) {
459 auto const taylor_0_bound = boost::math::tools::epsilon<double>();
460 if (abs(
x) >= taylor_0_bound) {
465 auto const taylor_2_bound = boost::math::tools::root_epsilon<double>();
466 if (abs(
x) >= taylor_2_bound) {
468 result += x2 * x2 * (3.0 / 40.0);
476 double Helix::calcATanXDividedByX(
const double& x)
481 BOOST_MATH_STD_USING;
483 auto const taylor_n_bound = boost::math::tools::forth_root_epsilon<double>();
485 if (abs(
x) >= taylor_n_bound) {
492 auto const taylor_0_bound = boost::math::tools::epsilon<double>();
493 if (abs(
x) >= taylor_0_bound) {
498 auto const taylor_2_bound = boost::math::tools::root_epsilon<double>();
499 if (abs(
x) >= taylor_2_bound) {
501 result += x2 * x2 * (1.0 / 5.0);
508 double Helix::calcDerivativeOfATanXDividedByX(
const double& x)
513 BOOST_MATH_STD_USING;
515 auto const taylor_n_bound = boost::math::tools::forth_root_epsilon<double>();
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);
526 auto const taylor_0_bound = boost::math::tools::epsilon<double>();
527 if (abs(
x) >= taylor_0_bound) {
529 result -= 2.0 *
x / 3.0;
531 auto const taylor_2_bound = boost::math::tools::root_epsilon<double>();
532 if (abs(
x) >= taylor_2_bound) {
534 result += x2 *
x * (4.0 / 5.0);
544 void Helix::calcArcLength2DAndDrAtXY(
const double& x,
const double& y,
double& arcLength2D,
double& dr)
const
547 const double omega = getOmega();
548 const double cosPhi0 = getCosPhi0();
549 const double sinPhi0 = getSinPhi0();
550 const double d0 = getD0();
554 const double deltaCylindricalR = hypot(deltaOrthogonal, deltaParallel);
555 const double deltaCylindricalRSquared = deltaCylindricalR * deltaCylindricalR;
557 const double A = 2 * deltaOrthogonal + omega * deltaCylindricalRSquared;
558 const double U = sqrt(1 + omega * A);
559 const double UOrthogonal = 1 + omega * deltaOrthogonal;
560 const double UParallel = omega * deltaParallel;
567 const double chi = -atan2(UParallel, UOrthogonal);
569 if (fabs(chi) < M_PI / 8) {
571 double principleArcLength2D = deltaParallel / UOrthogonal;
572 arcLength2D = principleArcLength2D * calcATanXDividedByX(principleArcLength2D * omega);
576 arcLength2D = -chi / omega;
580 double Helix::calcArcLength2DAtDeltaCylindricalRAndDr(
const double& deltaCylindricalR,
const double& dr)
const
582 const double omega = getOmega();
583 double secantLength = sqrt((deltaCylindricalR +
dr) * (deltaCylindricalR -
dr) / (1 +
dr * omega));
584 return calcArcLength2DFromSecantLength(secantLength);
587 void Helix::setCartesian(
const TVector3& position,
588 const TVector3& momentum,
589 const short int charge,
592 assert(abs(charge) <= 1);
593 const double alpha = getAlpha(bZ);
600 const double x = position.X();
601 const double y = position.Y();
602 const double z = position.Z();
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);
618 m_tanLambda = tanLambda;
621 passiveMoveBy(-
x, -y, 0);
634 <<
"d0=" << helix.getD0() <<
", "
635 <<
"phi0=" << helix.getPhi0() <<
", "
636 <<
"omega=" << helix.getOmega() <<
", "
637 <<
"z0=" << helix.getZ0() <<
", "
638 <<
"tanLambda=" << helix.getTanLambda() <<
")";
Helix(const HepPoint3D &pivot, const HepVector &a, const HepSymMatrix &Ea)
Constructor with pivot, helix parameter a, and its error matrix.
double phi0(void) const
Return helix parameter phi0.
const HepVector & a(void) const
Returns helix parameters.
double dr(void) const
Return helix parameter dr.
Hep3Vector momentum(double dPhi=0.) const
returns momentum vector after rotating angle dPhi in phi direction.
HepPoint3D x(double dPhi=0.) const
returns position after rotating angle dPhi in phi direction.
double cosPhi0(void) const
Return cos phi0.
double sinPhi0(void) const
Return sin phi0.
static const double speedOfLight
[cm/ns]
std::ostream & operator<<(std::ostream &output, const IntervalOfValidity &iov)
Abstract base class for different kinds of events.