Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 /** \file
0002  *
0003  *  \author S. Bolognesi - INFN TO
0004  */
0005 
0006 #include "FilterByLTC.h"
0007 
0008 #include "FWCore/Framework/interface/Event.h"
0009 // #include "DataFormats/Common/interface/Handle.h"
0010 // #include "FWCore/ParameterSet/interface/ParameterSet.h"
0011 
0012 #include "DataFormats/LTCDigi/interface/LTCDigi.h"
0013 
0014 #include "FWCore/Framework/interface/MakerMacros.h"
0015 
0016 #include <iostream>
0017 
0018 using namespace std;
0019 using namespace edm;
0020 
0021 FilterByLTC::FilterByLTC(const ParameterSet& pset)
0022     : nEventsProcessed(0), nEventsSelected(0), ltcTag_(pset.getParameter<edm::InputTag>("ltcTag")) {
0023   theTriggerSource = pset.getParameter<int>("triggerSource");
0024 }
0025 
0026 FilterByLTC::~FilterByLTC() {}
0027 
0028 bool FilterByLTC::filter(Event& event, const EventSetup& eventSetup) {
0029   bool selectThisEvent = false;
0030 
0031   nEventsProcessed++;
0032 
0033   edm::Handle<LTCDigiCollection> ltcdigis;
0034   event.getByLabel(ltcTag_, ltcdigis);
0035 
0036   bool DTtrig = false, CSCtrig = false, RPCtrig = false;
0037   for (std::vector<LTCDigi>::const_iterator ltc_it = ltcdigis->begin(); ltc_it != ltcdigis->end(); ltc_it++) {
0038     if ((*ltc_it).HasTriggered(0))
0039       DTtrig = true;
0040     if ((*ltc_it).HasTriggered(1))
0041       CSCtrig = true;
0042     if ((*ltc_it).HasTriggered(2))
0043       RPCtrig = true;
0044     if ((*ltc_it).HasTriggered(3))
0045       RPCtrig = true;
0046     if ((*ltc_it).HasTriggered(4))
0047       RPCtrig = true;
0048   }
0049 
0050   switch (theTriggerSource) {
0051     case 1:
0052       //only DT trigger
0053       if (DTtrig && !CSCtrig && !RPCtrig)
0054         selectThisEvent = true;
0055       break;
0056 
0057     case 2:
0058       //only CSC trigger
0059       if (!DTtrig && CSCtrig && !RPCtrig)
0060         selectThisEvent = true;
0061       break;
0062 
0063     case 3:
0064       //only RPC trigger
0065       if (!DTtrig && !CSCtrig && RPCtrig)
0066         selectThisEvent = true;
0067       break;
0068 
0069     case 4:
0070       //both DT && CSC trigger
0071       if (DTtrig && CSCtrig && !RPCtrig)
0072         selectThisEvent = true;
0073       break;
0074 
0075     case 5:
0076       // both DT && RPC trigger
0077       if (DTtrig && !CSCtrig && RPCtrig)
0078         selectThisEvent = true;
0079       break;
0080 
0081     case 6:
0082       // both CSC && RPC trigger
0083       if (!DTtrig && CSCtrig && RPCtrig)
0084         selectThisEvent = true;
0085       break;
0086 
0087     case 7:
0088       // all CSC && RPC && DT trigger
0089       if (DTtrig && CSCtrig && RPCtrig)
0090         selectThisEvent = true;
0091       break;
0092 
0093     case 8:
0094       // No DT
0095       if (!DTtrig)
0096         selectThisEvent = true;
0097       break;
0098 
0099     case 9:
0100       // No CSC
0101       if (!CSCtrig)
0102         selectThisEvent = true;
0103       break;
0104 
0105     case 10:
0106       // No RPC
0107       if (!RPCtrig)
0108         selectThisEvent = true;
0109       break;
0110 
0111     case 11:
0112       //DT at least
0113       if (DTtrig)
0114         selectThisEvent = true;
0115       break;
0116 
0117     case 12:
0118       //CSC at least
0119       if (CSCtrig)
0120         selectThisEvent = true;
0121       break;
0122 
0123     case 13:
0124       //RPC at least
0125       if (RPCtrig)
0126         selectThisEvent = true;
0127       break;
0128 
0129     default:
0130       cout << "[FilterByLTC] Wrong trigger source selected" << endl;
0131   }
0132 
0133   if (selectThisEvent)
0134     nEventsSelected++;
0135 
0136   //cout<<"DT "<<DTtrig<<"  CSC "<<CSCtrig<<"  RPC "<<RPCtrig<<endl;
0137   //cout<<"events selected "<<selectThisEvent<<endl;
0138   //cout<<"total "<<nEventsSelected<<" / "<<nEventsProcessed<<endl;
0139 
0140   return selectThisEvent;
0141 }
0142 
0143 DEFINE_FWK_MODULE(FilterByLTC);