CMS 3D CMS Logo

Functions | Variables
WeightedMeanFitter Namespace Reference

Functions

std::pair< GlobalPoint, double > nearestPoint (const GlobalPoint &vertex, reco::Track iclus)
 
TransientVertex weightedMeanOutlierRejection (const std::vector< std::pair< GlobalPoint, GlobalPoint >> &points, std::vector< reco::TransientTrack > iclus)
 
TransientVertex weightedMeanOutlierRejectionBeamSpot (const std::vector< std::pair< GlobalPoint, GlobalPoint >> &points, std::vector< reco::TransientTrack > iclus, const reco::BeamSpot &beamSpot)
 
TransientVertex weightedMeanOutlierRejectionVarianceAsError (const std::vector< std::pair< GlobalPoint, GlobalPoint >> &points, std::vector< std::vector< reco::TransientTrack >>::const_iterator iclus)
 

Variables

constexpr float corr_x = 1.2
 
constexpr float corr_x_bs = 1.0
 
constexpr float corr_z = 1.4
 
constexpr int maxIterations = 50
 
constexpr float muSquare = 9.
 
constexpr float precision = 1e-24
 
constexpr float startError = 20.0
 

Function Documentation

◆ nearestPoint()

std::pair<GlobalPoint, double> WeightedMeanFitter::nearestPoint ( const GlobalPoint vertex,
reco::Track  iclus 
)
inline

Definition at line 20 of file WeightedMeanFitter.h.

References AlCaHLTBitMon_ParallelJobs::p, conifer::pow(), reco::TrackBase::px(), reco::TrackBase::py(), reco::TrackBase::pz(), mathSSE::sqrt(), submitPVValidationJobs::t, bphysicsOniaDQM_cfi::vertex, run3scouting_cff::vx, reco::TrackBase::vx(), run3scouting_cff::vy, reco::TrackBase::vy(), run3scouting_cff::vz, and reco::TrackBase::vz().

Referenced by weightedMeanOutlierRejection(), and weightedMeanOutlierRejectionBeamSpot().

20  {
21  double ox = iclus.vx();
22  double oy = iclus.vy();
23  double oz = iclus.vz();
24 
25  double vx = iclus.px();
26  double vy = iclus.py();
27  double vz = iclus.pz();
28 
29  double opx = vertex.x() - ox;
30  double opy = vertex.y() - oy;
31  double opz = vertex.z() - oz;
32 
33  double vnorm2 = (vx * vx + vy * vy + vz * vz);
34  double t = (vx * opx + vy * opy + vz * opz) / (vnorm2);
35 
36  GlobalPoint p(ox + t * vx, oy + t * vy, oz + t * vz);
37  return std::pair<GlobalPoint, double>(
38  p,
39  std::sqrt(std::pow(p.x() - vertex.x(), 2) + std::pow(p.y() - vertex.y(), 2) + std::pow(p.z() - vertex.z(), 2)));
40  }
double vx() const
x coordinate of the reference point on track
Definition: TrackBase.h:655
double px() const
x coordinate of momentum vector
Definition: TrackBase.h:640
double py() const
y coordinate of momentum vector
Definition: TrackBase.h:643
constexpr int pow(int x)
Definition: conifer.h:24
double vz() const
z coordinate of the reference point on track
Definition: TrackBase.h:661
T sqrt(T t)
Definition: SSEVec.h:19
double pz() const
z coordinate of momentum vector
Definition: TrackBase.h:646
double vy() const
y coordinate of the reference point on track
Definition: TrackBase.h:658

◆ weightedMeanOutlierRejection()

TransientVertex WeightedMeanFitter::weightedMeanOutlierRejection ( const std::vector< std::pair< GlobalPoint, GlobalPoint >> &  points,
std::vector< reco::TransientTrack iclus 
)
inline

Definition at line 42 of file WeightedMeanFitter.h.

References funct::abs(), nano_mu_local_reco_cff::chi2, corr_x, corr_z, submitPVResolutionJobs::err, err_x, err_z, mps_fire::i, createfilelist::int, muSquare, nearestPoint(), AlCaHLTBitMon_ParallelJobs::p, hiPixelPairStep_cff::points, conifer::pow(), precision, edm::second(), startError, HLT_2023v12_cff::track, findQualityFiles::v, and x.

