CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
BrokenLine.h
Go to the documentation of this file.
1 #ifndef RecoPixelVertexing_PixelTrackFitting_interface_BrokenLine_h
2 #define RecoPixelVertexing_PixelTrackFitting_interface_BrokenLine_h
3 
4 #include <Eigen/Eigenvalues>
5 
7 
8 namespace brokenline {
9 
11 
19 
23  template <int n>
25  int qCharge;
28  // starting from the pre-fitted closest approach
32  };
33 
53  __host__ __device__ inline double multScatt(
54  const double& length, const double bField, const double radius, int layer, double slope) {
55  // limit R to 20GeV...
56  auto pt2 = std::min(20., bField * radius);
57  pt2 *= pt2;
58  constexpr double inv_X0 = 0.06 / 16.;
59  //if(Layer==1) XXI_0=0.06/16.;
60  // else XXI_0=0.06/16.;
61  //XX_0*=1;
62 
64  constexpr double geometry_factor = 0.7;
65  constexpr double fact = geometry_factor * riemannFit::sqr(13.6 / 1000.);
66  return fact / (pt2 * (1. + riemannFit::sqr(slope))) * (std::abs(length) * inv_X0) *
67  riemannFit::sqr(1. + 0.038 * log(std::abs(length) * inv_X0));
68  }
69 
79  rot(0, 0) = 1. / sqrt(1. + riemannFit::sqr(slope));
80  rot(0, 1) = slope * rot(0, 0);
81  rot(1, 0) = -rot(0, 1);
82  rot(1, 1) = rot(0, 0);
83  return rot;
84  }
85 
98  double x0,
99  double y0,
100  riemannFit::Matrix3d& jacobian) {
101  // Avoid multiple access to the circle.par vector.
102  using scalar = std::remove_reference<decltype(circle.par(0))>::type;
103  scalar phi = circle.par(0);
104  scalar dee = circle.par(1);
105  scalar rho = circle.par(2);
106 
107  // Avoid repeated trig. computations
108  scalar sinPhi = sin(phi);
109  scalar cosPhi = cos(phi);
110 
111  // Intermediate computations for the circle parameters
112  scalar deltaPara = x0 * cosPhi + y0 * sinPhi;
113  scalar deltaOrth = x0 * sinPhi - y0 * cosPhi + dee;
114  scalar tempSmallU = 1 + rho * dee;
115  scalar tempC = -rho * y0 + tempSmallU * cosPhi;
116  scalar tempB = rho * x0 + tempSmallU * sinPhi;
117  scalar tempA = 2. * deltaOrth + rho * (riemannFit::sqr(deltaOrth) + riemannFit::sqr(deltaPara));
118  scalar tempU = sqrt(1. + rho * tempA);
119 
120  // Intermediate computations for the error matrix transform
121  scalar xi = 1. / (riemannFit::sqr(tempB) + riemannFit::sqr(tempC));
122  scalar tempV = 1. + rho * deltaOrth;
123  scalar lambda = (0.5 * tempA) / (riemannFit::sqr(1. + tempU) * tempU);
124  scalar mu = 1. / (tempU * (1. + tempU)) + rho * lambda;
125  scalar zeta = riemannFit::sqr(deltaOrth) + riemannFit::sqr(deltaPara);
126  jacobian << xi * tempSmallU * tempV, -xi * riemannFit::sqr(rho) * deltaOrth, xi * deltaPara,
127  2. * mu * tempSmallU * deltaPara, 2. * mu * tempV, mu * zeta - lambda * tempA, 0, 0, 1.;
128 
129  // translated circle parameters
130  // phi
131  circle.par(0) = atan2(tempB, tempC);
132  // d
133  circle.par(1) = tempA / (1 + tempU);
134  // rho after translation. It is invariant, so noop
135  // circle.par(2)= rho;
136 
137  // translated error matrix
138  circle.cov = jacobian * circle.cov * jacobian.transpose();
139  }
140 
149  template <typename M3xN, typename V4, int n>
150  __host__ __device__ inline void prepareBrokenLineData(const M3xN& hits,
151  const V4& fast_fit,
152  const double bField,
156 
157  int mId = 1;
158 
159  if constexpr (n > 3) {
160  riemannFit::Vector2d middle = 0.5 * (hits.block(0, n - 1, 2, 1) + hits.block(0, 0, 2, 1));
161  auto d1 = (hits.block(0, n / 2, 2, 1) - middle).squaredNorm();
162  auto d2 = (hits.block(0, n / 2 - 1, 2, 1) - middle).squaredNorm();
163  mId = d1 < d2 ? n / 2 : n / 2 - 1;
164  }
165 
166  dVec = hits.block(0, mId, 2, 1) - hits.block(0, 0, 2, 1);
167  eVec = hits.block(0, n - 1, 2, 1) - hits.block(0, mId, 2, 1);
168  results.qCharge = riemannFit::cross2D(dVec, eVec) > 0 ? -1 : 1;
169 
170  const double slope = -results.qCharge / fast_fit(3);
171 
172  riemannFit::Matrix2d rotMat = rotationMatrix(slope);
173 
174  // calculate radii and s
175  results.radii = hits.block(0, 0, 2, n) - fast_fit.head(2) * riemannFit::MatrixXd::Constant(1, n, 1);
176  eVec = -fast_fit(2) * fast_fit.head(2) / fast_fit.head(2).norm();
177  for (u_int i = 0; i < n; i++) {
178  dVec = results.radii.block(0, i, 2, 1);
179  results.sTransverse(i) = results.qCharge * fast_fit(2) *
180  atan2(riemannFit::cross2D(dVec, eVec), dVec.dot(eVec)); // calculates the arc length
181  }
182  riemannFit::VectorNd<n> zVec = hits.block(2, 0, 1, n).transpose();
183 
184  //calculate sTotal and zVec
186  for (u_int i = 0; i < n; i++) {
187  pointsSZ(0, i) = results.sTransverse(i);
188  pointsSZ(1, i) = zVec(i);
189  pointsSZ.block(0, i, 2, 1) = rotMat * pointsSZ.block(0, i, 2, 1);
190  }
191  results.sTotal = pointsSZ.block(0, 0, 1, n).transpose();
192  results.zInSZplane = pointsSZ.block(1, 0, 1, n).transpose();
193 
194  //calculate varBeta
195  results.varBeta(0) = results.varBeta(n - 1) = 0;
196  for (u_int i = 1; i < n - 1; i++) {
197  results.varBeta(i) = multScatt(results.sTotal(i + 1) - results.sTotal(i), bField, fast_fit(2), i + 2, slope) +
198  multScatt(results.sTotal(i) - results.sTotal(i - 1), bField, fast_fit(2), i + 1, slope);
199  }
200  }
201 
214  template <int n>
216  const riemannFit::VectorNd<n>& sTotal,
217  const riemannFit::VectorNd<n>& varBeta) {
219  for (u_int i = 0; i < n; i++) {
220  c_uMat(i, i) = weights(i);
221  if (i > 1)
222  c_uMat(i, i) += 1. / (varBeta(i - 1) * riemannFit::sqr(sTotal(i) - sTotal(i - 1)));
223  if (i > 0 && i < n - 1)
224  c_uMat(i, i) +=
225  (1. / varBeta(i)) * riemannFit::sqr((sTotal(i + 1) - sTotal(i - 1)) /
226  ((sTotal(i + 1) - sTotal(i)) * (sTotal(i) - sTotal(i - 1))));
227  if (i < n - 2)
228  c_uMat(i, i) += 1. / (varBeta(i + 1) * riemannFit::sqr(sTotal(i + 1) - sTotal(i)));
229 
230  if (i > 0 && i < n - 1)
231  c_uMat(i, i + 1) =
232  1. / (varBeta(i) * (sTotal(i + 1) - sTotal(i))) *
233  (-(sTotal(i + 1) - sTotal(i - 1)) / ((sTotal(i + 1) - sTotal(i)) * (sTotal(i) - sTotal(i - 1))));
234  if (i < n - 2)
235  c_uMat(i, i + 1) +=
236  1. / (varBeta(i + 1) * (sTotal(i + 1) - sTotal(i))) *
237  (-(sTotal(i + 2) - sTotal(i)) / ((sTotal(i + 2) - sTotal(i + 1)) * (sTotal(i + 1) - sTotal(i))));
238 
239  if (i < n - 2)
240  c_uMat(i, i + 2) = 1. / (varBeta(i + 1) * (sTotal(i + 2) - sTotal(i + 1)) * (sTotal(i + 1) - sTotal(i)));
241 
242  c_uMat(i, i) *= 0.5;
243  }
244  return c_uMat + c_uMat.transpose();
245  }
246 
257  template <typename M3xN, typename V4>
258  __host__ __device__ inline void fastFit(const M3xN& hits, V4& result) {
259  constexpr uint32_t n = M3xN::ColsAtCompileTime;
260 
261  int mId = 1;
262 
263  if constexpr (n > 3) {
264  riemannFit::Vector2d middle = 0.5 * (hits.block(0, n - 1, 2, 1) + hits.block(0, 0, 2, 1));
265  auto d1 = (hits.block(0, n / 2, 2, 1) - middle).squaredNorm();
266  auto d2 = (hits.block(0, n / 2 - 1, 2, 1) - middle).squaredNorm();
267  mId = d1 < d2 ? n / 2 : n / 2 - 1;
268  }
269 
270  const riemannFit::Vector2d a = hits.block(0, mId, 2, 1) - hits.block(0, 0, 2, 1);
271  const riemannFit::Vector2d b = hits.block(0, n - 1, 2, 1) - hits.block(0, mId, 2, 1);
272  const riemannFit::Vector2d c = hits.block(0, 0, 2, 1) - hits.block(0, n - 1, 2, 1);
273 
274  auto tmp = 0.5 / riemannFit::cross2D(c, a);
275  result(0) = hits(0, 0) - (a(1) * c.squaredNorm() + c(1) * a.squaredNorm()) * tmp;
276  result(1) = hits(1, 0) + (a(0) * c.squaredNorm() + c(0) * a.squaredNorm()) * tmp;
277  // check Wikipedia for these formulas
278 
279  result(2) = sqrt(a.squaredNorm() * b.squaredNorm() * c.squaredNorm()) / (2. * std::abs(riemannFit::cross2D(b, a)));
280  // Using Math Olympiad's formula R=abc/(4A)
281 
282  const riemannFit::Vector2d d = hits.block(0, 0, 2, 1) - result.head(2);
283  const riemannFit::Vector2d e = hits.block(0, n - 1, 2, 1) - result.head(2);
284 
285  result(3) = result(2) * atan2(riemannFit::cross2D(d, e), d.dot(e)) / (hits(2, n - 1) - hits(2, 0));
286  // ds/dz slope between last and first point
287  }
288 
313  template <typename M3xN, typename M6xN, typename V4, int n>
314  __host__ __device__ inline void circleFit(const M3xN& hits,
315  const M6xN& hits_ge,
316  const V4& fast_fit,
317  const double bField,
319  karimaki_circle_fit& circle_results) {
320  circle_results.qCharge = data.qCharge;
321  auto& radii = data.radii;
322  const auto& sTransverse = data.sTransverse;
323  const auto& sTotal = data.sTotal;
324  auto& zInSZplane = data.zInSZplane;
325  auto& varBeta = data.varBeta;
326  const double slope = -circle_results.qCharge / fast_fit(3);
327  varBeta *= 1. + riemannFit::sqr(slope); // the kink angles are projected!
328 
329  for (u_int i = 0; i < n; i++) {
330  zInSZplane(i) = radii.block(0, i, 2, 1).norm() - fast_fit(2);
331  }
332 
333  riemannFit::Matrix2d vMat; // covariance matrix
334  riemannFit::VectorNd<n> weightsVec; // weights
335  riemannFit::Matrix2d rotMat; // rotation matrix point by point
336  for (u_int i = 0; i < n; i++) {
337  vMat(0, 0) = hits_ge.col(i)[0]; // x errors
338  vMat(0, 1) = vMat(1, 0) = hits_ge.col(i)[1]; // cov_xy
339  vMat(1, 1) = hits_ge.col(i)[2]; // y errors
340  rotMat = rotationMatrix(-radii(0, i) / radii(1, i));
341  weightsVec(i) =
342  1. / ((rotMat * vMat * rotMat.transpose())(1, 1)); // compute the orthogonal weight point by point
343  }
344 
346  r_uVec(n) = 0;
347  for (u_int i = 0; i < n; i++) {
348  r_uVec(i) = weightsVec(i) * zInSZplane(i);
349  }
350 
352  c_uMat.block(0, 0, n, n) = matrixC_u(weightsVec, sTransverse, varBeta);
353  c_uMat(n, n) = 0;
354  //add the border to the c_uMat matrix
355  for (u_int i = 0; i < n; i++) {
356  c_uMat(i, n) = 0;
357  if (i > 0 && i < n - 1) {
358  c_uMat(i, n) +=
359  -(sTransverse(i + 1) - sTransverse(i - 1)) * (sTransverse(i + 1) - sTransverse(i - 1)) /
360  (2. * varBeta(i) * (sTransverse(i + 1) - sTransverse(i)) * (sTransverse(i) - sTransverse(i - 1)));
361  }
362  if (i > 1) {
363  c_uMat(i, n) +=
364  (sTransverse(i) - sTransverse(i - 2)) / (2. * varBeta(i - 1) * (sTransverse(i) - sTransverse(i - 1)));
365  }
366  if (i < n - 2) {
367  c_uMat(i, n) +=
368  (sTransverse(i + 2) - sTransverse(i)) / (2. * varBeta(i + 1) * (sTransverse(i + 1) - sTransverse(i)));
369  }
370  c_uMat(n, i) = c_uMat(i, n);
371  if (i > 0 && i < n - 1)
372  c_uMat(n, n) += riemannFit::sqr(sTransverse(i + 1) - sTransverse(i - 1)) / (4. * varBeta(i));
373  }
374 
375 #ifdef CPP_DUMP
376  std::cout << "CU5\n" << c_uMat << std::endl;
377 #endif
379  math::cholesky::invert(c_uMat, iMat);
380 #ifdef CPP_DUMP
381  std::cout << "I5\n" << iMat << std::endl;
382 #endif
383 
384  riemannFit::VectorNplusONEd<n> uVec = iMat * r_uVec; // obtain the fitted parameters by solving the linear system
385 
386  // compute (phi, d_ca, k) in the system in which the midpoint of the first two corrected hits is the origin...
387 
388  radii.block(0, 0, 2, 1) /= radii.block(0, 0, 2, 1).norm();
389  radii.block(0, 1, 2, 1) /= radii.block(0, 1, 2, 1).norm();
390 
391  riemannFit::Vector2d dVec = hits.block(0, 0, 2, 1) + (-zInSZplane(0) + uVec(0)) * radii.block(0, 0, 2, 1);
392  riemannFit::Vector2d eVec = hits.block(0, 1, 2, 1) + (-zInSZplane(1) + uVec(1)) * radii.block(0, 1, 2, 1);
393  auto eMinusd = eVec - dVec;
394  auto eMinusd2 = eMinusd.squaredNorm();
395  auto tmp1 = 1. / eMinusd2;
396  auto tmp2 = sqrt(riemannFit::sqr(fast_fit(2)) - 0.25 * eMinusd2);
397 
398  circle_results.par << atan2(eMinusd(1), eMinusd(0)), circle_results.qCharge * (tmp2 - fast_fit(2)),
399  circle_results.qCharge * (1. / fast_fit(2) + uVec(n));
400 
401  tmp2 = 1. / tmp2;
402 
403  riemannFit::Matrix3d jacobian;
404  jacobian << (radii(1, 0) * eMinusd(0) - eMinusd(1) * radii(0, 0)) * tmp1,
405  (radii(1, 1) * eMinusd(0) - eMinusd(1) * radii(0, 1)) * tmp1, 0,
406  circle_results.qCharge * (eMinusd(0) * radii(0, 0) + eMinusd(1) * radii(1, 0)) * tmp2,
407  circle_results.qCharge * (eMinusd(0) * radii(0, 1) + eMinusd(1) * radii(1, 1)) * tmp2, 0, 0, 0,
408  circle_results.qCharge;
409 
410  circle_results.cov << iMat(0, 0), iMat(0, 1), iMat(0, n), iMat(1, 0), iMat(1, 1), iMat(1, n), iMat(n, 0),
411  iMat(n, 1), iMat(n, n);
412 
413  circle_results.cov = jacobian * circle_results.cov * jacobian.transpose();
414 
415  //...Translate in the system in which the first corrected hit is the origin, adding the m.s. correction...
416 
417  translateKarimaki(circle_results, 0.5 * eMinusd(0), 0.5 * eMinusd(1), jacobian);
418  circle_results.cov(0, 0) +=
419  (1 + riemannFit::sqr(slope)) * multScatt(sTotal(1) - sTotal(0), bField, fast_fit(2), 2, slope);
420 
421  //...And translate back to the original system
422 
423  translateKarimaki(circle_results, dVec(0), dVec(1), jacobian);
424 
425  // compute chi2
426  circle_results.chi2 = 0;
427  for (u_int i = 0; i < n; i++) {
428  circle_results.chi2 += weightsVec(i) * riemannFit::sqr(zInSZplane(i) - uVec(i));
429  if (i > 0 && i < n - 1)
430  circle_results.chi2 +=
431  riemannFit::sqr(uVec(i - 1) / (sTransverse(i) - sTransverse(i - 1)) -
432  uVec(i) * (sTransverse(i + 1) - sTransverse(i - 1)) /
433  ((sTransverse(i + 1) - sTransverse(i)) * (sTransverse(i) - sTransverse(i - 1))) +
434  uVec(i + 1) / (sTransverse(i + 1) - sTransverse(i)) +
435  (sTransverse(i + 1) - sTransverse(i - 1)) * uVec(n) / 2) /
436  varBeta(i);
437  }
438  }
439 
462  template <typename V4, typename M6xN, int n>
463  __host__ __device__ inline void lineFit(const M6xN& hits_ge,
464  const V4& fast_fit,
465  const double bField,
467  riemannFit::LineFit& line_results) {
468  const auto& radii = data.radii;
469  const auto& sTotal = data.sTotal;
470  const auto& zInSZplane = data.zInSZplane;
471  const auto& varBeta = data.varBeta;
472 
473  const double slope = -data.qCharge / fast_fit(3);
474  riemannFit::Matrix2d rotMat = rotationMatrix(slope);
475 
476  riemannFit::Matrix3d vMat = riemannFit::Matrix3d::Zero(); // covariance matrix XYZ
477  riemannFit::Matrix2x3d jacobXYZtosZ =
478  riemannFit::Matrix2x3d::Zero(); // jacobian for computation of the error on s (xyz -> sz)
480  for (u_int i = 0; i < n; i++) {
481  vMat(0, 0) = hits_ge.col(i)[0]; // x errors
482  vMat(0, 1) = vMat(1, 0) = hits_ge.col(i)[1]; // cov_xy
483  vMat(0, 2) = vMat(2, 0) = hits_ge.col(i)[3]; // cov_xz
484  vMat(1, 1) = hits_ge.col(i)[2]; // y errors
485  vMat(2, 1) = vMat(1, 2) = hits_ge.col(i)[4]; // cov_yz
486  vMat(2, 2) = hits_ge.col(i)[5]; // z errors
487  auto tmp = 1. / radii.block(0, i, 2, 1).norm();
488  jacobXYZtosZ(0, 0) = radii(1, i) * tmp;
489  jacobXYZtosZ(0, 1) = -radii(0, i) * tmp;
490  jacobXYZtosZ(1, 2) = 1.;
491  weights(i) = 1. / ((rotMat * jacobXYZtosZ * vMat * jacobXYZtosZ.transpose() * rotMat.transpose())(
492  1, 1)); // compute the orthogonal weight point by point
493  }
494 
496  for (u_int i = 0; i < n; i++) {
497  r_u(i) = weights(i) * zInSZplane(i);
498  }
499 #ifdef CPP_DUMP
500  std::cout << "CU4\n" << matrixC_u(w, sTotal, varBeta) << std::endl;
501 #endif
503  math::cholesky::invert(matrixC_u(weights, sTotal, varBeta), iMat);
504 #ifdef CPP_DUMP
505  std::cout << "I4\n" << iMat << std::endl;
506 #endif
507 
508  riemannFit::VectorNd<n> uVec = iMat * r_u; // obtain the fitted parameters by solving the linear system
509 
510  // line parameters in the system in which the first hit is the origin and with axis along SZ
511  line_results.par << (uVec(1) - uVec(0)) / (sTotal(1) - sTotal(0)), uVec(0);
512  auto idiff = 1. / (sTotal(1) - sTotal(0));
513  line_results.cov << (iMat(0, 0) - 2 * iMat(0, 1) + iMat(1, 1)) * riemannFit::sqr(idiff) +
514  multScatt(sTotal(1) - sTotal(0), bField, fast_fit(2), 2, slope),
515  (iMat(0, 1) - iMat(0, 0)) * idiff, (iMat(0, 1) - iMat(0, 0)) * idiff, iMat(0, 0);
516 
517  // translate to the original SZ system
518  riemannFit::Matrix2d jacobian;
519  jacobian(0, 0) = 1.;
520  jacobian(0, 1) = 0;
521  jacobian(1, 0) = -sTotal(0);
522  jacobian(1, 1) = 1.;
523  line_results.par(1) += -line_results.par(0) * sTotal(0);
524  line_results.cov = jacobian * line_results.cov * jacobian.transpose();
525 
526  // rotate to the original sz system
527  auto tmp = rotMat(0, 0) - line_results.par(0) * rotMat(0, 1);
528  jacobian(1, 1) = 1. / tmp;
529  jacobian(0, 0) = jacobian(1, 1) * jacobian(1, 1);
530  jacobian(0, 1) = 0;
531  jacobian(1, 0) = line_results.par(1) * rotMat(0, 1) * jacobian(0, 0);
532  line_results.par(1) = line_results.par(1) * jacobian(1, 1);
533  line_results.par(0) = (rotMat(0, 1) + line_results.par(0) * rotMat(0, 0)) * jacobian(1, 1);
534  line_results.cov = jacobian * line_results.cov * jacobian.transpose();
535 
536  // compute chi2
537  line_results.chi2 = 0;
538  for (u_int i = 0; i < n; i++) {
539  line_results.chi2 += weights(i) * riemannFit::sqr(zInSZplane(i) - uVec(i));
540  if (i > 0 && i < n - 1)
541  line_results.chi2 += riemannFit::sqr(uVec(i - 1) / (sTotal(i) - sTotal(i - 1)) -
542  uVec(i) * (sTotal(i + 1) - sTotal(i - 1)) /
543  ((sTotal(i + 1) - sTotal(i)) * (sTotal(i) - sTotal(i - 1))) +
544  uVec(i + 1) / (sTotal(i + 1) - sTotal(i))) /
545  varBeta(i);
546  }
547  }
548 
583  template <int n>
585  const Eigen::Matrix<float, 6, 4>& hits_ge,
586  const double bField) {
587  riemannFit::HelixFit helix;
588  riemannFit::Vector4d fast_fit;
589  fastFit(hits, fast_fit);
590 
592  karimaki_circle_fit circle;
594  riemannFit::Matrix3d jacobian;
595 
596  prepareBrokenLineData(hits, fast_fit, bField, data);
597  lineFit(hits_ge, fast_fit, bField, data, line);
598  circleFit(hits, hits_ge, fast_fit, bField, data, circle);
599 
600  // the circle fit gives k, but here we want p_t, so let's change the parameter and the covariance matrix
601  jacobian << 1., 0, 0, 0, 1., 0, 0, 0,
602  -std::abs(circle.par(2)) * bField / (riemannFit::sqr(circle.par(2)) * circle.par(2));
603  circle.par(2) = bField / std::abs(circle.par(2));
604  circle.cov = jacobian * circle.cov * jacobian.transpose();
605 
606  helix.par << circle.par, line.par;
607  helix.cov = riemannFit::MatrixXd::Zero(5, 5);
608  helix.cov.block(0, 0, 3, 3) = circle.cov;
609  helix.cov.block(3, 3, 2, 2) = line.cov;
610  helix.qCharge = circle.qCharge;
611  helix.chi2_circle = circle.chi2;
612  helix.chi2_line = line.chi2;
613 
614  return helix;
615  }
616 
617 } // namespace brokenline
618 
619 #endif // RecoPixelVertexing_PixelTrackFitting_interface_BrokenLine_h
__host__ __device__ void lineFit(const M6xN &hits_ge, const V4 &fast_fit, const double bField, const PreparedBrokenLineData< n > &data, riemannFit::LineFit &line_results)
Performs the Broken Line fit in the straight track case (that is, the fit parameters are only the int...
Definition: BrokenLine.h:463
int qCharge
particle charge
Definition: BrokenLine.h:25
Eigen::Matrix< double, 2, N > Matrix2xNd
Definition: FitUtils.h:25
static std::vector< std::string > checklist log
const edm::EventSetup & c
int32_t qCharge
particle charge
Definition: FitResult.h:34
#define __host__
static const double slope[3]
Vector5d par
(phi,Tip,pt,cotan(theta)),Zip)
Definition: FitResult.h:49
dictionary results
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Eigen::Matrix< double, N+1, 1 > VectorNplusONEd
Definition: FitUtils.h:35
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
Eigen::Matrix2d Matrix2d
Definition: FitResult.h:17
Eigen::Vector2d Vector2d
Definition: FitResult.h:13
constexpr std::array< uint8_t, layerIndexSize > layer
tuple result
Definition: mps_fire.py:311
__host__ __device__ void translateKarimaki(karimaki_circle_fit &circle, double x0, double y0, riemannFit::Matrix3d &jacobian)
Changes the Karimäki parameters (and consequently their covariance matrix) under a translation of the...
Definition: BrokenLine.h:97
tuple d
Definition: ztail.py:151
__host__ __device__ void prepareBrokenLineData(const M3xN &hits, const V4 &fast_fit, const double bField, PreparedBrokenLineData< n > &results)
Computes the data needed for the Broken Line fit procedure that are mainly common for the circle and ...
Definition: BrokenLine.h:150
Eigen::Matrix< double, 2, 3 > Matrix2x3d
Definition: FitUtils.h:45
__host__ __device__ double multScatt(const double &length, const double bField, const double radius, int layer, double slope)
Computes the Coulomb multiple scattering variance of the planar angle.
Definition: BrokenLine.h:53
Eigen::Matrix< double, N, 1 > VectorNd
Definition: FitUtils.h:33
riemannFit::VectorNd< n > sTransverse
total distance traveled in the transverse plane
Definition: BrokenLine.h:27
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
__host__ __device__ void fastFit(const M3xN &hits, V4 &result)
A very fast helix fit.
Definition: BrokenLine.h:258
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
riemannFit::HelixFit helixFit(const riemannFit::Matrix3xNd< n > &hits, const Eigen::Matrix< float, 6, 4 > &hits_ge, const double bField)
Helix fit by three step: -fast pre-fit (see Fast_fit() for further info); -circle fit of the hits pr...
Definition: BrokenLine.h:584
Eigen::Matrix< double, 3, N > Matrix3xNd
Definition: FitResult.h:24
const int mu
Definition: Constants.h:22
Polynomial< 0 > Constant
Definition: Constant.h:6
const double fact
__host__ __device__ riemannFit::Matrix2d rotationMatrix(double slope)
Computes the 2D rotation matrix that transforms the line y=slope*x into the line y=0.
Definition: BrokenLine.h:77
__host__ __device__ riemannFit::MatrixNd< n > matrixC_u(const riemannFit::VectorNd< n > &weights, const riemannFit::VectorNd< n > &sTotal, const riemannFit::VectorNd< n > &varBeta)
Computes the n-by-n band matrix obtained minimizing the Broken Line&#39;s cost function w...
Definition: BrokenLine.h:215
riemannFit::Matrix2xNd< n > radii
xy data in the system in which the pre-fitted center is the origin
Definition: BrokenLine.h:26
double b
Definition: hdecay.h:118
riemannFit::VectorNd< n > zInSZplane
orthogonal coordinate to the pre-fitted line in the sz plane
Definition: BrokenLine.h:30
Eigen::Matrix< double, N+1, N+1 > MatrixNplusONEd
Definition: FitUtils.h:17
Eigen::Vector4d Vector4d
Definition: FitResult.h:15
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:79
double a
Definition: hdecay.h:119
Vector2d par
(cotan(theta),Zip)
Definition: FitResult.h:39
riemannFit::VectorNd< n > varBeta
kink angles in the SZ plane
Definition: BrokenLine.h:31
tuple cout
Definition: gather_cfg.py:144
T w() const
Eigen::Matrix3d Matrix3d
Definition: FitResult.h:18
static constexpr float d1
__host__ __device__ void circleFit(const M3xN &hits, const M6xN &hits_ge, const V4 &fast_fit, const double bField, PreparedBrokenLineData< n > &data, karimaki_circle_fit &circle_results)
Performs the Broken Line fit in the curved track case (that is, the fit parameters are the intercepti...
Definition: BrokenLine.h:314
tmp
align.sh
Definition: createJobs.py:716
riemannFit::VectorNd< n > sTotal
total distance traveled (three-dimensional)
Definition: BrokenLine.h:29
data needed for the Broken Line fit procedure.
Definition: BrokenLine.h:24
#define __device__
double scalar(const CLHEP::HepGenMatrix &m)
Return the matrix as a scalar. Raise an assertion if the matris is not .
Definition: matutil.cc:166