Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:19:58

0001 /** \class CruzetL1DTfilter
0002 *
0003 *  
0004 *  This class is for use with Cosmic MC to simulate the L1 DT trigger in cruzet
0005 *
0006 *
0007 *  \author Ivan Mikulec
0008 *
0009 */
0010 #include "FWCore/Framework/interface/Event.h"
0011 #include "FWCore/Framework/interface/EventSetup.h"
0012 #include "FWCore/Framework/interface/global/EDFilter.h"
0013 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0014 #include "DataFormats/L1GlobalMuonTrigger/interface/L1MuRegionalCand.h"
0015 #include "DataFormats/L1GlobalMuonTrigger/interface/L1MuGMTReadoutCollection.h"
0016 #include "FWCore/Utilities/interface/InputTag.h"
0017 #include "FWCore/Framework/interface/Frameworkfwd.h"
0018 #include "FWCore/Framework/interface/MakerMacros.h"
0019 
0020 class CruzetL1DTfilter : public edm::global::EDFilter<> {
0021 public:
0022   explicit CruzetL1DTfilter(const edm::ParameterSet&);
0023   bool filter(edm::StreamID, edm::Event&, const edm::EventSetup&) const;
0024 
0025 private:
0026   // mode:
0027   // 1 - bottom only
0028   // 2 - top only
0029   // 3 - top or bottom
0030   // 4 - top and bottom
0031   int _mode;
0032   edm::InputTag _GMTInputTag;
0033 };
0034 
0035 CruzetL1DTfilter::CruzetL1DTfilter(const edm::ParameterSet& ps) {
0036   _mode = ps.getParameter<int>("mode");
0037   _GMTInputTag = ps.getParameter<edm::InputTag>("GMTInputTag");
0038 }
0039 
0040 bool CruzetL1DTfilter::filter(edm::StreamID, edm::Event& e, const edm::EventSetup& es) const {
0041   edm::Handle<L1MuGMTReadoutCollection> gmtrc;
0042   e.getByLabel(_GMTInputTag, gmtrc);
0043 
0044   bool result = false;
0045   L1MuGMTReadoutRecord gmtrr = gmtrc->getRecord();
0046 
0047   std::vector<L1MuRegionalCand>::const_iterator iter1;
0048   std::vector<L1MuRegionalCand> rmc = gmtrr.getDTBXCands();
0049   bool top = false;
0050   bool bot = false;
0051   for (iter1 = rmc.begin(); iter1 != rmc.end(); iter1++) {
0052     int phi = (*iter1).phi_packed();
0053     if (phi >= 18 && phi <= 53)
0054       top = true;
0055     if (phi >= 90 && phi <= 125)
0056       bot = true;
0057   }
0058   if (_mode == 1 && bot)
0059     result = true;
0060   if (_mode == 2 && top)
0061     result = true;
0062   if (_mode == 3 && (top || bot))
0063     result = true;
0064   if (_mode == 4 && (top && bot))
0065     result = true;
0066 
0067   return result;
0068 }
0069 
0070 DEFINE_FWK_MODULE(CruzetL1DTfilter);