Referenced by PrimaryVertexProducer::produce().

43  {
44  float x = 0., y = 0., z = 0.;
45  float s_wx = 0., s_wz = 0.;
46  float s2_wx = 0., s2_wz = 0.;
47  float wx = 0., wz = 0., chi2 = 0.;
48  float ndof_x = 0.;
49 
51  err(0, 0) = startError / 10 * startError / 10;
52  err(1, 1) = startError / 10 * startError / 10;
53  err(2, 2) = startError * startError; // error is 20 cm, so cov -> is 20 ^ 2
54  for (const auto& p : points) {
55  wx = p.second.x() <= precision ? 1. / std::pow(precision, 2) : 1. / std::pow(p.second.x(), 2);
56 
57  wz = p.second.z() <= precision ? 1. / std::pow(precision, 2) : 1. / std::pow(p.second.z(), 2);
58 
59  x += p.first.x() * wx;
60  y += p.first.y() * wx;
61  z += p.first.z() * wz;
62 
63  s_wx += wx;
64  s_wz += wz;
65  }
66 
67  if (s_wx == 0. || s_wz == 0.) {
68  edm::LogWarning("WeightedMeanFitter") << "Vertex fitting failed at beginning\n";
69  return TransientVertex(GlobalPoint(0, 0, 0), err, iclus, 0, 0);
70  }
71 
72  x /= s_wx;
73  y /= s_wx;
74  z /= s_wz;
75 
76  float old_x, old_y, old_z;
77 
78  float xpull;
79  int niter = 0;
80 
81  float err_x, err_z;
82 
83  err_x = 1. / s_wx;
84  err_z = 1. / s_wz;
85 
86  while ((niter++) < 2) {
87  old_x = x;
88  old_y = y;
89  old_z = z;
90  s_wx = 0;
91  s_wz = 0;
92  s2_wx = 0;
93  s2_wz = 0;
94  x = 0;
95  y = 0;
96  z = 0;
97  ndof_x = 0;
98 
99  for (unsigned int i = 0; i < (unsigned int)points.size(); i++) {
100  std::pair<GlobalPoint, double> p = nearestPoint(GlobalPoint(old_x, old_y, old_z), (iclus)[i].track());
101 
102  wx = points[i].second.x() <= precision ? std::pow(precision, 2) : std::pow(points[i].second.x(), 2);
103  wz = points[i].second.z() <= precision ? std::pow(precision, 2) : std::pow(points[i].second.z(), 2);
104 
105  xpull = 0.;
106 
107  if (std::pow(p.first.x() - old_x, 2) / (wx + err_x) < muSquare &&
108  std::pow(p.first.y() - old_y, 2) / (wx + err_x) < muSquare &&
109  std::pow(p.first.z() - old_z, 2) / (wz + err_z) < muSquare)
110  xpull = 1.;
111 
112  ndof_x += xpull;
113 
114  wx = xpull / wx;
115  wz = xpull / wz;
116 
117  x += wx * p.first.x();
118  y += wx * p.first.y();
119  z += wz * p.first.z();
120 
121  s_wx += wx;
122  s_wz += wz;
123 
124  s2_wx += wx * xpull;
125  s2_wz += wz * xpull;
126  }
127 
128  if (s_wx == 0. || s_wz == 0.) {
129  edm::LogWarning("WeightedMeanFitter")
130  << "Vertex fitting failed, either all tracks are outliers or they have a very large error\n";
131  return TransientVertex(GlobalPoint(0, 0, 0), err, iclus, 0, 0);
132  }
133  x /= s_wx;
134  y /= s_wx;
135  z /= s_wz;
136 
137  err_x = (s2_wx / std::pow(s_wx, 2));
138  err_z = (s2_wz / std::pow(s_wz, 2));
139 
140  if (std::abs(x - old_x) < (precision / 1.) && std::abs(y - old_y) < (precision / 1.) &&
141  std::abs(z - old_z) < (precision / 1.)) {
142  break;
143  }
144  }
145 
146  err(0, 0) = err_x * corr_x * corr_x;
147  err(1, 1) = err_x * corr_x * corr_x;
148  err(2, 2) = err_z * corr_z * corr_z;
149 
150  float dist = 0;
151  for (const auto& p : points) {
152  wx = p.second.x();
153  wx = wx <= precision ? precision : wx;
154 
155  wz = p.second.z();
156  wz = wz <= precision ? precision : wz;
157 
158  dist = std::pow(p.first.x() - x, 2) / (std::pow(wx, 2) + err(0, 0));
159  dist += std::pow(p.first.y() - y, 2) / (std::pow(wx, 2) + err(1, 1));
160  dist += std::pow(p.first.z() - z, 2) / (std::pow(wz, 2) + err(2, 2));
161  chi2 += dist;
162  }
163  TransientVertex v(GlobalPoint(x, y, z), err, iclus, chi2, (int)ndof_x);
164  return v;
165  }
constexpr float corr_x
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
constexpr int pow(int x)
Definition: conifer.h:24
float float float z
U second(std::pair< T, U > const &p)
constexpr float precision
float err_x
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
float err_z
constexpr float corr_z
std::pair< GlobalPoint, double > nearestPoint(const GlobalPoint &vertex, reco::Track iclus)
constexpr float muSquare
float x
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33
Log< level::Warning, false > LogWarning
constexpr float startError

