CMS 3D CMS Logo

Classes | Typedefs | Functions
brokenline Namespace Reference

Classes

struct  PreparedBrokenLineData
 data needed for the Broken Line fit procedure. More...
 

Typedefs

using karimaki_circle_fit = riemannFit::CircleFit
 < Karimäki's parameters: (phi, d, k=1/R) More...
 

Functions

template<typename M3xN , typename M6xN , typename V4 , int n>
__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 interceptions u and the curvature correction \Delta\kappa). More...
 
template<typename M3xN , typename V4 >
__host__ __device__ void fastFit (const M3xN &hits, V4 &result)
 A very fast helix fit. More...
 
template<int n>
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 projected in the transverse plane by Broken Line algorithm (see BL_Circle_fit() for further info);
-line fit of the hits projected on the (pre-fitted) cilinder surface by Broken Line algorithm (see BL_Line_fit() for further info);
Points must be passed ordered (from inner to outer layer). More...
 
template<typename V4 , typename M6xN , int n>
__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 interceptions u). More...
 
template<int n>
__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's cost function w.r.t u. This is the whole matrix in the case of the line fit and the main n-by-n block in the case of the circle fit. More...
 
__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. More...
 
template<typename M3xN , typename V4 , int n>
__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 the line fit. More...
 
__host__ __device__ riemannFit::Matrix2d rotationMatrix (double slope)
 Computes the 2D rotation matrix that transforms the line y=slope*x into the line y=0. More...
 
__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 coordinate system, such that the old origin has coordinates (x0,y0) in the new coordinate system. The formulas are taken from Karimäki V., 1990, Effective circle fitting for particle trajectories, Nucl. Instr. and Meth. A305 (1991) 187. More...
 

Typedef Documentation

◆ karimaki_circle_fit

< Karimäki's parameters: (phi, d, k=1/R)

< covariance matrix:
|cov(phi,phi)|cov( d ,phi)|cov( k ,phi)|
|cov(phi, d )|cov( d , d )|cov( k , d )|
|cov(phi, k )|cov( d , k )|cov( k , k )|
as defined in Karimäki V., 1990, Effective circle fitting for particle trajectories, Nucl. Instr. and Meth. A305 (1991) 187.

Definition at line 18 of file BrokenLine.h.

Function Documentation

◆ circleFit()

template<typename M3xN , typename M6xN , typename V4 , int n>
__host__ __device__ void brokenline::circleFit ( const M3xN &  hits,
const M6xN &  hits_ge,
const V4 &  fast_fit,
const double  bField,
PreparedBrokenLineData< n > &  data,
karimaki_circle_fit circle_results 
)
inline

Performs the Broken Line fit in the curved track case (that is, the fit parameters are the interceptions u and the curvature correction \Delta\kappa).

Parameters
hitshits coordinates.
hits_covhits covariance matrix.
fast_fitpre-fit result in the form (X0,Y0,R,tan(theta)).
bFieldmagnetic field in Gev/cm/c.
dataPreparedBrokenLineData.
circle_resultsstruct to be filled with the results in this form: -par parameter of the line in this form: (phi, d, k);
-cov covariance matrix of the fitted parameter;
-chi2 value of the cost function in the minimum.

The function implements the steps 2 and 3 of the Broken Line fit with the curvature correction.
The step 2 is the least square fit, done by imposing the minimum constraint on the cost function and solving the consequent linear system. It determines the fitted parameters u and \Delta\kappa and their covariance matrix. The step 3 is the correction of the fast pre-fitted parameters for the innermost part of the track. It is first done in a comfortable coordinate system (the one in which the first hit is the origin) and then the parameters and their covariance matrix are transformed to the original coordinate system.

Definition at line 296 of file BrokenLine.h.

