CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MatcherUsingTracksAlgorithm.cc
Go to the documentation of this file.
2 
5 
8 
15 
18 
19 // a few template-related workarounds
20 namespace reco {
21 template<>
22 inline double deltaR<GlobalVector>(const GlobalVector &v1, const GlobalVector &v2) {
23  return deltaR<float>(v1.eta(),v1.phi(),v2.eta(),v2.phi());
24 }
25 template<>
26 inline double deltaR2<GlobalVector>(const GlobalVector &v1, const GlobalVector &v2) {
27  return deltaR2<float>(v1.eta(),v1.phi(),v2.eta(),v2.phi());
28 }
29 }
30 
32  whichTrack1_(None), whichTrack2_(None),
33  whichState1_(AtVertex), whichState2_(AtVertex),
34  srcCut_(iConfig.existsAs<std::string>("srcPreselection") ? iConfig.getParameter<std::string>("srcPreselection") : ""),
35  matchedCut_(iConfig.existsAs<std::string>("matchedPreselection") ? iConfig.getParameter<std::string>("matchedPreselection") : ""),
36  requireSameCharge_(iConfig.existsAs<bool>("requireSameCharge") ? iConfig.getParameter<bool>("requireSameCharge") : false)
37 {
38  std::string algo = iConfig.getParameter<std::string>("algorithm");
39  if (algo == "byTrackRef") { algo_ = ByTrackRef; }
40  else if (algo == "byPropagatingSrc") { algo_ = ByPropagatingSrc; }
41  else if (algo == "byPropagatingMatched") { algo_ = ByPropagatingMatched; }
42  else if (algo == "byDirectComparison") { algo_ = ByDirectComparison; }
43  else throw cms::Exception("Configuration") << "Value '" << algo << "' for algorithm not yet implemented.\n";
44 
45  getConf(iConfig, "src", whichTrack1_, whichState1_);
46  getConf(iConfig, "matched", whichTrack2_, whichState2_);
47 
48  if (algo_ == ByTrackRef) {
49  // validate the config
50  if (whichTrack1_ == None || whichTrack2_ == None) throw cms::Exception("Configuration") << "Algorithm 'byTrackRef' needs tracks not to be 'none'.\n";
52  // read matching cuts
53  maxLocalPosDiff_ = iConfig.getParameter<double>("maxDeltaLocalPos");
54  maxGlobalMomDeltaR_ = iConfig.getParameter<double>("maxDeltaR");
55  maxGlobalMomDeltaEta_ = iConfig.existsAs<double>("maxDeltaEta") ? iConfig.getParameter<double>("maxDeltaEta") : maxGlobalMomDeltaR_;
56  maxGlobalMomDeltaPhi_ = iConfig.existsAs<double>("maxDeltaPhi") ? iConfig.getParameter<double>("maxDeltaPhi") : maxGlobalMomDeltaR_;
57  maxGlobalDPtRel_ = iConfig.getParameter<double>("maxDeltaPtRel");
58 
59  // choice of sorting variable
60  std::string sortBy = iConfig.getParameter<std::string>("sortBy");
61  if (sortBy == "deltaLocalPos") sortBy_ = LocalPosDiff;
62  else if (sortBy == "deltaPtRel") sortBy_ = GlobalDPtRel;
63  else if (sortBy == "deltaR") sortBy_ = GlobalMomDeltaR;
64  else if (sortBy == "deltaEta") sortBy_ = GlobalMomDeltaEta;
65  else if (sortBy == "deltaPhi") sortBy_ = GlobalMomDeltaPhi;
66  else if (sortBy == "chi2") sortBy_ = Chi2;
67  else throw cms::Exception("Configuration") << "Parameter 'sortBy' must be one of: deltaLocalPos, deltaPtRel, deltaR, chi2.\n";
68  // validate the config
69  if (algo_ == ByPropagatingSrc) {
70  if (whichTrack2_ == None || whichState2_ == AtVertex) {
72  //throw cms::Exception("Configuration") << "Destination track must be non-null, and state must not be 'AtVertex' (not yet).\n";
73  }
74  } else if (algo_ == ByPropagatingMatched) {
75  if (whichTrack1_ == None || whichState1_ == AtVertex) {
77  //throw cms::Exception("Configuration") << "Destination track must be non-null, and state must not be 'AtVertex' (not yet).\n";
78  }
79  } else if (algo_ == ByDirectComparison) {
80  bool firstAtVertex = (whichTrack1_ == None || whichState1_ == AtVertex);
81  bool secAtVertex = (whichTrack2_ == None || whichState2_ == AtVertex);
82  if (firstAtVertex) {
83  if (!secAtVertex) throw cms::Exception("Configuration") << "When using 'byDirectComparison' with 'src' at vertex (or None), 'matched' must be at vertex too.\n";
84  } else {
85  if (secAtVertex) throw cms::Exception("Configuration") << "When using 'byDirectComparison' with 'src' not at vertex, 'matched' can't be at vertex or None.\n";
86  if (whichState1_ != whichState2_) throw cms::Exception("Configuration") << "You can't use 'byDirectComparison' with non-matching states.\n";
87  }
88  }
89 
90  useChi2_ = iConfig.existsAs<bool>("computeChi2") ? iConfig.getParameter<bool>("computeChi2") : false;
91  if (useChi2_) {
92  if (whichTrack1_ == None && whichTrack2_ == None) throw cms::Exception("Configuration") << "Can't compute chi2s if both tracks are set to 'none'.\n";
93  maxChi2_ = iConfig.getParameter<double>("maxChi2");
94  chi2DiagonalOnly_ = iConfig.getParameter<bool>("chi2DiagonalOnly");
96  chi2UseVertex_ = iConfig.getParameter<bool>("chi2UsePosition");
97  } else {
98  chi2UseVertex_ = iConfig.getParameter<bool>("chi2UseVertex");
99  if (algo_ == ByDirectComparison) {
100  std::string choice = iConfig.getParameter<std::string>("chi2MomentumForDxy");
101  if (choice == "src") chi2FirstMomentum_ = true;
102  else if (choice != "matched") throw cms::Exception("Configuration") << "chi2MomentumForDxy must be 'src' or 'matched'\n";
103  }
104  }
105  } else maxChi2_ = 1;
106 
107  if (sortBy_ == Chi2 && !useChi2_) throw cms::Exception("Configuration") << "Can't sort by chi2s if 'computeChi2s' is not set to true.\n";
108  }
109 }
110 
111 void
112 MatcherUsingTracksAlgorithm::getConf(const edm::ParameterSet & iConfig, const std::string &whatFor, WhichTrack &whichTrack, WhichState &whichState) {
113  std::string s_whichTrack = iConfig.getParameter<std::string>(whatFor+"Track");
114  if (s_whichTrack == "none") { whichTrack = None; }
115  else if (s_whichTrack == "tracker") { whichTrack = TrackerTk; }
116  else if (s_whichTrack == "muon") { whichTrack = MuonTk; }
117  else if (s_whichTrack == "global") { whichTrack = GlobalTk; }
118  else throw cms::Exception("Configuration") << "Parameter 'useTrack' must be 'none', 'tracker', 'muon', 'global'\n";
119  if ((whichTrack != None) && (algo_ != ByTrackRef)) {
120  std::string s_whichState = iConfig.getParameter<std::string>(whatFor+"State");
121  if (s_whichState == "atVertex") { whichState = AtVertex; }
122  else if (s_whichState == "innermost") { whichState = Innermost; }
123  else if (s_whichState == "outermost") { whichState = Outermost; }
124  else throw cms::Exception("Configuration") << "Parameter 'useState' must be 'atVertex', 'innermost', 'outermost'\n";
125  }
126 }
127 
130 bool
131 MatcherUsingTracksAlgorithm::match(const reco::Candidate &c1, const reco::Candidate &c2, float &deltR, float &deltEta, float &deltPhi, float &deltaLocalPos, float &deltaPtRel, float &chi2) const {
132  if (!(srcCut_(c1) && matchedCut_(c2))) return false;
133  if (requireSameCharge_ && (c1.charge() != c2.charge())) return false;
134  switch (algo_) {
135  case ByTrackRef: {
138  if (t1.isNonnull()) {
139  if (t1 == t2) return true;
140  if (t1.id() != t2.id()) { edm::LogWarning("MatcherUsingTracksAlgorithm") << "Trying to match by reference tracks coming from different collections.\n"; }
141  }
142  }
143  case ByPropagatingSrc: {
146  return matchWithPropagation(start, target, deltR, deltEta, deltPhi, deltaLocalPos, deltaPtRel, chi2);
147  }
148  case ByPropagatingMatched: {
151  return matchWithPropagation(start, target, deltR, deltEta, deltPhi, deltaLocalPos, deltaPtRel, chi2);
152  }
153  case ByPropagatingSrcTSCP: {
156  return matchWithPropagation(start, target, deltR, deltEta, deltPhi, deltaLocalPos, deltaPtRel, chi2);
157  }
158  case ByPropagatingMatchedTSCP: {
161  return matchWithPropagation(start, target, deltR, deltEta, deltPhi, deltaLocalPos, deltaPtRel, chi2);
162  }
163 
164  case ByDirectComparison: {
167  return matchByDirectComparison(start, otherstart, deltR, deltEta, deltPhi, deltaLocalPos, deltaPtRel, chi2);
168  }
169  }
170  return false;
171 }
172 
176 int
177 MatcherUsingTracksAlgorithm::match(const reco::Candidate &c1, const edm::View<reco::Candidate> &c2s, float &deltR, float &deltEta, float &deltPhi, float &deltaLocalPos, float &deltaPtRel, float &chi2) const {
178  if (!srcCut_(c1)) return -1;
179 
180  // working and output variables
182  int match = -1;
183 
184  // pre-fetch some states if needed
187  } else if (algo_ == ByPropagatingMatched) target = targetState(c1, whichTrack1_, whichState1_);
188 
189  // loop on the collection
191  for (it = c2s.begin(), ed = c2s.end(), i = 0; it != ed; ++it, ++i) {
192  if (!matchedCut_(*it)) continue;
193  if (requireSameCharge_ && (c1.charge() != it->charge()) ) continue;
194  bool exit = false;
195  switch (algo_) {
196  case ByTrackRef: {
199  if (t1.isNonnull()) {
200  if (t1 == t2) { match = i; exit = true; }
201  if (t1.id() != t2.id()) { edm::LogWarning("MatcherUsingTracksAlgorithm") << "Trying to match by reference tracks coming from different collections.\n"; }
202  }
203  } break;
204  case ByPropagatingSrc:
205  case ByPropagatingMatched: {
207  else if (algo_ == ByPropagatingSrc) target = targetState(*it, whichTrack2_, whichState2_);
208  if (matchWithPropagation(start, target, deltR, deltEta, deltPhi, deltaLocalPos, deltaPtRel, chi2)) {
209  match = i;
210  }
211  } break;
212  case ByDirectComparison: {
214  if (matchByDirectComparison(start, otherstart, deltR, deltEta, deltPhi, deltaLocalPos, deltaPtRel, chi2)) {
215  match = i;
216  }
217  } break;
218  case ByPropagatingSrcTSCP: {
220  if (matchWithPropagation(start, otherstart, deltR, deltEta, deltPhi, deltaLocalPos, deltaPtRel, chi2)) {
221  match = i;
222  }
223  } break;
226  if (matchWithPropagation(otherstart, start, deltR, deltEta, deltPhi, deltaLocalPos, deltaPtRel, chi2)) {
227  match = i;
228  }
229  } break;
230  }
231  if (exit) break;
232  }
233 
234  return match;
235 }
236 
237 void
239  iSetup.get<IdealMagneticFieldRecord>().get(magfield_);
240  iSetup.get<TrackingComponentsRecord>().get("SteppingHelixPropagatorAny", propagator_);
242 }
243 
244 
247  reco::TrackRef tk;
248  const reco::RecoCandidate *rc = dynamic_cast<const reco::RecoCandidate *>(&reco);
249  if (rc == 0) throw cms::Exception("Invalid Data") << "Input object is not a RecoCandidate.\n";
250  switch (whichTrack) {
251  case TrackerTk: tk = rc->track(); break;
252  case MuonTk : tk = rc->standAloneMuon(); break;
253  case GlobalTk : tk = rc->combinedMuon(); break;
254  default: break; // just to make gcc happy
255  }
256  return tk;
257 }
258 
262  if (whichTrack != None) {
263  reco::TrackRef tk = getTrack(reco, whichTrack);
264  if (tk.isNull()) {
265  ret = FreeTrajectoryState();
266  } else {
267  switch (whichState) {
271  }
272  }
273  } else {
274  ret = FreeTrajectoryState( GlobalPoint( reco.vx(), reco.vy(), reco.vz()),
275  GlobalVector(reco.px(), reco.py(), reco.pz()),
276  reco.charge(),
277  magfield_.product());
278  }
279  return ret;
280 }
281 
285  reco::TrackRef tk = getTrack(reco, whichTrack);
286  if (tk.isNonnull()) {
287  switch (whichState) {
290  default: break; // just to make gcc happy
291  }
292  }
293  return ret;
294 }
295 
296 bool
299  float &lastDeltaR,
300  float &lastDeltaEta,
301  float &lastDeltaPhi,
302  float &lastDeltaLocalPos,
303  float &lastGlobalDPtRel,
304  float &lastChi2) const
305 {
306  if ((start.momentum().mag() == 0) || !target.isValid()) return false;
307 
308  TrajectoryStateOnSurface tsos = propagator_->propagate(start, target.surface());
309 
310  bool isBest = false;
311  if (tsos.isValid()) {
312  float thisLocalPosDiff = (tsos.localPosition()-target.localPosition()).mag();
313  float thisGlobalMomDeltaR = deltaR(tsos.globalMomentum(), target.globalMomentum());
314  float thisGlobalMomDeltaPhi = fabs(deltaPhi(tsos.globalMomentum().phi(), target.globalMomentum().phi()));
315  float thisGlobalMomDeltaEta = fabs(tsos.globalMomentum().eta() - target.globalMomentum().eta());
316  float thisGlobalDPtRel = (tsos.globalMomentum().perp() - target.globalMomentum().perp())/target.globalMomentum().perp();
317 
318  if ((thisLocalPosDiff < maxLocalPosDiff_) &&
319  (thisGlobalMomDeltaR < maxGlobalMomDeltaR_) &&
320  (thisGlobalMomDeltaEta < maxGlobalMomDeltaEta_) &&
321  (thisGlobalMomDeltaPhi < maxGlobalMomDeltaPhi_) &&
322  (fabs(thisGlobalDPtRel) < maxGlobalDPtRel_)) {
323  float thisChi2 = useChi2_ ? getChi2(target, tsos, chi2DiagonalOnly_, chi2UseVertex_) : 0;
324  if (thisChi2 >= maxChi2_) return false;
325  switch (sortBy_) {
326  case LocalPosDiff: isBest = (thisLocalPosDiff < lastDeltaLocalPos); break;
327  case GlobalMomDeltaR: isBest = (thisGlobalMomDeltaR < lastDeltaR); break;
328  case GlobalMomDeltaEta: isBest = (thisGlobalMomDeltaEta < lastDeltaEta); break;
329  case GlobalMomDeltaPhi: isBest = (thisGlobalMomDeltaPhi < lastDeltaPhi); break;
330  case GlobalDPtRel: isBest = (thisGlobalDPtRel < lastGlobalDPtRel); break;
331  case Chi2: isBest = (thisChi2 < lastChi2); break;
332  }
333  if (isBest) {
334  lastDeltaLocalPos = thisLocalPosDiff;
335  lastDeltaR = thisGlobalMomDeltaR;
336  lastDeltaEta = thisGlobalMomDeltaEta;
337  lastDeltaPhi = thisGlobalMomDeltaPhi;
338  lastGlobalDPtRel = thisGlobalDPtRel;
339  lastChi2 = thisChi2;
340  }
341  } // if match
342  }
343 
344  return isBest;
345 }
346 
347 bool
349  const FreeTrajectoryState &target,
350  float &lastDeltaR,
351  float &lastDeltaEta,
352  float &lastDeltaPhi,
353  float &lastDeltaLocalPos,
354  float &lastGlobalDPtRel,
355  float &lastChi2) const
356 {
357  if ((start.momentum().mag() == 0) || (target.momentum().mag() == 0)) return false;
359  /*2.2.X*/ try {
360  TrajectoryStateClosestToPoint tscp = propagator(start, target.position());
361  // if (!tscp.isValid()) return false; // in 3.1.X
362 
363  bool isBest = false;
364  float thisLocalPosDiff = (tscp.position()-target.position()).mag();
365  float thisGlobalMomDeltaR = deltaR(tscp.momentum(), target.momentum());
366  float thisGlobalMomDeltaPhi = fabs(deltaPhi(tscp.momentum().phi(), target.momentum().phi()));
367  float thisGlobalMomDeltaEta = fabs(tscp.momentum().eta() - target.momentum().eta());
368  float thisGlobalDPtRel = (tscp.momentum().perp() - target.momentum().perp())/target.momentum().perp();
369 
370  if ((thisLocalPosDiff < maxLocalPosDiff_) &&
371  (thisGlobalMomDeltaR < maxGlobalMomDeltaR_) &&
372  (thisGlobalMomDeltaEta < maxGlobalMomDeltaEta_) &&
373  (thisGlobalMomDeltaPhi < maxGlobalMomDeltaPhi_) &&
374  (fabs(thisGlobalDPtRel) < maxGlobalDPtRel_)) {
375  float thisChi2 = useChi2_ ? getChi2(target, tscp, chi2DiagonalOnly_, chi2UseVertex_) : 0;
376  if (thisChi2 >= maxChi2_) return false;
377  switch (sortBy_) {
378  case LocalPosDiff: isBest = (thisLocalPosDiff < lastDeltaLocalPos); break;
379  case GlobalMomDeltaR: isBest = (thisGlobalMomDeltaR < lastDeltaR); break;
380  case GlobalMomDeltaEta: isBest = (thisGlobalMomDeltaEta < lastDeltaEta); break;
381  case GlobalMomDeltaPhi: isBest = (thisGlobalMomDeltaPhi < lastDeltaPhi); break;
382  case GlobalDPtRel: isBest = (thisGlobalDPtRel < lastGlobalDPtRel); break;
383  case Chi2: isBest = (thisChi2 < lastChi2); break;
384  }
385  if (isBest) {
386  lastDeltaLocalPos = thisLocalPosDiff;
387  lastDeltaR = thisGlobalMomDeltaR;
388  lastDeltaEta = thisGlobalMomDeltaEta;
389  lastDeltaPhi = thisGlobalMomDeltaPhi;
390  lastGlobalDPtRel = thisGlobalDPtRel;
391  lastChi2 = thisChi2;
392  }
393  } // if match
394 
395  return isBest;
396  /*2.2.X*/ } catch (const TrajectoryStateException &err) { return false; }
397 }
398 
399 
400 bool
402  const FreeTrajectoryState &target,
403  float &lastDeltaR,
404  float &lastDeltaEta,
405  float &lastDeltaPhi,
406  float &lastDeltaLocalPos,
407  float &lastGlobalDPtRel,
408  float &lastChi2) const
409 {
410  if ((start.momentum().mag() == 0) || target.momentum().mag() == 0) return false;
411 
412  bool isBest = false;
413  float thisLocalPosDiff = (start.position()-target.position()).mag();
414  float thisGlobalMomDeltaR = deltaR(start.momentum(), target.momentum());
415  float thisGlobalMomDeltaPhi = fabs(deltaPhi(start.momentum().phi(), target.momentum().phi()));
416  float thisGlobalMomDeltaEta = fabs(start.momentum().eta() - target.momentum().eta());
417  float thisGlobalDPtRel = (start.momentum().perp() - target.momentum().perp())/target.momentum().perp();
418 
419  if ((thisLocalPosDiff < maxLocalPosDiff_) &&
420  (thisGlobalMomDeltaR < maxGlobalMomDeltaR_) &&
421  (thisGlobalMomDeltaEta < maxGlobalMomDeltaEta_) &&
422  (thisGlobalMomDeltaPhi < maxGlobalMomDeltaPhi_) &&
423  (fabs(thisGlobalDPtRel) < maxGlobalDPtRel_)) {
424  float thisChi2 = useChi2_ ? getChi2(start, target, chi2DiagonalOnly_, chi2UseVertex_, chi2FirstMomentum_) : 0;
425  if (thisChi2 >= maxChi2_) return false;
426  switch (sortBy_) {
427  case LocalPosDiff: isBest = (thisLocalPosDiff < lastDeltaLocalPos); break;
428  case GlobalMomDeltaR: isBest = (thisGlobalMomDeltaR < lastDeltaR); break;
429  case GlobalMomDeltaEta: isBest = (thisGlobalMomDeltaEta < lastDeltaEta); break;
430  case GlobalMomDeltaPhi: isBest = (thisGlobalMomDeltaPhi < lastDeltaPhi); break;
431  case GlobalDPtRel: isBest = (thisGlobalDPtRel < lastGlobalDPtRel); break;
432  case Chi2: isBest = (thisChi2 < lastChi2); break;
433  }
434  if (isBest) {
435  lastDeltaLocalPos = thisLocalPosDiff;
436  lastDeltaR = thisGlobalMomDeltaR;
437  lastDeltaEta = thisGlobalMomDeltaEta;
438  lastDeltaPhi = thisGlobalMomDeltaPhi;
439  lastGlobalDPtRel = thisGlobalDPtRel;
440  lastChi2 = thisChi2;
441  }
442  } // if match
443 
444  return isBest;
445 }
446 
447 double
448 MatcherUsingTracksAlgorithm::getChi2(const FreeTrajectoryState &start, const FreeTrajectoryState &other, bool diagonalOnly, bool useVertex, bool useFirstMomentum) {
449  if (!start.hasError() && !other.hasError()) throw cms::Exception("LogicError") << "At least one of the two states must have errors to make chi2s.\n";
451  if (start.hasError()) cov += start.curvilinearError().matrix();
452  if (other.hasError()) cov += other.curvilinearError().matrix();
453  cropAndInvert(cov, diagonalOnly, !useVertex);
454  GlobalVector p1 = start.momentum(), p2 = other.momentum();
455  GlobalPoint x1 = start.position(), x2 = other.position();
456  GlobalVector p = useFirstMomentum ? p1 : p2; double pt = p.perp(), pm = p.mag();
457  double dsz = (x1.z()-x2.z())*pt/pm - ( (x1.x()-x2.x())*p.x() + (x1.y()-x2.y())*p.y() ) / pt * p.z()/pm;
458  double dxy = ( -(x1.x()-x2.x())*p.y() + (x1.y()-x2.y())*p.x() )/pt;
459  AlgebraicVector5 diff(start.charge()/p1.mag()-other.charge()/p2.mag(),
460  p1.theta()-p2.theta(),
461  (p1.phi()-p2.phi()).value(),
462  dxy,
463  dsz);
464  return ROOT::Math::Similarity(diff, cov);
465 }
466 
467 double
468 MatcherUsingTracksAlgorithm::getChi2(const FreeTrajectoryState &start, const TrajectoryStateClosestToPoint &other, bool diagonalOnly, bool useVertex) {
469  if (!start.hasError() && !other.hasError()) throw cms::Exception("LogicError") << "At least one of the two states must have errors to make chi2s.\n";
470  PerigeeConversions pgconvert; double pt; // needed by pgconvert
472  if (start.hasError()) cov += pgconvert.ftsToPerigeeError(start).covarianceMatrix();
473  if (other.hasError()) cov += other.perigeeError().covarianceMatrix();
474  cropAndInvert(cov, diagonalOnly, !useVertex);
475  AlgebraicVector5 pgpar1 = pgconvert.ftsToPerigeeParameters(start,other.referencePoint(),pt).vector();
476  AlgebraicVector5 pgpar2 = other.perigeeParameters().vector();
477  AlgebraicVector5 diff(pgpar1-pgpar2);
478  return ROOT::Math::Similarity(diff, cov);
479 }
480 
481 double
482 MatcherUsingTracksAlgorithm::getChi2(const TrajectoryStateOnSurface &start, const TrajectoryStateOnSurface &other, bool diagonalOnly, bool usePosition) {
483  if (!start.hasError() && !other.hasError()) throw cms::Exception("LogicError") << "At least one of the two states must have errors to make chi2s.\n";
485  if (start.hasError()) cov += start.localError().matrix();
486  if (other.hasError()) cov += other.localError().matrix();
487  cropAndInvert(cov, diagonalOnly, !usePosition);
489  return ROOT::Math::Similarity(diff, cov);
490 }
491 
492 void
493 MatcherUsingTracksAlgorithm::cropAndInvert(AlgebraicSymMatrix55 &cov, bool diagonalOnly, bool top3by3only) {
494  if (!top3by3only) {
495  if (diagonalOnly) {
496  for (size_t i = 0; i < 5; ++i) { for (size_t j = i+1; j < 5; ++j) {
497  cov(i,j) = 0;
498  } }
499  }
500  cov.Invert();
501  } else {
502  // get 3x3 covariance
503  AlgebraicSymMatrix33 momCov = cov.Sub<AlgebraicSymMatrix33>(0,0); // get 3x3 matrix
504  if (diagonalOnly) { momCov(0,1) = 0; momCov(0,2) = 0; momCov(1,2) = 0; }
505  // invert
506  momCov.Invert();
507  // place it
508  cov.Place_at(momCov,0,0);
509  // zero the rest
510  for (size_t i = 3; i < 5; ++i) { for (size_t j = i; j < 5; ++j) {
511  cov(i,j) = 0;
512  } }
513  }
514 }
515 
516 
T getParameter(std::string const &) const
const AlgebraicVector5 & vector() const
int i
Definition: DBlmapReader.cc:9
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:81
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:187
T perp() const
Definition: PV3DBase.h:71
const PerigeeTrajectoryError & perigeeError() const
TrajectoryStateOnSurface innerStateOnSurface(const reco::Track &tk, const TrackingGeometry &geom, const MagneticField *field)
const LocalTrajectoryParameters & localParameters() const
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
void getConf(const edm::ParameterSet &iConfig, const std::string &whatFor, WhichTrack &whichTrack, WhichState &whichState)
Parse some configuration.
FreeTrajectoryState innerFreeState(const reco::Track &tk, const MagneticField *field)
PerigeeTrajectoryParameters ftsToPerigeeParameters(const FTS &originalFTS, const GlobalPoint &referencePoint, double &pt) const
virtual double pz() const =0
z coordinate of momentum vector
Geom::Phi< T > phi() const
Definition: PV3DBase.h:68
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
virtual double vx() const =0
x coordinate of vertex position
T y() const
Definition: PV3DBase.h:62
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
TrackCharge charge() const
double deltaR2< GlobalVector >(const GlobalVector &v1, const GlobalVector &v2)
virtual reco::TrackRef standAloneMuon() const
reference to a stand-alone muon Track
static double getChi2(const FreeTrajectoryState &start, const FreeTrajectoryState &other, bool diagonalOnly, bool useVertex, bool useFirstMomentum)
const CurvilinearTrajectoryError & curvilinearError() const
bool matchWithPropagation(const FreeTrajectoryState &start, const FreeTrajectoryState &target, float &lastDeltaR, float &lastDeltaEta, float &lastDeltaPhi, float &lastDeltaLocalPos, float &lastGlobalDPtRel, float &lastChi2) const
virtual double vy() const =0
y coordinate of vertex position
virtual reco::TrackRef track() const
reference to a Track
Geom::Theta< T > theta() const
Definition: PV3DBase.h:74
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:250
bool isNull() const
Checks for null.
Definition: Ref.h:247
T mag() const
Definition: PV3DBase.h:66
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33
const PerigeeTrajectoryParameters & perigeeParameters() const
edm::ESHandle< GlobalTrackingGeometry > geometry_
T z() const
Definition: PV3DBase.h:63
edm::ESHandle< MagneticField > magfield_
int j
Definition: DBlmapReader.cc:9
const AlgebraicSymMatrix55 & matrix() const
virtual int charge() const =0
electric charge
MatcherUsingTracksAlgorithm(const edm::ParameterSet &iConfig)
const LocalTrajectoryError & localError() const
virtual double py() const =0
y coordinate of momentum vector
double p2[4]
Definition: TauolaWrapper.h:90
GlobalVector momentum() const
edm::ESHandle< Propagator > propagator_
const AlgebraicSymMatrix55 & covarianceMatrix() const
virtual double px() const =0
x coordinate of momentum vector
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< float >, ROOT::Math::GlobalCoordinateSystemTag > GlobalVector
vector in glovbal coordinate system
Definition: Vector3D.h:28
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
TrajectoryStateOnSurface targetState(const reco::Candidate &reco, WhichTrack whichTrack, WhichState whichState) const
End state for the propagation.
GlobalPoint position() const
StringCutObjectSelector< reco::Candidate, true > srcCut_
ROOT::Math::SVector< double, 5 > AlgebraicVector5
FreeTrajectoryState initialFreeState(const reco::Track &tk, const MagneticField *field)
const GlobalPoint & referencePoint() const
FreeTrajectoryState startingState(const reco::Candidate &reco, WhichTrack whichTrack, WhichState whichState) const
Starting state for the propagation.
const T & get() const
Definition: EventSetup.h:55
TrajectoryStateOnSurface outerStateOnSurface(const reco::Track &tk, const TrackingGeometry &geom, const MagneticField *field)
T const * product() const
Definition: ESHandle.h:62
static void cropAndInvert(AlgebraicSymMatrix55 &cov, bool diagonalOnly, bool top3by3only)
Possibly crop the 3x3 part of the matrix or remove off-diagonal terms, then invert.
AlgebraicVector5 mixedFormatVector() const
T eta() const
Definition: PV3DBase.h:75
virtual double vz() const =0
z coordinate of vertex position
double p1[4]
Definition: TauolaWrapper.h:89
const AlgebraicSymMatrix55 & matrix() const
GlobalVector globalMomentum() const
const Surface & surface() const
PerigeeTrajectoryError ftsToPerigeeError(const FTS &originalFTS) const
bool match(const reco::Candidate &c1, const reco::Candidate &c2, float &deltaR, float &deltaEta, float &deltaPhi, float &deltaLocalPos, float &deltaPtRel, float &chi2) const
void init(const edm::EventSetup &iSetup)
Call this method at the beginning of each run, to initialize geometry, magnetic field and propagators...
const_iterator begin() const
const_iterator end() const
FreeTrajectoryState outerFreeState(const reco::Track &tk, const MagneticField *field)
ProductID id() const
Accessor for product ID.
Definition: Ref.h:256
reco::TrackRef getTrack(const reco::Candidate &reco, WhichTrack which) const
Get track reference out of a Candidate (via dynamic_cast to reco::RecoCandidate)
double deltaR< GlobalVector >(const GlobalVector &v1, const GlobalVector &v2)
StringCutObjectSelector< reco::Candidate, true > matchedCut_
LimitAlgo * algo
Definition: Combine.cc:60
Definition: Chi2.h:17
bool matchByDirectComparison(const FreeTrajectoryState &start, const FreeTrajectoryState &other, float &lastDeltaR, float &lastDeltaEta, float &lastDeltaPhi, float &lastDeltaLocalPos, float &lastGlobalDPtRel, float &lastChi2) const
Compare directly two states. return true if current pair is the new best match (in that case...
T x() const
Definition: PV3DBase.h:61
Global3DVector GlobalVector
Definition: GlobalVector.h:10
virtual reco::TrackRef combinedMuon() const
reference to a stand-alone muon Track