Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-09-07 04:36:38

0001 #include "HLTCSCOverlapFilter.h"
0002 
0003 #include "FWCore/Framework/interface/ESHandle.h"
0004 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0005 #include "FWCore/ServiceRegistry/interface/Service.h"
0006 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0007 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0008 #include "CommonTools/UtilAlgos/interface/TFileService.h"
0009 #include "DataFormats/Common/interface/Handle.h"
0010 #include "DataFormats/MuonDetId/interface/CSCDetId.h"
0011 #include "Geometry/CSCGeometry/interface/CSCGeometry.h"
0012 #include "Geometry/Records/interface/MuonGeometryRecord.h"
0013 
0014 HLTCSCOverlapFilter::HLTCSCOverlapFilter(const edm::ParameterSet& iConfig)
0015     : HLTFilter(iConfig),
0016       muonGeometryRecordToken_(esConsumes()),
0017       m_input(iConfig.getParameter<edm::InputTag>("input")),
0018       m_minHits(iConfig.getParameter<unsigned int>("minHits")),
0019       m_xWindow(iConfig.getParameter<double>("xWindow")),
0020       m_yWindow(iConfig.getParameter<double>("yWindow")),
0021       m_ring1(iConfig.getParameter<bool>("ring1")),
0022       m_ring2(iConfig.getParameter<bool>("ring2")),
0023       m_fillHists(iConfig.getParameter<bool>("fillHists")),
0024       m_nhitsNoWindowCut(nullptr),
0025       m_xdiff(nullptr),
0026       m_ydiff(nullptr),
0027       m_pairsWithWindowCut(nullptr) {
0028   cscrechitsToken = consumes<CSCRecHit2DCollection>(m_input);
0029   if (m_fillHists) {
0030     edm::Service<TFileService> tfile;
0031     m_nhitsNoWindowCut = tfile->make<TH1F>("nhitsNoWindowCut", "nhitsNoWindowCut", 16, -0.5, 15.5);
0032     m_xdiff = tfile->make<TH1F>("xdiff", "xdiff", 200, 0., 10.);
0033     m_ydiff = tfile->make<TH1F>("ydiff", "ydiff", 200, 0., 10.);
0034     m_pairsWithWindowCut = tfile->make<TH1F>("pairsWithWindowCut", "pairsWithWindowCut", 226, -0.5, 225.5);
0035   }
0036 }
0037 
0038 HLTCSCOverlapFilter::~HLTCSCOverlapFilter() = default;
0039 
0040 void HLTCSCOverlapFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0041   edm::ParameterSetDescription desc;
0042   makeHLTFilterDescription(desc);
0043   desc.add<edm::InputTag>("input", edm::InputTag("hltCsc2DRecHits"));
0044   desc.add<unsigned int>("minHits", 4);
0045   desc.add<double>("xWindow", 1000.);
0046   desc.add<double>("yWindow", 1000.);
0047   desc.add<bool>("ring1", true);
0048   desc.add<bool>("ring2", false);
0049   desc.add<bool>("fillHists", false);
0050   descriptions.add("hltCSCOverlapFilter", desc);
0051 }
0052 
0053 bool HLTCSCOverlapFilter::hltFilter(edm::Event& iEvent,
0054                                     const edm::EventSetup& iSetup,
0055                                     trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0056   edm::Handle<CSCRecHit2DCollection> hits;
0057   iEvent.getByToken(cscrechitsToken, hits);
0058 
0059   edm::ESHandle<CSCGeometry> cscGeometry;
0060   bool got_cscGeometry = false;
0061 
0062   std::map<int, std::vector<const CSCRecHit2D*> > chamber_tohit;
0063 
0064   for (auto const& hit : *hits) {
0065     CSCDetId id(hit.geographicalId());
0066     int chamber_id = CSCDetId(id.endcap(), id.station(), id.ring(), id.chamber(), 0).rawId();
0067 
0068     if ((m_ring1 && (id.ring() == 1 || id.ring() == 4)) || (m_ring2 && id.ring() == 2)) {
0069       std::map<int, std::vector<const CSCRecHit2D*> >::const_iterator chamber_iter = chamber_tohit.find(chamber_id);
0070       if (chamber_iter == chamber_tohit.end()) {
0071         std::vector<const CSCRecHit2D*> newlist;
0072         newlist.push_back(&hit);
0073       }
0074       chamber_tohit[chamber_id].push_back(&hit);
0075     }  // end if this ring is selected
0076   }  // end loop over hits
0077 
0078   bool keep = false;
0079   unsigned int minHitsSquared = m_minHits * m_minHits;
0080   for (std::map<int, std::vector<const CSCRecHit2D*> >::const_iterator chamber_iter = chamber_tohit.begin();
0081        chamber_iter != chamber_tohit.end();
0082        ++chamber_iter) {
0083     if (m_fillHists) {
0084       m_nhitsNoWindowCut->Fill(chamber_iter->second.size());
0085     }
0086 
0087     if (chamber_iter->second.size() >= m_minHits) {
0088       CSCDetId chamber_id(chamber_iter->first);
0089       int chamber = chamber_id.chamber();
0090       int next = chamber + 1;
0091 
0092       // Some rings have 36 chambers, others have 18.  This will still be valid when ME4/2 is added.
0093       if (next == 37 && (std::abs(chamber_id.station()) == 1 || chamber_id.ring() == 2))
0094         next = 1;
0095       if (next == 19 && (std::abs(chamber_id.station()) != 1 && chamber_id.ring() == 1))
0096         next = 1;
0097 
0098       int next_id = CSCDetId(chamber_id.endcap(), chamber_id.station(), chamber_id.ring(), next, 0).rawId();
0099 
0100       std::map<int, std::vector<const CSCRecHit2D*> >::const_iterator chamber_next = chamber_tohit.find(next_id);
0101       if (chamber_next != chamber_tohit.end() && chamber_next->second.size() >= m_minHits) {
0102         if (!got_cscGeometry) {
0103           cscGeometry = iSetup.getHandle(muonGeometryRecordToken_);
0104           got_cscGeometry = true;
0105         }
0106 
0107         unsigned int pairs_in_window = 0;
0108         for (auto hit1 = chamber_iter->second.begin(); hit1 != chamber_iter->second.end(); ++hit1) {
0109           for (auto hit2 : chamber_next->second) {
0110             GlobalPoint pos1 =
0111                 cscGeometry->idToDet((*hit1)->geographicalId())->surface().toGlobal((*hit1)->localPosition());
0112             GlobalPoint pos2 = cscGeometry->idToDet(hit2->geographicalId())->surface().toGlobal(hit2->localPosition());
0113 
0114             if (m_fillHists) {
0115               m_xdiff->Fill(pos1.x() - pos2.x());
0116               m_ydiff->Fill(pos1.y() - pos2.y());
0117             }
0118 
0119             if (fabs(pos1.x() - pos2.x()) < m_xWindow && fabs(pos1.y() - pos2.y()) < m_yWindow)
0120               pairs_in_window++;
0121 
0122             if (pairs_in_window >= minHitsSquared) {
0123               keep = true;
0124               if (!m_fillHists)
0125                 return true;
0126             }
0127           }  // end loop over hits in chamber 2
0128         }  // end loop over hits in chamber 1
0129 
0130         if (m_fillHists) {
0131           m_pairsWithWindowCut->Fill(pairs_in_window);
0132         }
0133 
0134       }  // end if chamber 2 has enough hits
0135     }  // end if chamber 1 has enough hits
0136   }  // end loop over chambers
0137 
0138   return keep;
0139 }
0140 
0141 // declare this class as a framework plugin
0142 #include "FWCore/Framework/interface/MakerMacros.h"
0143 DEFINE_FWK_MODULE(HLTCSCOverlapFilter);