301  {
302  circle_results.qCharge = data.qCharge;
303  auto& radii = data.radii;
304  const auto& sTransverse = data.sTransverse;
305  const auto& sTotal = data.sTotal;
306  auto& zInSZplane = data.zInSZplane;
307  auto& varBeta = data.varBeta;
308  const double slope = -circle_results.qCharge / fast_fit(3);
309  varBeta *= 1. + riemannFit::sqr(slope); // the kink angles are projected!
310 
311  for (u_int i = 0; i < n; i++) {
312  zInSZplane(i) = radii.block(0, i, 2, 1).norm() - fast_fit(2);
313  }
314 
315  riemannFit::Matrix2d vMat; // covariance matrix
316  riemannFit::VectorNd<n> weightsVec; // weights
317  riemannFit::Matrix2d rotMat; // rotation matrix point by point
318  for (u_int i = 0; i < n; i++) {
319  vMat(0, 0) = hits_ge.col(i)[0]; // x errors
320  vMat(0, 1) = vMat(1, 0) = hits_ge.col(i)[1]; // cov_xy
321  vMat(1, 1) = hits_ge.col(i)[2]; // y errors
322  rotMat = rotationMatrix(-radii(0, i) / radii(1, i));
323  weightsVec(i) =
324  1. / ((rotMat * vMat * rotMat.transpose())(1, 1)); // compute the orthogonal weight point by point
325  }
326 
328  r_uVec(n) = 0;
329  for (u_int i = 0; i < n; i++) {
330  r_uVec(i) = weightsVec(i) * zInSZplane(i);
331  }
332 
334  c_uMat.block(0, 0, n, n) = matrixC_u(weightsVec, sTransverse, varBeta);
335  c_uMat(n, n) = 0;
336  //add the border to the c_uMat matrix
337  for (u_int i = 0; i < n; i++) {
338  c_uMat(i, n) = 0;
339  if (i > 0 && i < n - 1) {
340  c_uMat(i, n) +=
341  -(sTransverse(i + 1) - sTransverse(i - 1)) * (sTransverse(i + 1) - sTransverse(i - 1)) /
342  (2. * varBeta(i) * (sTransverse(i + 1) - sTransverse(i)) * (sTransverse(i) - sTransverse(i - 1)));
343  }
344  if (i > 1) {
345  c_uMat(i, n) +=
346  (sTransverse(i) - sTransverse(i - 2)) / (2. * varBeta(i - 1) * (sTransverse(i) - sTransverse(i - 1)));
347  }
348  if (i < n - 2) {
349  c_uMat(i, n) +=
350  (sTransverse(i + 2) - sTransverse(i)) / (2. * varBeta(i + 1) * (sTransverse(i + 1) - sTransverse(i)));
351  }
352  c_uMat(n, i) = c_uMat(i, n);
353  if (i > 0 && i < n - 1)
354  c_uMat(n, n) += riemannFit::sqr(sTransverse(i + 1) - sTransverse(i - 1)) / (4. * varBeta(i));
355  }
356 
357 #ifdef CPP_DUMP
358  std::cout << "CU5\n" << c_uMat << std::endl;
359 #endif
361  math::cholesky::invert(c_uMat, iMat);
362 #ifdef CPP_DUMP
363  std::cout << "I5\n" << iMat << std::endl;
364 #endif
365 
366  riemannFit::VectorNplusONEd<n> uVec = iMat * r_uVec; // obtain the fitted parameters by solving the linear system
367 
368  // compute (phi, d_ca, k) in the system in which the midpoint of the first two corrected hits is the origin...
369 
370  radii.block(0, 0, 2, 1) /= radii.block(0, 0, 2, 1).norm();
371  radii.block(0, 1, 2, 1) /= radii.block(0, 1, 2, 1).norm();
372 
373  riemannFit::Vector2d dVec = hits.block(0, 0, 2, 1) + (-zInSZplane(0) + uVec(0)) * radii.block(0, 0, 2, 1);
374  riemannFit::Vector2d eVec = hits.block(0, 1, 2, 1) + (-zInSZplane(1) + uVec(1)) * radii.block(0, 1, 2, 1);
375  auto eMinusd = eVec - dVec;
376  auto eMinusd2 = eMinusd.squaredNorm();
377  auto tmp1 = 1. / eMinusd2;
378  auto tmp2 = sqrt(riemannFit::sqr(fast_fit(2)) - 0.25 * eMinusd2);
379 
380  circle_results.par << atan2(eMinusd(1), eMinusd(0)), circle_results.qCharge * (tmp2 - fast_fit(2)),
381  circle_results.qCharge * (1. / fast_fit(2) + uVec(n));
382 
383  tmp2 = 1. / tmp2;
384 
385  riemannFit::Matrix3d jacobian;
386  jacobian << (radii(1, 0) * eMinusd(0) - eMinusd(1) * radii(0, 0)) * tmp1,
387  (radii(1, 1) * eMinusd(0) - eMinusd(1) * radii(0, 1)) * tmp1, 0,
388  circle_results.qCharge * (eMinusd(0) * radii(0, 0) + eMinusd(1) * radii(1, 0)) * tmp2,
389  circle_results.qCharge * (eMinusd(0) * radii(0, 1) + eMinusd(1) * radii(1, 1)) * tmp2, 0, 0, 0,
390  circle_results.qCharge;
391 
392  circle_results.cov << iMat(0, 0), iMat(0, 1), iMat(0, n), iMat(1, 0), iMat(1, 1), iMat(1, n), iMat(n, 0),
393  iMat(n, 1), iMat(n, n);
394 
395  circle_results.cov = jacobian * circle_results.cov * jacobian.transpose();
396 
397  //...Translate in the system in which the first corrected hit is the origin, adding the m.s. correction...
398 
399  translateKarimaki(circle_results, 0.5 * eMinusd(0), 0.5 * eMinusd(1), jacobian);
400  circle_results.cov(0, 0) +=
401  (1 + riemannFit::sqr(slope)) * multScatt(sTotal(1) - sTotal(0), bField, fast_fit(2), 2, slope);
402 
403  //...And translate back to the original system
404 
405  translateKarimaki(circle_results, dVec(0), dVec(1), jacobian);
406 
407  // compute chi2
408  circle_results.chi2 = 0;
409  for (u_int i = 0; i < n; i++) {
410  circle_results.chi2 += weightsVec(i) * riemannFit::sqr(zInSZplane(i) - uVec(i));
411  if (i > 0 && i < n - 1)
412  circle_results.chi2 +=
413  riemannFit::sqr(uVec(i - 1) / (sTransverse(i) - sTransverse(i - 1)) -
414  uVec(i) * (sTransverse(i + 1) - sTransverse(i - 1)) /
415  ((sTransverse(i + 1) - sTransverse(i)) * (sTransverse(i) - sTransverse(i - 1))) +
416  uVec(i + 1) / (sTransverse(i + 1) - sTransverse(i)) +
417  (sTransverse(i + 1) - sTransverse(i - 1)) * uVec(n) / 2) /
418  varBeta(i);
419  }
420  }

