CMS 3D CMS Logo

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