![]() |
![]() |
00001 // 00002 // Original Author: Marco Zanetti 00003 // Created: Tue Sep 9 15:56:24 CEST 2008 00004 00005 00006 // system include files 00007 #include <memory> 00008 00009 // user include files 00010 #include "FWCore/Framework/interface/Frameworkfwd.h" 00011 #include "FWCore/Framework/interface/EDFilter.h" 00012 00013 #include "FWCore/Framework/interface/Event.h" 00014 #include "FWCore/Framework/interface/MakerMacros.h" 00015 00016 #include "FWCore/ParameterSet/interface/ParameterSet.h" 00017 00018 #include "EventFilter/FEDInterface/interface/GlobalEventNumber.h" 00019 #include "DataFormats/FEDRawData/interface/FEDRawDataCollection.h" 00020 00021 #include <vector> 00022 00023 class BxNumberFilter : public edm::EDFilter { 00024 public: 00025 explicit BxNumberFilter(const edm::ParameterSet&); 00026 ~BxNumberFilter(); 00027 00028 private: 00029 virtual void beginJob() ; 00030 virtual bool filter(edm::Event&, const edm::EventSetup&); 00031 virtual void endJob() ; 00032 00033 edm::InputTag inputLabel; 00034 std::vector<int> goldenBXIds; 00035 unsigned int range; 00036 bool debug; 00037 00038 00039 }; 00040 00041 BxNumberFilter::BxNumberFilter(const edm::ParameterSet& iConfig) { 00042 00043 inputLabel = iConfig.getUntrackedParameter<edm::InputTag>("inputLabel",edm::InputTag("source")); 00044 goldenBXIds = iConfig.getParameter<std::vector<int> >("goldenBXIds"); 00045 range = iConfig.getUntrackedParameter<unsigned int>("range", 1); 00046 debug = iConfig.getUntrackedParameter<unsigned int>("debug", false); 00047 } 00048 00049 00050 BxNumberFilter::~BxNumberFilter() { } 00051 00052 00053 bool BxNumberFilter::filter(edm::Event& iEvent, const edm::EventSetup& iSetup) { 00054 using namespace edm; 00055 using namespace std; 00056 00057 bool result = false; 00058 00059 unsigned int GTEVMId= 812; 00060 00061 Handle<FEDRawDataCollection> rawdata; 00062 iEvent.getByLabel(inputLabel, rawdata); 00063 const FEDRawData& data = rawdata->FEDData(GTEVMId); 00064 evf::evtn::evm_board_setformat(data.size()); 00065 // loop over the predefined BX's 00066 for (vector<int>::const_iterator i = goldenBXIds.begin(); i != goldenBXIds.end(); i++) { 00067 00068 // Select the BX 00069 if ( evf::evtn::getfdlbx(data.data()) <= (*i) + range 00070 && 00071 evf::evtn::getfdlbx(data.data()) >= (*i) - range ) { 00072 result = true; 00073 00074 if (debug) { 00075 cout << "Event # " << evf::evtn::get(data.data(),true) << endl; 00076 cout << "LS # " << evf::evtn::getlbn(data.data()) << endl; 00077 cout << "ORBIT # " << evf::evtn::getorbit(data.data()) << endl; 00078 cout << "GPS LOW # " << evf::evtn::getgpslow(data.data()) << endl; 00079 cout << "GPS HI # " << evf::evtn::getgpshigh(data.data()) << endl; 00080 cout << "BX FROM FDL 0-xing # " << evf::evtn::getfdlbx(data.data()) << endl; 00081 } 00082 00083 } 00084 } 00085 return result; 00086 } 00087 00088 // ------------ method called once each job just before starting event loop ------------ 00089 void BxNumberFilter::beginJob() { 00090 } 00091 00092 // ------------ method called once each job just after ending the event loop ------------ 00093 void BxNumberFilter::endJob() { 00094 } 00095 00096 //define this as a plug-in 00097 DEFINE_FWK_MODULE(BxNumberFilter);