References Calorimetry_cff::bField, riemannFit::CircleFit::chi2, gather_cfg::cout, riemannFit::CircleFit::cov, data, hfClusterShapes_cfi::hits, mps_fire::i, l1tstage2_dqm_sourceclient-live_cfg::invert, matrixC_u(), multScatt(), dqmiodumpmetadata::n, riemannFit::CircleFit::par, riemannFit::CircleFit::qCharge, rotationMatrix(), slope, riemannFit::sqr(), mathSSE::sqrt(), and translateKarimaki().

Referenced by __attribute__(), and helixFit().

◆ fastFit()

template<typename M3xN , typename V4 >
__host__ __device__ void brokenline::fastFit ( const M3xN &  hits,
V4 &  result 
)
inline

A very fast helix fit.

Parameters
hitsthe measured hits.
Returns
(X0,Y0,R,tan(theta)).
Warning
sign of theta is (intentionally, for now) mistaken for negative charges.

Definition at line 249 of file BrokenLine.h.

249  {
250  constexpr uint32_t n = M3xN::ColsAtCompileTime;
251 
252  const riemannFit::Vector2d a = hits.block(0, n / 2, 2, 1) - hits.block(0, 0, 2, 1);
253  const riemannFit::Vector2d b = hits.block(0, n - 1, 2, 1) - hits.block(0, n / 2, 2, 1);
254  const riemannFit::Vector2d c = hits.block(0, 0, 2, 1) - hits.block(0, n - 1, 2, 1);
255 
256  auto tmp = 0.5 / riemannFit::cross2D(c, a);
257  result(0) = hits(0, 0) - (a(1) * c.squaredNorm() + c(1) * a.squaredNorm()) * tmp;
258  result(1) = hits(1, 0) + (a(0) * c.squaredNorm() + c(0) * a.squaredNorm()) * tmp;
259  // check Wikipedia for these formulas
260 
261  result(2) = sqrt(a.squaredNorm() * b.squaredNorm() * c.squaredNorm()) / (2. * std::abs(riemannFit::cross2D(b, a)));
262  // Using Math Olympiad's formula R=abc/(4A)
263 
264  const riemannFit::Vector2d d = hits.block(0, 0, 2, 1) - result.head(2);
265  const riemannFit::Vector2d e = hits.block(0, n - 1, 2, 1) - result.head(2);
266 
267  result(3) = result(2) * atan2(riemannFit::cross2D(d, e), d.dot(e)) / (hits(2, n - 1) - hits(2, 0));
268  // ds/dz slope between last and first point
269  }

References a, funct::abs(), b, c, riemannFit::cross2D(), ztail::d, MillePedeFileConverter_cfg::e, hfClusterShapes_cfi::hits, dqmiodumpmetadata::n, mps_fire::result, mathSSE::sqrt(), and createJobs::tmp.

Referenced by for(), and helixFit().

◆ helixFit()

template<int n>
riemannFit::HelixFit brokenline::helixFit ( const riemannFit::Matrix3xNd< n > &  hits,
const Eigen::Matrix< float, 6, 4 > &  hits_ge,
const double  bField 
)
inline

Helix fit by three step: -fast pre-fit (see Fast_fit() for further info);
-circle fit of the hits projected in the transverse plane by Broken Line algorithm (see BL_Circle_fit() for further info);
-line fit of the hits projected on the (pre-fitted) cilinder surface by Broken Line algorithm (see BL_Line_fit() for further info);
Points must be passed ordered (from inner to outer layer).

