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