Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:32:48

0001 #include "Validation/MuonCSCDigis/interface/CSCStubMatcher.h"
0002 #include "DataFormats/CSCDigi/interface/CSCConstants.h"
0003 #include <algorithm>
0004 
0005 using namespace std;
0006 
0007 CSCStubMatcher::CSCStubMatcher(const edm::ParameterSet& pSet, edm::ConsumesCollector&& iC) {
0008   useGEMs_ = pSet.getParameter<bool>("useGEMs");
0009 
0010   const auto& cscCLCT = pSet.getParameter<edm::ParameterSet>("cscCLCT");
0011   minBXCLCT_ = cscCLCT.getParameter<int>("minBX");
0012   maxBXCLCT_ = cscCLCT.getParameter<int>("maxBX");
0013   verboseCLCT_ = cscCLCT.getParameter<int>("verbose");
0014   minNHitsChamberCLCT_ = cscCLCT.getParameter<int>("minNHitsChamber");
0015 
0016   const auto& cscALCT = pSet.getParameter<edm::ParameterSet>("cscALCT");
0017   minBXALCT_ = cscALCT.getParameter<int>("minBX");
0018   maxBXALCT_ = cscALCT.getParameter<int>("maxBX");
0019   verboseALCT_ = cscALCT.getParameter<int>("verbose");
0020   minNHitsChamberALCT_ = cscALCT.getParameter<int>("minNHitsChamber");
0021 
0022   const auto& cscLCT = pSet.getParameter<edm::ParameterSet>("cscLCT");
0023   minBXLCT_ = cscLCT.getParameter<int>("minBX");
0024   maxBXLCT_ = cscLCT.getParameter<int>("maxBX");
0025   matchTypeTightLCT_ = cscLCT.getParameter<bool>("matchTypeTight");
0026   verboseLCT_ = cscLCT.getParameter<int>("verbose");
0027   minNHitsChamberLCT_ = cscLCT.getParameter<int>("minNHitsChamber");
0028   addGhostLCTs_ = cscLCT.getParameter<bool>("addGhosts");
0029 
0030   const auto& cscMPLCT = pSet.getParameter<edm::ParameterSet>("cscMPLCT");
0031   minBXMPLCT_ = cscMPLCT.getParameter<int>("minBX");
0032   maxBXMPLCT_ = cscMPLCT.getParameter<int>("maxBX");
0033   verboseMPLCT_ = cscMPLCT.getParameter<int>("verbose");
0034   minNHitsChamberMPLCT_ = cscMPLCT.getParameter<int>("minNHitsChamber");
0035 
0036   if (useGEMs_)
0037     gemDigiMatcher_.reset(new GEMDigiMatcher(pSet, std::move(iC)));
0038   cscDigiMatcher_.reset(new CSCDigiMatcher(pSet, std::move(iC)));
0039 
0040   clctInputTag_ = cscCLCT.getParameter<edm::InputTag>("inputTag");
0041   alctInputTag_ = cscALCT.getParameter<edm::InputTag>("inputTag");
0042   lctInputTag_ = cscLCT.getParameter<edm::InputTag>("inputTag");
0043   mplctInputTag_ = cscMPLCT.getParameter<edm::InputTag>("inputTag");
0044 
0045   clctToken_ = iC.consumes<CSCCLCTDigiCollection>(clctInputTag_);
0046   alctToken_ = iC.consumes<CSCALCTDigiCollection>(alctInputTag_);
0047   lctToken_ = iC.consumes<CSCCorrelatedLCTDigiCollection>(lctInputTag_);
0048   mplctToken_ = iC.consumes<CSCCorrelatedLCTDigiCollection>(mplctInputTag_);
0049 
0050   geomToken_ = iC.esConsumes<CSCGeometry, MuonGeometryRecord>();
0051 }
0052 
0053 void CSCStubMatcher::init(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0054   if (useGEMs_)
0055     gemDigiMatcher_->init(iEvent, iSetup);
0056   cscDigiMatcher_->init(iEvent, iSetup);
0057 
0058   iEvent.getByToken(clctToken_, clctsH_);
0059   iEvent.getByToken(alctToken_, alctsH_);
0060   iEvent.getByToken(lctToken_, lctsH_);
0061   iEvent.getByToken(mplctToken_, mplctsH_);
0062 
0063   cscGeometry_ = &iSetup.getData(geomToken_);
0064 }
0065 
0066 // do the matching
0067 void CSCStubMatcher::match(const SimTrack& t, const SimVertex& v) {
0068   // match simhits first
0069   if (useGEMs_)
0070     gemDigiMatcher_->match(t, v);
0071   cscDigiMatcher_->match(t, v);
0072 
0073   const CSCCLCTDigiCollection& clcts = *clctsH_.product();
0074   const CSCALCTDigiCollection& alcts = *alctsH_.product();
0075   const CSCCorrelatedLCTDigiCollection& lcts = *lctsH_.product();
0076   const CSCCorrelatedLCTDigiCollection& mplcts = *mplctsH_.product();
0077 
0078   // clear collections
0079   clear();
0080 
0081   if (!alctsH_.isValid()) {
0082     edm::LogError("CSCStubMatcher") << "Cannot get ALCTs with label " << alctInputTag_.encode();
0083   } else {
0084     matchALCTsToSimTrack(alcts);
0085   }
0086 
0087   if (!clctsH_.isValid()) {
0088     edm::LogError("CSCStubMatcher") << "Cannot get CLCTs with label " << clctInputTag_.encode();
0089   } else {
0090     matchCLCTsToSimTrack(clcts);
0091   }
0092 
0093   if (!lctsH_.isValid()) {
0094     edm::LogError("CSCStubMatcher") << "Cannot get LCTs with label " << lctInputTag_.encode();
0095   } else {
0096     matchLCTsToSimTrack(lcts);
0097   }
0098 
0099   if (!mplctsH_.isValid()) {
0100     edm::LogError("CSCStubMatcher") << "Cannot get MPLCTs with label " << mplctInputTag_.encode();
0101   } else {
0102     matchMPLCTsToSimTrack(mplcts);
0103   }
0104 }
0105 
0106 void CSCStubMatcher::matchCLCTsToSimTrack(const CSCCLCTDigiCollection& clcts) {
0107   const auto& cathode_ids = cscDigiMatcher_->chamberIdsStrip(0);
0108 
0109   for (const auto& id : cathode_ids) {
0110     CSCDetId ch_id(id);
0111     if (verboseCLCT_) {
0112       edm::LogInfo("CSCStubMatcher") << "To check CSC chamber " << ch_id;
0113     }
0114 
0115     int ring = ch_id.ring();
0116 
0117     // do not consider CSCs with too few hits
0118     if (cscDigiMatcher_->nLayersWithStripInChamber(ch_id) < minNHitsChamberCLCT_)
0119       continue;
0120 
0121     // get the comparator digis in this chamber
0122     std::vector<CSCComparatorDigiContainer> comps;
0123     for (int ilayer = CSCDetId::minLayerId(); ilayer <= CSCDetId::maxLayerId(); ilayer++) {
0124       CSCDetId layerid(ch_id.endcap(), ch_id.station(), ring, ch_id.chamber(), ilayer);
0125       comps.push_back(cscDigiMatcher_->comparatorDigisInDetId(layerid));
0126     }
0127 
0128     // print out the digis
0129     if (verboseCLCT_) {
0130       edm::LogInfo("CSCStubMatcher") << "clct: comparators " << ch_id;
0131       int layer = 0;
0132       for (const auto& p : comps) {
0133         layer++;
0134         for (const auto& q : p) {
0135           edm::LogInfo("CSCStubMatcher") << "L" << layer << " " << q << " " << q.getHalfStrip() << " ";
0136         }
0137       }
0138     }
0139 
0140     //use ME1b id to get CLCTs
0141     const bool isME1a(ch_id.station() == 1 and ch_id.ring() == 4);
0142     if (isME1a)
0143       ring = 1;
0144     CSCDetId ch_id2(ch_id.endcap(), ch_id.station(), ring, ch_id.chamber(), 0);
0145     auto id2 = ch_id2.rawId();  // CLCTs should be sorted into the det of the CLCTs.
0146 
0147     const auto& clcts_in_det = clcts.get(ch_id2);
0148 
0149     for (auto c = clcts_in_det.first; c != clcts_in_det.second; ++c) {
0150       if (verboseCLCT_)
0151         edm::LogInfo("CSCStubMatcher") << "clct " << ch_id2 << " " << *c;
0152 
0153       if (!c->isValid())
0154         continue;
0155 
0156       // check that the BX for this stub wasn't too early or too late
0157       if (c->getBX() < minBXCLCT_ || c->getBX() > maxBXCLCT_)
0158         continue;
0159 
0160       // store all CLCTs in this chamber
0161       chamber_to_clcts_all_[id2].push_back(*c);
0162 
0163       // check that at least 3 comparator digis were matched!
0164       int nMatches = 0;
0165       int layer = 0;
0166       for (const auto& p : comps) {
0167         layer++;
0168         for (const auto& q : p) {
0169           if (verboseCLCT_)
0170             edm::LogInfo("CSCStubMatcher") << "L" << layer << " " << q << " " << q.getHalfStrip() << " " << std::endl;
0171           for (const auto& clctComp : (*c).getHits()[layer - 1]) {
0172             if (clctComp == 65535)
0173               continue;
0174             if (verboseCLCT_) {
0175               edm::LogInfo("CSCStubMatcher") << "\t" << clctComp << " ";
0176             }
0177             if (q.getHalfStrip() == clctComp or (isME1a and q.getHalfStrip() + 128 == clctComp)) {
0178               nMatches++;
0179               if (verboseCLCT_) {
0180                 edm::LogInfo("CSCStubMatcher") << "\t\tnMatches " << nMatches << std::endl;
0181               }
0182             }
0183           }
0184         }
0185       }
0186 
0187       // require at least 3 good matches
0188       if (nMatches < 3)
0189         continue;
0190 
0191       if (verboseCLCT_)
0192         edm::LogInfo("CSCStubMatcher") << "clctGOOD";
0193 
0194       // store matching CLCTs in this chamber
0195       if (std::find(chamber_to_clcts_[id2].begin(), chamber_to_clcts_[id2].end(), *c) == chamber_to_clcts_[id2].end()) {
0196         chamber_to_clcts_[id2].push_back(*c);
0197       }
0198     }
0199     if (chamber_to_clcts_[id2].size() > 2) {
0200       edm::LogInfo("CSCStubMatcher") << "WARNING!!! too many CLCTs " << chamber_to_clcts_[id2].size() << " in "
0201                                      << ch_id2;
0202       for (auto& c : chamber_to_clcts_[id2])
0203         edm::LogInfo("CSCStubMatcher") << "  " << c;
0204     }
0205   }
0206 }
0207 
0208 void CSCStubMatcher::matchALCTsToSimTrack(const CSCALCTDigiCollection& alcts) {
0209   const auto& anode_ids = cscDigiMatcher_->chamberIdsWire(0);
0210   for (const auto& id : anode_ids) {
0211     CSCDetId ch_id(id);
0212 
0213     // fill 1 WG wide gaps
0214     const auto& digi_wgs = cscDigiMatcher_->wiregroupsInChamber(id, 1);
0215     if (verboseALCT_) {
0216       cout << "alct: digi_wgs " << ch_id << " ";
0217       copy(digi_wgs.begin(), digi_wgs.end(), ostream_iterator<int>(cout, " "));
0218       cout << endl;
0219     }
0220 
0221     int ring = ch_id.ring();
0222     if (ring == 4)
0223       ring = 1;  //use ME1b id to get ALCTs
0224     CSCDetId ch_id2(ch_id.endcap(), ch_id.station(), ring, ch_id.chamber(), 0);
0225     auto id2 = ch_id2.rawId();  // ALCTs should be sorted into the det of the ALCTs.
0226 
0227     const auto& alcts_in_det = alcts.get(ch_id2);
0228     for (auto a = alcts_in_det.first; a != alcts_in_det.second; ++a) {
0229       if (!a->isValid())
0230         continue;
0231 
0232       if (verboseALCT_)
0233         edm::LogInfo("CSCStubMatcher") << "alct " << ch_id << " " << *a;
0234 
0235       // check that the BX for stub wasn't too early or too late
0236       if (a->getBX() < minBXALCT_ || a->getBX() > maxBXALCT_)
0237         continue;
0238 
0239       int wg = a->getKeyWG() + 1;  // as ALCT wiregroups numbers start from 0
0240 
0241       // store all ALCTs in this chamber
0242       chamber_to_alcts_all_[id2].push_back(*a);
0243 
0244       // match by wiregroup with the digis
0245       if (digi_wgs.find(wg) == digi_wgs.end()) {
0246         continue;
0247       }
0248       if (verboseALCT_)
0249         edm::LogInfo("CSCStubMatcher") << "alctGOOD";
0250 
0251       // store matching ALCTs in this chamber
0252       if (std::find(chamber_to_alcts_[id2].begin(), chamber_to_alcts_[id2].end(), *a) == chamber_to_alcts_[id2].end()) {
0253         chamber_to_alcts_[id2].push_back(*a);
0254       }
0255     }
0256     if (chamber_to_alcts_[id2].size() > 2) {
0257       edm::LogInfo("CSCStubMatcher") << "WARNING!!! too many ALCTs " << chamber_to_alcts_[id2].size() << " in "
0258                                      << ch_id;
0259       for (auto& a : chamber_to_alcts_[id2])
0260         edm::LogInfo("CSCStubMatcher") << "  " << a;
0261     }
0262   }
0263 }
0264 
0265 void CSCStubMatcher::matchLCTsToSimTrack(const CSCCorrelatedLCTDigiCollection& lcts) {
0266   // only look for stubs in chambers that already have CLCT and ALCT
0267   const auto& cathode_ids = chamberIdsAllCLCT(0);
0268   const auto& anode_ids = chamberIdsAllALCT(0);
0269 
0270   std::set<int> cathode_and_anode_ids;
0271   std::set_union(cathode_ids.begin(),
0272                  cathode_ids.end(),
0273                  anode_ids.begin(),
0274                  anode_ids.end(),
0275                  std::inserter(cathode_and_anode_ids, cathode_and_anode_ids.end()));
0276 
0277   for (const auto& id : cathode_and_anode_ids) {
0278     CSCDetId ch_id(id);
0279 
0280     //use ME1b id to get LCTs
0281     int ring = ch_id.ring();
0282     if (ring == 4)
0283       ring = 1;
0284     CSCDetId ch_id2(ch_id.endcap(), ch_id.station(), ring, ch_id.chamber(), 0);
0285     auto id2 = ch_id2.rawId();  // LCTs should be sorted into the det of the LCTs.
0286 
0287     const auto& lcts_in_det = lcts.get(ch_id2);
0288 
0289     std::map<int, CSCCorrelatedLCTDigiContainer> bx_to_lcts;
0290 
0291     // collect all valid LCTs in a handy container
0292     CSCCorrelatedLCTDigiContainer lcts_tmp;
0293     for (auto lct = lcts_in_det.first; lct != lcts_in_det.second; ++lct) {
0294       if (!lct->isValid())
0295         continue;
0296       lcts_tmp.push_back(*lct);
0297       int bx = lct->getBX();
0298       bx_to_lcts[bx].push_back(*lct);
0299 
0300       // Add ghost LCTs when there are two in bx
0301       // and the two don't share half-strip or wiregroup
0302       if (bx_to_lcts[bx].size() == 2 and addGhostLCTs_) {
0303         auto lct11 = bx_to_lcts[bx][0];
0304         auto lct22 = bx_to_lcts[bx][1];
0305         addGhostLCTs(lct11, lct22, lcts_tmp);
0306       }
0307     }
0308 
0309     for (const auto& lct : lcts_tmp) {
0310       bool lct_clct_match(false);
0311       bool lct_alct_match(false);
0312       bool lct_gem1_match(false);
0313       bool lct_gem2_match(false);
0314 
0315       if (verboseLCT_) {
0316         edm::LogInfo("CSCStubMatcher") << ch_id << " " << ch_id2;
0317         edm::LogInfo("CSCStubMatcher") << lct;
0318         edm::LogInfo("CSCStubMatcher") << "getCLCT " << lct.getCLCT() << "\ngetALCT " << lct.getALCT() << "\ngetGEM1 "
0319                                        << lct.getGEM1() << "\ngetGEM2 " << lct.getGEM2();
0320       }
0321       // Check if matched to an CLCT
0322       for (const auto& p : clctsInChamber(id)) {
0323         if (p == lct.getCLCT()) {
0324           lct_clct_match = true;
0325           if (verboseLCT_)
0326             edm::LogInfo("CSCStubMatcher") << "\t...lct_clct_match";
0327           break;
0328         }
0329       }
0330 
0331       // Check if matched to an ALCT
0332       for (const auto& p : alctsInChamber(id)) {
0333         if (p == lct.getALCT()) {
0334           lct_alct_match = true;
0335           if (verboseLCT_)
0336             edm::LogInfo("CSCStubMatcher") << "\t...lct_alct_match";
0337           break;
0338         }
0339       }
0340 
0341       if (useGEMs_) {
0342         // fixME here: double check the timing of GEMPad
0343         if (ch_id.ring() == 1 and (ch_id.station() == 1 or ch_id.station() == 2)) {
0344           // Check if matched to an GEM pad L1
0345           const GEMDetId gemDetIdL1(ch_id.zendcap(), 1, ch_id.station(), 1, ch_id.chamber(), 0);
0346           for (const auto& p : gemDigiMatcher_->padsInChamber(gemDetIdL1.rawId())) {
0347             if (p == lct.getGEM1()) {
0348               lct_gem1_match = true;
0349               if (verboseLCT_)
0350                 edm::LogInfo("CSCStubMatcher") << "\t...lct_gem1_match";
0351               break;
0352             }
0353           }
0354 
0355           // Check if matched to an GEM pad L2
0356           const GEMDetId gemDetIdL2(ch_id.zendcap(), 1, ch_id.station(), 2, ch_id.chamber(), 0);
0357           for (const auto& p : gemDigiMatcher_->padsInChamber(gemDetIdL2.rawId())) {
0358             if (p == lct.getGEM2()) {
0359               lct_gem2_match = true;
0360               if (verboseLCT_)
0361                 edm::LogInfo("CSCStubMatcher") << "\t...lct_gem2_match";
0362               break;
0363             }
0364           }
0365         }
0366       }
0367 
0368       const bool alct_clct = lct_clct_match and lct_alct_match;
0369       const bool alct_gem = lct_alct_match and lct_gem1_match and lct_gem2_match;
0370       const bool clct_gem = lct_clct_match and lct_gem1_match and lct_gem2_match;
0371 
0372       bool lct_tight_matched = alct_clct or alct_gem or clct_gem;
0373       bool lct_loose_matched = lct_clct_match or lct_alct_match;
0374       bool lct_matched = matchTypeTightLCT_ ? lct_tight_matched : lct_loose_matched;
0375 
0376       if (lct_matched) {
0377         if (verboseLCT_)
0378           edm::LogInfo("CSCStubMatcher") << "...was matched";
0379         if (std::find(chamber_to_lcts_[id2].begin(), chamber_to_lcts_[id2].end(), lct) == chamber_to_lcts_[id2].end()) {
0380           chamber_to_lcts_[id2].emplace_back(lct);
0381         }
0382       }
0383     }  // lct loop over
0384   }
0385 }
0386 
0387 void CSCStubMatcher::matchMPLCTsToSimTrack(const CSCCorrelatedLCTDigiCollection& mplcts) {
0388   // match simtrack to MPC LCT by looking only in chambers
0389   // that already have LCTs matched to this simtrack
0390   const auto& lcts_ids = chamberIdsLCT(0);
0391 
0392   // loop on the detids
0393   for (const auto& id : lcts_ids) {
0394     const auto& mplcts_in_det = mplcts.get(id);
0395 
0396     // loop on the MPC LCTs in this detid
0397     for (auto lct = mplcts_in_det.first; lct != mplcts_in_det.second; ++lct) {
0398       if (!lct->isValid())
0399         continue;
0400 
0401       chamber_to_mplcts_all_[id].emplace_back(*lct);
0402 
0403       // check if this stub corresponds with a previously matched stub
0404       for (const auto& sim_stub : lctsInChamber(id)) {
0405         if (sim_stub == *lct) {
0406           if (std::find(chamber_to_mplcts_[id].begin(), chamber_to_mplcts_[id].end(), *lct) ==
0407               chamber_to_mplcts_[id].end()) {
0408             chamber_to_mplcts_[id].emplace_back(*lct);
0409           }
0410         }
0411       }
0412     }
0413   }
0414 }
0415 
0416 std::set<unsigned int> CSCStubMatcher::chamberIdsAllCLCT(int csc_type) const {
0417   return selectDetIds(chamber_to_clcts_all_, csc_type);
0418 }
0419 
0420 std::set<unsigned int> CSCStubMatcher::chamberIdsAllALCT(int csc_type) const {
0421   return selectDetIds(chamber_to_alcts_all_, csc_type);
0422 }
0423 
0424 std::set<unsigned int> CSCStubMatcher::chamberIdsAllLCT(int csc_type) const {
0425   return selectDetIds(chamber_to_lcts_all_, csc_type);
0426 }
0427 
0428 std::set<unsigned int> CSCStubMatcher::chamberIdsAllMPLCT(int csc_type) const {
0429   return selectDetIds(chamber_to_mplcts_all_, csc_type);
0430 }
0431 
0432 std::set<unsigned int> CSCStubMatcher::chamberIdsCLCT(int csc_type) const {
0433   return selectDetIds(chamber_to_clcts_, csc_type);
0434 }
0435 
0436 std::set<unsigned int> CSCStubMatcher::chamberIdsALCT(int csc_type) const {
0437   return selectDetIds(chamber_to_alcts_, csc_type);
0438 }
0439 
0440 std::set<unsigned int> CSCStubMatcher::chamberIdsLCT(int csc_type) const {
0441   return selectDetIds(chamber_to_lcts_, csc_type);
0442 }
0443 
0444 std::set<unsigned int> CSCStubMatcher::chamberIdsMPLCT(int csc_type) const {
0445   return selectDetIds(chamber_to_mplcts_, csc_type);
0446 }
0447 
0448 const CSCCLCTDigiContainer& CSCStubMatcher::allCLCTsInChamber(unsigned int detid) const {
0449   if (chamber_to_clcts_all_.find(detid) == chamber_to_clcts_all_.end())
0450     return no_clcts_;
0451   return chamber_to_clcts_all_.at(detid);
0452 }
0453 
0454 const CSCALCTDigiContainer& CSCStubMatcher::allALCTsInChamber(unsigned int detid) const {
0455   if (chamber_to_alcts_all_.find(detid) == chamber_to_alcts_all_.end())
0456     return no_alcts_;
0457   return chamber_to_alcts_all_.at(detid);
0458 }
0459 
0460 const CSCCorrelatedLCTDigiContainer& CSCStubMatcher::allLCTsInChamber(unsigned int detid) const {
0461   if (chamber_to_lcts_all_.find(detid) == chamber_to_lcts_all_.end())
0462     return no_lcts_;
0463   return chamber_to_lcts_all_.at(detid);
0464 }
0465 
0466 const CSCCorrelatedLCTDigiContainer& CSCStubMatcher::allMPLCTsInChamber(unsigned int detid) const {
0467   if (chamber_to_mplcts_all_.find(detid) == chamber_to_mplcts_all_.end())
0468     return no_mplcts_;
0469   return chamber_to_mplcts_all_.at(detid);
0470 }
0471 
0472 const CSCCLCTDigiContainer& CSCStubMatcher::clctsInChamber(unsigned int detid) const {
0473   if (chamber_to_clcts_.find(detid) == chamber_to_clcts_.end())
0474     return no_clcts_;
0475   return chamber_to_clcts_.at(detid);
0476 }
0477 
0478 const CSCALCTDigiContainer& CSCStubMatcher::alctsInChamber(unsigned int detid) const {
0479   if (chamber_to_alcts_.find(detid) == chamber_to_alcts_.end())
0480     return no_alcts_;
0481   return chamber_to_alcts_.at(detid);
0482 }
0483 
0484 const CSCCorrelatedLCTDigiContainer& CSCStubMatcher::lctsInChamber(unsigned int detid) const {
0485   if (chamber_to_lcts_.find(detid) == chamber_to_lcts_.end())
0486     return no_lcts_;
0487   return chamber_to_lcts_.at(detid);
0488 }
0489 
0490 const CSCCorrelatedLCTDigiContainer& CSCStubMatcher::mplctsInChamber(unsigned int detid) const {
0491   if (chamber_to_mplcts_.find(detid) == chamber_to_mplcts_.end())
0492     return no_mplcts_;
0493   return chamber_to_mplcts_.at(detid);
0494 }
0495 
0496 CSCCLCTDigi CSCStubMatcher::bestClctInChamber(unsigned int detid) const {
0497   //sort stubs based on quality
0498   const auto& input(clctsInChamber(detid));
0499   int bestQ = 0;
0500   int index = -1;
0501   for (unsigned int i = 0; i < input.size(); ++i) {
0502     int quality = input[i].getQuality();
0503     if (quality > bestQ) {
0504       bestQ = quality;
0505       index = i;
0506     }
0507   }
0508   if (index != -1)
0509     return input[index];
0510   return CSCCLCTDigi();
0511 }
0512 
0513 CSCALCTDigi CSCStubMatcher::bestAlctInChamber(unsigned int detid) const {
0514   //sort stubs based on quality
0515   const auto& input(alctsInChamber(detid));
0516   int bestQ = 0;
0517   int index = -1;
0518   for (unsigned int i = 0; i < input.size(); ++i) {
0519     int quality = input[i].getQuality();
0520     if (quality > bestQ) {
0521       bestQ = quality;
0522       index = i;
0523     }
0524   }
0525   if (index != -1)
0526     return input[index];
0527   return CSCALCTDigi();
0528 }
0529 
0530 CSCCorrelatedLCTDigi CSCStubMatcher::bestLctInChamber(unsigned int detid) const {
0531   //sort stubs based on quality
0532   const auto& input(lctsInChamber(detid));
0533   int bestQ = 0;
0534   int index = -1;
0535   for (unsigned int i = 0; i < input.size(); ++i) {
0536     int quality = input[i].getQuality();
0537     if (quality > bestQ) {
0538       bestQ = quality;
0539       index = i;
0540     }
0541   }
0542   if (index != -1)
0543     return input[index];
0544   return CSCCorrelatedLCTDigi();
0545 }
0546 
0547 float CSCStubMatcher::zpositionOfLayer(unsigned int detid, int layer) const {
0548   const auto& id = CSCDetId(detid);
0549   const auto& chamber(cscGeometry_->chamber(id));
0550   return fabs(chamber->layer(layer)->centerOfStrip(20).z());
0551 }
0552 
0553 int CSCStubMatcher::nChambersWithCLCT(int min_quality) const {
0554   int result = 0;
0555   const auto& chamber_ids = chamberIdsCLCT();
0556   for (const auto& id : chamber_ids) {
0557     int nStubChamber = 0;
0558     const auto& clcts = clctsInChamber(id);
0559     for (const auto& clct : clcts) {
0560       if (!clct.isValid())
0561         continue;
0562       if (clct.getQuality() >= min_quality) {
0563         nStubChamber++;
0564       }
0565     }
0566     if (nStubChamber > 0) {
0567       ++result;
0568     }
0569   }
0570   return result;
0571 }
0572 
0573 int CSCStubMatcher::nChambersWithALCT(int min_quality) const {
0574   int result = 0;
0575   const auto& chamber_ids = chamberIdsALCT();
0576   for (const auto& id : chamber_ids) {
0577     int nStubChamber = 0;
0578     const auto& alcts = alctsInChamber(id);
0579     for (const auto& alct : alcts) {
0580       if (!alct.isValid())
0581         continue;
0582       if (alct.getQuality() >= min_quality) {
0583         nStubChamber++;
0584       }
0585     }
0586     if (nStubChamber > 0) {
0587       ++result;
0588     }
0589   }
0590   return result;
0591 }
0592 
0593 int CSCStubMatcher::nChambersWithLCT(int min_quality) const {
0594   int result = 0;
0595   const auto& chamber_ids = chamberIdsLCT();
0596   for (const auto& id : chamber_ids) {
0597     int nStubChamber = 0;
0598     const auto& lcts = lctsInChamber(id);
0599     for (const auto& lct : lcts) {
0600       if (!lct.isValid())
0601         continue;
0602       if (lct.getQuality() >= min_quality) {
0603         nStubChamber++;
0604       }
0605     }
0606     if (nStubChamber > 0) {
0607       ++result;
0608     }
0609   }
0610   return result;
0611 }
0612 
0613 int CSCStubMatcher::nChambersWithMPLCT(int min_quality) const {
0614   int result = 0;
0615   const auto& chamber_ids = chamberIdsMPLCT();
0616   for (const auto& id : chamber_ids) {
0617     int nStubChamber = 0;
0618     const auto& mplcts = mplctsInChamber(id);
0619     for (const auto& mplct : mplcts) {
0620       if (!mplct.isValid())
0621         continue;
0622       if (mplct.getQuality() >= min_quality) {
0623         nStubChamber++;
0624       }
0625     }
0626     if (nStubChamber > 0) {
0627       ++result;
0628     }
0629   }
0630   return result;
0631 }
0632 
0633 bool CSCStubMatcher::lctInChamber(const CSCDetId& id, const CSCCorrelatedLCTDigi& lct) const {
0634   for (const auto& stub : lctsInChamber(id.rawId())) {
0635     if (stub == lct)
0636       return true;
0637   }
0638   return false;
0639 }
0640 
0641 GlobalPoint CSCStubMatcher::getGlobalPosition(unsigned int rawId, const CSCCorrelatedLCTDigi& lct) const {
0642   CSCDetId cscId(rawId);
0643   CSCDetId keyId(cscId.endcap(), cscId.station(), cscId.ring(), cscId.chamber(), CSCConstants::KEY_CLCT_LAYER);
0644   float fractional_strip = lct.getFractionalStrip();
0645   // case ME1/1
0646   if (cscId.station() == 1 and (cscId.ring() == 4 || cscId.ring() == 1)) {
0647     int ring = 1;  // Default to ME1/b
0648     if (lct.getStrip() > CSCConstants::MAX_HALF_STRIP_ME1B) {
0649       ring = 4;  // Change to ME1/a if the HalfStrip Number exceeds the range of ME1/b
0650       fractional_strip -= CSCConstants::NUM_STRIPS_ME1B;
0651     }
0652     CSCDetId cscId_(cscId.endcap(), cscId.station(), ring, cscId.chamber(), cscId.layer());
0653     cscId = cscId_;
0654   }
0655   // regular cases
0656   const auto& chamber = cscGeometry_->chamber(cscId);
0657   const auto& layer_geo = chamber->layer(CSCConstants::KEY_CLCT_LAYER)->geometry();
0658   // LCT::getKeyWG() also starts from 0
0659   float wire = layer_geo->middleWireOfGroup(lct.getKeyWG() + 1);
0660   const LocalPoint& csc_intersect = layer_geo->intersectionOfStripAndWire(fractional_strip, wire);
0661   const GlobalPoint& csc_gp = cscGeometry_->idToDet(keyId)->surface().toGlobal(csc_intersect);
0662   return csc_gp;
0663 }
0664 
0665 void CSCStubMatcher::clear() {
0666   chamber_to_clcts_all_.clear();
0667   chamber_to_alcts_all_.clear();
0668   chamber_to_lcts_all_.clear();
0669   chamber_to_mplcts_all_.clear();
0670 
0671   chamber_to_clcts_.clear();
0672   chamber_to_alcts_.clear();
0673   chamber_to_lcts_.clear();
0674   chamber_to_mplcts_.clear();
0675 }
0676 
0677 void CSCStubMatcher::addGhostLCTs(const CSCCorrelatedLCTDigi& lct11,
0678                                   const CSCCorrelatedLCTDigi& lct22,
0679                                   CSCCorrelatedLCTDigiContainer& lcts_tmp) const {
0680   int wg1 = lct11.getKeyWG();
0681   int wg2 = lct22.getKeyWG();
0682   int hs1 = lct11.getStrip();
0683   int hs2 = lct22.getStrip();
0684 
0685   if (!(wg1 == wg2 || hs1 == hs2)) {
0686     // flip the ALCTs
0687     CSCCorrelatedLCTDigi lct12 = lct11;
0688     lct12.setWireGroup(wg2);
0689     lct12.setALCT(lct22.getALCT());
0690     lct12.setCLCT(lct11.getCLCT());
0691     lcts_tmp.push_back(lct12);
0692 
0693     CSCCorrelatedLCTDigi lct21 = lct22;
0694     lct21.setWireGroup(wg1);
0695     lct21.setALCT(lct11.getALCT());
0696     lct21.setCLCT(lct22.getCLCT());
0697     lcts_tmp.push_back(lct21);
0698   }
0699 }