Parameters
hitsMatrix3xNd hits coordinates in this form:
|x1|x2|x3|...|xn|
|y1|y2|y3|...|yn|
|z1|z2|z3|...|zn|
hits_covMatrix3Nd covariance matrix in this form (()->cov()):
|(x1,x1)|(x2,x1)|(x3,x1)|(x4,x1)|.|(y1,x1)|(y2,x1)|(y3,x1)|(y4,x1)|.|(z1,x1)|(z2,x1)|(z3,x1)|(z4,x1)|
|(x1,x2)|(x2,x2)|(x3,x2)|(x4,x2)|.|(y1,x2)|(y2,x2)|(y3,x2)|(y4,x2)|.|(z1,x2)|(z2,x2)|(z3,x2)|(z4,x2)|
|(x1,x3)|(x2,x3)|(x3,x3)|(x4,x3)|.|(y1,x3)|(y2,x3)|(y3,x3)|(y4,x3)|.|(z1,x3)|(z2,x3)|(z3,x3)|(z4,x3)|
|(x1,x4)|(x2,x4)|(x3,x4)|(x4,x4)|.|(y1,x4)|(y2,x4)|(y3,x4)|(y4,x4)|.|(z1,x4)|(z2,x4)|(z3,x4)|(z4,x4)|
. . . . . . . . . . . . . . .
|(x1,y1)|(x2,y1)|(x3,y1)|(x4,y1)|.|(y1,y1)|(y2,y1)|(y3,x1)|(y4,y1)|.|(z1,y1)|(z2,y1)|(z3,y1)|(z4,y1)|
|(x1,y2)|(x2,y2)|(x3,y2)|(x4,y2)|.|(y1,y2)|(y2,y2)|(y3,x2)|(y4,y2)|.|(z1,y2)|(z2,y2)|(z3,y2)|(z4,y2)|
|(x1,y3)|(x2,y3)|(x3,y3)|(x4,y3)|.|(y1,y3)|(y2,y3)|(y3,x3)|(y4,y3)|.|(z1,y3)|(z2,y3)|(z3,y3)|(z4,y3)|
|(x1,y4)|(x2,y4)|(x3,y4)|(x4,y4)|.|(y1,y4)|(y2,y4)|(y3,x4)|(y4,y4)|.|(z1,y4)|(z2,y4)|(z3,y4)|(z4,y4)|
. . . . . . . . . . . . . . .
|(x1,z1)|(x2,z1)|(x3,z1)|(x4,z1)|.|(y1,z1)|(y2,z1)|(y3,z1)|(y4,z1)|.|(z1,z1)|(z2,z1)|(z3,z1)|(z4,z1)|
|(x1,z2)|(x2,z2)|(x3,z2)|(x4,z2)|.|(y1,z2)|(y2,z2)|(y3,z2)|(y4,z2)|.|(z1,z2)|(z2,z2)|(z3,z2)|(z4,z2)|
|(x1,z3)|(x2,z3)|(x3,z3)|(x4,z3)|.|(y1,z3)|(y2,z3)|(y3,z3)|(y4,z3)|.|(z1,z3)|(z2,z3)|(z3,z3)|(z4,z3)|
|(x1,z4)|(x2,z4)|(x3,z4)|(x4,z4)|.|(y1,z4)|(y2,z4)|(y3,z4)|(y4,z4)|.|(z1,z4)|(z2,z4)|(z3,z4)|(z4,z4)|
bFieldmagnetic field in the center of the detector in Gev/cm/c, in order to perform the p_t calculation.
Warning
see BL_Circle_fit(), BL_Line_fit() and Fast_fit() warnings.
Returns
(phi,Tip,p_t,cot(theta)),Zip), their covariance matrix and the chi2's of the circle and line fits.

Definition at line 566 of file BrokenLine.h.

568  {
569  riemannFit::HelixFit helix;
570  riemannFit::Vector4d fast_fit;
571  fastFit(hits, fast_fit);
572 
573  PreparedBrokenLineData<n> data;
574  karimaki_circle_fit circle;
576  riemannFit::Matrix3d jacobian;
577 
578  prepareBrokenLineData(hits, fast_fit, bField, data);
579  lineFit(hits_ge, fast_fit, bField, data, line);
580  circleFit(hits, hits_ge, fast_fit, bField, data, circle);
581 
582  // the circle fit gives k, but here we want p_t, so let's change the parameter and the covariance matrix
583  jacobian << 1., 0, 0, 0, 1., 0, 0, 0,
584  -std::abs(circle.par(2)) * bField / (riemannFit::sqr(circle.par(2)) * circle.par(2));
585  circle.par(2) = bField / std::abs(circle.par(2));
586  circle.cov = jacobian * circle.cov * jacobian.transpose();
587 
588  helix.par << circle.par, line.par;
589  helix.cov = riemannFit::MatrixXd::Zero(5, 5);
590  helix.cov.block(0, 0, 3, 3) = circle.cov;
591  helix.cov.block(3, 3, 2, 2) = line.cov;
592  helix.qCharge = circle.qCharge;
593  helix.chi2_circle = circle.chi2;
594  helix.chi2_line = line.chi2;
595 
596  return helix;
597  }

References funct::abs(), Calorimetry_cff::bField, riemannFit::CircleFit::chi2, riemannFit::HelixFit::chi2_circle, riemannFit::HelixFit::chi2_line, circleFit(), riemannFit::CircleFit::cov, riemannFit::HelixFit::cov, data, fastFit(), hfClusterShapes_cfi::hits, mps_splice::line, lineFit(), riemannFit::CircleFit::par, riemannFit::HelixFit::par, prepareBrokenLineData(), riemannFit::CircleFit::qCharge, riemannFit::HelixFit::qCharge, and riemannFit::sqr().

Referenced by PixelNtupletsFitter::run().

◆ lineFit()

template<typename V4 , typename M6xN , int n>
__host__ __device__ void brokenline::lineFit ( const M6xN &  hits_ge,
const V4 &  fast_fit,
const double  bField,
const PreparedBrokenLineData< n > &  data,
riemannFit::LineFit line_results 
)
inline

Performs the Broken Line fit in the straight track case (that is, the fit parameters are only the interceptions u).

Parameters
hitshits coordinates.
fast_fitpre-fit result in the form (X0,Y0,R,tan(theta)).
bFieldmagnetic field in Gev/cm/c.
dataPreparedBrokenLineData.
line_resultsstruct to be filled with the results in this form: -par parameter of the line in this form: (cot(theta), Zip);
-cov covariance matrix of the fitted parameter;
-chi2 value of the cost function in the minimum.

