Go to the documentation of this file.00001
00018 #include <vector>
00019
00020 #include "FWCore/Framework/interface/Frameworkfwd.h"
00021 #include "FWCore/Framework/interface/Event.h"
00022 #include "FWCore/Framework/interface/EDFilter.h"
00023 #include "FWCore/Framework/interface/ESWatcher.h"
00024 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00025 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
00026 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
00027 #include "FWCore/Utilities/interface/InputTag.h"
00028 #include "CondFormats/DataRecord/interface/L1GtTriggerMaskTechTrigRcd.h"
00029 #include "CondFormats/DataRecord/interface/L1GtTriggerMaskAlgoTrigRcd.h"
00030
00031
00032 #define PHYSICS_BITS_SIZE 128
00033 #define TECHNICAL_BITS_SIZE 64
00034
00035
00036
00037
00038
00039 class HLTLevel1Activity : public edm::EDFilter {
00040 public:
00041 explicit HLTLevel1Activity(const edm::ParameterSet&);
00042 ~HLTLevel1Activity();
00043 static void fillDescriptions(edm::ConfigurationDescriptions & descriptions);
00044 virtual bool filter(edm::Event&, const edm::EventSetup&);
00045
00046 private:
00047 edm::InputTag m_gtReadoutRecord;
00048 std::vector<int> m_bunchCrossings;
00049 std::vector<bool> m_selectPhysics;
00050 std::vector<bool> m_selectTechnical;
00051 std::vector<bool> m_maskedPhysics;
00052 std::vector<bool> m_maskedTechnical;
00053 unsigned int m_daqPartitions;
00054 bool m_ignoreL1Mask;
00055 bool m_invert;
00056
00057 edm::ESWatcher<L1GtTriggerMaskAlgoTrigRcd> m_watchPhysicsMask;
00058 edm::ESWatcher<L1GtTriggerMaskTechTrigRcd> m_watchTechnicalMask;
00059 };
00060
00061 #include <boost/foreach.hpp>
00062
00063 #include "FWCore/Framework/interface/EventSetup.h"
00064 #include "FWCore/Framework/interface/ESHandle.h"
00065 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00066 #include "CondFormats/L1TObjects/interface/L1GtTriggerMask.h"
00067 #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerReadoutRecord.h"
00068
00069
00070
00071
00072 HLTLevel1Activity::HLTLevel1Activity(const edm::ParameterSet & config) :
00073 m_gtReadoutRecord( config.getParameter<edm::InputTag> ("L1GtReadoutRecordTag") ),
00074 m_bunchCrossings( config.getParameter<std::vector<int> > ("bunchCrossings") ),
00075 m_selectPhysics( PHYSICS_BITS_SIZE ),
00076 m_selectTechnical( TECHNICAL_BITS_SIZE ),
00077 m_maskedPhysics( PHYSICS_BITS_SIZE ),
00078 m_maskedTechnical( TECHNICAL_BITS_SIZE ),
00079 m_daqPartitions( config.getParameter<unsigned int> ("daqPartitions") ),
00080 m_ignoreL1Mask( config.getParameter<bool> ("ignoreL1Mask") ),
00081 m_invert( config.getParameter<bool> ("invert") )
00082 {
00083 unsigned long long low = config.getParameter<unsigned long long>("physicsLoBits");
00084 unsigned long long high = config.getParameter<unsigned long long>("physicsHiBits");
00085 unsigned long long tech = config.getParameter<unsigned long long>("technicalBits");
00086 for (unsigned int i = 0; i < 64; i++) {
00087 m_selectPhysics[i] = low & (0x01ULL << (unsigned long long) i);
00088 m_maskedPhysics[i] = low & (0x01ULL << (unsigned long long) i);
00089 }
00090 for (unsigned int i = 0; i < 64; i++) {
00091 m_selectPhysics[i+64] = high & (0x01ULL << (unsigned long long) i);
00092 m_maskedPhysics[i+64] = high & (0x01ULL << (unsigned long long) i);
00093 }
00094 for (unsigned int i = 0; i < 64; i++) {
00095 m_selectTechnical[i] = tech & (0x01ULL << (unsigned long long) i);
00096 m_maskedTechnical[i] = tech & (0x01ULL << (unsigned long long) i);
00097 }
00098 }
00099
00100 HLTLevel1Activity::~HLTLevel1Activity()
00101 {
00102 }
00103
00104 void
00105 HLTLevel1Activity::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
00106 edm::ParameterSetDescription desc;
00107 desc.add<edm::InputTag>("L1GtReadoutRecordTag",edm::InputTag("hltGtDigis"));
00108 {
00109 std::vector<int> temp1;
00110 temp1.reserve(3);
00111 temp1.push_back(0);
00112 temp1.push_back(-1);
00113 temp1.push_back(1);
00114 desc.add<std::vector<int> >("bunchCrossings",temp1);
00115 }
00116 desc.add<unsigned int>("daqPartitions",1);
00117 desc.add<bool>("ignoreL1Mask",false);
00118 desc.add<bool>("invert",false);
00119 desc.add<unsigned long long int>("physicsLoBits",1);
00120 desc.add<unsigned long long int>("physicsHiBits",262144);
00121 desc.add<unsigned long long int>("technicalBits",1);
00122 descriptions.add("hltLevel1Activity",desc);
00123 }
00124
00125
00126
00127
00128
00129
00130 bool
00131 HLTLevel1Activity::filter(edm::Event& event, const edm::EventSetup& setup)
00132 {
00133
00134
00135
00136
00137 if (not m_ignoreL1Mask and m_watchPhysicsMask.check(setup)) {
00138 edm::ESHandle<L1GtTriggerMask> h_mask;
00139 setup.get<L1GtTriggerMaskAlgoTrigRcd>().get(h_mask);
00140 const std::vector<unsigned int> & mask = h_mask->gtTriggerMask();
00141 for (unsigned int i = 0; i < PHYSICS_BITS_SIZE; ++i)
00142 m_maskedPhysics[i] = m_selectPhysics[i] and ((mask[i] & m_daqPartitions) != m_daqPartitions);
00143 }
00144
00145
00146
00147
00148
00149 if (not m_ignoreL1Mask and m_watchTechnicalMask.check(setup)) {
00150 edm::ESHandle<L1GtTriggerMask> h_mask;
00151 setup.get<L1GtTriggerMaskTechTrigRcd>().get(h_mask);
00152 const std::vector<unsigned int> & mask = h_mask->gtTriggerMask();
00153 for (unsigned int i = 0; i < TECHNICAL_BITS_SIZE; ++i)
00154 m_maskedTechnical[i] = m_selectTechnical[i] and ((mask[i] & m_daqPartitions) != m_daqPartitions);
00155 }
00156
00157
00158 edm::Handle<L1GlobalTriggerReadoutRecord> h_gtReadoutRecord;
00159 event.getByLabel(m_gtReadoutRecord, h_gtReadoutRecord);
00160
00161
00162 BOOST_FOREACH(int bx, m_bunchCrossings) {
00163 const std::vector<bool> & physics = h_gtReadoutRecord->decisionWord(bx);
00164 if (physics.size() != PHYSICS_BITS_SIZE)
00165
00166 return m_invert;
00167 for (unsigned int i = 0; i < PHYSICS_BITS_SIZE; ++i)
00168 if (m_maskedPhysics[i] and physics[i])
00169 return not m_invert;
00170 const std::vector<bool> & technical = h_gtReadoutRecord->technicalTriggerWord(bx);
00171 if (technical.size() != TECHNICAL_BITS_SIZE)
00172
00173 return m_invert;
00174 for (unsigned int i = 0; i < TECHNICAL_BITS_SIZE; ++i)
00175 if (m_maskedTechnical[i] and technical[i])
00176 return not m_invert;
00177 }
00178
00179 return m_invert;
00180 }
00181
00182
00183 #include "FWCore/Framework/interface/MakerMacros.h"
00184 DEFINE_FWK_MODULE(HLTLevel1Activity);