Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:27:17

0001 #include "RecoMuon/TrackerSeedGenerator/interface/L1MuonSeedsMerger.h"
0002 
0003 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0004 #include "DataFormats/TrackReco/interface/Track.h"
0005 
0006 L1MuonSeedsMerger::L1MuonSeedsMerger(const edm::ParameterSet& cfg) {
0007   theDeltaEtaCut = cfg.getParameter<double>("deltaEtaCut");
0008   theDiffRelPtCut = cfg.getParameter<double>("diffRelPtCut");
0009 }
0010 
0011 void L1MuonSeedsMerger::resolve(std::vector<TrackAndHits>& tracks) const {
0012   sort(tracks.begin(), tracks.end(), Less());
0013   typedef std::vector<TrackAndHits>::iterator Tracks_Itr;
0014   Tracks_Itr it1 = tracks.begin();
0015   while (it1 != tracks.end()) {
0016     for (Tracks_Itr it2 = it1 + 1; it1->first && it2 < tracks.end(); it2++) {
0017       if (!it2->first)
0018         continue;
0019       if (it2->first->eta() - it1->first->eta() > theDeltaEtaCut)
0020         break;
0021       switch (compare(&(*it1), &(*it2))) {
0022         case killFirst:
0023           delete it1->first;
0024           it1->first = nullptr;
0025           break;
0026         case killSecond:
0027           delete it2->first;
0028           it2->first = nullptr;
0029           break;
0030         case mergeTwo:
0031           *it2 = *(merge(&(*it1), &(*it2)));
0032           it1->first = nullptr;
0033           break;
0034         case goAhead:
0035         default:
0036           break;
0037       }
0038     }
0039     if (nullptr == it1->first)
0040       tracks.erase(it1);
0041     else
0042       it1++;
0043   }
0044 }
0045 
0046 bool L1MuonSeedsMerger::Less::operator()(const TrackAndHits& a, const TrackAndHits& b) const {
0047   return (a.first->eta() < b.first->eta());
0048 }
0049 
0050 const L1MuonSeedsMerger::TrackAndHits* L1MuonSeedsMerger::merge(const TrackAndHits* a, const TrackAndHits* b) const {
0051   // temporary algorith, takes track with bigger pt, to be reimplemented
0052   if (a->first->pt() > b->first->pt()) {
0053     delete b->first;
0054     return a;
0055   } else {
0056     delete a->first;
0057     return b;
0058   }
0059 }
0060 
0061 L1MuonSeedsMerger::Action L1MuonSeedsMerger::compare(const TrackAndHits* a, const TrackAndHits* b) const {
0062   int nshared = 0;
0063   const SeedingHitSet& hitsA = a->second;
0064   const SeedingHitSet& hitsB = b->second;
0065   for (unsigned int iHitA = 0, nHitsA = hitsA.size(); iHitA < nHitsA; ++iHitA) {
0066     const TrackingRecHit* trhA = hitsA[iHitA]->hit();
0067     for (unsigned int iHitB = 0, nHitsB = hitsB.size(); iHitB < nHitsB; ++iHitB) {
0068       const TrackingRecHit* trhB = hitsB[iHitB]->hit();
0069       if (trhA == trhB)
0070         nshared++;
0071     }
0072   }
0073 
0074   if (nshared >= 2) {
0075     if (hitsA.size() >= 3 && hitsB.size() >= 3)
0076       return (a->first->chi2() > b->first->chi2()) ? killFirst : killSecond;
0077     else if (hitsB.size() >= 3)
0078       return killFirst;
0079     else
0080       return killSecond;
0081   } else if (nshared >= 1) {
0082     if (hitsA.size() != hitsB.size())
0083       return (hitsA.size() < hitsB.size()) ? killFirst : killSecond;
0084     else if (hitsA.size() >= 3 && a->first->charge() == b->first->charge() &&
0085              fabs(a->first->pt() - b->first->pt()) / b->first->pt() < theDiffRelPtCut)
0086       return (a->first->chi2() > b->first->chi2()) ? killFirst : killSecond;
0087     else if (hitsA.size() == 2)
0088       return (a->first->pt() < b->first->pt()) ? killFirst : killSecond;
0089     else
0090       return goAhead;
0091   } else
0092     return goAhead;
0093 }