The function implements the steps 2 and 3 of the Broken Line fit without the curvature correction.
The step 2 is the least square fit, done by imposing the minimum constraint on the cost function and solving the consequent linear system. It determines the fitted parameters u and their covariance matrix. The step 3 is the correction of the fast pre-fitted parameters for the innermost part of the track. It is first done in a comfortable coordinate system (the one in which the first hit is the origin) and then the parameters and their covariance matrix are transformed to the original coordinate system.

Definition at line 445 of file BrokenLine.h.

449  {
450  const auto& radii = data.radii;
451  const auto& sTotal = data.sTotal;
452  const auto& zInSZplane = data.zInSZplane;
453  const auto& varBeta = data.varBeta;
454 
455  const double slope = -data.qCharge / fast_fit(3);
457 
458  riemannFit::Matrix3d vMat = riemannFit::Matrix3d::Zero(); // covariance matrix XYZ
459  riemannFit::Matrix2x3d jacobXYZtosZ =
460  riemannFit::Matrix2x3d::Zero(); // jacobian for computation of the error on s (xyz -> sz)
462  for (u_int i = 0; i < n; i++) {
463  vMat(0, 0) = hits_ge.col(i)[0]; // x errors
464  vMat(0, 1) = vMat(1, 0) = hits_ge.col(i)[1]; // cov_xy
465  vMat(0, 2) = vMat(2, 0) = hits_ge.col(i)[3]; // cov_xz
466  vMat(1, 1) = hits_ge.col(i)[2]; // y errors
467  vMat(2, 1) = vMat(1, 2) = hits_ge.col(i)[4]; // cov_yz
468  vMat(2, 2) = hits_ge.col(i)[5]; // z errors
469  auto tmp = 1. / radii.block(0, i, 2, 1).norm();
470  jacobXYZtosZ(0, 0) = radii(1, i) * tmp;
471  jacobXYZtosZ(0, 1) = -radii(0, i) * tmp;
472  jacobXYZtosZ(1, 2) = 1.;
473  weights(i) = 1. / ((rotMat * jacobXYZtosZ * vMat * jacobXYZtosZ.transpose() * rotMat.transpose())(
474  1, 1)); // compute the orthogonal weight point by point
475  }
476 
478  for (u_int i = 0; i < n; i++) {
479  r_u(i) = weights(i) * zInSZplane(i);
480  }
481 #ifdef CPP_DUMP
482  std::cout << "CU4\n" << matrixC_u(w, sTotal, varBeta) << std::endl;
483 #endif
485  math::cholesky::invert(matrixC_u(weights, sTotal, varBeta), iMat);
486 #ifdef CPP_DUMP
487  std::cout << "I4\n" << iMat << std::endl;
488 #endif
489 
490  riemannFit::VectorNd<n> uVec = iMat * r_u; // obtain the fitted parameters by solving the linear system
491 
492  // line parameters in the system in which the first hit is the origin and with axis along SZ
493  line_results.par << (uVec(1) - uVec(0)) / (sTotal(1) - sTotal(0)), uVec(0);
494  auto idiff = 1. / (sTotal(1) - sTotal(0));
495  line_results.cov << (iMat(0, 0) - 2 * iMat(0, 1) + iMat(1, 1)) * riemannFit::sqr(idiff) +
496  multScatt(sTotal(1) - sTotal(0), bField, fast_fit(2), 2, slope),
497  (iMat(0, 1) - iMat(0, 0)) * idiff, (iMat(0, 1) - iMat(0, 0)) * idiff, iMat(0, 0);
498 
499  // translate to the original SZ system
500  riemannFit::Matrix2d jacobian;
501  jacobian(0, 0) = 1.;
502  jacobian(0, 1) = 0;
503  jacobian(1, 0) = -sTotal(0);
504  jacobian(1, 1) = 1.;
505  line_results.par(1) += -line_results.par(0) * sTotal(0);
506  line_results.cov = jacobian * line_results.cov * jacobian.transpose();
507 
508  // rotate to the original sz system
509  auto tmp = rotMat(0, 0) - line_results.par(0) * rotMat(0, 1);
510  jacobian(1, 1) = 1. / tmp;
511  jacobian(0, 0) = jacobian(1, 1) * jacobian(1, 1);
512  jacobian(0, 1) = 0;
513  jacobian(1, 0) = line_results.par(1) * rotMat(0, 1) * jacobian(0, 0);
514  line_results.par(1) = line_results.par(1) * jacobian(1, 1);
515  line_results.par(0) = (rotMat(0, 1) + line_results.par(0) * rotMat(0, 0)) * jacobian(1, 1);
516  line_results.cov = jacobian * line_results.cov * jacobian.transpose();
517 
518  // compute chi2
519  line_results.chi2 = 0;
520  for (u_int i = 0; i < n; i++) {
521  line_results.chi2 += weights(i) * riemannFit::sqr(zInSZplane(i) - uVec(i));
522  if (i > 0 && i < n - 1)
523  line_results.chi2 += riemannFit::sqr(uVec(i - 1) / (sTotal(i) - sTotal(i - 1)) -
524  uVec(i) * (sTotal(i + 1) - sTotal(i - 1)) /
525  ((sTotal(i + 1) - sTotal(i)) * (sTotal(i) - sTotal(i - 1))) +
526  uVec(i + 1) / (sTotal(i + 1) - sTotal(i))) /
527  varBeta(i);
528  }
529  }

