CMS 3D CMS Logo

TrigObjTnPSource.cc
Go to the documentation of this file.
1 
11 
12 #include "TrigObjTnPHistColl.h"
13 
14 class TrigObjTnPSource : public DQMGlobalEDAnalyzer<std::vector<TrigObjTnPHistColl>> {
15 public:
16  explicit TrigObjTnPSource(const edm::ParameterSet&);
17  ~TrigObjTnPSource() override = default;
18  TrigObjTnPSource(const TrigObjTnPSource&) = delete;
20 
21  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
22 
23  void dqmAnalyze(const edm::Event&, const edm::EventSetup&, const std::vector<TrigObjTnPHistColl>&) const override;
25  const edm::Run&,
26  const edm::EventSetup&,
27  std::vector<TrigObjTnPHistColl>&) const override;
28  void dqmBeginRun(const edm::Run&, const edm::EventSetup&, std::vector<TrigObjTnPHistColl>&) const override;
29 
30 private:
34  //it would be more memory efficient to save this as a vector of TrigTnPHistColls and then use this to instance
35  //the TrigTnPHistColls for each run, something to do for the future
36  std::vector<edm::ParameterSet> histCollConfigs_;
37 };
38 
40  : trigEvtToken_(consumes<trigger::TriggerEvent>(config.getParameter<edm::InputTag>("triggerEvent"))),
41  trigResultsToken_(consumes<edm::TriggerResults>(config.getParameter<edm::InputTag>("triggerResults"))),
42  hltProcess_(config.getParameter<edm::InputTag>("triggerResults").process()),
43  histCollConfigs_(config.getParameter<std::vector<edm::ParameterSet>>("histColls")) {}
44 
47  desc.add<edm::InputTag>("triggerEvent", edm::InputTag("hltTriggerSummaryAOD::HLT"));
48  desc.add<edm::InputTag>("triggerResults", edm::InputTag("TriggerResults::HLT"));
49  desc.addVPSet("histColls", TrigObjTnPHistColl::makePSetDescription(), std::vector<edm::ParameterSet>());
50  descriptions.add("trigObjTnPSource", desc);
51 }
52 
54  const edm::Run& run,
55  const edm::EventSetup& setup,
56  std::vector<TrigObjTnPHistColl>& tnpHistColls) const {
57  tnpHistColls.clear();
59  bool hltChanged = false;
60  hltConfig.init(run, setup, hltProcess_, hltChanged);
61  for (auto& histCollConfig : histCollConfigs_)
62  tnpHistColls.emplace_back(TrigObjTnPHistColl(histCollConfig));
63  for (auto& histColl : tnpHistColls) {
64  histColl.init(hltConfig);
65  histColl.bookHists(iBooker);
66  }
67 }
68 
70  const edm::EventSetup& setup,
71  std::vector<TrigObjTnPHistColl>& tnpHistColls) const {}
72 
74  const edm::EventSetup& setup,
75  const std::vector<TrigObjTnPHistColl>& tnpHistColls) const {
77  event.getByToken(trigEvtToken_, trigEvtHandle);
78  edm::Handle<edm::TriggerResults> trigResultsHandle;
79  event.getByToken(trigResultsToken_, trigResultsHandle);
80  //DQM should never crash the HLT under any circumstances (except at configure)
81  if (trigEvtHandle.isValid() && trigResultsHandle.isValid()) {
82  for (auto& histColl : tnpHistColls) {
83  histColl.fill(*trigEvtHandle, *trigResultsHandle, event.triggerNames(*trigResultsHandle));
84  }
85  }
86 }
87 
TrigObjTnPSource(const edm::ParameterSet &)
std::vector< edm::ParameterSet > histCollConfigs_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: config.py:1
TrigObjTnPSource & operator=(const TrigObjTnPSource &)=delete
static edm::ParameterSetDescription makePSetDescription()
void dqmAnalyze(const edm::Event &, const edm::EventSetup &, const std::vector< TrigObjTnPHistColl > &) const override
void bookHistograms(DQMStore::IBooker &, const edm::Run &, const edm::EventSetup &, std::vector< TrigObjTnPHistColl > &) const override
edm::EDGetTokenT< edm::TriggerResults > trigResultsToken_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
~TrigObjTnPSource() override=default
void dqmBeginRun(const edm::Run &, const edm::EventSetup &, std::vector< TrigObjTnPHistColl > &) const override
void add(std::string const &label, ParameterSetDescription const &psetDescription)
bool isValid() const
Definition: HandleBase.h:70
std::string hltProcess_
HLT enums.
edm::EDGetTokenT< trigger::TriggerEvent > trigEvtToken_
Definition: event.py:1
Definition: Run.h:45