CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
FastCircleFit Class Reference

#include <FastCircleFit.h>

Public Member Functions

float chi2 () const
 
template<typename P , typename E >
 FastCircleFit (const P &points, const E &errors)
 
template<size_t N>
 FastCircleFit (const std::array< GlobalPoint, N > &points, const std::array< GlobalError, N > &errors)
 
float rho () const
 
float x0 () const
 
float y0 () const
 
 ~FastCircleFit ()=default
 

Private Member Functions

template<typename P , typename E , typename C >
void calculate (const P &points, const E &errors, C &x, C &y, C &z, C &weight)
 
template<typename T >
T sqr (T t)
 

Private Attributes

float chi2_
 
float rho_
 
float x0_
 
float y0_
 

Detailed Description

Same (almost) as FastCircle but with arbitrary number of hits.

Strandlie, Wroldsen, Frühwirth NIM A 488 (2002) 332-341. Frühwirth, Strandlie, Waltenberger NIM A 490 (2002) 366-378.

Definition at line 27 of file FastCircleFit.h.

Constructor & Destructor Documentation

◆ FastCircleFit() [1/2]

template<typename P , typename E >
FastCircleFit::FastCircleFit ( const P points,
const E &  errors 
)
inline

Constructor for containers of GlobalPoint and GlobalError

Template Parameters
PContainer of GlobalPoint
EContainer of GlobalError

Container can be e.g. std::vector or DynArray.

Definition at line 38 of file FastCircleFit.h.

References calculate(), declareDynArray, N, hiPixelPairStep_cff::points, x, y, and z.

38  {
39  const auto N = points.size();
40  declareDynArray(float, N, x);
41  declareDynArray(float, N, y);
42  declareDynArray(float, N, z);
43  declareDynArray(float, N, weight); // 1/sigma^2
45  }
Definition: weight.py:1
void calculate(const P &points, const E &errors, C &x, C &y, C &z, C &weight)
Definition: FastCircleFit.h:84
#define N
Definition: blowfish.cc:9
Definition: errors.py:1
#define declareDynArray(T, n, x)
Definition: DynArray.h:91

◆ FastCircleFit() [2/2]

template<size_t N>
FastCircleFit::FastCircleFit ( const std::array< GlobalPoint, N > &  points,
const std::array< GlobalError, N > &  errors 
)
inline

Constructor for std::array of GlobalPoint and GlobalError

Definition at line 51 of file FastCircleFit.h.

References calculate(), hiPixelPairStep_cff::points, mps_merge::weight, x, y, and z.

51  {
52  std::array<float, N> x;
53  std::array<float, N> y;
54  std::array<float, N> z;
55  std::array<float, N> weight; // 1/sigma^2
57  }
Definition: weight.py:1
void calculate(const P &points, const E &errors, C &x, C &y, C &z, C &weight)
Definition: FastCircleFit.h:84
Definition: errors.py:1

◆ ~FastCircleFit()

FastCircleFit::~FastCircleFit ( )
default

Member Function Documentation

◆ calculate()

template<typename P , typename E , typename C >
void FastCircleFit::calculate ( const P points,
const E &  errors,
C &  x,
C &  y,
C &  z,
C &  weight 
)
private

Definition at line 84 of file FastCircleFit.h.

References A, HltBtagPostValidation_cff::c, chi2_, HLT_2022v15_cff::distance, mps_fire::i, SiStripPI::mean, N, dqmiodumpmetadata::n, AlCaHLTBitMon_ParallelJobs::p, point, hiPixelPairStep_cff::points, PRINT, rho_, sqr(), mathSSE::sqrt(), createJobs::tmp, w(), x, x0_, y, y0_, and z.

Referenced by FastCircleFit().

