CMS 3D CMS Logo

TriggerRatesMonitorClient.cc
Go to the documentation of this file.
2 
4 
5 //
6 // -------------------------------------- Constructor --------------------------------------------
7 //
9  : m_dqm_path(iConfig.getUntrackedParameter<std::string>("dqmPath")) {
10  edm::LogInfo("TriggerRatesMonitorClient")
11  << "Constructor TriggerRatesMonitorClient::TriggerRatesMonitorClient " << std::endl;
12 }
13 
14 //
15 // -------------------------------------- beginJob --------------------------------------------
16 //
18  edm::LogInfo("TriggerRatesMonitorClient") << "TriggerRatesMonitorClient::beginJob " << std::endl;
19 }
20 
21 //
22 // -------------------------------------- get and book in the endJob --------------------------------------------
23 //
25  // create and cd into new folder
26  ibooker_.setCurrentFolder(m_dqm_path);
27 
28  //get available histograms
29  std::vector<std::string> directories = igetter_.getSubdirs();
30  m_hltXpd_counts.resize(directories.size());
31 
32  int i = 0;
33  for (auto const& dir : directories) {
34  // std::cout << "dir: " << dir << std::endl;
35  ibooker_.setCurrentFolder(m_dqm_path + "/" + dir);
36 
37  std::vector<std::string> const& all_mes = igetter_.getMEs();
38  std::vector<std::string> mes;
39  for (auto const& me : all_mes)
40  if (me.find("accept") != std::string::npos)
41  mes.push_back(me);
42 
43  int nbinsY = mes.size();
44  double xminY = 0.;
45  double xmaxY = xminY + 1. * nbinsY;
46  int nbinsX = 0;
47  int ibinY = 1;
48  for (auto const& me : mes) {
49  // std::cout << "me: " << me << std::endl;
50  ibooker_.setCurrentFolder(m_dqm_path + "/" + dir);
51  TH1F* histo = igetter_.get(me)->getTH1F();
52 
53  if (m_hltXpd_counts[i] == nullptr) {
54  // get range and binning for new MEs
55  nbinsX = histo->GetNbinsX();
56  double xminX = histo->GetXaxis()->GetXmin();
57  double xmaxX = histo->GetXaxis()->GetXmax();
58 
59  //book new histogram
60  std::string hname = dir + "_summary";
61  ibooker_.setCurrentFolder(m_dqm_path);
62  m_hltXpd_counts[i] = ibooker_.book2D(hname, hname, nbinsX, xminX, xmaxX, nbinsY, xminY, xmaxY)->getTH2F();
63  } else {
64  m_hltXpd_counts[i]->GetYaxis()->SetBinLabel(ibinY, me.c_str());
65  }
66 
67  // handle mes
68  for (int ibinX = 1; ibinX <= nbinsX; ++ibinX) {
69  float rate = histo->GetBinContent(ibinX);
70  m_hltXpd_counts[i]->SetBinContent(ibinX, ibinY, rate);
71  }
72  ibinY++;
73  }
74 
75  i++;
76  }
77 }
78 
79 //
80 // -------------------------------------- get in the endLumi if needed --------------------------------------------
81 //
83  DQMStore::IGetter& igetter_,
84  edm::LuminosityBlock const& iLumi,
85  edm::EventSetup const& iSetup) {
86  edm::LogInfo("TriggerRatesMonitorClient") << "TriggerRatesMonitorClient::endLumi " << std::endl;
87 }
88 
91  desc.addUntracked<std::string>("dqmPath", "HLT/Datasets");
92  descriptions.add("triggerRatesMonitorClient", desc);
93 }
94 
95 // Define this as a plug-in
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
virtual std::vector< std::string > getMEs() const
Definition: DQMStore.cc:759
void dqmEndLuminosityBlock(DQMStore::IBooker &, DQMStore::IGetter &, edm::LuminosityBlock const &, edm::EventSetup const &) override
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
virtual TH2F * getTH2F() const
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
Log< level::Info, false > LogInfo
TriggerRatesMonitorClient(const edm::ParameterSet &ps)
virtual TH1F * getTH1F() const
MonitorElement * book2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:221
void add(std::string const &label, ParameterSetDescription const &psetDescription)
virtual MonitorElement * get(std::string const &fullpath) const
Definition: DQMStore.cc:712
std::vector< TH2F * > m_hltXpd_counts
double rate(double x)
Definition: Constants.cc:3
void dqmEndJob(DQMStore::IBooker &, DQMStore::IGetter &) override
virtual DQM_DEPRECATED std::vector< std::string > getSubdirs() const
Definition: DQMStore.cc:739