CMS 3D CMS Logo

DuplicateTrackMerger.cc
Go to the documentation of this file.
1 
9 
14 
18 
20 
24 #include "TrackMerger.h"
25 
27 #include <vector>
28 #include <algorithm>
29 #include <string>
30 #include <iostream>
31 #include <atomic>
32 
34 
35 using namespace reco;
36 namespace {
37  class DuplicateTrackMerger final : public edm::stream::EDProducer<> {
38  public:
40  explicit DuplicateTrackMerger(const edm::ParameterSet &iPara);
42  ~DuplicateTrackMerger() override;
43 
44  using CandidateToDuplicate = std::vector<std::pair<int, int>>;
45 
46  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
47 
48  private:
50  void produce(edm::Event &, const edm::EventSetup &) override;
51 
52  bool checkForDisjointTracks(const reco::Track *t1, const reco::Track *t2, TSCPBuilderNoMaterial &tscpBuilder) const;
53  bool checkForOverlappingTracks(
54  const reco::Track *t1, const reco::Track *t2, unsigned int nvh1, unsigned int nvh2, double cosT) const;
55 
56  private:
58  const GBRForest *forest_;
59 
61  std::string dbFileName_;
62  bool useForestFromDB_;
63  std::string forestLabel_;
64 
65  std::string propagatorName_;
66  std::string chi2EstimatorName_;
67 
71  double minDeltaR3d2_;
73  double minBDTG_;
75  double minpT2_;
77  double minP_;
79  float maxDCA2_;
81  float maxDPhi_;
83  float maxDLambda_;
85  float maxDdxy_;
87  float maxDdsz_;
89  float maxDQoP_;
91  unsigned int overlapCheckMaxHits_;
93  unsigned int overlapCheckMaxMissingLayers_;
95  double overlapCheckMinCosT_;
96 
97  const MagneticField *magfield_;
98  const TrackerTopology *ttopo_;
99  const TrackerGeometry *geom_;
100  const Propagator *propagator_;
101  const Chi2MeasurementEstimatorBase *chi2Estimator_;
102 
104  TrackMerger merger_;
105 
106 #ifdef EDM_ML_DEBUG
107  bool debug_;
108 #endif
109  };
110 } // namespace
111 
120 #include "TFile.h"
121 
124 
125 #include "DuplicateTrackType.h"
126 
127 namespace {
128 
131  desc.add<edm::InputTag>("source", edm::InputTag());
132  desc.add<double>("minDeltaR3d", -4.0);
133  desc.add<double>("minBDTG", -0.1);
134  desc.add<double>("minpT", 0.2);
135  desc.add<double>("minP", 0.4);
136  desc.add<double>("maxDCA", 30.0);
137  desc.add<double>("maxDPhi", 0.30);
138  desc.add<double>("maxDLambda", 0.30);
139  desc.add<double>("maxDdsz", 10.0);
140  desc.add<double>("maxDdxy", 10.0);
141  desc.add<double>("maxDQoP", 0.25);
142  desc.add<unsigned>("overlapCheckMaxHits", 4);
143  desc.add<unsigned>("overlapCheckMaxMissingLayers", 1);
144  desc.add<double>("overlapCheckMinCosT", 0.99);
145  desc.add<std::string>("forestLabel", "MVADuplicate");
146  desc.add<std::string>("GBRForestFileName", "");
147  desc.add<bool>("useInnermostState", true);
148  desc.add<std::string>("ttrhBuilderName", "WithAngleAndTemplate");
149  desc.add<std::string>("propagatorName", "PropagatorWithMaterial");
150  desc.add<std::string>("chi2EstimatorName", "DuplicateTrackMergerChi2Est");
151  descriptions.add("DuplicateTrackMerger", desc);
152  }
153 
154  DuplicateTrackMerger::DuplicateTrackMerger(const edm::ParameterSet &iPara) : forest_(nullptr), merger_(iPara) {
155  trackSource_ = consumes<reco::TrackCollection>(iPara.getParameter<edm::InputTag>("source"));
156  minDeltaR3d2_ = iPara.getParameter<double>("minDeltaR3d");
157  minDeltaR3d2_ *= std::abs(minDeltaR3d2_);
158  minBDTG_ = iPara.getParameter<double>("minBDTG");
159  minpT2_ = iPara.getParameter<double>("minpT");
160  minpT2_ *= minpT2_;
161  minP_ = iPara.getParameter<double>("minP");
162  maxDCA2_ = iPara.getParameter<double>("maxDCA");
163  maxDCA2_ *= maxDCA2_;
164  maxDPhi_ = iPara.getParameter<double>("maxDPhi");
165  maxDLambda_ = iPara.getParameter<double>("maxDLambda");
166  maxDdsz_ = iPara.getParameter<double>("maxDdsz");
167  maxDdxy_ = iPara.getParameter<double>("maxDdxy");
168  maxDQoP_ = iPara.getParameter<double>("maxDQoP");
169  overlapCheckMaxHits_ = iPara.getParameter<unsigned>("overlapCheckMaxHits");
170  overlapCheckMaxMissingLayers_ = iPara.getParameter<unsigned>("overlapCheckMaxMissingLayers");
171  overlapCheckMinCosT_ = iPara.getParameter<double>("overlapCheckMinCosT");
172 
173  produces<std::vector<TrackCandidate>>("candidates");
174  produces<CandidateToDuplicate>("candidateMap");
175 
176  forestLabel_ = iPara.getParameter<std::string>("forestLabel");
177 
178  dbFileName_ = iPara.getParameter<std::string>("GBRForestFileName");
179  useForestFromDB_ = dbFileName_.empty();
180 
181  propagatorName_ = iPara.getParameter<std::string>("propagatorName");
182  chi2EstimatorName_ = iPara.getParameter<std::string>("chi2EstimatorName");
183 
184  /*
185  tmvaReader_ = new TMVA::Reader("!Color:Silent");
186  tmvaReader_->AddVariable("ddsz",&tmva_ddsz_);
187  tmvaReader_->AddVariable("ddxy",&tmva_ddxy_);
188  tmvaReader_->AddVariable("dphi",&tmva_dphi_);
189  tmvaReader_->AddVariable("dlambda",&tmva_dlambda_);
190  tmvaReader_->AddVariable("dqoverp",&tmva_dqoverp_);
191  tmvaReader_->AddVariable("delta3d_r",&tmva_d3dr_);
192  tmvaReader_->AddVariable("delta3d_z",&tmva_d3dz_);
193  tmvaReader_->AddVariable("outer_nMissingInner",&tmva_outer_nMissingInner_);
194  tmvaReader_->AddVariable("inner_nMissingOuter",&tmva_inner_nMissingOuter_);
195  tmvaReader_->BookMVA("BDTG",mvaFilePath);
196  */
197  }
198 
199  DuplicateTrackMerger::~DuplicateTrackMerger() {
200  if (!useForestFromDB_)
201  delete forest_;
202  }
203 
204 #ifdef VI_STAT
205  struct Stat {
206  Stat() : maxCos(1.1), nCand(0), nLoop0(0) {}
207  ~Stat() { std::cout << "Stats " << nCand << ' ' << nLoop0 << ' ' << maxCos << std::endl; }
208  std::atomic<float> maxCos;
209  std::atomic<int> nCand, nLoop0;
210  };
211  Stat stat;
212 #endif
213 
214  template <typename T>
215  void update_maximum(std::atomic<T> &maximum_value, T const &value) noexcept {
216  T prev_value = maximum_value;
217  while (prev_value < value && !maximum_value.compare_exchange_weak(prev_value, value))
218  ;
219  }
220 
221  template <typename T>
222  void update_minimum(std::atomic<T> &minimum_value, T const &value) noexcept {
223  T prev_value = minimum_value;
224  while (prev_value > value && !minimum_value.compare_exchange_weak(prev_value, value))
225  ;
226  }
227 
228  void DuplicateTrackMerger::produce(edm::Event &iEvent, const edm::EventSetup &iSetup) {
229  merger_.init(iSetup);
230 
231  if (!forest_) {
232  if (useForestFromDB_) {
233  edm::ESHandle<GBRForest> forestHandle;
234  iSetup.get<GBRWrapperRcd>().get(forestLabel_, forestHandle);
235  forest_ = forestHandle.product();
236  } else {
237  TFile gbrfile(dbFileName_.c_str());
238  forest_ = dynamic_cast<const GBRForest *>(gbrfile.Get(forestLabel_.c_str()));
239  }
240  }
241 
242  //edm::Handle<edm::View<reco::Track> >handle;
244  iEvent.getByToken(trackSource_, handle);
245  auto const &tracks = *handle;
246 
248  iSetup.get<IdealMagneticFieldRecord>().get(hmagfield);
249  magfield_ = hmagfield.product();
250 
252  iSetup.get<TrackerTopologyRcd>().get(httopo);
253  ttopo_ = httopo.product();
254 
256  iSetup.get<TrackerDigiGeometryRecord>().get(hgeom);
257  geom_ = hgeom.product();
258 
259  edm::ESHandle<Propagator> hpropagator;
260  iSetup.get<TrackingComponentsRecord>().get(propagatorName_, hpropagator);
261  propagator_ = hpropagator.product();
262 
264  iSetup.get<TrackingComponentsRecord>().get(chi2EstimatorName_, hestimator);
265  chi2Estimator_ = hestimator.product();
266 
267  TSCPBuilderNoMaterial tscpBuilder;
268  auto out_duplicateCandidates = std::make_unique<std::vector<TrackCandidate>>();
269 
270  auto out_candidateMap = std::make_unique<CandidateToDuplicate>();
271  LogDebug("DuplicateTrackMerger") << "Number of tracks to be checked for merging: " << tracks.size();
272 
273 #ifdef EDM_ML_DEBUG
274  auto test = [&](const reco::Track *a, const reco::Track *b) {
275  const auto ev = iEvent.id().event();
276  const auto aOriAlgo = a->originalAlgo();
277  const auto bOriAlgo = b->originalAlgo();
278  const auto aSeed = a->seedRef().key();
279  const auto bSeed = b->seedRef().key();
280  return ((ev == 6903 && ((aOriAlgo == 23 && aSeed == 695 && bOriAlgo == 5 && bSeed == 652) ||
281  (aOriAlgo == 23 && aSeed == 400 && bOriAlgo == 7 && bSeed == 156) ||
282  (aOriAlgo == 4 && aSeed == 914 && bOriAlgo == 22 && bSeed == 503) ||
283  (aOriAlgo == 5 && aSeed == 809 && bOriAlgo == 4 && bSeed == 1030) ||
284  (aOriAlgo == 23 && aSeed == 749 && bOriAlgo == 5 && bSeed == 659) ||
285  (aOriAlgo == 4 && aSeed == 1053 && bOriAlgo == 23 && bSeed == 1035) ||
286  (aOriAlgo == 4 && aSeed == 810 && bOriAlgo == 5 && bSeed == 666) ||
287  (aOriAlgo == 4 && aSeed == 974 && bOriAlgo == 5 && bSeed == 778))) ||
288  (ev == 6904 && ((aOriAlgo == 23 && aSeed == 526 && bOriAlgo == 5 && bSeed == 307) ||
289  (aOriAlgo == 4 && aSeed == 559 && bOriAlgo == 22 && bSeed == 472))) ||
290  (ev == 6902 && ((aOriAlgo == 4 && aSeed == 750 && bOriAlgo == 22 && bSeed == 340) ||
291  (aOriAlgo == 4 && aSeed == 906 && bOriAlgo == 5 && bSeed == 609) ||
292  (aOriAlgo == 4 && aSeed == 724 && bOriAlgo == 5 && bSeed == 528) ||
293  (aOriAlgo == 4 && aSeed == 943 && bOriAlgo == 23 && bSeed == 739) ||
294  (aOriAlgo == 8 && aSeed == 2 && bOriAlgo == 9 && bSeed == 2282) ||
295  (aOriAlgo == 23 && aSeed == 827 && bOriAlgo == 5 && bSeed == 656) ||
296  (aOriAlgo == 22 && aSeed == 667 && bOriAlgo == 7 && bSeed == 516))));
297  };
298 #endif
299 
300  // cache few "heavy to compute quantities
301  int nTracks = 0;
302  declareDynArray(const reco::Track *, tracks.size(), selTracks);
303  declareDynArray(unsigned int, tracks.size(), nValidHits);
304  declareDynArray(unsigned int, tracks.size(), oriIndex);
305  for (auto i = 0U; i < tracks.size(); i++) {
306  const reco::Track *rt1 = &tracks[i];
307  if (rt1->innerMomentum().perp2() < minpT2_)
308  continue;
309  selTracks[nTracks] = rt1;
310  nValidHits[nTracks] = rt1->numberOfValidHits(); // yes it is extremely heavy!
311  oriIndex[nTracks] = i;
312  ++nTracks;
313  }
314 
315  for (int i = 0; i < nTracks; i++) {
316  const reco::Track *rt1 = selTracks[i];
317  for (int j = i + 1; j < nTracks; j++) {
318  const reco::Track *rt2 = selTracks[j];
319 
320 #ifdef EDM_ML_DEBUG
321  debug_ = false;
322  if (test(rt1, rt2) || test(rt2, rt1)) {
323  debug_ = true;
324  LogTrace("DuplicateTrackMerger")
325  << "Track1 " << i << " originalAlgo " << rt1->originalAlgo() << " seed " << rt1->seedRef().key() << " pT "
326  << std::sqrt(rt1->innerMomentum().perp2()) << " charge " << rt1->charge() << " outerPosition2 "
327  << rt1->outerPosition().perp2() << "\n"
328  << "Track2 " << j << " originalAlgo " << rt2->originalAlgo() << " seed " << rt2->seedRef().key() << " pT "
329  << std::sqrt(rt2->innerMomentum().perp2()) << " charge " << rt2->charge() << " outerPosition2 "
330  << rt2->outerPosition().perp2();
331  }
332 #endif
333 
334  if (rt1->charge() != rt2->charge())
335  continue;
336  auto cosT = (*rt1).momentum().Dot((*rt2).momentum()); // not normalized!
337  IfLogTrace(debug_, "DuplicateTrackMerger") << " cosT " << cosT;
338  if (cosT < 0.)
339  continue;
340  cosT /= std::sqrt((*rt1).momentum().Mag2() * (*rt2).momentum().Mag2());
341 
342  const reco::Track *t1, *t2;
343  unsigned int nhv1, nhv2;
344  if (rt1->outerPosition().perp2() < rt2->outerPosition().perp2()) {
345  t1 = rt1;
346  nhv1 = nValidHits[i];
347  t2 = rt2;
348  nhv2 = nValidHits[j];
349  } else {
350  t1 = rt2;
351  nhv1 = nValidHits[j];
352  t2 = rt1;
353  nhv2 = nValidHits[i];
354  }
355  auto deltaR3d2 = (t1->outerPosition() - t2->innerPosition()).mag2();
356 
357  if (t1->outerPosition().perp2() > t2->innerPosition().perp2())
358  deltaR3d2 *= -1.0;
359  IfLogTrace(debug_, "DuplicateTrackMerger")
360  << " deltaR3d2 " << deltaR3d2 << " t1.outerPos2 " << t1->outerPosition().perp2() << " t2.innerPos2 "
361  << t2->innerPosition().perp2();
362 
363  bool compatible = false;
364  DuplicateTrackType duplType;
365  if (deltaR3d2 >= minDeltaR3d2_) {
366  compatible = checkForDisjointTracks(t1, t2, tscpBuilder);
367  duplType = DuplicateTrackType::Disjoint;
368  } else {
369  compatible = checkForOverlappingTracks(t1, t2, nhv1, nhv2, cosT);
371  }
372  if (!compatible)
373  continue;
374 
375  IfLogTrace(debug_, "DuplicateTrackMerger") << " marking as duplicates" << oriIndex[i] << ',' << oriIndex[j];
376  out_duplicateCandidates->push_back(merger_.merge(*t1, *t2, duplType));
377  out_candidateMap->emplace_back(oriIndex[i], oriIndex[j]);
378 
379 #ifdef VI_STAT
380  ++stat.nCand;
381  // auto cosT = float((*t1).momentum().unit().Dot((*t2).momentum().unit()));
382  if (cosT > 0)
383  update_minimum(stat.maxCos, float(cosT));
384  else
385  ++stat.nLoop0;
386 #endif
387  }
388  }
389  iEvent.put(std::move(out_duplicateCandidates), "candidates");
390  iEvent.put(std::move(out_candidateMap), "candidateMap");
391  }
392 
393  bool DuplicateTrackMerger::checkForDisjointTracks(const reco::Track *t1,
394  const reco::Track *t2,
395  TSCPBuilderNoMaterial &tscpBuilder) const {
396  IfLogTrace(debug_, "DuplicateTrackMerger") << " Checking for disjoint duplicates";
397 
400  GlobalPoint avgPoint((t1->outerPosition().x() + t2->innerPosition().x()) * 0.5,
401  (t1->outerPosition().y() + t2->innerPosition().y()) * 0.5,
402  (t1->outerPosition().z() + t2->innerPosition().z()) * 0.5);
403  TrajectoryStateClosestToPoint TSCP1 = tscpBuilder(fts1, avgPoint);
404  TrajectoryStateClosestToPoint TSCP2 = tscpBuilder(fts2, avgPoint);
405  IfLogTrace(debug_, "DuplicateTrackMerger")
406  << " TSCP1.isValid " << TSCP1.isValid() << " TSCP2.isValid " << TSCP2.isValid();
407  if (!TSCP1.isValid())
408  return false;
409  if (!TSCP2.isValid())
410  return false;
411 
412  const FreeTrajectoryState &ftsn1 = TSCP1.theState();
413  const FreeTrajectoryState &ftsn2 = TSCP2.theState();
414 
415  IfLogTrace(debug_, "DuplicateTrackMerger") << " DCA2 " << (ftsn2.position() - ftsn1.position()).mag2();
416  if ((ftsn2.position() - ftsn1.position()).mag2() > maxDCA2_)
417  return false;
418 
419  auto qoverp1 = ftsn1.signedInverseMomentum();
420  auto qoverp2 = ftsn2.signedInverseMomentum();
421  float tmva_dqoverp_ = qoverp1 - qoverp2;
422  IfLogTrace(debug_, "DuplicateTrackMerger") << " dqoverp " << tmva_dqoverp_;
423  if (std::abs(tmva_dqoverp_) > maxDQoP_)
424  return false;
425 
426  //auto pp = [&](TrajectoryStateClosestToPoint const & ts) { std::cout << ' ' << ts.perigeeParameters().vector()[0] << '/' << ts.perigeeError().transverseCurvatureError();};
427  //if(qoverp1*qoverp2 <0) { std::cout << "charge different " << qoverp1 <<',' << qoverp2; pp(TSCP1); pp(TSCP2); std::cout << std::endl;}
428 
429  // lambda = pi/2 - theta.... so l1-l2 == t2-t1
430  float tmva_dlambda_ = ftsn2.momentum().theta() - ftsn1.momentum().theta();
431  IfLogTrace(debug_, "DuplicateTrackMerger") << " dlambda " << tmva_dlambda_;
432  if (std::abs(tmva_dlambda_) > maxDLambda_)
433  return false;
434 
435  auto phi1 = ftsn1.momentum().phi();
436  auto phi2 = ftsn2.momentum().phi();
437  float tmva_dphi_ = phi1 - phi2;
438  if (std::abs(tmva_dphi_) > float(M_PI))
439  tmva_dphi_ = 2.f * float(M_PI) - std::abs(tmva_dphi_);
440  IfLogTrace(debug_, "DuplicateTrackMerger") << " dphi " << tmva_dphi_;
441  if (std::abs(tmva_dphi_) > maxDPhi_)
442  return false;
443 
444  auto dxy1 =
445  (-ftsn1.position().x() * ftsn1.momentum().y() + ftsn1.position().y() * ftsn1.momentum().x()) / TSCP1.pt();
446  auto dxy2 =
447  (-ftsn2.position().x() * ftsn2.momentum().y() + ftsn2.position().y() * ftsn2.momentum().x()) / TSCP2.pt();
448  float tmva_ddxy_ = dxy1 - dxy2;
449  IfLogTrace(debug_, "DuplicateTrackMerger") << " ddxy " << tmva_ddxy_;
450  if (std::abs(tmva_ddxy_) > maxDdxy_)
451  return false;
452 
453  auto dsz1 = ftsn1.position().z() * TSCP1.pt() / TSCP1.momentum().mag() -
454  (ftsn1.position().x() * ftsn1.momentum().y() + ftsn1.position().y() * ftsn1.momentum().x()) /
455  TSCP1.pt() * ftsn1.momentum().z() / ftsn1.momentum().mag();
456  auto dsz2 = ftsn2.position().z() * TSCP2.pt() / TSCP2.momentum().mag() -
457  (ftsn2.position().x() * ftsn2.momentum().y() + ftsn2.position().y() * ftsn2.momentum().x()) /
458  TSCP2.pt() * ftsn2.momentum().z() / ftsn2.momentum().mag();
459  float tmva_ddsz_ = dsz1 - dsz2;
460  IfLogTrace(debug_, "DuplicateTrackMerger") << " ddsz " << tmva_ddsz_;
461  if (std::abs(tmva_ddsz_) > maxDdsz_)
462  return false;
463 
464  float tmva_d3dr_ = avgPoint.perp();
465  float tmva_d3dz_ = avgPoint.z();
466  float tmva_outer_nMissingInner_ = t2->hitPattern().numberOfLostHits(reco::HitPattern::MISSING_INNER_HITS);
467  float tmva_inner_nMissingOuter_ = t1->hitPattern().numberOfLostHits(reco::HitPattern::MISSING_OUTER_HITS);
468 
469  float gbrVals_[9];
470  gbrVals_[0] = tmva_ddsz_;
471  gbrVals_[1] = tmva_ddxy_;
472  gbrVals_[2] = tmva_dphi_;
473  gbrVals_[3] = tmva_dlambda_;
474  gbrVals_[4] = tmva_dqoverp_;
475  gbrVals_[5] = tmva_d3dr_;
476  gbrVals_[6] = tmva_d3dz_;
477  gbrVals_[7] = tmva_outer_nMissingInner_;
478  gbrVals_[8] = tmva_inner_nMissingOuter_;
479 
480  auto mvaBDTG = forest_->GetClassifier(gbrVals_);
481  IfLogTrace(debug_, "DuplicateTrackMerger") << " mvaBDTG " << mvaBDTG;
482  if (mvaBDTG < minBDTG_)
483  return false;
484 
485  // std::cout << "to merge " << mvaBDTG << ' ' << std::copysign(std::sqrt(std::abs(deltaR3d2)),deltaR3d2) << ' ' << tmva_dphi_ << ' ' << TSCP1.pt() <<'/'<<TSCP2.pt() << std::endl;
486  return true;
487  }
488 
489  bool DuplicateTrackMerger::checkForOverlappingTracks(
490  const reco::Track *t1, const reco::Track *t2, unsigned int nvh1, unsigned int nvh2, double cosT) const {
491  // ensure t1 is the shorter track
492  if (nvh2 < nvh1) {
493  std::swap(t1, t2);
494  std::swap(nvh1, nvh2);
495  }
496 
497  IfLogTrace(debug_, "DuplicateTrackMerger")
498  << " Checking for overlapping duplicates, cosT " << cosT << " t1 hits " << nvh1;
499  if (cosT < overlapCheckMinCosT_)
500  return false;
501  if (nvh1 > overlapCheckMaxHits_)
502  return false;
503 
504  // find the hit on the longer track on layer of the first hit of the shorter track
505  auto findHitOnT2 = [&](const TrackingRecHit *hit1) {
506  const auto hitDet = hit1->geographicalId().det();
507  const auto hitSubdet = hit1->geographicalId().subdetId();
508  const auto hitLayer = ttopo_->layer(hit1->geographicalId());
509  return std::find_if(t2->recHitsBegin(), t2->recHitsEnd(), [&](const TrackingRecHit *hit2) {
510  const auto &detId = hit2->geographicalId();
511  return (detId.det() == hitDet && detId.subdetId() == hitSubdet && ttopo_->layer(detId) == hitLayer);
512  });
513  };
514 
515  auto t1HitIter = t1->recHitsBegin();
516  if (!(*t1HitIter)->isValid()) {
517  IfLogTrace(debug_, "DuplicateTrackMerger") << " first t1 hit invalid";
518  return false;
519  }
520  auto t2HitIter = findHitOnT2(*t1HitIter);
521  if (t2HitIter == t2->recHitsEnd()) {
522  // if first hit not found, try with second
523  // if that fails, then reject
524  ++t1HitIter;
525  assert(t1HitIter != t1->recHitsEnd());
526 
527  if (!(*t1HitIter)->isValid()) {
528  IfLogTrace(debug_, "DuplicateTrackMerger") << " second t1 hit invalid";
529  return false;
530  }
531  t2HitIter = findHitOnT2(*t1HitIter);
532  if (t2HitIter == t2->recHitsEnd())
533  return false;
534  }
535  IfLogTrace(debug_, "DuplicateTrackMerger")
536  << " starting overlap check from t1 hit " << std::distance(t1->recHitsBegin(), t1HitIter) << " t2 hit "
537  << std::distance(t2->recHitsBegin(), t2HitIter);
538 
539  auto prevSubdet = (*t1HitIter)->geographicalId().subdetId();
540  const TrajectoryStateOnSurface tsosInner = trajectoryStateTransform::innerStateOnSurface(*t2, *geom_, magfield_);
541 
542  ++t1HitIter;
543  ++t2HitIter;
544  unsigned int missedLayers = 0;
545  while (t1HitIter != t1->recHitsEnd() && t2HitIter != t2->recHitsEnd()) {
546  // in case of invalid hits, reject immediately
547  if ((*t1HitIter)->getType() != TrackingRecHit::valid || trackerHitRTTI::isUndef(**t1HitIter) ||
548  (*t2HitIter)->getType() != TrackingRecHit::valid || trackerHitRTTI::isUndef(**t2HitIter)) {
549  IfLogTrace(debug_, "DuplicateTrackMerger")
550  << " t1 hit " << std::distance(t1->recHitsBegin(), t1HitIter) << " t2 hit "
551  << std::distance(t2->recHitsBegin(), t2HitIter) << " either is invalid, types t1 "
552  << (*t1HitIter)->getType() << " t2 " << (*t2HitIter)->getType();
553  return false;
554  }
555 
556  const auto &t1DetId = (*t1HitIter)->geographicalId();
557  const auto &t2DetId = (*t2HitIter)->geographicalId();
558 
559  const auto t1Det = t1DetId.det();
560  const auto t2Det = t2DetId.det();
561  if (t1Det != DetId::Tracker || t2Det != DetId::Tracker) {
562  IfLogTrace(debug_, "DuplicateTrackMerger") << " t1 hit " << std::distance(t1->recHitsBegin(), t1HitIter)
563  << " t2 hit " << std::distance(t2->recHitsBegin(), t2HitIter)
564  << " either not from Tracker, dets t1 " << t1Det << " t2 " << t2Det;
565  return false;
566  }
567 
568  const auto t1Subdet = t1DetId.subdetId();
569  const auto t1Layer = ttopo_->layer(t1DetId);
570 
571  // reject if hits have the same DetId but are different
572  if (t1DetId == t2DetId) {
573  if (!(*t1HitIter)->sharesInput(*t2HitIter, TrackingRecHit::all)) {
574  IfLogTrace(debug_, "DuplicateTrackMerger")
575  << " t1 hit " << std::distance(t1->recHitsBegin(), t1HitIter) << " t2 hit "
576  << std::distance(t2->recHitsBegin(), t2HitIter) << " same DetId (" << t1DetId.rawId()
577  << ") but do not share all input";
578  return false;
579  }
580  } else {
581  const auto t2Subdet = t2DetId.subdetId();
582  const auto t2Layer = ttopo_->layer(t2DetId);
583 
584  // reject if hits are on different layers
585  if (t1Subdet != t2Subdet || t1Layer != t2Layer) {
586  bool recovered = false;
587  // but try to recover first by checking if either one has skipped over a layer
588  if (t1Subdet == prevSubdet && t2Subdet != prevSubdet) {
589  // t1 has a layer t2 doesn't
590  ++t1HitIter;
591  recovered = true;
592  } else if (t1Subdet != prevSubdet && t2Subdet == prevSubdet) {
593  // t2 has a layer t1 doesn't
594  ++t2HitIter;
595  recovered = true;
596  } else if (t1Subdet == t2Subdet) {
597  prevSubdet = t1Subdet;
598  // same subdet, so layer must be different
599  if (t2Layer > t1Layer) {
600  // t1 has a layer t2 doesn't
601  ++t1HitIter;
602  recovered = true;
603  } else if (t1Layer > t2Layer) {
604  // t2 has a layer t1 doesn't
605  ++t2HitIter;
606  recovered = true;
607  }
608  }
609  if (recovered) {
610  ++missedLayers;
611  if (missedLayers > overlapCheckMaxMissingLayers_) {
612  IfLogTrace(debug_, "DuplicateTrackMerger") << " max number of missed layers exceeded";
613  return false;
614  }
615  continue;
616  }
617 
618  IfLogTrace(debug_, "DuplicateTrackMerger")
619  << " t1 hit " << std::distance(t1->recHitsBegin(), t1HitIter) << " t2 hit "
620  << std::distance(t2->recHitsBegin(), t2HitIter) << " are on different layers (subdet, layer) t1 "
621  << t1Subdet << "," << t1Layer << " t2 " << t2Subdet << "," << t2Layer;
622  return false;
623  }
624  // reject if same layer (but not same hit) in non-pixel detector
625  else if (t1Subdet != PixelSubdetector::PixelBarrel && t1Subdet != PixelSubdetector::PixelEndcap) {
626  IfLogTrace(debug_, "DuplicateTrackMerger")
627  << " t1 hit " << std::distance(t1->recHitsBegin(), t1HitIter) << " t2 hit "
628  << std::distance(t2->recHitsBegin(), t2HitIter) << " are on same layer, but in non-pixel detector (det "
629  << t1Det << " subdet " << t1Subdet << " layer " << t1Layer << ")";
630  return false;
631  }
632  }
633 
634  // Propagate longer track to the shorter track hit surface, check compatibility
635  TrajectoryStateOnSurface tsosPropagated = propagator_->propagate(tsosInner, (*t1HitIter)->det()->surface());
636  if (!tsosPropagated.isValid()) { // reject immediately if TSOS is not valid
637  IfLogTrace(debug_, "DuplicateTrackMerger")
638  << " t1 hit " << std::distance(t1->recHitsBegin(), t1HitIter) << " t2 hit "
639  << std::distance(t2->recHitsBegin(), t2HitIter) << " TSOS not valid";
640  return false;
641  }
642  auto passChi2Pair = chi2Estimator_->estimate(tsosPropagated, **t1HitIter);
643  if (!passChi2Pair.first) {
644  IfLogTrace(debug_, "DuplicateTrackMerger")
645  << " t1 hit " << std::distance(t1->recHitsBegin(), t1HitIter) << " t2 hit "
646  << std::distance(t2->recHitsBegin(), t2HitIter) << " hit chi2 compatibility failed with chi2 "
647  << passChi2Pair.second;
648  return false;
649  }
650 
651  prevSubdet = t1Subdet;
652  ++t1HitIter;
653  ++t2HitIter;
654  }
655  if (t1HitIter != t1->recHitsEnd()) {
656  IfLogTrace(debug_, "DuplicateTrackMerger") << " hits on t2 ended before hits on t1";
657  return false;
658  }
659 
660  IfLogTrace(debug_, "DuplicateTrackMerger") << " all hits on t2 are on layers whre t1 has also a hit";
661  return true;
662  }
663 } // namespace
664 
667 
reco::Track::outerPosition
const math::XYZPoint & outerPosition() const
position of the outermost hit
Definition: Track.h:62
ConfigurationDescriptions.h
Chi2MeasurementEstimatorBase.h
Propagator.h
edm::ESHandle::product
T const * product() const
Definition: ESHandle.h:86
RandomServiceHelper.t2
t2
Definition: RandomServiceHelper.py:257
reco::HitPattern::MISSING_OUTER_HITS
Definition: HitPattern.h:155
FreeTrajectoryState::momentum
GlobalVector momentum() const
Definition: FreeTrajectoryState.h:68
TrajectoryStateOnSurface.h
nTracks
const unsigned int nTracks(const reco::Vertex &sv)
Definition: TemplatedVertexArbitrator.h:44
PDWG_EXOHSCP_cff.tracks
tracks
Definition: PDWG_EXOHSCP_cff.py:28
mps_fire.i
i
Definition: mps_fire.py:355
edm::ParameterSetDescription::add
ParameterDescriptionBase * add(U const &iLabel, T const &value)
Definition: ParameterSetDescription.h:95
PixelSubdetector.h
FreeTrajectoryState.h
MessageLogger.h
dqmMemoryStats.float
float
Definition: dqmMemoryStats.py:127
PixelSubdetector::PixelEndcap
Definition: PixelSubdetector.h:11
trajectoryStateTransform::outerFreeState
FreeTrajectoryState outerFreeState(const reco::Track &tk, const MagneticField *field, bool withErr=true)
Definition: TrajectoryStateTransform.cc:98
PixelSubdetector::PixelBarrel
Definition: PixelSubdetector.h:11
ESHandle.h
edm::RefToBase::key
size_t key() const
Definition: RefToBase.h:219
PV3DBase::x
T x() const
Definition: PV3DBase.h:59
DuplicateTrackType.h
patZpeak.handle
handle
Definition: patZpeak.py:23
edm::EDGetTokenT< reco::TrackCollection >
GBRWrapperRcd.h
TrackerTopology
Definition: TrackerTopology.h:16
PV3DBase::theta
Geom::Theta< T > theta() const
Definition: PV3DBase.h:72
GBRForest
Definition: GBRForest.h:25
gather_cfg.cout
cout
Definition: gather_cfg.py:144
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
cms::cuda::assert
assert(be >=bs)
HLT_2018_cff.distance
distance
Definition: HLT_2018_cff.py:6417
EDProducer.h
trackerHitRTTI::isUndef
bool isUndef(TrackingRecHit const &hit)
Definition: trackerHitRTTI.h:23
FreeTrajectoryState::signedInverseMomentum
double signedInverseMomentum() const
Definition: FreeTrajectoryState.h:70
reco::TrackBase::originalAlgo
TrackAlgorithm originalAlgo() const
Definition: TrackBase.h:533
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
FreeTrajectoryState::position
GlobalPoint position() const
Definition: FreeTrajectoryState.h:67
GBRForest.h
edm::Handle< reco::TrackCollection >
TrackMerger
Definition: TrackMerger.h:14
reco::TrackBase::numberOfValidHits
unsigned short numberOfValidHits() const
number of valid hits found
Definition: TrackBase.h:751
PV3DBase::z
T z() const
Definition: PV3DBase.h:61
IdealMagneticFieldRecord
Definition: IdealMagneticFieldRecord.h:11
Propagator
Definition: Propagator.h:44
TrajectoryStateOnSurface
Definition: TrajectoryStateOnSurface.h:16
MakerMacros.h
reco::Track::innerMomentum
const math::XYZVector & innerMomentum() const
momentum vector at the innermost hit position
Definition: Track.h:59
test
Definition: SmallWORMDict.h:13
RandomServiceHelper.t1
t1
Definition: RandomServiceHelper.py:256
ctpps_dqm_sourceclient-live_cfg.test
test
Definition: ctpps_dqm_sourceclient-live_cfg.py:7
Track.h
edm::EventSetup::get
T get() const
Definition: EventSetup.h:73
TrackFwd.h
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
hgcalPlots.stat
stat
Definition: hgcalPlots.py:1111
TrackMerger.h
std::swap
void swap(edm::DataFrameContainer &lhs, edm::DataFrameContainer &rhs)
Definition: DataFrameContainer.h:209
TrackerDigiGeometryRecord
Definition: TrackerDigiGeometryRecord.h:15
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
reco::Track
Definition: Track.h:27
edm::ESHandle< GBRForest >
declareDynArray
#define declareDynArray(T, n, x)
Definition: DynArray.h:91
reco::TrackBase::charge
int charge() const
track electric charge
Definition: TrackBase.h:581
IfLogTrace
#define IfLogTrace(cond, cat)
Definition: MessageLogger.h:699
Point3DBase< float, GlobalTag >
ParameterSetDescription.h
b
double b
Definition: hdecay.h:118
mitigatedMETSequence_cff.U
U
Definition: mitigatedMETSequence_cff.py:36
edm::ConfigurationDescriptions
Definition: ConfigurationDescriptions.h:28
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
HLT_2018_cff.InputTag
InputTag
Definition: HLT_2018_cff.py:79016
LogDebug
#define LogDebug(id)
Definition: MessageLogger.h:670
edm::ParameterSet
Definition: ParameterSet.h:36
a
double a
Definition: hdecay.h:119
DetId::Tracker
Definition: DetId.h:25
TrackingRecHit::all
Definition: TrackingRecHit.h:59
Event.h
PV3DBase::y
T y() const
Definition: PV3DBase.h:60
TrajectoryStateClosestToPoint::pt
double pt() const
Definition: TrajectoryStateClosestToPoint.h:78
ModuleDef.h
DuplicateTrackType
DuplicateTrackType
Definition: DuplicateTrackType.h:4
iEvent
int iEvent
Definition: GenABIO.cc:224
value
Definition: value.py:1
mag2
T mag2() const
The vector magnitude squared. Equivalent to vec.dot(vec)
Definition: Basic3DVectorLD.h:124
M_PI
#define M_PI
Definition: BXVectorInputProducer.cc:50
edm::stream::EDProducer
Definition: EDProducer.h:38
edm::EventSetup
Definition: EventSetup.h:57
TrajectoryStateClosestToPoint
Definition: TrajectoryStateClosestToPoint.h:18
TSCPBuilderNoMaterial.h
TrackingRecHit::valid
Definition: TrackingRecHit.h:46
DynArray.h
Chi2MeasurementEstimatorBase
Definition: Chi2MeasurementEstimatorBase.h:14
reco::get
T get(const Candidate &c)
Definition: component.h:60
TrackingRecHit
Definition: TrackingRecHit.h:21
BaseTrackerRecHit.h
PV3DBase::mag
T mag() const
Definition: PV3DBase.h:64
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
trajectoryStateTransform::innerFreeState
FreeTrajectoryState innerFreeState(const reco::Track &tk, const MagneticField *field, bool withErr=true)
Definition: TrajectoryStateTransform.cc:86
eostools.move
def move(src, dest)
Definition: eostools.py:511
TSCPBuilderNoMaterial
Definition: TSCPBuilderNoMaterial.h:17
FreeTrajectoryState
Definition: FreeTrajectoryState.h:27
reco::Track::seedRef
const edm::RefToBase< TrajectorySeed > & seedRef() const
Definition: Track.h:155
T
long double T
Definition: Basic3DVectorLD.h:48
DuplicateTrackType::Overlapping
TrackingComponentsRecord.h
ev
bool ev
Definition: Hydjet2Hadronizer.cc:95
reco::HitPattern::MISSING_INNER_HITS
Definition: HitPattern.h:155
GBRWrapperRcd
Definition: GBRWrapperRcd.h:24
EventSetup.h
TrajectoryStateTransform.h
DuplicateTrackType::Disjoint
TrajectoryStateClosestToPoint::isValid
bool isValid() const
Definition: TrajectoryStateClosestToPoint.h:111
DuplicateTrackMerger
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
LogTrace
#define LogTrace(id)
Definition: MessageLogger.h:671
TrackerTopologyRcd
Definition: TrackerTopologyRcd.h:10
ParameterSet.h
TrackCandidate.h
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
trajectoryStateTransform::innerStateOnSurface
TrajectoryStateOnSurface innerStateOnSurface(const reco::Track &tk, const TrackingGeometry &geom, const MagneticField *field, bool withErr=true)
Definition: TrajectoryStateTransform.cc:110
TrajectoryStateClosestToPoint::theState
const FreeTrajectoryState & theState() const
Definition: TrajectoryStateClosestToPoint.h:96
TrajectoryStateClosestToPoint::momentum
GlobalVector momentum() const
Definition: TrajectoryStateClosestToPoint.h:92
edm::Event
Definition: Event.h:73
MagneticField
Definition: MagneticField.h:19
edm::InputTag
Definition: InputTag.h:15
PV3DBase::phi
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
TrajectoryStateOnSurface::isValid
bool isValid() const
Definition: TrajectoryStateOnSurface.h:54
TrackerGeometry
Definition: TrackerGeometry.h:14
TrackingComponentsRecord
Definition: TrackingComponentsRecord.h:12