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