References Calorimetry_cff::bField, riemannFit::LineFit::chi2, gather_cfg::cout, riemannFit::LineFit::cov, data, mps_fire::i, l1tstage2_dqm_sourceclient-live_cfg::invert, matrixC_u(), multScatt(), dqmiodumpmetadata::n, riemannFit::LineFit::par, rotationMatrix(), slope, riemannFit::sqr(), createJobs::tmp, w, and HLT_FULL_cff::weights.

Referenced by __attribute__(), and helixFit().

◆ matrixC_u()

template<int n>
__host__ __device__ riemannFit::MatrixNd<n> brokenline::matrixC_u ( const riemannFit::VectorNd< n > &  weights,
const riemannFit::VectorNd< n > &  sTotal,
const riemannFit::VectorNd< n > &  varBeta 
)
inline

Computes the n-by-n band matrix obtained minimizing the Broken Line's cost function w.r.t u. This is the whole matrix in the case of the line fit and the main n-by-n block in the case of the circle fit.

Parameters
weightsweights of the first part of the cost function, the one with the measurements and not the angles (\sum_{i=1}^n w*(y_i-u_i)^2).
sTotaltotal distance traveled by the particle from the pre-fitted closest approach.
varBetakink angles' variance.
Returns
the n-by-n matrix of the linear system

Definition at line 206 of file BrokenLine.h.

208  {
210  for (u_int i = 0; i < n; i++) {
211  c_uMat(i, i) = weights(i);
212  if (i > 1)
213  c_uMat(i, i) += 1. / (varBeta(i - 1) * riemannFit::sqr(sTotal(i) - sTotal(i - 1)));
214  if (i > 0 && i < n - 1)
215  c_uMat(i, i) +=
216  (1. / varBeta(i)) * riemannFit::sqr((sTotal(i + 1) - sTotal(i - 1)) /
217  ((sTotal(i + 1) - sTotal(i)) * (sTotal(i) - sTotal(i - 1))));
218  if (i < n - 2)
219  c_uMat(i, i) += 1. / (varBeta(i + 1) * riemannFit::sqr(sTotal(i + 1) - sTotal(i)));
220 
221  if (i > 0 && i < n - 1)
222  c_uMat(i, i + 1) =
223  1. / (varBeta(i) * (sTotal(i + 1) - sTotal(i))) *
224  (-(sTotal(i + 1) - sTotal(i - 1)) / ((sTotal(i + 1) - sTotal(i)) * (sTotal(i) - sTotal(i - 1))));
225  if (i < n - 2)
226  c_uMat(i, i + 1) +=
227  1. / (varBeta(i + 1) * (sTotal(i + 1) - sTotal(i))) *
228  (-(sTotal(i + 2) - sTotal(i)) / ((sTotal(i + 2) - sTotal(i + 1)) * (sTotal(i + 1) - sTotal(i))));
229 
230  if (i < n - 2)
231  c_uMat(i, i + 2) = 1. / (varBeta(i + 1) * (sTotal(i + 2) - sTotal(i + 1)) * (sTotal(i + 1) - sTotal(i)));
232 
233  c_uMat(i, i) *= 0.5;
234  }
235  return c_uMat + c_uMat.transpose();
236  }

References mps_fire::i, dqmiodumpmetadata::n, riemannFit::sqr(), and HLT_FULL_cff::weights.

Referenced by circleFit(), and lineFit().

◆ multScatt()

__host__ __device__ double brokenline::multScatt ( const double &  length,
const double  bField,
const double  radius,
int  layer,
double  slope 
)
inline

Computes the Coulomb multiple scattering variance of the planar angle.

Parameters
lengthlength of the track in the material.
bFieldmagnetic field in Gev/cm/c.
radiusradius of curvature (needed to evaluate p).
layerdenotes which of the four layers of the detector is the endpoint of the multiple scattered track. For example, if Layer=3, then the particle has just gone through the material between the second and the third layer.
Warning
the formula used here assumes beta=1, and so neglects the dependence of theta_0 on the mass of the particle at fixed momentum.
Returns
the variance of the planar angle ((theta_0)^2 /3).

< inverse of radiation length of the material in cm

number between 1/3 (uniform material) and 1 (thin scatterer) to be manually tuned

Definition at line 53 of file BrokenLine.h.

54  {
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  }

References funct::abs(), Calorimetry_cff::bField, fact, dqm-mbProfile::log, min(), HLT_FULL_cff::pt2, CosmicsPD_Skims::radius, slope, and riemannFit::sqr().

Referenced by circleFit(), lineFit(), and prepareBrokenLineData().

◆ prepareBrokenLineData()

template<typename M3xN , typename V4 , int n>
__host__ __device__ void brokenline::prepareBrokenLineData ( const M3xN &  hits,
const V4 &  fast_fit,
const double  bField,
PreparedBrokenLineData< n > &  results 
)
inline

Computes the data needed for the Broken Line fit procedure that are mainly common for the circle and the line fit.