◆ weightedMeanOutlierRejectionBeamSpot()

TransientVertex WeightedMeanFitter::weightedMeanOutlierRejectionBeamSpot ( const std::vector< std::pair< GlobalPoint, GlobalPoint >> &  points,
std::vector< reco::TransientTrack iclus,
const reco::BeamSpot beamSpot 
)
inline

Definition at line 167 of file WeightedMeanFitter.h.

References funct::abs(), pwdgSkimBPark_cfi::beamSpot, nano_mu_local_reco_cff::chi2, corr_x_bs, corr_z, submitPVResolutionJobs::err, err_x, err_y, err_z, mps_fire::i, createfilelist::int, muSquare, nearestPoint(), AlCaHLTBitMon_ParallelJobs::p, hiPixelPairStep_cff::points, conifer::pow(), precision, edm::second(), startError, HLT_2023v12_cff::track, findQualityFiles::v, and x.

Referenced by PrimaryVertexProducer::produce().

170  {
171  float x = 0., y = 0., z = 0.;
172  float s_wx = 0., s_wz = 0.;
173  float s2_wx = 0., s2_wz = 0.;
174  float wx = 0., wz = 0., chi2 = 0.;
175  float wy = 0., s_wy = 0., s2_wy = 0.;
176  float ndof_x = 0.;
177 
179  err(0, 0) = startError / 10 * startError / 10;
180  err(1, 1) = startError / 10 * startError / 10;
181  err(2, 2) = startError * startError; // error is 20 cm, so cov -> is 20 ^ 2
182 
183  GlobalError bse(beamSpot.rotatedCovariance3D());
184  GlobalPoint bsp(Basic3DVector<float>(beamSpot.position()));
185 
186  for (const auto& p : points) {
187  wx = p.second.x() <= precision ? 1. / std::pow(precision, 2) : 1. / std::pow(p.second.x(), 2);
188  wy = p.second.y() <= precision ? 1. / std::pow(precision, 2) : 1. / std::pow(p.second.y(), 2);
189 
190  wz = p.second.z() <= precision ? 1. / std::pow(precision, 2) : 1. / std::pow(p.second.z(), 2);
191 
192  x += p.first.x() * wx;
193  y += p.first.y() * wy;
194  z += p.first.z() * wz;
195 
196  s_wx += wx;
197  s_wy += wy;
198  s_wz += wz;
199  }
200 
201  if (s_wx == 0. || s_wy == 0. || s_wz == 0.) {
202  edm::LogWarning("WeightedMeanFitter") << "Vertex fitting failed at beginning\n";
203  return TransientVertex(GlobalPoint(0, 0, 0), err, iclus, 0, 0);
204  }
205  // use the square of covariance element to increase it's weight: it will be the most important
206  wx = bse.cxx() <= precision ? 1. / std::pow(precision, 2) : 1. / std::pow(bse.cxx(), 2);
207  wy = bse.cyy() <= precision ? 1. / std::pow(precision, 2) : 1. / std::pow(bse.cyy(), 2);
208 
209  x += bsp.x() * wx;
210  y += bsp.y() * wy;
211 
212  x /= (s_wx + wx);
213  y /= (s_wy + wy);
214  z /= s_wz;
215 
216  float old_x, old_y, old_z;
217 
218  float xpull;
219  int niter = 0;
220 
221  float err_x, err_y, err_z;
222 
223  err_x = 1. / s_wx;
224  err_y = 1. / s_wy;
225  err_z = 1. / s_wz;
226 
227  while ((niter++) < 2) {
228  old_x = x;
229  old_y = y;
230  old_z = z;
231  s_wx = 0;
232  s_wz = 0;
233  s2_wx = 0;
234  s2_wz = 0;
235 
236  s_wy = 0;
237  s2_wy = 0;
238 
239  x = 0;
240  y = 0;
241  z = 0;
242  ndof_x = 0;
243 
244  for (unsigned int i = 0; i < (unsigned int)points.size(); i++) {
245  std::pair<GlobalPoint, double> p = nearestPoint(GlobalPoint(old_x, old_y, old_z), (iclus)[i].track());
246 
247  wx = points[i].second.x() <= precision ? std::pow(precision, 2) : std::pow(points[i].second.x(), 2);
248  wy = points[i].second.y() <= precision ? std::pow(precision, 2) : std::pow(points[i].second.y(), 2);
249 
250  wz = points[i].second.z() <= precision ? std::pow(precision, 2) : std::pow(points[i].second.z(), 2);
251 
252  xpull = 0.;
253  if (std::pow(p.first.x() - old_x, 2) / (wx + err_x) < muSquare &&
254  std::pow(p.first.y() - old_y, 2) / (wy + err_y) < muSquare &&
255  std::pow(p.first.z() - old_z, 2) / (wz + err_z) < muSquare)
256  xpull = 1.;
257 
258  ndof_x += xpull;
259 
260  wx = xpull / wx;
261  wy = xpull / wy;
262  wz = xpull / wz;
263 
264  x += wx * p.first.x();
265  y += wy * p.first.y();
266  z += wz * p.first.z();
267 
268  s_wx += wx;
269  s_wy += wy;
270  s_wz += wz;
271 
272  s2_wx += wx * xpull;
273  s2_wy += wy * xpull;
274  s2_wz += wz * xpull;
275  }
276 
277  if (s_wx == 0. || s_wy == 0. || s_wz == 0.) {
278  edm::LogWarning("WeightedMeanFitter")
279  << "Vertex fitting failed, either all tracks are outliers or they have a very large error\n";
280  return TransientVertex(GlobalPoint(0, 0, 0), err, iclus, 0, 0);
281  }
282  wx = bse.cxx() <= std::pow(precision, 2) ? 1. / std::pow(precision, 2) : 1. / bse.cxx();
283  wy = bse.cyy() <= std::pow(precision, 2) ? 1. / std::pow(precision, 2) : 1. / bse.cyy();
284 
285  x += bsp.x() * wx;
286  y += bsp.y() * wy;
287  s_wx += wx;
288  s2_wx += wx;
289  s_wy += wy;
290  s2_wy += wy;
291 
292  x /= s_wx;
293  y /= s_wy;
294  z /= s_wz;
295 
296  err_x = (s2_wx / std::pow(s_wx, 2));
297  err_y = (s2_wy / std::pow(s_wy, 2));
298  err_z = (s2_wz / std::pow(s_wz, 2));
299 
300  if (std::abs(x - old_x) < (precision) && std::abs(y - old_y) < (precision) && std::abs(z - old_z) < (precision)) {
301  break;
302  }
303  }
304  err(0, 0) = err_x * corr_x_bs * corr_x_bs;
305  err(1, 1) = err_y * corr_x_bs * corr_x_bs;
306  err(2, 2) = err_z * corr_z * corr_z;
307 
308  float dist = 0;
309  for (const auto& p : points) {
310  wx = p.second.x();
311  wx = wx <= precision ? precision : wx;
312 
313  wz = p.second.z();
314  wz = wz <= precision ? precision : wz;
315 
316  dist = std::pow(p.first.x() - x, 2) / (std::pow(wx, 2) + err(0, 0));
317  dist += std::pow(p.first.y() - y, 2) / (std::pow(wx, 2) + err(1, 1));
318  dist += std::pow(p.first.z() - z, 2) / (std::pow(wz, 2) + err(2, 2));
319  chi2 += dist;
320  }
321  TransientVertex v(GlobalPoint(x, y, z), err, iclus, chi2, (int)ndof_x);
322  return v;
323  }
constexpr float corr_x_bs
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
constexpr int pow(int x)
Definition: conifer.h:24
float float float z
U second(std::pair< T, U > const &p)
constexpr float precision
float err_x
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
float err_z
constexpr float corr_z
std::pair< GlobalPoint, double > nearestPoint(const GlobalPoint &vertex, reco::Track iclus)
constexpr float muSquare
float x
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33
float err_y
Log< level::Warning, false > LogWarning
constexpr float startError