84  {
85  const auto N = points.size();
86 
87  // transform
88  for (size_t i = 0; i < N; ++i) {
89  const auto& point = points[i];
90  const auto& p = point.basicVector();
91  x[i] = p.x();
92  y[i] = p.y();
93  z[i] = sqr(p.x()) + sqr(p.y());
94 
95  // TODO: The following accounts only the hit uncertainty in phi.
96  // Albeit it "works nicely", should we account also the
97  // uncertainty in r in FPix (and the correlation between phi and r)?
98  const float phiErr2 = errors[i].phierr(point);
99  const float w = 1.f / (point.perp2() * phiErr2);
100  weight[i] = w;
101 
102  PRINT << " point " << p.x() << " " << p.y() << " transformed " << x[i] << " " << y[i] << " " << z[i] << " weight "
103  << w << std::endl;
104  }
105  const float invTotWeight = 1.f / std::accumulate(weight.begin(), weight.end(), 0.f);
106  PRINT << " invTotWeight " << invTotWeight;
107 
108  Eigen::Vector3f mean = Eigen::Vector3f::Zero();
109  for (size_t i = 0; i < N; ++i) {
110  const float w = weight[i];
111  mean[0] += w * x[i];
112  mean[1] += w * y[i];
113  mean[2] += w * z[i];
114  }
115  mean *= invTotWeight;
116  PRINT << " mean " << mean[0] << " " << mean[1] << " " << mean[2] << std::endl;
117 
118  Eigen::Matrix3f A = Eigen::Matrix3f::Zero();
119  for (size_t i = 0; i < N; ++i) {
120  const auto diff_x = x[i] - mean[0];
121  const auto diff_y = y[i] - mean[1];
122  const auto diff_z = z[i] - mean[2];
123  const auto w = weight[i];
124 
125  // exploit the fact that only lower triangular part of the matrix
126  // is used in the eigenvalue calculation
127  A(0, 0) += w * diff_x * diff_x;
128  A(1, 0) += w * diff_y * diff_x;
129  A(1, 1) += w * diff_y * diff_y;
130  A(2, 0) += w * diff_z * diff_x;
131  A(2, 1) += w * diff_z * diff_y;
132  A(2, 2) += w * diff_z * diff_z;
133  PRINT << " i " << i << " A " << A(0, 0) << " " << A(0, 1) << " " << A(0, 2) << std::endl
134  << " " << A(1, 0) << " " << A(1, 1) << " " << A(1, 2) << std::endl
135  << " " << A(2, 0) << " " << A(2, 1) << " " << A(2, 2) << std::endl;
136  }
137  A *= 1. / static_cast<float>(N);
138 
139  PRINT << " A " << A(0, 0) << " " << A(0, 1) << " " << A(0, 2) << std::endl
140  << " " << A(1, 0) << " " << A(1, 1) << " " << A(1, 2) << std::endl
141  << " " << A(2, 0) << " " << A(2, 1) << " " << A(2, 2) << std::endl;
142 
143  // find eigenvalues and vectors
144  Eigen::SelfAdjointEigenSolver<Eigen::Matrix3f> eigen(A);
145  const auto& eigenValues = eigen.eigenvalues();
146  const auto& eigenVectors = eigen.eigenvectors();
147 
148  // in Eigen the first eigenvalue is the smallest
149  PRINT << " eigenvalues " << eigenValues[0] << " " << eigenValues[1] << " " << eigenValues[2] << std::endl;
150 
151  PRINT << " eigen " << eigenVectors(0, 0) << " " << eigenVectors(0, 1) << " " << eigenVectors(0, 2) << std::endl
152  << " " << eigenVectors(1, 0) << " " << eigenVectors(1, 1) << " " << eigenVectors(1, 2) << std::endl
153  << " " << eigenVectors(2, 0) << " " << eigenVectors(2, 1) << " " << eigenVectors(2, 2) << std::endl;
154 
155  // eivenvector corresponding smallest eigenvalue
156  auto n = eigenVectors.col(0);
157  PRINT << " n (unit) " << n[0] << " " << n[1] << " " << n[2] << std::endl;
158 
159  const float c = -1.f * (n[0] * mean[0] + n[1] * mean[1] + n[2] * mean[2]);
160 
161  PRINT << " c " << c << std::endl;
162 
163  // calculate fit parameters
164  const auto tmp = 0.5f * 1.f / n[2];
165  x0_ = -n[0] * tmp;
166  y0_ = -n[1] * tmp;
167  const float rho2 = (1 - sqr(n[2]) - 4 * c * n[2]) * sqr(tmp);
168  rho_ = std::sqrt(rho2);
169 
170  // calculate square sum
171  float S = 0;
172  for (size_t i = 0; i < N; ++i) {
173  const float incr = sqr(c + n[0] * x[i] + n[1] * y[i] + n[2] * z[i]) * weight[i];
174 #if defined(MK_DEBUG) || defined(EDM_ML_DEBUG)
175  const float distance = std::sqrt(sqr(x[i] - x0_) + sqr(y[i] - y0_)) - rho_;
176  PRINT << " i " << i << " chi2 incr " << incr << " d(hit, circle) " << distance << " sigma "
177  << 1.f / std::sqrt(weight[i]) << std::endl;
178 #endif
179  S += incr;
180  }
181  chi2_ = S;
182 }
Eigen::Vector3f Vector3f
Definition: FitUtils.h:48
#define PRINT
Definition: FastCircleFit.h:15
Eigen::Matrix3f Matrix3f
Definition: FitUtils.h:47
T w() const
Definition: weight.py:1
T sqrt(T t)
Definition: SSEVec.h:19
#define N
Definition: blowfish.cc:9
Definition: errors.py:1
Definition: APVGainStruct.h:7
tmp
align.sh
Definition: createJobs.py:716
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5

◆ chi2()

float FastCircleFit::chi2 ( void  ) const
inline

Definition at line 66 of file FastCircleFit.h.

References chi2_.

66 { return chi2_; }

◆ rho()

float FastCircleFit::rho ( ) const
inline

Definition at line 63 of file FastCircleFit.h.

References rho_.

Referenced by Lepton.Lepton::absIsoFromEA(), and Muon.Muon::absIsoWithFSR().

63 { return rho_; }

◆ sqr()

template<typename T >
T FastCircleFit::sqr ( T  t)
inlineprivate

Definition at line 73 of file FastCircleFit.h.

References submitPVValidationJobs::t.

Referenced by calculate().

73  {
74  return t * t;
75  }

◆ x0()

float FastCircleFit::x0 ( ) const
inline

Definition at line 61 of file FastCircleFit.h.

References x0_.

61 { return x0_; }

◆ y0()

float FastCircleFit::y0 ( ) const
inline

Definition at line 62 of file FastCircleFit.h.

References y0_.

62 { return y0_; }

Member Data Documentation

◆ chi2_

float FastCircleFit::chi2_
private

Definition at line 80 of file FastCircleFit.h.

Referenced by calculate(), and chi2().

◆ rho_

float FastCircleFit::rho_
private

Definition at line 79 of file FastCircleFit.h.

Referenced by calculate(), and rho().

◆ x0_

float FastCircleFit::x0_
private

Definition at line 77 of file FastCircleFit.h.

Referenced by calculate(), and x0().

◆ y0_

float FastCircleFit::y0_
private

Definition at line 78 of file FastCircleFit.h.

Referenced by calculate(), and y0().