Parameters
hitshits coordinates.
fast_fitpre-fit result in the form (X0,Y0,R,tan(theta)).
bFieldmagnetic field in Gev/cm/c.
resultsPreparedBrokenLineData to be filled (see description of PreparedBrokenLineData).

Definition at line 150 of file BrokenLine.h.

153  {
156 
157  dVec = hits.block(0, 1, 2, 1) - hits.block(0, 0, 2, 1);
158  eVec = hits.block(0, n - 1, 2, 1) - hits.block(0, n - 2, 2, 1);
159  results.qCharge = riemannFit::cross2D(dVec, eVec) > 0 ? -1 : 1;
160 
161  const double slope = -results.qCharge / fast_fit(3);
162 
164 
165  // calculate radii and s
166  results.radii = hits.block(0, 0, 2, n) - fast_fit.head(2) * riemannFit::MatrixXd::Constant(1, n, 1);
167  eVec = -fast_fit(2) * fast_fit.head(2) / fast_fit.head(2).norm();
168  for (u_int i = 0; i < n; i++) {
169  dVec = results.radii.block(0, i, 2, 1);
170  results.sTransverse(i) = results.qCharge * fast_fit(2) *
171  atan2(riemannFit::cross2D(dVec, eVec), dVec.dot(eVec)); // calculates the arc length
172  }
173  riemannFit::VectorNd<n> zVec = hits.block(2, 0, 1, n).transpose();
174 
175  //calculate sTotal and zVec
177  for (u_int i = 0; i < n; i++) {
178  pointsSZ(0, i) = results.sTransverse(i);
179  pointsSZ(1, i) = zVec(i);
180  pointsSZ.block(0, i, 2, 1) = rotMat * pointsSZ.block(0, i, 2, 1);
181  }
182  results.sTotal = pointsSZ.block(0, 0, 1, n).transpose();
183  results.zInSZplane = pointsSZ.block(1, 0, 1, n).transpose();
184 
185  //calculate varBeta
186  results.varBeta(0) = results.varBeta(n - 1) = 0;
187  for (u_int i = 1; i < n - 1; i++) {
188  results.varBeta(i) = multScatt(results.sTotal(i + 1) - results.sTotal(i), bField, fast_fit(2), i + 2, slope) +
189  multScatt(results.sTotal(i) - results.sTotal(i - 1), bField, fast_fit(2), i + 1, slope);
190  }
191  }

References Calorimetry_cff::bField, riemannFit::cross2D(), hfClusterShapes_cfi::hits, mps_fire::i, multScatt(), dqmiodumpmetadata::n, bookConverter::results, rotationMatrix(), and slope.

Referenced by __attribute__(), and helixFit().

◆ rotationMatrix()

__host__ __device__ riemannFit::Matrix2d brokenline::rotationMatrix ( double  slope)
inline

Computes the 2D rotation matrix that transforms the line y=slope*x into the line y=0.

Parameters
slopetangent of the angle of rotation.
Returns
2D rotation matrix.

Definition at line 77 of file BrokenLine.h.

77  {
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  }

References makeMuonMisalignmentScenario::rot, slope, riemannFit::sqr(), and mathSSE::sqrt().

Referenced by circleFit(), lineFit(), prepareBrokenLineData(), and reco::BeamSpot::rotatedCovariance3D().

◆ translateKarimaki()

__host__ __device__ void brokenline::translateKarimaki ( karimaki_circle_fit circle,
double  x0,
double  y0,
riemannFit::Matrix3d jacobian 
)
inline

Changes the Karimäki parameters (and consequently their covariance matrix) under a translation of the coordinate system, such that the old origin has coordinates (x0,y0) in the new coordinate system. The formulas are taken from Karimäki V., 1990, Effective circle fitting for particle trajectories, Nucl. Instr. and Meth. A305 (1991) 187.

Parameters
circlecircle fit in the old coordinate system. circle.par(0) is phi, circle.par(1) is d and circle.par(2) is rho.
x0x coordinate of the translation vector.
y0y coordinate of the translation vector.
jacobianpassed by reference in order to save stack.

Definition at line 97 of file BrokenLine.h.

100  {
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  }

References funct::cos(), Phase1L1TJetProducer_cfi::cosPhi, riemannFit::CircleFit::cov, amptDefaultParameters_cff::mu, riemannFit::CircleFit::par, hitfit::scalar(), funct::sin(), Phase1L1TJetProducer_cfi::sinPhi, riemannFit::sqr(), mathSSE::sqrt(), and protons_cff::xi.

Referenced by circleFit().

