1 #ifndef RecoTracker_PixelTrackFitting_RiemannFit_h 2 #define RecoTracker_PixelTrackFitting_RiemannFit_h 5 #include <Eigen/Eigenvalues> 33 template <
typename VNd1,
typename VNd2>
38 uint n = length_values.rows();
39 rad_lengths(0) = length_values(0) * xx_0_inv;
41 rad_lengths(
j) =
std::abs(length_values(
j) - length_values(
j - 1)) * xx_0_inv;
64 template <
typename V4,
typename VNd1,
typename VNd2,
int N>
77 double p_2 = p_t * p_t * (1. + 1. /
sqr(
fast_fit(3)));
82 VectorNd<N> s_values = s_arcs.array() * s_arcs.array() + z_values.array() * z_values.array();
83 s_values = s_values.array().sqrt();
86 sig2_S = .000225 / p_2 * (1. + 0.038 * rad_lengths_S.array().log()).abs2() * rad_lengths_S.array();
92 tmp(
k,
k) = cov_sz[
k](0, 0);
125 template <
typename M2xN,
typename V4,
int N>
132 double p_2 = p_t * p_t * (1. + 1. /
sqr(
fast_fit(3)));
141 Vector2d pVec = p2D.block(0,
i, 2, 1) - oVec;
143 const double dot = (-oVec).
dot(pVec);
144 const double tempAtan2 = atan2(
cross,
dot);
149 VectorNd<N> sig2 = (1. + 0.038 * rad_lengths.array().log()).abs2() * rad_lengths.array();
173 template <
typename M2xN,
int N>
178 printf(
"Address of p2D: %p\n", &p2D);
180 printIt(&p2D,
"cov_radtocart - p2D:");
184 printIt(&rad_inv,
"cov_radtocart - rad_inv:");
187 cov_cart(
i,
j) = cov_rad(
i,
j) * p2D(1,
i) * rad_inv(
i) * p2D(1,
j) * rad_inv(
j);
188 cov_cart(
i +
n,
j +
n) = cov_rad(
i,
j) * p2D(0,
i) * rad_inv(
i) * p2D(0,
j) * rad_inv(
j);
189 cov_cart(
i,
j +
n) = -cov_rad(
i,
j) * p2D(1,
i) * rad_inv(
i) * p2D(0,
j) * rad_inv(
j);
190 cov_cart(
i +
n,
j) = -cov_rad(
i,
j) * p2D(0,
i) * rad_inv(
i) * p2D(1,
j) * rad_inv(
j);
191 cov_cart(
j,
i) = cov_cart(
i,
j);
192 cov_cart(
j +
n,
i +
n) = cov_cart(
i +
n,
j +
n);
193 cov_cart(
j +
n,
i) = cov_cart(
i,
j +
n);
194 cov_cart(
j,
i +
n) = cov_cart(
i +
n,
j);
210 template <
typename M2xN,
int N>
216 const VectorNd<N> rad_inv2 = rad.cwiseInverse().array().square();
220 cov_rad(
i) = cov_cart(
i,
i);
222 cov_rad(
i) = rad_inv2(
i) * (cov_cart(
i,
i) *
sqr(p2D(1,
i)) + cov_cart(
i +
n,
i +
n) *
sqr(p2D(0,
i)) -
223 2. * cov_cart(
i,
i +
n) * p2D(0,
i) * p2D(1,
i));
241 template <
typename M2xN,
typename V4,
int N>
251 cov_rad(
i) = cov_cart(
i,
i);
255 const double x2 =
a.dot(
b);
257 const double tan_c = -
y2 / x2;
258 const double tan_c2 =
sqr(tan_c);
260 1. / (1. + tan_c2) * (cov_cart(
i,
i) + cov_cart(
i +
n,
i +
n) * tan_c2 + 2 * cov_cart(
i,
i +
n) * tan_c);
279 return cov_rad_inv.colwise().sum().transpose();
290 template <
typename M2xN>
292 return ((p2D(0, 1) - p2D(0, 0)) * (par_uvr.y() - p2D(1, 0)) - (p2D(1, 1) - p2D(1, 0)) * (par_uvr.x() - p2D(0, 0)) >
315 printf(
"min_eigen3D - enter\n");
317 Eigen::SelfAdjointEigenSolver<Matrix3d> solver(3);
318 solver.computeDirect(
A);
320 chi2 = solver.eigenvalues().minCoeff(&min_index);
322 printf(
"min_eigen3D - exit\n");
324 return solver.eigenvectors().col(min_index);
339 Eigen::SelfAdjointEigenSolver<Matrix3f> solver(3);
340 solver.computeDirect(
A.cast<
float>());
342 solver.eigenvalues().minCoeff(&min_index);
343 return solver.eigenvectors().col(min_index).cast<
double>();
357 Eigen::SelfAdjointEigenSolver<Matrix2d> solver(2);
358 solver.computeDirect(aMat);
360 chi2 = solver.eigenvalues().minCoeff(&min_index);
361 return solver.eigenvectors().col(min_index);
377 template <
typename M3xN,
typename V4>
379 constexpr uint32_t
N = M3xN::ColsAtCompileTime;
387 printIt(&bVec,
"Fast_fit - b: ");
388 printIt(&cVec,
"Fast_fit - c: ");
390 auto b2 = bVec.squaredNorm();
391 auto c2 = cVec.squaredNorm();
399 auto bx =
flip ? bVec.y() : bVec.x();
400 auto by =
flip ? bVec.x() : bVec.y();
401 auto cx =
flip ? cVec.y() : cVec.x();
402 auto cy =
flip ? cVec.x() : cVec.y();
404 auto div = 2. * (cx * by -
bx * cy);
406 auto y0 = (cx *
b2 -
bx * c2) / div;
407 auto x0 = (0.5 *
b2 - y0 * by) /
bx;
457 template <
typename M2xN,
typename V4,
int N>
465 printf(
"circle_fit - enter\n");
470 printIt(&hits2D,
"circle_fit - hits2D:");
471 printIt(&hits_cov2D,
"circle_fit - hits_cov2D:");
474 printf(
"circle_fit - WEIGHT COMPUTATION\n");
483 printIt(&scatterCovRadMat,
"circle_fit - scatter_cov_rad:");
484 printIt(&hits2D,
"circle_fit - hits2D bis:");
486 printf(
"Address of hits2D: a) %p\n", &hits2D);
489 printIt(&vMat,
"circle_fit - V:");
490 cov_rad += scatterCovRadMat;
491 printIt(&cov_rad,
"circle_fit - cov_rad:");
502 printf(
"circle_fit - SPACE TRANSFORMATION\n");
507 printf(
"Address of hits2D: b) %p\n", &hits2D);
509 const Vector2d hCentroid = hits2D.rowwise().mean();
510 printIt(&hCentroid,
"circle_fit - h_:");
512 p3D.block(0, 0, 2,
n) = hits2D.colwise() - hCentroid;
513 printIt(&p3D,
"circle_fit - p3D: a)");
515 mc << p3D.row(0).transpose(), p3D.row(1).transpose();
516 printIt(&
mc,
"circle_fit - mc(centered hits):");
519 const double tempQ =
mc.squaredNorm();
520 const double tempS =
sqrt(
n * 1. / tempQ);
521 p3D.block(0, 0, 2,
n) *= tempS;
524 p3D.row(2) = p3D.block(0, 0, 2,
n).colwise().squaredNorm();
525 printIt(&p3D,
"circle_fit - p3D: b)");
528 printf(
"circle_fit - COST FUNCTION\n");
534 r0.noalias() = p3D *
weight;
536 Matrix3d aMat = xMat * gMat * xMat.transpose();
537 printIt(&aMat,
"circle_fit - A:");
540 printf(
"circle_fit - MINIMIZE\n");
546 printf(
"circle_fit - AFTER MIN_EIGEN\n");
548 printIt(&vVec,
"v BEFORE INVERSION");
549 vVec *= (vVec(2) > 0) ? 1 : -1;
550 printIt(&vVec,
"v AFTER INVERSION");
554 printf(
"circle_fit - AFTER MIN_EIGEN 1\n");
556 Eigen::Matrix<double, 1, 1> cm;
558 printf(
"circle_fit - AFTER MIN_EIGEN 2\n");
560 cm = -vVec.transpose() * r0;
562 printf(
"circle_fit - AFTER MIN_EIGEN 3\n");
564 const double tempC = cm(0, 0);
567 printf(
"circle_fit - COMPUTE CIRCLE PARAMETER\n");
572 const double tempH =
sqrt(1. -
sqr(vVec(2)) - 4. * tempC * vVec(2));
573 const double v2x2_inv = 1. / (2. * vVec(2));
574 const double s_inv = 1. / tempS;
576 par_uvr << -vVec(0) * v2x2_inv, -vVec(1) * v2x2_inv, tempH * v2x2_inv;
579 circle.
par << par_uvr(0) * s_inv + hCentroid(0), par_uvr(1) * s_inv + hCentroid(1), par_uvr(2) * s_inv;
581 circle.
chi2 =
abs(
chi2) * renorm /
sqr(2 * vVec(2) * par_uvr(2) * tempS);
582 printIt(&circle.
par,
"circle_fit - CIRCLE PARAMETERS:");
583 printIt(&circle.
cov,
"circle_fit - CIRCLE COVARIANCE:");
585 printf(
"circle_fit - CIRCLE CHARGE: %d\n", circle.
qCharge);
589 printf(
"circle_fit - ERROR PROPAGATION\n");
594 printf(
"circle_fit - ERROR PRPAGATION ACTIVATED\n");
599 printf(
"circle_fit - ERROR PRPAGATION ACTIVATED 2\n");
602 Eigen::Matrix<double, 1, 1> cm;
603 Eigen::Matrix<double, 1, 1> cm2;
604 cm =
mc.transpose() * vMat *
mc;
605 const double tempC2 = cm(0, 0);
607 tempVcsMat.template triangularView<Eigen::Upper>() =
608 (
sqr(tempS) * vMat +
sqr(
sqr(tempS)) * 1. / (4. * tempQ *
n) *
609 (2. * vMat.squaredNorm() + 4. * tempC2) *
610 (
mc *
mc.transpose()));
612 printIt(&tempVcsMat,
"circle_fit - Vcs:");
613 cMat[0][0] = tempVcsMat.block(0, 0,
n,
n).template selfadjointView<Eigen::Upper>();
614 vcsMat[0][1] = tempVcsMat.block(0,
n,
n,
n);
615 cMat[1][1] = tempVcsMat.block(
n,
n,
n,
n).template selfadjointView<Eigen::Upper>();
616 vcsMat[1][0] = vcsMat[0][1].transpose();
617 printIt(&tempVcsMat,
"circle_fit - Vcs:");
623 const ArrayNd<N> t00 = p3D.row(0).transpose() * p3D.row(0);
624 const ArrayNd<N> t01 = p3D.row(0).transpose() * p3D.row(1);
625 const ArrayNd<N> t11 = p3D.row(1).transpose() * p3D.row(1);
627 vcsMat[0][0] = cMat[0][0];
628 cMat[0][1] = vcsMat[0][1];
629 cMat[0][2] = 2. * (vcsMat[0][0] *
t0 + vcsMat[0][1] *
t1);
630 vcsMat[1][1] = cMat[1][1];
631 cMat[1][2] = 2. * (vcsMat[1][0] *
t0 + vcsMat[1][1] *
t1);
633 tmp.template triangularView<Eigen::Upper>() =
634 (2. * (vcsMat[0][0] * vcsMat[0][0] + vcsMat[0][0] * vcsMat[0][1] + vcsMat[1][1] * vcsMat[1][0] +
635 vcsMat[1][1] * vcsMat[1][1]) +
636 4. * (vcsMat[0][0] * t00 + vcsMat[0][1] * t01 + vcsMat[1][0] * t10 + vcsMat[1][1] * t11))
638 cMat[2][2] =
tmp.template selfadjointView<Eigen::Upper>();
640 printIt(&cMat[0][0],
"circle_fit - C[0][0]:");
645 Eigen::Matrix<double, 1, 1>
tmp;
648 const double tempC =
tmp(0, 0);
650 c0Mat(
j,
i) = c0Mat(
i,
j);
653 printIt(&c0Mat,
"circle_fit - C0:");
658 printIt(&wMat,
"circle_fit - W:");
659 printIt(&hMat,
"circle_fit - H:");
660 printIt(&s_v,
"circle_fit - s_v:");
663 dMat[0][0] = (hMat * cMat[0][0] * hMat.transpose()).cwiseProduct(wMat);
664 dMat[0][1] = (hMat * cMat[0][1] * hMat.transpose()).cwiseProduct(wMat);
665 dMat[0][2] = (hMat * cMat[0][2] * hMat.transpose()).cwiseProduct(wMat);
666 dMat[1][1] = (hMat * cMat[1][1] * hMat.transpose()).cwiseProduct(wMat);
667 dMat[1][2] = (hMat * cMat[1][2] * hMat.transpose()).cwiseProduct(wMat);
668 dMat[2][2] = (hMat * cMat[2][2] * hMat.transpose()).cwiseProduct(wMat);
669 dMat[1][0] = dMat[0][1].transpose();
670 dMat[2][0] = dMat[0][2].transpose();
671 dMat[2][1] = dMat[1][2].transpose();
672 printIt(&dMat[0][0],
"circle_fit - D_[0][0]:");
674 constexpr uint nu[6][2] = {{0, 0}, {0, 1}, {0, 2}, {1, 1}, {1, 2}, {2, 2}};
678 const uint i = nu[
a][0],
j = nu[
a][1];
680 const uint k = nu[
b][0],
l = nu[
b][1];
684 t0 = 2. * dMat[
j][
l] * s_v.col(
l);
688 t1 = 2. * dMat[
i][
l] * s_v.col(
l);
690 t0 = dMat[
j][
l] * s_v.col(
k) + dMat[
j][
k] * s_v.col(
l);
694 t1 = dMat[
i][
l] * s_v.col(
k) + dMat[
i][
k] * s_v.col(
l);
698 Eigen::Matrix<double, 1, 1> cm;
699 cm = s_v.col(
i).transpose() * (
t0 +
t1);
701 const double tempC = cm(0, 0);
702 eMat(
a,
b) = 0. + tempC;
704 Eigen::Matrix<double, 1, 1> cm;
705 cm = (s_v.col(
i).transpose() *
t0) + (s_v.col(
j).transpose() *
t1);
707 const double tempC = cm(0, 0);
708 eMat(
a,
b) = 0. + tempC;
711 eMat(
b,
a) = eMat(
a,
b);
714 printIt(&eMat,
"circle_fit - E:");
716 Eigen::Matrix<double, 3, 6> j2Mat;
718 const uint i = nu[
a][0],
j = nu[
a][1];
722 const int sign = (j2Mat.col(
a)(2) > 0) ? 1 : -1;
725 printIt(&j2Mat,
"circle_fit - J2:");
729 Matrix3d t0 = j2Mat * eMat * j2Mat.transpose();
731 cvcMat.block(0, 0, 3, 3) =
t0;
732 cvcMat.block(0, 3, 3, 1) =
t1;
733 cvcMat.block(3, 0, 1, 3) =
t1.transpose();
734 Eigen::Matrix<double, 1, 1> cm1;
735 Eigen::Matrix<double, 1, 1> cm3;
736 cm1 = (vVec.transpose() * c0Mat * vVec);
738 cm3 = (r0.transpose() *
t0 * r0);
740 const double tempC = cm1(0, 0) + (c0Mat.cwiseProduct(
t0)).sum() + cm3(0, 0);
741 cvcMat(3, 3) = tempC;
744 printIt(&cvcMat,
"circle_fit - Cvc:");
746 Eigen::Matrix<double, 3, 4> j3Mat;
748 const double t = 1. / tempH;
749 j3Mat << -v2x2_inv, 0, vVec(0) *
sqr(v2x2_inv) * 2., 0, 0, -v2x2_inv, vVec(1) *
sqr(v2x2_inv) * 2., 0,
750 vVec(0) * v2x2_inv *
t, vVec(1) * v2x2_inv *
t,
751 -tempH *
sqr(v2x2_inv) * 2. - (2. * tempC + vVec(2)) * v2x2_inv *
t, -
t;
753 printIt(&j3Mat,
"circle_fit - J3:");
756 printIt(&Jq,
"circle_fit - Jq:");
758 Matrix3d cov_uvr = j3Mat * cvcMat * j3Mat.transpose() *
sqr(s_inv)
759 + (par_uvr * par_uvr.transpose()) * (Jq * vMat * Jq.transpose());
761 circle.
cov = cov_uvr;
766 printf(
"circle_fit - exit\n");
787 template <
typename M3xN,
typename M6xN,
typename V4>
794 constexpr uint32_t
N = M3xN::ColsAtCompileTime;
800 Eigen::Matrix<double, 2, 2>
rot;
812 Eigen::Matrix<double, 2, 6> jxMat;
815 printf(
"Line_fit - B: %g\n",
bField);
817 printIt(&hits_ge,
"Line_fit covs: ");
833 const double dot = (-oVec).
dot(pVec);
838 p2D(0,
i) = tempQAtan2 * circle.
par(2);
842 double d_X0 = 0., d_Y0 = 0., d_R = 0.;
844 d_X0 = -temp0 * ((pVec(1) + oVec(1)) *
dot - (pVec(0) - oVec(0)) *
cross);
845 d_Y0 = temp0 * ((pVec(0) + oVec(0)) *
dot - (oVec(1) - pVec(1)) *
cross);
848 const double d_x = temp0 * (oVec(1) *
dot + oVec(0) *
cross);
849 const double d_y = temp0 * (-oVec(0) *
dot + oVec(1) *
cross);
850 jxMat << d_X0, d_Y0, d_R,
d_x,
d_y, 0., 0., 0., 0., 0., 0., 1.;
852 covMat.block(0, 0, 3, 3) = circle.
cov;
853 covMat(3, 3) = hits_ge.col(
i)[0];
854 covMat(4, 4) = hits_ge.col(
i)[2];
855 covMat(5, 5) = hits_ge.col(
i)[5];
856 covMat(3, 4) = covMat(4, 3) = hits_ge.col(
i)[1];
857 covMat(3, 5) = covMat(5, 3) = hits_ge.col(
i)[3];
858 covMat(4, 5) = covMat(5, 4) = hits_ge.col(
i)[4];
860 cov_sz[
i].noalias() =
rot *
tmp *
rot.transpose();
863 p2D.row(1) =
hits.row(2);
870 printIt(cov_sz,
"line_fit - cov_sz:");
871 printIt(&cov_with_ms,
"line_fit - cov_with_ms: ");
878 printf(
"Fast fit Tan(theta): %g\n",
fast_fit(3));
879 printf(
"Rotation angle: %g\n",
theta);
881 printIt(&p2D,
"Original Hits(s,z):");
882 printIt(&p2D_rot,
"Rotated hits(S3D, Z'):");
888 aMat << MatrixXd::Ones(1,
n), p2D_rot.row(0);
898 Eigen::Matrix<double, 2, 2> covParamsMat = aMat * vyInvMat * aMat.transpose();
905 Eigen::Matrix<double, 2, 1>
sol = covParamsMat * aMat * vyInvMat * p2D_rot.row(1).transpose();
914 auto common_factor = 1. / (sinTheta -
sol(1, 0) * cosTheta);
915 Eigen::Matrix<double, 2, 2> jMat;
916 jMat << 0., common_factor * common_factor, common_factor,
sol(0, 0) * cosTheta * common_factor * common_factor;
918 double tempM = common_factor * (
sol(1, 0) * sinTheta + cosTheta);
919 double tempQ = common_factor *
sol(0, 0);
920 auto cov_mq = jMat * covParamsMat * jMat.transpose();
923 double chi2 =
res.transpose() * vyInvMat *
res;
926 line.par << tempM, tempQ;
931 printf(
"Common_factor: %g\n", common_factor);
934 printIt(&covParamsMat,
"Cov_params:");
935 printIt(&cov_mq,
"Rotated Covariance Matrix:");
938 printf(
"Chi2: %g\n",
chi2);
979 const Eigen::Matrix<float, 6, N>& hits_ge,
998 helix.
cov = MatrixXd::Zero(5, 5);
999 helix.
cov.block(0, 0, 3, 3) = circle.
cov;
1000 helix.
cov.block(3, 3, 2, 2) =
line.cov;
1011 #endif // RecoTracker_PixelTrackFitting_RiemannFit_h
Eigen::Matrix< double, 2, N > Matrix2xNd
Basic3DVector cross(const Basic3DVector &v) const
Vector product, or "cross" product, with a vector of same type.
int32_t qCharge
particle charge
Eigen::Matrix< double, 1, 2 *N > RowVector2Nd
Eigen::Matrix< double, 2 *N, 1 > Vector2Nd
ret
prodAgent to be discontinued
Vector5d par
(phi,Tip,pt,cotan(theta)),Zip)
Sin< T >::type sin(const T &t)
__host__ __device__ VectorNd< N > cov_carttorad_prefit(const M2xN &p2D, const Matrix2Nd< N > &cov_cart, V4 &fast_fit, const VectorNd< N > &rad)
Transform covariance matrix from Cartesian coordinates (only transverse plane component) to coordinat...
Eigen::Array< double, N, N > ArrayNd
__host__ __device__ Vector2d min_eigen2D(const Matrix2d &aMat, double &chi2)
2D version of min_eigen3D().
Vector3d par
parameter: (X0,Y0,R)
__host__ __device__ double cross2D(const Vector2d &a, const Vector2d &b)
Compute cross product of two 2D vector (assuming z component 0), returning z component of the result...
__host__ __device__ LineFit lineFit(const M3xN &hits, const M6xN &hits_ge, const CircleFit &circle, const V4 &fast_fit, const double bField, const bool error)
Perform an ordinary least square fit in the s-z plane to compute the parameters cotTheta and Zip...
__host__ __device__ CircleFit circleFit(const M2xN &hits2D, const Matrix2Nd< N > &hits_cov2D, const V4 &fast_fit, const VectorNd< N > &rad, const double bField, const bool error)
Fit a generic number of 2D points with a circle using Riemann-Chernov algorithm. Covariance matrix of...
T dot(const Basic3DVector &v) const
Scalar product, or "dot" product, with a vector of same type.
Eigen::Matrix< double, 6, 6 > Matrix6d
__host__ __device__ void par_uvrtopak(CircleFit &circle, const double B, const bool error)
Transform circle parameter from (X0,Y0,R) to (phi,Tip,p_t) and consequently covariance matrix...
Eigen::Matrix< double, N, 1 > VectorNd
int32_t qCharge
particle charge
Cos< T >::type cos(const T &t)
Eigen::Matrix< double, N, N > MatrixNd
constexpr T sqr(const T a)
raise to square.
constexpr double epsilon
used in numerical derivative (J2 in Circle_fit())
Abs< T >::type abs(const T &t)
__host__ __device__ Vector3d min_eigen3D(const Matrix3d &A, double &chi2)
Compute the eigenvector associated to the minimum eigenvalue.
Eigen::Matrix< double, 3, N > Matrix3xNd
__host__ __device__ void computeRadLenUniformMaterial(const VNd1 &length_values, VNd2 &rad_lengths)
ALPAKA_FN_ACC ALPAKA_FN_INLINE void const M3xN const V4 & fast_fit
Eigen::Matrix< double, 2 *N, 2 *N > Matrix2Nd
__host__ __device__ VectorNd< N > weightCircle(const MatrixNd< N > &cov_rad_inv)
Compute the points' weights' vector for the circle fit when multiple scattering is managed...
__host__ __device__ Vector3d min_eigen3D_fast(const Matrix3d &A)
A faster version of min_eigen3D() where double precision is not needed.
Eigen::Matrix< double, N, 3 > MatrixNx3d
riemannFit::Vector2d eVec
__host__ __device__ VectorNd< N > cov_carttorad(const M2xN &p2D, const Matrix2Nd< N > &cov_cart, const VectorNd< N > &rad)
Transform covariance matrix from Cartesian coordinates (only transverse plane component) to radial co...
__host__ __device__ void loadCovariance2D(M6xNf const &ge, M2Nd &hits_cov)
__host__ __device__ MatrixNd< N > scatter_cov_rad(const M2xN &p2D, const V4 &fast_fit, VectorNd< N > const &rad, double B)
Compute the covariance matrix (in radial coordinates) of points in the transverse plane due to multip...
__host__ __device__ auto scatterCovLine(Matrix2d const *cov_sz, const V4 &fast_fit, VNd1 const &s_arcs, VNd2 const &z_values, const double theta, const double bField, MatrixNd< N > &ret)
Compute the covariance matrix along cartesian S-Z of points due to multiple Coulomb scattering to be ...
__host__ __device__ Matrix2Nd< N > cov_radtocart(const M2xN &p2D, const MatrixNd< N > &cov_rad, const VectorNd< N > &rad)
Transform covariance matrix from radial (only tangential component) to Cartesian coordinates (only tr...
__host__ __device__ int32_t charge(const M2xN &p2D, const Vector3d &par_uvr)
Find particle q considering the sign of cross product between particles velocity (estimated by the fi...
Geom::Theta< T > theta() const
HelixFit helixFit(const Matrix3xNd< N > &hits, const Eigen::Matrix< float, 6, N > &hits_ge, const double bField, const bool error)
Helix fit by three step: -fast pre-fit (see Fast_fit() for further info); -circle fit of hits proje...
__host__ __device__ void fastFit(const M3xN &hits, V4 &result)
A very fast helix fit: it fits a circle by three points (first, middle and last point) and a line by ...
__host__ __device__ void printIt(C *m, const char *prefix="")