1 #ifndef RecoPixelVertexing_PixelTrackFitting_interface_RiemannFit_h 2 #define RecoPixelVertexing_PixelTrackFitting_interface_RiemannFit_h 30 template <
typename VNd1,
typename VNd2>
34 constexpr
double xx_0_inv = 0.06 / 16.;
35 uint n = length_values.rows();
36 rad_lengths(0) = length_values(0) * xx_0_inv;
38 rad_lengths(
j) =
std::abs(length_values(
j) - length_values(
j - 1)) * xx_0_inv;
61 template <
typename V4,
typename VNd1,
typename VNd2,
int N>
74 double p_2 = p_t * p_t * (1. + 1. /
sqr(fast_fit(3)));
79 VectorNd<N> s_values = s_arcs.array() * s_arcs.array() + z_values.array() * z_values.array();
80 s_values = s_values.array().sqrt();
83 sig2_S = .000225 / p_2 * (1. + 0.038 * rad_lengths_S.array().log()).abs2() * rad_lengths_S.array();
89 tmp(
k,
k) = cov_sz[
k](0, 0);
122 template <
typename M2xN,
typename V4,
int N>
128 double p_t =
std::min(20., fast_fit(2) *
B);
129 double p_2 = p_t * p_t * (1. + 1. /
sqr(fast_fit(3)));
130 double theta = atan(fast_fit(3));
134 const Vector2d oVec(fast_fit(0), fast_fit(1));
138 Vector2d pVec = p2D.block(0,
i, 2, 1) - oVec;
140 const double dot = (-oVec).
dot(pVec);
141 const double tempAtan2 = atan2(
cross,
dot);
142 s_values(
i) =
std::abs(tempAtan2 * fast_fit(2));
146 VectorNd<N> sig2 = (1. + 0.038 * rad_lengths.array().log()).abs2() * rad_lengths.array();
170 template <
typename M2xN,
int N>
175 printf(
"Address of p2D: %p\n", &p2D);
177 printIt(&p2D,
"cov_radtocart - p2D:");
181 printIt(&rad_inv,
"cov_radtocart - rad_inv:");
184 cov_cart(
i,
j) = cov_rad(
i,
j) * p2D(1,
i) * rad_inv(
i) * p2D(1,
j) * rad_inv(
j);
185 cov_cart(
i +
n,
j +
n) = cov_rad(
i,
j) * p2D(0,
i) * rad_inv(
i) * p2D(0,
j) * rad_inv(
j);
186 cov_cart(
i,
j +
n) = -cov_rad(
i,
j) * p2D(1,
i) * rad_inv(
i) * p2D(0,
j) * rad_inv(
j);
187 cov_cart(
i +
n,
j) = -cov_rad(
i,
j) * p2D(0,
i) * rad_inv(
i) * p2D(1,
j) * rad_inv(
j);
188 cov_cart(
j,
i) = cov_cart(
i,
j);
189 cov_cart(
j +
n,
i +
n) = cov_cart(
i +
n,
j +
n);
190 cov_cart(
j +
n,
i) = cov_cart(
i,
j +
n);
191 cov_cart(
j,
i +
n) = cov_cart(
i +
n,
j);
207 template <
typename M2xN,
int N>
213 const VectorNd<N> rad_inv2 = rad.cwiseInverse().array().square();
217 cov_rad(
i) = cov_cart(
i,
i);
219 cov_rad(
i) = rad_inv2(
i) * (cov_cart(
i,
i) *
sqr(p2D(1,
i)) + cov_cart(
i +
n,
i +
n) *
sqr(p2D(0,
i)) -
220 2. * cov_cart(
i,
i +
n) * p2D(0,
i) * p2D(1,
i));
238 template <
typename M2xN,
typename V4,
int N>
248 cov_rad(
i) = cov_cart(
i,
i);
252 const double x2 =
a.dot(
b);
254 const double tan_c = -
y2 / x2;
255 const double tan_c2 =
sqr(tan_c);
257 1. / (1. + tan_c2) * (cov_cart(
i,
i) + cov_cart(
i +
n,
i +
n) * tan_c2 + 2 * cov_cart(
i,
i +
n) * tan_c);
276 return cov_rad_inv.colwise().sum().transpose();
287 template <
typename M2xN>
289 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)) >
312 printf(
"min_eigen3D - enter\n");
314 Eigen::SelfAdjointEigenSolver<Matrix3d> solver(3);
315 solver.computeDirect(
A);
317 chi2 = solver.eigenvalues().minCoeff(&min_index);
319 printf(
"min_eigen3D - exit\n");
321 return solver.eigenvectors().col(min_index);
336 Eigen::SelfAdjointEigenSolver<Matrix3f> solver(3);
337 solver.computeDirect(
A.cast<
float>());
339 solver.eigenvalues().minCoeff(&min_index);
340 return solver.eigenvectors().col(min_index).cast<
double>();
354 Eigen::SelfAdjointEigenSolver<Matrix2d> solver(2);
355 solver.computeDirect(aMat);
357 chi2 = solver.eigenvalues().minCoeff(&min_index);
358 return solver.eigenvectors().col(min_index);
374 template <
typename M3xN,
typename V4>
376 constexpr uint32_t
N = M3xN::ColsAtCompileTime;
377 constexpr
auto n =
N;
384 printIt(&bVec,
"Fast_fit - b: ");
385 printIt(&cVec,
"Fast_fit - c: ");
387 auto b2 = bVec.squaredNorm();
388 auto c2 = cVec.squaredNorm();
396 auto bx =
flip ? bVec.y() : bVec.x();
397 auto by =
flip ? bVec.x() : bVec.y();
398 auto cx =
flip ? cVec.y() : cVec.x();
399 auto cy =
flip ? cVec.x() : cVec.y();
401 auto div = 2. * (cx * by -
bx * cy);
403 auto y0 = (cx *
b2 -
bx * c2) / div;
404 auto x0 = (0.5 *
b2 - y0 * by) /
bx;
413 printIt(&eVec,
"Fast_fit - e: ");
414 printIt(&dVec,
"Fast_fit - d: ");
416 auto dr =
result(2) * atan2(
cross2D(dVec, eVec), dVec.dot(eVec));
454 template <
typename M2xN,
typename V4,
int N>
462 printf(
"circle_fit - enter\n");
467 printIt(&hits2D,
"circle_fit - hits2D:");
468 printIt(&hits_cov2D,
"circle_fit - hits_cov2D:");
471 printf(
"circle_fit - WEIGHT COMPUTATION\n");
480 printIt(&scatterCovRadMat,
"circle_fit - scatter_cov_rad:");
481 printIt(&hits2D,
"circle_fit - hits2D bis:");
483 printf(
"Address of hits2D: a) %p\n", &hits2D);
486 printIt(&vMat,
"circle_fit - V:");
487 cov_rad += scatterCovRadMat;
488 printIt(&cov_rad,
"circle_fit - cov_rad:");
499 printf(
"circle_fit - SPACE TRANSFORMATION\n");
504 printf(
"Address of hits2D: b) %p\n", &hits2D);
506 const Vector2d hCentroid = hits2D.rowwise().mean();
507 printIt(&hCentroid,
"circle_fit - h_:");
509 p3D.block(0, 0, 2,
n) = hits2D.colwise() - hCentroid;
510 printIt(&p3D,
"circle_fit - p3D: a)");
512 mc << p3D.row(0).transpose(), p3D.row(1).transpose();
513 printIt(&
mc,
"circle_fit - mc(centered hits):");
516 const double tempQ =
mc.squaredNorm();
517 const double tempS =
sqrt(
n * 1. / tempQ);
518 p3D.block(0, 0, 2,
n) *= tempS;
521 p3D.row(2) = p3D.block(0, 0, 2,
n).colwise().squaredNorm();
522 printIt(&p3D,
"circle_fit - p3D: b)");
525 printf(
"circle_fit - COST FUNCTION\n");
531 r0.noalias() = p3D *
weight;
533 Matrix3d aMat = xMat * gMat * xMat.transpose();
534 printIt(&aMat,
"circle_fit - A:");
537 printf(
"circle_fit - MINIMIZE\n");
543 printf(
"circle_fit - AFTER MIN_EIGEN\n");
545 printIt(&vVec,
"v BEFORE INVERSION");
546 vVec *= (vVec(2) > 0) ? 1 : -1;
547 printIt(&vVec,
"v AFTER INVERSION");
551 printf(
"circle_fit - AFTER MIN_EIGEN 1\n");
553 Eigen::Matrix<double, 1, 1> cm;
555 printf(
"circle_fit - AFTER MIN_EIGEN 2\n");
557 cm = -vVec.transpose() * r0;
559 printf(
"circle_fit - AFTER MIN_EIGEN 3\n");
561 const double tempC = cm(0, 0);
564 printf(
"circle_fit - COMPUTE CIRCLE PARAMETER\n");
569 const double tempH =
sqrt(1. -
sqr(vVec(2)) - 4. * tempC * vVec(2));
570 const double v2x2_inv = 1. / (2. * vVec(2));
571 const double s_inv = 1. / tempS;
573 par_uvr << -vVec(0) * v2x2_inv, -vVec(1) * v2x2_inv, tempH * v2x2_inv;
576 circle.
par << par_uvr(0) * s_inv + hCentroid(0), par_uvr(1) * s_inv + hCentroid(1), par_uvr(2) * s_inv;
578 circle.
chi2 =
abs(
chi2) * renorm /
sqr(2 * vVec(2) * par_uvr(2) * tempS);
579 printIt(&circle.
par,
"circle_fit - CIRCLE PARAMETERS:");
580 printIt(&circle.
cov,
"circle_fit - CIRCLE COVARIANCE:");
582 printf(
"circle_fit - CIRCLE CHARGE: %d\n", circle.
qCharge);
586 printf(
"circle_fit - ERROR PROPAGATION\n");
591 printf(
"circle_fit - ERROR PRPAGATION ACTIVATED\n");
596 printf(
"circle_fit - ERROR PRPAGATION ACTIVATED 2\n");
599 Eigen::Matrix<double, 1, 1> cm;
600 Eigen::Matrix<double, 1, 1> cm2;
601 cm =
mc.transpose() * vMat *
mc;
602 const double tempC2 = cm(0, 0);
604 tempVcsMat.template triangularView<Eigen::Upper>() =
605 (
sqr(tempS) * vMat +
sqr(
sqr(tempS)) * 1. / (4. * tempQ *
n) *
606 (2. * vMat.squaredNorm() + 4. * tempC2) *
607 (
mc *
mc.transpose()));
609 printIt(&tempVcsMat,
"circle_fit - Vcs:");
610 cMat[0][0] = tempVcsMat.block(0, 0,
n,
n).template selfadjointView<Eigen::Upper>();
611 vcsMat[0][1] = tempVcsMat.block(0,
n,
n,
n);
612 cMat[1][1] = tempVcsMat.block(
n,
n,
n,
n).template selfadjointView<Eigen::Upper>();
613 vcsMat[1][0] = vcsMat[0][1].transpose();
614 printIt(&tempVcsMat,
"circle_fit - Vcs:");
620 const ArrayNd<N> t00 = p3D.row(0).transpose() * p3D.row(0);
621 const ArrayNd<N> t01 = p3D.row(0).transpose() * p3D.row(1);
622 const ArrayNd<N> t11 = p3D.row(1).transpose() * p3D.row(1);
624 vcsMat[0][0] = cMat[0][0];
625 cMat[0][1] = vcsMat[0][1];
626 cMat[0][2] = 2. * (vcsMat[0][0] *
t0 + vcsMat[0][1] *
t1);
627 vcsMat[1][1] = cMat[1][1];
628 cMat[1][2] = 2. * (vcsMat[1][0] *
t0 + vcsMat[1][1] *
t1);
630 tmp.template triangularView<Eigen::Upper>() =
631 (2. * (vcsMat[0][0] * vcsMat[0][0] + vcsMat[0][0] * vcsMat[0][1] + vcsMat[1][1] * vcsMat[1][0] +
632 vcsMat[1][1] * vcsMat[1][1]) +
633 4. * (vcsMat[0][0] * t00 + vcsMat[0][1] * t01 + vcsMat[1][0] * t10 + vcsMat[1][1] * t11))
635 cMat[2][2] =
tmp.template selfadjointView<Eigen::Upper>();
637 printIt(&cMat[0][0],
"circle_fit - C[0][0]:");
642 Eigen::Matrix<double, 1, 1>
tmp;
645 const double tempC =
tmp(0, 0);
647 c0Mat(
j,
i) = c0Mat(
i,
j);
650 printIt(&c0Mat,
"circle_fit - C0:");
655 printIt(&wMat,
"circle_fit - W:");
656 printIt(&hMat,
"circle_fit - H:");
657 printIt(&s_v,
"circle_fit - s_v:");
660 dMat[0][0] = (hMat * cMat[0][0] * hMat.transpose()).cwiseProduct(wMat);
661 dMat[0][1] = (hMat * cMat[0][1] * hMat.transpose()).cwiseProduct(wMat);
662 dMat[0][2] = (hMat * cMat[0][2] * hMat.transpose()).cwiseProduct(wMat);
663 dMat[1][1] = (hMat * cMat[1][1] * hMat.transpose()).cwiseProduct(wMat);
664 dMat[1][2] = (hMat * cMat[1][2] * hMat.transpose()).cwiseProduct(wMat);
665 dMat[2][2] = (hMat * cMat[2][2] * hMat.transpose()).cwiseProduct(wMat);
666 dMat[1][0] = dMat[0][1].transpose();
667 dMat[2][0] = dMat[0][2].transpose();
668 dMat[2][1] = dMat[1][2].transpose();
669 printIt(&dMat[0][0],
"circle_fit - D_[0][0]:");
671 constexpr
uint nu[6][2] = {{0, 0}, {0, 1}, {0, 2}, {1, 1}, {1, 2}, {2, 2}};
675 const uint i = nu[
a][0],
j = nu[
a][1];
677 const uint k = nu[
b][0],
l = nu[
b][1];
681 t0 = 2. * dMat[
j][
l] * s_v.col(
l);
685 t1 = 2. * dMat[
i][
l] * s_v.col(
l);
687 t0 = dMat[
j][
l] * s_v.col(
k) + dMat[
j][
k] * s_v.col(
l);
691 t1 = dMat[
i][
l] * s_v.col(
k) + dMat[
i][
k] * s_v.col(
l);
695 Eigen::Matrix<double, 1, 1> cm;
696 cm = s_v.col(
i).transpose() * (
t0 +
t1);
698 const double tempC = cm(0, 0);
699 eMat(
a,
b) = 0. + tempC;
701 Eigen::Matrix<double, 1, 1> cm;
702 cm = (s_v.col(
i).transpose() *
t0) + (s_v.col(
j).transpose() *
t1);
704 const double tempC = cm(0, 0);
705 eMat(
a,
b) = 0. + tempC;
708 eMat(
b,
a) = eMat(
a,
b);
711 printIt(&eMat,
"circle_fit - E:");
713 Eigen::Matrix<double, 3, 6> j2Mat;
715 const uint i = nu[
a][0],
j = nu[
a][1];
719 const int sign = (j2Mat.col(
a)(2) > 0) ? 1 : -1;
722 printIt(&j2Mat,
"circle_fit - J2:");
726 Matrix3d t0 = j2Mat * eMat * j2Mat.transpose();
728 cvcMat.block(0, 0, 3, 3) =
t0;
729 cvcMat.block(0, 3, 3, 1) =
t1;
730 cvcMat.block(3, 0, 1, 3) =
t1.transpose();
731 Eigen::Matrix<double, 1, 1> cm1;
732 Eigen::Matrix<double, 1, 1> cm3;
733 cm1 = (vVec.transpose() * c0Mat * vVec);
735 cm3 = (r0.transpose() *
t0 * r0);
737 const double tempC = cm1(0, 0) + (c0Mat.cwiseProduct(
t0)).sum() + cm3(0, 0);
738 cvcMat(3, 3) = tempC;
741 printIt(&cvcMat,
"circle_fit - Cvc:");
743 Eigen::Matrix<double, 3, 4> j3Mat;
745 const double t = 1. / tempH;
746 j3Mat << -v2x2_inv, 0, vVec(0) *
sqr(v2x2_inv) * 2., 0, 0, -v2x2_inv, vVec(1) *
sqr(v2x2_inv) * 2., 0,
747 vVec(0) * v2x2_inv *
t, vVec(1) * v2x2_inv *
t,
748 -tempH *
sqr(v2x2_inv) * 2. - (2. * tempC + vVec(2)) * v2x2_inv *
t, -
t;
750 printIt(&j3Mat,
"circle_fit - J3:");
753 printIt(&Jq,
"circle_fit - Jq:");
755 Matrix3d cov_uvr = j3Mat * cvcMat * j3Mat.transpose() *
sqr(s_inv)
756 + (par_uvr * par_uvr.transpose()) * (Jq * vMat * Jq.transpose());
758 circle.
cov = cov_uvr;
763 printf(
"circle_fit - exit\n");
784 template <
typename M3xN,
typename M6xN,
typename V4>
791 constexpr uint32_t
N = M3xN::ColsAtCompileTime;
792 constexpr
auto n =
N;
797 Eigen::Matrix<double, 2, 2>
rot;
809 Eigen::Matrix<double, 2, 6> jxMat;
812 printf(
"Line_fit - B: %g\n",
bField);
814 printIt(&hits_ge,
"Line_fit covs: ");
830 const double dot = (-oVec).
dot(pVec);
835 p2D(0,
i) = tempQAtan2 * circle.
par(2);
839 double d_X0 = 0., d_Y0 = 0., d_R = 0.;
841 d_X0 = -temp0 * ((pVec(1) + oVec(1)) *
dot - (pVec(0) - oVec(0)) *
cross);
842 d_Y0 = temp0 * ((pVec(0) + oVec(0)) *
dot - (oVec(1) - pVec(1)) *
cross);
845 const double d_x = temp0 * (oVec(1) *
dot + oVec(0) *
cross);
846 const double d_y = temp0 * (-oVec(0) *
dot + oVec(1) *
cross);
847 jxMat << d_X0, d_Y0, d_R,
d_x,
d_y, 0., 0., 0., 0., 0., 0., 1.;
849 covMat.block(0, 0, 3, 3) = circle.
cov;
850 covMat(3, 3) = hits_ge.col(
i)[0];
851 covMat(4, 4) = hits_ge.col(
i)[2];
852 covMat(5, 5) = hits_ge.col(
i)[5];
853 covMat(3, 4) = covMat(4, 3) = hits_ge.col(
i)[1];
854 covMat(3, 5) = covMat(5, 3) = hits_ge.col(
i)[3];
855 covMat(4, 5) = covMat(5, 4) = hits_ge.col(
i)[4];
857 cov_sz[
i].noalias() =
rot *
tmp *
rot.transpose();
860 p2D.row(1) =
hits.row(2);
867 printIt(cov_sz,
"line_fit - cov_sz:");
868 printIt(&cov_with_ms,
"line_fit - cov_with_ms: ");
875 printf(
"Fast fit Tan(theta): %g\n", fast_fit(3));
876 printf(
"Rotation angle: %g\n",
theta);
878 printIt(&p2D,
"Original Hits(s,z):");
879 printIt(&p2D_rot,
"Rotated hits(S3D, Z'):");
885 aMat << MatrixXd::Ones(1,
n), p2D_rot.row(0);
895 Eigen::Matrix<double, 2, 2> covParamsMat = aMat * vyInvMat * aMat.transpose();
902 Eigen::Matrix<double, 2, 1>
sol = covParamsMat * aMat * vyInvMat * p2D_rot.row(1).transpose();
911 auto common_factor = 1. / (sinTheta -
sol(1, 0) * cosTheta);
912 Eigen::Matrix<double, 2, 2> jMat;
913 jMat << 0., common_factor * common_factor, common_factor,
sol(0, 0) * cosTheta * common_factor * common_factor;
915 double tempM = common_factor * (
sol(1, 0) * sinTheta + cosTheta);
916 double tempQ = common_factor *
sol(0, 0);
917 auto cov_mq = jMat * covParamsMat * jMat.transpose();
920 double chi2 =
res.transpose() * vyInvMat *
res;
923 line.par << tempM, tempQ;
928 printf(
"Common_factor: %g\n", common_factor);
931 printIt(&covParamsMat,
"Cov_params:");
932 printIt(&cov_mq,
"Rotated Covariance Matrix:");
935 printf(
"Chi2: %g\n",
chi2);
976 const Eigen::Matrix<float, 6, N>& hits_ge,
995 helix.
cov = MatrixXd::Zero(5, 5);
996 helix.
cov.block(0, 0, 3, 3) = circle.
cov;
997 helix.
cov.block(3, 3, 2, 2) =
line.cov;
1008 #endif // RecoPixelVertexing_PixelTrackFitting_interface_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, 6, 6 > Matrix6d
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.
__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)
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
__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="")