CMS 3D CMS Logo

RiemannFit.h
Go to the documentation of this file.
1 #ifndef RecoTracker_PixelTrackFitting_RiemannFit_h
2 #define RecoTracker_PixelTrackFitting_RiemannFit_h
3 
4 #include <Eigen/Core>
5 #include <Eigen/Eigenvalues>
6 
8 
9 namespace riemannFit {
10 
33  template <typename VNd1, typename VNd2>
34  __host__ __device__ inline void computeRadLenUniformMaterial(const VNd1& length_values, VNd2& rad_lengths) {
35  // Radiation length of the pixel detector in the uniform assumption, with
36  // 0.06 rad_len at 16 cm
37  constexpr double xx_0_inv = 0.06 / 16.;
38  uint n = length_values.rows();
39  rad_lengths(0) = length_values(0) * xx_0_inv;
40  for (uint j = 1; j < n; ++j) {
41  rad_lengths(j) = std::abs(length_values(j) - length_values(j - 1)) * xx_0_inv;
42  }
43  }
44 
64  template <typename V4, typename VNd1, typename VNd2, int N>
65  __host__ __device__ inline auto scatterCovLine(Matrix2d const* cov_sz,
66  const V4& fast_fit,
67  VNd1 const& s_arcs,
68  VNd2 const& z_values,
69  const double theta,
70  const double bField,
71  MatrixNd<N>& ret) {
72 #ifdef RFIT_DEBUG
73  riemannFit::printIt(&s_arcs, "Scatter_cov_line - s_arcs: ");
74 #endif
75  constexpr uint n = N;
76  double p_t = std::min(20., fast_fit(2) * bField); // limit pt to avoid too small error!!!
77  double p_2 = p_t * p_t * (1. + 1. / sqr(fast_fit(3)));
78  VectorNd<N> rad_lengths_S;
79  // See documentation at http://eigen.tuxfamily.org/dox/group__TutorialArrayClass.html
80  // Basically, to perform cwise operations on Matrices and Vectors, you need
81  // to transform them into Array-like objects.
82  VectorNd<N> s_values = s_arcs.array() * s_arcs.array() + z_values.array() * z_values.array();
83  s_values = s_values.array().sqrt();
84  computeRadLenUniformMaterial(s_values, rad_lengths_S);
85  VectorNd<N> sig2_S;
86  sig2_S = .000225 / p_2 * (1. + 0.038 * rad_lengths_S.array().log()).abs2() * rad_lengths_S.array();
87 #ifdef RFIT_DEBUG
88  riemannFit::printIt(cov_sz, "Scatter_cov_line - cov_sz: ");
89 #endif
91  for (uint k = 0; k < n; ++k) {
92  tmp(k, k) = cov_sz[k](0, 0);
93  tmp(k + n, k + n) = cov_sz[k](1, 1);
94  tmp(k, k + n) = tmp(k + n, k) = cov_sz[k](0, 1);
95  }
96  for (uint k = 0; k < n; ++k) {
97  for (uint l = k; l < n; ++l) {
98  for (uint i = 0; i < std::min(k, l); ++i) {
99  tmp(k + n, l + n) += std::abs(s_values(k) - s_values(i)) * std::abs(s_values(l) - s_values(i)) * sig2_S(i);
100  }
101  tmp(l + n, k + n) = tmp(k + n, l + n);
102  }
103  }
104  // We are interested only in the errors orthogonal to the rotated s-axis
105  // which, in our formalism, are in the lower square matrix.
106 #ifdef RFIT_DEBUG
107  riemannFit::printIt(&tmp, "Scatter_cov_line - tmp: ");
108 #endif
109  ret = tmp.block(n, n, n, n);
110  }
111 
125  template <typename M2xN, typename V4, int N>
127  const V4& fast_fit,
128  VectorNd<N> const& rad,
129  double B) {
130  constexpr uint n = N;
131  double p_t = std::min(20., fast_fit(2) * B); // limit pt to avoid too small error!!!
132  double p_2 = p_t * p_t * (1. + 1. / sqr(fast_fit(3)));
133  double theta = atan(fast_fit(3));
134  theta = theta < 0. ? theta + M_PI : theta;
135  VectorNd<N> s_values;
136  VectorNd<N> rad_lengths;
137  const Vector2d oVec(fast_fit(0), fast_fit(1));
138 
139  // associated Jacobian, used in weights and errors computation
140  for (uint i = 0; i < n; ++i) { // x
141  Vector2d pVec = p2D.block(0, i, 2, 1) - oVec;
142  const double cross = cross2D(-oVec, pVec);
143  const double dot = (-oVec).dot(pVec);
144  const double tempAtan2 = atan2(cross, dot);
145  s_values(i) = std::abs(tempAtan2 * fast_fit(2));
146  }
147  computeRadLenUniformMaterial(s_values * sqrt(1. + 1. / sqr(fast_fit(3))), rad_lengths);
149  VectorNd<N> sig2 = (1. + 0.038 * rad_lengths.array().log()).abs2() * rad_lengths.array();
150  sig2 *= 0.000225 / (p_2 * sqr(sin(theta)));
151  for (uint k = 0; k < n; ++k) {
152  for (uint l = k; l < n; ++l) {
153  for (uint i = 0; i < std::min(k, l); ++i) {
154  scatter_cov_rad(k, l) += (rad(k) - rad(i)) * (rad(l) - rad(i)) * sig2(i);
155  }
157  }
158  }
159 #ifdef RFIT_DEBUG
160  riemannFit::printIt(&scatter_cov_rad, "Scatter_cov_rad - scatter_cov_rad: ");
161 #endif
162  return scatter_cov_rad;
163  }
164 
173  template <typename M2xN, int N>
175  const MatrixNd<N>& cov_rad,
176  const VectorNd<N>& rad) {
177 #ifdef RFIT_DEBUG
178  printf("Address of p2D: %p\n", &p2D);
179 #endif
180  printIt(&p2D, "cov_radtocart - p2D:");
181  constexpr uint n = N;
182  Matrix2Nd<N> cov_cart = Matrix2Nd<N>::Zero();
183  VectorNd<N> rad_inv = rad.cwiseInverse();
184  printIt(&rad_inv, "cov_radtocart - rad_inv:");
185  for (uint i = 0; i < n; ++i) {
186  for (uint j = i; j < n; ++j) {
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);
195  }
196  }
197  return cov_cart;
198  }
199 
210  template <typename M2xN, int N>
212  const Matrix2Nd<N>& cov_cart,
213  const VectorNd<N>& rad) {
214  constexpr uint n = N;
215  VectorNd<N> cov_rad;
216  const VectorNd<N> rad_inv2 = rad.cwiseInverse().array().square();
217  for (uint i = 0; i < n; ++i) {
219  if (rad(i) < 1.e-4)
220  cov_rad(i) = cov_cart(i, i);
221  else {
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));
224  }
225  }
226  return cov_rad;
227  }
228 
241  template <typename M2xN, typename V4, int N>
243  const Matrix2Nd<N>& cov_cart,
244  V4& fast_fit,
245  const VectorNd<N>& rad) {
246  constexpr uint n = N;
247  VectorNd<N> cov_rad;
248  for (uint i = 0; i < n; ++i) {
250  if (rad(i) < 1.e-4)
251  cov_rad(i) = cov_cart(i, i); // TO FIX
252  else {
253  Vector2d a = p2D.col(i);
254  Vector2d b = p2D.col(i) - fast_fit.head(2);
255  const double x2 = a.dot(b);
256  const double y2 = cross2D(a, b);
257  const double tan_c = -y2 / x2;
258  const double tan_c2 = sqr(tan_c);
259  cov_rad(i) =
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);
261  }
262  }
263  return cov_rad;
264  }
265 
277  template <int N>
279  return cov_rad_inv.colwise().sum().transpose();
280  }
281 
290  template <typename M2xN>
291  __host__ __device__ inline int32_t charge(const M2xN& p2D, const Vector3d& par_uvr) {
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)) >
293  0)
294  ? -1
295  : 1;
296  }
297 
314 #ifdef RFIT_DEBUG
315  printf("min_eigen3D - enter\n");
316 #endif
317  Eigen::SelfAdjointEigenSolver<Matrix3d> solver(3);
318  solver.computeDirect(A);
319  int min_index;
320  chi2 = solver.eigenvalues().minCoeff(&min_index);
321 #ifdef RFIT_DEBUG
322  printf("min_eigen3D - exit\n");
323 #endif
324  return solver.eigenvectors().col(min_index);
325  }
326 
339  Eigen::SelfAdjointEigenSolver<Matrix3f> solver(3);
340  solver.computeDirect(A.cast<float>());
341  int min_index;
342  solver.eigenvalues().minCoeff(&min_index);
343  return solver.eigenvectors().col(min_index).cast<double>();
344  }
345 
356  __host__ __device__ inline Vector2d min_eigen2D(const Matrix2d& aMat, double& chi2) {
357  Eigen::SelfAdjointEigenSolver<Matrix2d> solver(2);
358  solver.computeDirect(aMat);
359  int min_index;
360  chi2 = solver.eigenvalues().minCoeff(&min_index);
361  return solver.eigenvectors().col(min_index);
362  }
363 
377  template <typename M3xN, typename V4>
378  __host__ __device__ inline void fastFit(const M3xN& hits, V4& result) {
379  constexpr uint32_t N = M3xN::ColsAtCompileTime;
380  constexpr auto n = N; // get the number of hits
381  printIt(&hits, "Fast_fit - hits: ");
382 
383  // CIRCLE FIT
384  // Make segments between middle-to-first(b) and last-to-first(c) hits
385  const Vector2d bVec = hits.block(0, n / 2, 2, 1) - hits.block(0, 0, 2, 1);
386  const Vector2d cVec = hits.block(0, n - 1, 2, 1) - hits.block(0, 0, 2, 1);
387  printIt(&bVec, "Fast_fit - b: ");
388  printIt(&cVec, "Fast_fit - c: ");
389  // Compute their lengths
390  auto b2 = bVec.squaredNorm();
391  auto c2 = cVec.squaredNorm();
392  // The algebra has been verified (MR). The usual approach has been followed:
393  // * use an orthogonal reference frame passing from the first point.
394  // * build the segments (chords)
395  // * build orthogonal lines through mid points
396  // * make a system and solve for X0 and Y0.
397  // * add the initial point
398  bool flip = abs(bVec.x()) < abs(bVec.y());
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);
405  // if aligned TO FIX
406  auto y0 = (cx * b2 - bx * c2) / div;
407  auto x0 = (0.5 * b2 - y0 * by) / bx;
408  result(0) = hits(0, 0) + (flip ? y0 : x0);
409  result(1) = hits(1, 0) + (flip ? x0 : y0);
410  result(2) = sqrt(sqr(x0) + sqr(y0));
411  printIt(&result, "Fast_fit - result: ");
412 
413  // LINE FIT
414  const Vector2d dVec = hits.block(0, 0, 2, 1) - result.head(2);
415  const Vector2d eVec = hits.block(0, n - 1, 2, 1) - result.head(2);
416  printIt(&eVec, "Fast_fit - e: ");
417  printIt(&dVec, "Fast_fit - d: ");
418  // Compute the arc-length between first and last point: L = R * theta = R * atan (tan (Theta) )
419  auto dr = result(2) * atan2(cross2D(dVec, eVec), dVec.dot(eVec));
420  // Simple difference in Z between last and first hit
421  auto dz = hits(2, n - 1) - hits(2, 0);
422 
423  result(3) = (dr / dz);
424 
425 #ifdef RFIT_DEBUG
426  printf("Fast_fit: [%f, %f, %f, %f]\n", result(0), result(1), result(2), result(3));
427 #endif
428  }
429 
457  template <typename M2xN, typename V4, int N>
458  __host__ __device__ inline CircleFit circleFit(const M2xN& hits2D,
459  const Matrix2Nd<N>& hits_cov2D,
460  const V4& fast_fit,
461  const VectorNd<N>& rad,
462  const double bField,
463  const bool error) {
464 #ifdef RFIT_DEBUG
465  printf("circle_fit - enter\n");
466 #endif
467  // INITIALIZATION
468  Matrix2Nd<N> vMat = hits_cov2D;
469  constexpr uint n = N;
470  printIt(&hits2D, "circle_fit - hits2D:");
471  printIt(&hits_cov2D, "circle_fit - hits_cov2D:");
472 
473 #ifdef RFIT_DEBUG
474  printf("circle_fit - WEIGHT COMPUTATION\n");
475 #endif
476  // WEIGHT COMPUTATION
478  MatrixNd<N> gMat;
479  double renorm;
480  {
481  MatrixNd<N> cov_rad = cov_carttorad_prefit(hits2D, vMat, fast_fit, rad).asDiagonal();
482  MatrixNd<N> scatterCovRadMat = scatter_cov_rad(hits2D, fast_fit, rad, bField);
483  printIt(&scatterCovRadMat, "circle_fit - scatter_cov_rad:");
484  printIt(&hits2D, "circle_fit - hits2D bis:");
485 #ifdef RFIT_DEBUG
486  printf("Address of hits2D: a) %p\n", &hits2D);
487 #endif
488  vMat += cov_radtocart(hits2D, scatterCovRadMat, rad);
489  printIt(&vMat, "circle_fit - V:");
490  cov_rad += scatterCovRadMat;
491  printIt(&cov_rad, "circle_fit - cov_rad:");
492  math::cholesky::invert(cov_rad, gMat);
493  // gMat = cov_rad.inverse();
494  renorm = gMat.sum();
495  gMat *= 1. / renorm;
496  weight = weightCircle(gMat);
497  }
498  printIt(&weight, "circle_fit - weight:");
499 
500  // SPACE TRANSFORMATION
501 #ifdef RFIT_DEBUG
502  printf("circle_fit - SPACE TRANSFORMATION\n");
503 #endif
504 
505  // center
506 #ifdef RFIT_DEBUG
507  printf("Address of hits2D: b) %p\n", &hits2D);
508 #endif
509  const Vector2d hCentroid = hits2D.rowwise().mean(); // centroid
510  printIt(&hCentroid, "circle_fit - h_:");
511  Matrix3xNd<N> p3D;
512  p3D.block(0, 0, 2, n) = hits2D.colwise() - hCentroid;
513  printIt(&p3D, "circle_fit - p3D: a)");
514  Vector2Nd<N> mc; // centered hits, used in error computation
515  mc << p3D.row(0).transpose(), p3D.row(1).transpose();
516  printIt(&mc, "circle_fit - mc(centered hits):");
517 
518  // scale
519  const double tempQ = mc.squaredNorm();
520  const double tempS = sqrt(n * 1. / tempQ); // scaling factor
521  p3D.block(0, 0, 2, n) *= tempS;
522 
523  // project on paraboloid
524  p3D.row(2) = p3D.block(0, 0, 2, n).colwise().squaredNorm();
525  printIt(&p3D, "circle_fit - p3D: b)");
526 
527 #ifdef RFIT_DEBUG
528  printf("circle_fit - COST FUNCTION\n");
529 #endif
530  // COST FUNCTION
531 
532  // compute
533  Vector3d r0;
534  r0.noalias() = p3D * weight; // center of gravity
535  const Matrix3xNd<N> xMat = p3D.colwise() - r0;
536  Matrix3d aMat = xMat * gMat * xMat.transpose();
537  printIt(&aMat, "circle_fit - A:");
538 
539 #ifdef RFIT_DEBUG
540  printf("circle_fit - MINIMIZE\n");
541 #endif
542  // minimize
543  double chi2;
544  Vector3d vVec = min_eigen3D(aMat, chi2);
545 #ifdef RFIT_DEBUG
546  printf("circle_fit - AFTER MIN_EIGEN\n");
547 #endif
548  printIt(&vVec, "v BEFORE INVERSION");
549  vVec *= (vVec(2) > 0) ? 1 : -1; // TO FIX dovrebbe essere N(3)>0
550  printIt(&vVec, "v AFTER INVERSION");
551  // This hack to be able to run on GPU where the automatic assignment to a
552  // double from the vector multiplication is not working.
553 #ifdef RFIT_DEBUG
554  printf("circle_fit - AFTER MIN_EIGEN 1\n");
555 #endif
556  Eigen::Matrix<double, 1, 1> cm;
557 #ifdef RFIT_DEBUG
558  printf("circle_fit - AFTER MIN_EIGEN 2\n");
559 #endif
560  cm = -vVec.transpose() * r0;
561 #ifdef RFIT_DEBUG
562  printf("circle_fit - AFTER MIN_EIGEN 3\n");
563 #endif
564  const double tempC = cm(0, 0);
565 
566 #ifdef RFIT_DEBUG
567  printf("circle_fit - COMPUTE CIRCLE PARAMETER\n");
568 #endif
569  // COMPUTE CIRCLE PARAMETER
570 
571  // auxiliary quantities
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;
575  Vector3d par_uvr; // used in error propagation
576  par_uvr << -vVec(0) * v2x2_inv, -vVec(1) * v2x2_inv, tempH * v2x2_inv;
577 
578  CircleFit circle;
579  circle.par << par_uvr(0) * s_inv + hCentroid(0), par_uvr(1) * s_inv + hCentroid(1), par_uvr(2) * s_inv;
580  circle.qCharge = charge(hits2D, circle.par);
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:");
584 #ifdef RFIT_DEBUG
585  printf("circle_fit - CIRCLE CHARGE: %d\n", circle.qCharge);
586 #endif
587 
588 #ifdef RFIT_DEBUG
589  printf("circle_fit - ERROR PROPAGATION\n");
590 #endif
591  // ERROR PROPAGATION
592  if (error) {
593 #ifdef RFIT_DEBUG
594  printf("circle_fit - ERROR PRPAGATION ACTIVATED\n");
595 #endif
596  ArrayNd<N> vcsMat[2][2]; // cov matrix of center & scaled points
597  MatrixNd<N> cMat[3][3]; // cov matrix of 3D transformed points
598 #ifdef RFIT_DEBUG
599  printf("circle_fit - ERROR PRPAGATION ACTIVATED 2\n");
600 #endif
601  {
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);
606  Matrix2Nd<N> tempVcsMat;
607  tempVcsMat.template triangularView<Eigen::Upper>() =
608  (sqr(tempS) * vMat + sqr(sqr(tempS)) * 1. / (4. * tempQ * n) *
609  (2. * vMat.squaredNorm() + 4. * tempC2) * // mc.transpose() * V * mc) *
610  (mc * mc.transpose()));
611 
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:");
618  }
619 
620  {
621  const ArrayNd<N> t0 = (VectorXd::Constant(n, 1.) * p3D.row(0));
622  const ArrayNd<N> t1 = (VectorXd::Constant(n, 1.) * p3D.row(1));
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);
626  const ArrayNd<N> t10 = t01.transpose();
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))
637  .matrix();
638  cMat[2][2] = tmp.template selfadjointView<Eigen::Upper>();
639  }
640  printIt(&cMat[0][0], "circle_fit - C[0][0]:");
641 
642  Matrix3d c0Mat; // cov matrix of center of gravity (r0.x,r0.y,r0.z)
643  for (uint i = 0; i < 3; ++i) {
644  for (uint j = i; j < 3; ++j) {
645  Eigen::Matrix<double, 1, 1> tmp;
646  tmp = weight.transpose() * cMat[i][j] * weight;
647  // Workaround to get things working in GPU
648  const double tempC = tmp(0, 0);
649  c0Mat(i, j) = tempC; //weight.transpose() * C[i][j] * weight;
650  c0Mat(j, i) = c0Mat(i, j);
651  }
652  }
653  printIt(&c0Mat, "circle_fit - C0:");
654 
655  const MatrixNd<N> wMat = weight * weight.transpose();
656  const MatrixNd<N> hMat = MatrixNd<N>::Identity().rowwise() - weight.transpose();
657  const MatrixNx3d<N> s_v = hMat * p3D.transpose();
658  printIt(&wMat, "circle_fit - W:");
659  printIt(&hMat, "circle_fit - H:");
660  printIt(&s_v, "circle_fit - s_v:");
661 
662  MatrixNd<N> dMat[3][3]; // cov(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]:");
673 
674  constexpr uint nu[6][2] = {{0, 0}, {0, 1}, {0, 2}, {1, 1}, {1, 2}, {2, 2}};
675 
676  Matrix6d eMat; // cov matrix of the 6 independent elements of A
677  for (uint a = 0; a < 6; ++a) {
678  const uint i = nu[a][0], j = nu[a][1];
679  for (uint b = a; b < 6; ++b) {
680  const uint k = nu[b][0], l = nu[b][1];
681  VectorNd<N> t0(n);
682  VectorNd<N> t1(n);
683  if (l == k) {
684  t0 = 2. * dMat[j][l] * s_v.col(l);
685  if (i == j)
686  t1 = t0;
687  else
688  t1 = 2. * dMat[i][l] * s_v.col(l);
689  } else {
690  t0 = dMat[j][l] * s_v.col(k) + dMat[j][k] * s_v.col(l);
691  if (i == j)
692  t1 = t0;
693  else
694  t1 = dMat[i][l] * s_v.col(k) + dMat[i][k] * s_v.col(l);
695  }
696 
697  if (i == j) {
698  Eigen::Matrix<double, 1, 1> cm;
699  cm = s_v.col(i).transpose() * (t0 + t1);
700  // Workaround to get things working in GPU
701  const double tempC = cm(0, 0);
702  eMat(a, b) = 0. + tempC;
703  } else {
704  Eigen::Matrix<double, 1, 1> cm;
705  cm = (s_v.col(i).transpose() * t0) + (s_v.col(j).transpose() * t1);
706  // Workaround to get things working in GPU
707  const double tempC = cm(0, 0);
708  eMat(a, b) = 0. + tempC; //(s_v.col(i).transpose() * t0) + (s_v.col(j).transpose() * t1);
709  }
710  if (b != a)
711  eMat(b, a) = eMat(a, b);
712  }
713  }
714  printIt(&eMat, "circle_fit - E:");
715 
716  Eigen::Matrix<double, 3, 6> j2Mat; // Jacobian of min_eigen() (numerically computed)
717  for (uint a = 0; a < 6; ++a) {
718  const uint i = nu[a][0], j = nu[a][1];
719  Matrix3d delta = Matrix3d::Zero();
720  delta(i, j) = delta(j, i) = abs(aMat(i, j) * epsilon);
721  j2Mat.col(a) = min_eigen3D_fast(aMat + delta);
722  const int sign = (j2Mat.col(a)(2) > 0) ? 1 : -1;
723  j2Mat.col(a) = (j2Mat.col(a) * sign - vVec) / delta(i, j);
724  }
725  printIt(&j2Mat, "circle_fit - J2:");
726 
727  Matrix4d cvcMat; // joint cov matrix of (v0,v1,v2,c)
728  {
729  Matrix3d t0 = j2Mat * eMat * j2Mat.transpose();
730  Vector3d t1 = -t0 * r0;
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);
737  // cm2 = (c0Mat.cwiseProduct(t0)).sum();
738  cm3 = (r0.transpose() * t0 * r0);
739  // Workaround to get things working in GPU
740  const double tempC = cm1(0, 0) + (c0Mat.cwiseProduct(t0)).sum() + cm3(0, 0);
741  cvcMat(3, 3) = tempC;
742  // (v.transpose() * c0Mat * v) + (c0Mat.cwiseProduct(t0)).sum() + (r0.transpose() * t0 * r0);
743  }
744  printIt(&cvcMat, "circle_fit - Cvc:");
745 
746  Eigen::Matrix<double, 3, 4> j3Mat; // Jacobian (v0,v1,v2,c)->(X0,Y0,R)
747  {
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;
752  }
753  printIt(&j3Mat, "circle_fit - J3:");
754 
755  const RowVector2Nd<N> Jq = mc.transpose() * tempS * 1. / n; // var(q)
756  printIt(&Jq, "circle_fit - Jq:");
757 
758  Matrix3d cov_uvr = j3Mat * cvcMat * j3Mat.transpose() * sqr(s_inv) // cov(X0,Y0,R)
759  + (par_uvr * par_uvr.transpose()) * (Jq * vMat * Jq.transpose());
760 
761  circle.cov = cov_uvr;
762  }
763 
764  printIt(&circle.cov, "Circle cov:");
765 #ifdef RFIT_DEBUG
766  printf("circle_fit - exit\n");
767 #endif
768  return circle;
769  }
770 
787  template <typename M3xN, typename M6xN, typename V4>
788  __host__ __device__ inline LineFit lineFit(const M3xN& hits,
789  const M6xN& hits_ge,
790  const CircleFit& circle,
791  const V4& fast_fit,
792  const double bField,
793  const bool error) {
794  constexpr uint32_t N = M3xN::ColsAtCompileTime;
795  constexpr auto n = N;
796  double theta = -circle.qCharge * atan(fast_fit(3));
797  theta = theta < 0. ? theta + M_PI : theta;
798 
799  // Prepare the Rotation Matrix to rotate the points
800  Eigen::Matrix<double, 2, 2> rot;
801  rot << sin(theta), cos(theta), -cos(theta), sin(theta);
802 
803  // PROJECTION ON THE CILINDER
804  //
805  // p2D will be:
806  // [s1, s2, s3, ..., sn]
807  // [z1, z2, z3, ..., zn]
808  // s values will be ordinary x-values
809  // z values will be ordinary y-values
810 
812  Eigen::Matrix<double, 2, 6> jxMat;
813 
814 #ifdef RFIT_DEBUG
815  printf("Line_fit - B: %g\n", bField);
816  printIt(&hits, "Line_fit points: ");
817  printIt(&hits_ge, "Line_fit covs: ");
818  printIt(&rot, "Line_fit rot: ");
819 #endif
820  // x & associated Jacobian
821  // cfr https://indico.cern.ch/event/663159/contributions/2707659/attachments/1517175/2368189/Riemann_fit.pdf
822  // Slide 11
823  // a ==> -o i.e. the origin of the circle in XY plane, negative
824  // b ==> p i.e. distances of the points wrt the origin of the circle.
825  const Vector2d oVec(circle.par(0), circle.par(1));
826 
827  // associated Jacobian, used in weights and errors computation
828  Matrix6d covMat = Matrix6d::Zero();
829  Matrix2d cov_sz[N];
830  for (uint i = 0; i < n; ++i) {
831  Vector2d pVec = hits.block(0, i, 2, 1) - oVec;
832  const double cross = cross2D(-oVec, pVec);
833  const double dot = (-oVec).dot(pVec);
834  // atan2(cross, dot) give back the angle in the transverse plane so tha the
835  // final equation reads: x_i = -q*R*theta (theta = angle returned by atan2)
836  const double tempQAtan2 = -circle.qCharge * atan2(cross, dot);
837  // p2D.coeffRef(1, i) = atan2_ * circle.par(2);
838  p2D(0, i) = tempQAtan2 * circle.par(2);
839 
840  // associated Jacobian, used in weights and errors- computation
841  const double temp0 = -circle.qCharge * circle.par(2) * 1. / (sqr(dot) + sqr(cross));
842  double d_X0 = 0., d_Y0 = 0., d_R = 0.; // good approximation for big pt and eta
843  if (error) {
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);
846  d_R = tempQAtan2;
847  }
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.;
851 
852  covMat.block(0, 0, 3, 3) = circle.cov;
853  covMat(3, 3) = hits_ge.col(i)[0]; // x errors
854  covMat(4, 4) = hits_ge.col(i)[2]; // y errors
855  covMat(5, 5) = hits_ge.col(i)[5]; // z errors
856  covMat(3, 4) = covMat(4, 3) = hits_ge.col(i)[1]; // cov_xy
857  covMat(3, 5) = covMat(5, 3) = hits_ge.col(i)[3]; // cov_xz
858  covMat(4, 5) = covMat(5, 4) = hits_ge.col(i)[4]; // cov_yz
859  Matrix2d tmp = jxMat * covMat * jxMat.transpose();
860  cov_sz[i].noalias() = rot * tmp * rot.transpose();
861  }
862  // Math of d_{X0,Y0,R,x,y} all verified by hand
863  p2D.row(1) = hits.row(2);
864 
865  // The following matrix will contain errors orthogonal to the rotated S
866  // component only, with the Multiple Scattering properly treated!!
867  MatrixNd<N> cov_with_ms;
868  scatterCovLine(cov_sz, fast_fit, p2D.row(0), p2D.row(1), theta, bField, cov_with_ms);
869 #ifdef RFIT_DEBUG
870  printIt(cov_sz, "line_fit - cov_sz:");
871  printIt(&cov_with_ms, "line_fit - cov_with_ms: ");
872 #endif
873 
874  // Rotate Points with the shape [2, n]
875  Matrix2xNd<N> p2D_rot = rot * p2D;
876 
877 #ifdef RFIT_DEBUG
878  printf("Fast fit Tan(theta): %g\n", fast_fit(3));
879  printf("Rotation angle: %g\n", theta);
880  printIt(&rot, "Rotation Matrix:");
881  printIt(&p2D, "Original Hits(s,z):");
882  printIt(&p2D_rot, "Rotated hits(S3D, Z'):");
883  printIt(&rot, "Rotation Matrix:");
884 #endif
885 
886  // Build the A Matrix
887  Matrix2xNd<N> aMat;
888  aMat << MatrixXd::Ones(1, n), p2D_rot.row(0); // rotated s values
889 
890 #ifdef RFIT_DEBUG
891  printIt(&aMat, "A Matrix:");
892 #endif
893 
894  // Build A^T V-1 A, where V-1 is the covariance of only the Y components.
895  MatrixNd<N> vyInvMat;
896  math::cholesky::invert(cov_with_ms, vyInvMat);
897  // MatrixNd<N> vyInvMat = cov_with_ms.inverse();
898  Eigen::Matrix<double, 2, 2> covParamsMat = aMat * vyInvMat * aMat.transpose();
899  // Compute the Covariance Matrix of the fit parameters
900  math::cholesky::invert(covParamsMat, covParamsMat);
901 
902  // Now Compute the Parameters in the form [2,1]
903  // The first component is q.
904  // The second component is m.
905  Eigen::Matrix<double, 2, 1> sol = covParamsMat * aMat * vyInvMat * p2D_rot.row(1).transpose();
906 
907 #ifdef RFIT_DEBUG
908  printIt(&sol, "Rotated solutions:");
909 #endif
910 
911  // We need now to transfer back the results in the original s-z plane
912  const auto sinTheta = sin(theta);
913  const auto cosTheta = cos(theta);
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;
917 
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();
921 
922  VectorNd<N> res = p2D_rot.row(1).transpose() - aMat.transpose() * sol;
923  double chi2 = res.transpose() * vyInvMat * res;
924 
925  LineFit line;
926  line.par << tempM, tempQ;
927  line.cov << cov_mq;
928  line.chi2 = chi2;
929 
930 #ifdef RFIT_DEBUG
931  printf("Common_factor: %g\n", common_factor);
932  printIt(&jMat, "Jacobian:");
933  printIt(&sol, "Rotated solutions:");
934  printIt(&covParamsMat, "Cov_params:");
935  printIt(&cov_mq, "Rotated Covariance Matrix:");
936  printIt(&(line.par), "Real Parameters:");
937  printIt(&(line.cov), "Real Covariance Matrix:");
938  printf("Chi2: %g\n", chi2);
939 #endif
940 
941  return line;
942  }
943 
977  template <int N>
979  const Eigen::Matrix<float, 6, N>& hits_ge,
980  const double bField,
981  const bool error) {
982  constexpr uint n = N;
983  VectorNd<4> rad = (hits.block(0, 0, 2, n).colwise().norm());
984 
985  // Fast_fit gives back (X0, Y0, R, theta) w/o errors, using only 3 points.
988  riemannFit::Matrix2Nd<N> hits_cov = MatrixXd::Zero(2 * n, 2 * n);
989  riemannFit::loadCovariance2D(hits_ge, hits_cov);
990  CircleFit circle = circleFit(hits.block(0, 0, 2, n), hits_cov, fast_fit, rad, bField, error);
991  LineFit line = lineFit(hits, hits_ge, circle, fast_fit, bField, error);
992 
993  par_uvrtopak(circle, bField, error);
994 
995  HelixFit helix;
996  helix.par << circle.par, line.par;
997  if (error) {
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;
1001  }
1002  helix.qCharge = circle.qCharge;
1003  helix.chi2_circle = circle.chi2;
1004  helix.chi2_line = line.chi2;
1005 
1006  return helix;
1007  }
1008 
1009 } // namespace riemannFit
1010 
1011 #endif // RecoTracker_PixelTrackFitting_RiemannFit_h
Eigen::Vector4d Vector4d
Definition: FitResult.h:12
Eigen::Matrix< double, 2, N > Matrix2xNd
Definition: FitUtils.h:29
Eigen::Vector3d Vector3d
Definition: FitResult.h:11
Basic3DVector cross(const Basic3DVector &v) const
Vector product, or "cross" product, with a vector of same type.
Eigen::Matrix3d Matrix3d
Definition: FitResult.h:15
int32_t qCharge
particle charge
Definition: FitResult.h:31
Eigen::Matrix< double, 1, 2 *N > RowVector2Nd
Definition: FitUtils.h:47
Definition: APVGainStruct.h:7
Eigen::Matrix< double, 2 *N, 1 > Vector2Nd
Definition: FitUtils.h:41
#define __host__
ret
prodAgent to be discontinued
Vector5d par
(phi,Tip,pt,cotan(theta)),Zip)
Definition: FitResult.h:46
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
__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...
Definition: RiemannFit.h:242
Eigen::Array< double, N, N > ArrayNd
Definition: FitUtils.h:23
Definition: weight.py:1
__host__ __device__ Vector2d min_eigen2D(const Matrix2d &aMat, double &chi2)
2D version of min_eigen3D().
Definition: RiemannFit.h:356
Vector3d par
parameter: (X0,Y0,R)
Definition: FitResult.h:24
__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...
Definition: FitUtils.h:79
__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...
Definition: RiemannFit.h:788
__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...
Definition: RiemannFit.h:458
Definition: Electron.h:6
T dot(const Basic3DVector &v) const
Scalar product, or "dot" product, with a vector of same type.
Eigen::Matrix< double, 6, 6 > Matrix6d
Definition: FitResult.h:18
__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...
Definition: FitUtils.h:173
Eigen::Matrix< double, N, 1 > VectorNd
Definition: FitUtils.h:37
int32_t qCharge
particle charge
Definition: FitResult.h:58
T sqrt(T t)
Definition: SSEVec.h:19
Eigen::Matrix2d Matrix2d
Definition: FitResult.h:14
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Eigen::Matrix< double, N, N > MatrixNd
Definition: FitUtils.h:19
constexpr T sqr(const T a)
raise to square.
Definition: FitUtils.h:67
Eigen::Matrix4d Matrix4d
Definition: FitResult.h:16
constexpr double epsilon
used in numerical derivative (J2 in Circle_fit())
Definition: FitUtils.h:14
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
__host__ __device__ Vector3d min_eigen3D(const Matrix3d &A, double &chi2)
Compute the eigenvector associated to the minimum eigenvalue.
Definition: RiemannFit.h:313
Eigen::Matrix< double, 3, N > Matrix3xNd
Definition: HelixFit.h:30
__host__ __device__ void computeRadLenUniformMaterial(const VNd1 &length_values, VNd2 &rad_lengths)
Definition: RiemannFit.h:34
#define M_PI
bias2_t b2[25]
Definition: b2.h:9
Polynomial< 0 > Constant
Definition: Constant.h:6
ALPAKA_FN_ACC ALPAKA_FN_INLINE void const M3xN const V4 & fast_fit
Definition: BrokenLine.h:158
Eigen::Matrix< double, 2 *N, 2 *N > Matrix2Nd
Definition: FitUtils.h:25
__host__ __device__ VectorNd< N > weightCircle(const MatrixNd< N > &cov_rad_inv)
Compute the points&#39; weights&#39; vector for the circle fit when multiple scattering is managed...
Definition: RiemannFit.h:278
constexpr float sol
Definition: Config.h:13
__host__ __device__ Vector3d min_eigen3D_fast(const Matrix3d &A)
A faster version of min_eigen3D() where double precision is not needed.
Definition: RiemannFit.h:338
#define N
Definition: blowfish.cc:9
Eigen::Vector2d Vector2d
Definition: FitResult.h:10
double b
Definition: hdecay.h:120
Eigen::Matrix< double, N, 3 > MatrixNx3d
Definition: FitUtils.h:33
__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...
Definition: RiemannFit.h:211
double a
Definition: hdecay.h:121
__host__ __device__ void loadCovariance2D(M6xNf const &ge, M2Nd &hits_cov)
Definition: FitUtils.h:88
__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...
Definition: RiemannFit.h:126
__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 ...
Definition: RiemannFit.h:65
Definition: APVGainStruct.h:7
__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...
Definition: RiemannFit.h:174
__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...
Definition: RiemannFit.h:291
tmp
align.sh
Definition: createJobs.py:716
Geom::Theta< T > theta() const
#define __device__
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...
Definition: RiemannFit.h:978
__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 ...
Definition: RiemannFit.h:378
__host__ __device__ void printIt(C *m, const char *prefix="")
Definition: FitUtils.h:53