CMS 3D CMS Logo

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