riemannFit::Vector4d
Eigen::Vector4d Vector4d
Definition: FitResult.h:15
riemannFit::LineFit::par
Vector2d par
(cotan(theta),Zip)
Definition: FitResult.h:39
riemannFit::sqr
constexpr T sqr(const T a)
raise to square.
Definition: FitUtils.h:67
mps_fire.i
i
Definition: mps_fire.py:428
funct::Constant
Polynomial< 0 > Constant
Definition: Constant.h:6
dqmiodumpmetadata.n
n
Definition: dqmiodumpmetadata.py:28
riemannFit::LineFit::cov
Matrix2d cov
Definition: FitResult.h:40
riemannFit::LineFit::chi2
double chi2
Definition: FitResult.h:45
hfClusterShapes_cfi.hits
hits
Definition: hfClusterShapes_cfi.py:5
brokenline::prepareBrokenLineData
__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
amptDefaultParameters_cff.mu
mu
Definition: amptDefaultParameters_cff.py:16
min
T min(T a, T b)
Definition: MathUtil.h:58
riemannFit::HelixFit::qCharge
int32_t qCharge
particle charge
Definition: FitResult.h:61
gather_cfg.cout
cout
Definition: gather_cfg.py:144
ZElectronSkim_cff.rho
rho
Definition: ZElectronSkim_cff.py:38
riemannFit::VectorNplusONEd
Eigen::Matrix< double, N+1, 1 > VectorNplusONEd
Definition: FitUtils.h:35
bookConverter.results
results
Definition: bookConverter.py:144
brokenline::matrixC_u
__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's cost function w....
Definition: BrokenLine.h:206
brokenline::circleFit
__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:296
createJobs.tmp
tmp
align.sh
Definition: createJobs.py:716
brokenline::karimaki_circle_fit
riemannFit::CircleFit karimaki_circle_fit
< Karimäki's parameters: (phi, d, k=1/R)
Definition: BrokenLine.h:18
riemannFit::Matrix2d
Eigen::Matrix2d Matrix2d
Definition: FitResult.h:17
Phase1L1TJetProducer_cfi.cosPhi
cosPhi
Definition: Phase1L1TJetProducer_cfi.py:40
funct::sin
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
riemannFit::HelixFit::cov
Matrix5d cov
Definition: FitResult.h:50
brokenline::rotationMatrix
__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
funct::cos
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
riemannFit::MatrixNd
Eigen::Matrix< double, N, N > MatrixNd
Definition: FitUtils.h:15
brokenline::lineFit
__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:445
w
const double w
Definition: UKUtility.cc:23
HLT_FULL_cff.weights
weights
Definition: HLT_FULL_cff.py:99207
Phase1L1TJetProducer_cfi.sinPhi
sinPhi
Definition: Phase1L1TJetProducer_cfi.py:39
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
hitfit::scalar
double scalar(const CLHEP::HepGenMatrix &m)
Return the matrix as a scalar. Raise an assertion if the matris is not .
Definition: matutil.cc:166
riemannFit::Matrix2xNd
Eigen::Matrix< double, 2, N > Matrix2xNd
Definition: FitUtils.h:25
riemannFit::HelixFit::chi2_line
float chi2_line
Definition: FitResult.h:59
b
double b
Definition: hdecay.h:118
riemannFit::MatrixNplusONEd
Eigen::Matrix< double, N+1, N+1 > MatrixNplusONEd
Definition: FitUtils.h:17
a
double a
Definition: hdecay.h:119
PVValHelper::phi
Definition: PVValidationHelpers.h:69
type
type
Definition: SiPixelVCal_PayloadInspector.cc:39
riemannFit::HelixFit::chi2_circle
float chi2_circle
Definition: FitResult.h:58
fact
const double fact
Definition: NuclearInteractionFTFSimulator.cc:74
riemannFit::LineFit
Definition: FitResult.h:38
HLT_FULL_cff.pt2
pt2
Definition: HLT_FULL_cff.py:9895
riemannFit::HelixFit::par
Vector5d par
(phi,Tip,pt,cotan(theta)),Zip)
Definition: FitResult.h:49
riemannFit::VectorNd
Eigen::Matrix< double, N, 1 > VectorNd
Definition: FitUtils.h:33
protons_cff.xi
xi
Definition: protons_cff.py:31
brokenline::translateKarimaki
__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
Calorimetry_cff.bField
bField
Definition: Calorimetry_cff.py:284
l1tstage2_dqm_sourceclient-live_cfg.invert
invert
Definition: l1tstage2_dqm_sourceclient-live_cfg.py:86
riemannFit::HelixFit
Definition: FitResult.h:48
riemannFit::Vector2d
Eigen::Vector2d Vector2d
Definition: FitResult.h:13
makeMuonMisalignmentScenario.rot
rot
Definition: makeMuonMisalignmentScenario.py:322
CosmicsPD_Skims.radius
radius
Definition: CosmicsPD_Skims.py:135
dqm-mbProfile.log
log
Definition: dqm-mbProfile.py:17
data
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:79
brokenline::multScatt
__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
ztail.d
d
Definition: ztail.py:151
riemannFit::Matrix2x3d
Eigen::Matrix< double, 2, 3 > Matrix2x3d
Definition: FitUtils.h:45
mps_fire.result
result
Definition: mps_fire.py:311
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
brokenline::fastFit
__host__ __device__ void fastFit(const M3xN &hits, V4 &result)
A very fast helix fit.
Definition: BrokenLine.h:249
riemannFit::cross2D
__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
riemannFit::Matrix3d
Eigen::Matrix3d Matrix3d
Definition: FitResult.h:18
c
auto & c
Definition: CAHitNtupletGeneratorKernelsImpl.h:56
slope
static const double slope[3]
Definition: CastorTimeSlew.cc:6
mps_splice.line
line
Definition: mps_splice.py:76
MillePedeFileConverter_cfg.e
e
Definition: MillePedeFileConverter_cfg.py:37