◆ weightedMeanOutlierRejectionVarianceAsError()

TransientVertex WeightedMeanFitter::weightedMeanOutlierRejectionVarianceAsError ( const std::vector< std::pair< GlobalPoint, GlobalPoint >> &  points,
std::vector< std::vector< reco::TransientTrack >>::const_iterator  iclus 
)
inline

Definition at line 325 of file WeightedMeanFitter.h.

References funct::abs(), nano_mu_local_reco_cff::chi2, submitPVResolutionJobs::err, err_x, err_y, err_z, JetChargeProducer_cfi::exp, maxIterations, muSquare, AlCaHLTBitMon_ParallelJobs::p, hiPixelPairStep_cff::points, conifer::pow(), precision, startError, findQualityFiles::v, and x.

327  {
328  float x = 0, y = 0, z = 0, s_wx = 0, s_wy = 0, s_wz = 0, wx = 0, wy = 0, wz = 0, chi2 = 0;
329  float ndof_x = 0;
331  err(2, 2) = startError * startError; // error is 20 cm, so cov -> is 20 ^ 2
332  err(0, 0) = err(1, 1) = err(2, 2) / 100.;
333 
334  for (const auto& p : points) {
335  wx = p.second.x();
336  wx = wx <= precision ? 1. / std::pow(precision, 2) : 1. / std::pow(wx, 2);
337 
338  wz = p.second.z();
339  wz = wz <= precision ? 1. / std::pow(precision, 2) : 1. / std::pow(wz, 2);
340 
341  x += p.first.x() * wx;
342  y += p.first.y() * wx;
343  z += p.first.z() * wz;
344 
345  s_wx += wx;
346  s_wz += wz;
347  }
348 
349  if (s_wx == 0. || s_wz == 0.) {
350  edm::LogWarning("WeightedMeanFitter") << "Vertex fitting failed at beginning\n";
351  return TransientVertex(GlobalPoint(0, 0, 0), err, *iclus, 0, 0);
352  }
353 
354  x /= s_wx;
355  y /= s_wx;
356  z /= s_wz;
357 
358  float old_x, old_y, old_z;
359  float xpull;
360  int niter = 0;
361  float err_x, err_y, err_z;
362  err_x = 1. / s_wx;
363  err_y = 1. / s_wx;
364  err_z = 1. / s_wz;
365  float s_err_x = 0., s_err_y = 0., s_err_z = 0.;
366  while ((niter++) < maxIterations) {
367  old_x = x;
368  old_y = y;
369  old_z = z;
370  s_wx = 0;
371  s_wy = 0;
372  s_wz = 0;
373  x = 0;
374  y = 0;
375  z = 0;
376  s_err_x = 0.;
377  s_err_y = 0.;
378  s_err_z = 0.;
379 
380  for (const auto& p : points) {
381  wx = p.second.x();
382  wx = wx <= precision ? precision : wx;
383 
384  wy = wx * wx + err_y;
385  wx = wx * wx + err_x;
386 
387  wz = p.second.z();
388  wz = wz <= precision ? precision : wz;
389  wz = wz * wz + err_z;
390 
391  xpull = std::pow((p.first.x() - old_x), 2) / wx;
392  xpull += std::pow((p.first.y() - old_y), 2) / wy;
393  xpull += std::pow((p.first.z() - old_z), 2) / wz;
394  xpull = 1. / (1. + std::exp(-0.5 * ((muSquare)-xpull)));
395  ndof_x += xpull;
396 
397  wx = 1. / wx;
398  wy = 1. / wy;
399  wz = 1. / wz;
400 
401  wx *= xpull;
402  wy *= xpull;
403  wz *= xpull;
404 
405  x += wx * p.first.x();
406  y += wy * p.first.y();
407  z += wz * p.first.z();
408 
409  s_wx += wx;
410  s_wy += wy;
411  s_wz += wz;
412 
413  s_err_x += wx * pow(p.first.x() - old_x, 2);
414  s_err_y += wy * pow(p.first.y() - old_y, 2);
415  s_err_z += wz * pow(p.first.z() - old_z, 2);
416  }
417  if (s_wx == 0. || s_wy == 0. || s_wz == 0.) {
418  edm::LogWarning("WeightedMeanFitter")
419  << "Vertex fitting failed, either all tracks are outliers or they have a very large error\n";
420  return TransientVertex(GlobalPoint(0, 0, 0), err, *iclus, 0, 0);
421  }
422  x /= s_wx;
423  y /= s_wy;
424  z /= s_wz;
425 
426  err_x = s_err_x / s_wx;
427  err_y = s_err_y / s_wy;
428  err_z = s_err_z / s_wz;
429 
430  if (std::abs(x - old_x) < (precision / 1.) && std::abs(y - old_y) < (precision / 1.) &&
431  std::abs(z - old_z) < (precision / 1.))
432  break;
433  }
434 
435  err(0, 0) = err_x;
436  err(1, 1) = err_y;
437  err(2, 2) = err_z;
438 
439  float dist = 0.f;
440  for (const auto& p : points) {
441  wx = p.second.x();
442  wx = wx <= precision ? precision : wx;
443 
444  wz = p.second.z();
445  wz = wz <= precision ? precision : wz;
446 
447  dist = std::pow(p.first.x() - x, 2) / (std::pow(wx, 2) + std::pow(err(0, 0), 2));
448  dist += std::pow(p.first.y() - y, 2) / (std::pow(wx, 2) + std::pow(err(1, 1), 2));
449  dist += std::pow(p.first.z() - z, 2) / (std::pow(wz, 2) + std::pow(err(2, 2), 2));
450  chi2 += dist;
451  }
452  TransientVertex v(GlobalPoint(x, y, z), err, *iclus, chi2, (int)ndof_x);
453  return v;
454  }
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
constexpr int pow(int x)
Definition: conifer.h:24
float float float z
constexpr float precision
float err_x
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
float err_z
constexpr int maxIterations
constexpr float muSquare
float x
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33
float err_y
Log< level::Warning, false > LogWarning
constexpr float startError

Variable Documentation

◆ corr_x

constexpr float WeightedMeanFitter::corr_x = 1.2

Definition at line 14 of file WeightedMeanFitter.h.

Referenced by weightedMeanOutlierRejection().

◆ corr_x_bs

constexpr float WeightedMeanFitter::corr_x_bs = 1.0

Definition at line 15 of file WeightedMeanFitter.h.

Referenced by weightedMeanOutlierRejectionBeamSpot().

◆ corr_z

constexpr float WeightedMeanFitter::corr_z = 1.4

◆ maxIterations

constexpr int WeightedMeanFitter::maxIterations = 50

Definition at line 17 of file WeightedMeanFitter.h.

Referenced by weightedMeanOutlierRejectionVarianceAsError().

◆ muSquare

constexpr float WeightedMeanFitter::muSquare = 9.

◆ precision

constexpr float WeightedMeanFitter::precision = 1e-24

◆ startError

constexpr float WeightedMeanFitter::startError = 20.0