CMS 3D CMS Logo

HLTPhysicsDeclared.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: HLTPhysicsDeclared
4 // Class: HLTPhysicsDeclared
5 //
6 // Original Author: Luca Malgeri
7 // Adapted for HLT by: Andrea Bocci
8 
9 // user include files
17 
20 //
21 // class declaration
22 //
23 
25 public:
26  explicit HLTPhysicsDeclared(const edm::ParameterSet&);
27  ~HLTPhysicsDeclared() override;
28  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
29 
30 private:
31  bool filter(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
32 
33  bool m_invert;
36 };
37 
38 // system include files
39 #include <memory>
40 
41 // user include files
43 
44 using namespace edm;
45 
47  : m_invert(config.getParameter<bool>("invert")),
48  m_gtDigis(config.getParameter<edm::InputTag>("L1GtReadoutRecordTag")) {
49  m_gtDigisToken = consumes<L1GlobalTriggerReadoutRecord>(m_gtDigis);
50 }
51 
53 
56  desc.add<edm::InputTag>("L1GtReadoutRecordTag", edm::InputTag("hltGtDigis"));
57  desc.add<bool>("invert", false);
58  descriptions.add("hltPhysicsDeclared", desc);
59 }
60 
62  bool accept = false;
63 
64  if (event.isRealData()) {
65  // for real data, access the "physics enabled" bit in the L1 GT data
67  if (not event.getByToken(m_gtDigisToken, h_gtDigis)) {
68  edm::LogWarning(h_gtDigis.whyFailed()->category()) << h_gtDigis.whyFailed()->what();
69  // not enough informations to make a decision - reject the event
70  return false;
71  } else {
72  L1GtFdlWord fdlWord = h_gtDigis->gtFdlWord();
73  if (fdlWord.physicsDeclared() == 1)
74  accept = true;
75  }
76  } else {
77  // for MC, assume the "physics enabled" bit to be always set
78  accept = true;
79  }
80 
81  // if requested, invert the filter decision
82  if (m_invert)
83  accept = not accept;
84 
85  return accept;
86 }
87 
88 // declare this class as a framework plugin
~HLTPhysicsDeclared() override
const cms_uint16_t physicsDeclared() const
get/set "physics declared" bit
Definition: L1GtFdlWord.h:169
Definition: config.py:1
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:31
const L1GtFdlWord gtFdlWord(int bxInEventValue) const
get / set FDL word (record) in the GT readout record
HLTPhysicsDeclared(const edm::ParameterSet &)
std::shared_ptr< cms::Exception > whyFailed() const
Definition: HandleBase.h:91
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
HLT enums.
bool filter(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
Log< level::Warning, false > LogWarning
edm::EDGetTokenT< L1GlobalTriggerReadoutRecord > m_gtDigisToken
Definition: event.py:1