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 {
11  edm::LogInfo("TriggerRatesMonitorClient") << "Constructor TriggerRatesMonitorClient::TriggerRatesMonitorClient " << std::endl;
12 
13 }
14 
15 //
16 // -------------------------------------- beginJob --------------------------------------------
17 //
19 {
20  edm::LogInfo("TriggerRatesMonitorClient") << "TriggerRatesMonitorClient::beginJob " << std::endl;
21 }
22 
23 //
24 // -------------------------------------- get and book in the endJob --------------------------------------------
25 //
27 {
28  // create and cd into new folder
29  ibooker_.setCurrentFolder(m_dqm_path);
30 
31  //get available histograms
32  std::vector<std::string> directories = igetter_.getSubdirs();
33  m_hltXpd_counts.resize( directories.size() );
34 
35  int i = 0;
36  for ( auto const & dir : directories) {
37  // std::cout << "dir: " << dir << std::endl;
38  ibooker_.setCurrentFolder(m_dqm_path + "/" + dir);
39 
40  std::vector<std::string> const & all_mes = igetter_.getMEs();
41  std::vector<std::string> mes;
42  for ( auto const & me : all_mes )
43  if ( me.find("accept") != std::string::npos ) mes.push_back(me);
44 
45  int nbinsY = mes.size();
46  double xminY = 0.;
47  double xmaxY = xminY + 1.*nbinsY;
48  int nbinsX = 0;
49  int ibinY = 1;
50  for ( auto const & me : mes ) {
51  // std::cout << "me: " << me << std::endl;
52  ibooker_.setCurrentFolder(m_dqm_path + "/" + dir);
53  TH1F* histo = igetter_.get( me )->getTH1F();
54 
55  if (m_hltXpd_counts[i] == nullptr ) {
56  // get range and binning for new MEs
57  nbinsX = histo->GetNbinsX();
58  double xminX = histo->GetXaxis()->GetXmin();
59  double xmaxX = histo->GetXaxis()->GetXmax();
60 
61  //book new histogram
62  std::string hname = dir + "_summary";
63  ibooker_.setCurrentFolder(m_dqm_path);
64  m_hltXpd_counts[i] = ibooker_.book2D(hname,hname,nbinsX,xminX,xmaxX,nbinsY,xminY,xmaxY)->getTH2F();
65  } else {
66  m_hltXpd_counts[i]->GetYaxis()->SetBinLabel(ibinY,me.c_str());
67  }
68 
69  // handle mes
70  for ( int ibinX=1; ibinX <= nbinsX; ++ibinX ) {
71  float rate = histo->GetBinContent(ibinX);
72  m_hltXpd_counts[i]->SetBinContent(ibinX,ibinY,rate);
73  }
74  ibinY++;
75  }
76 
77  i++;
78  }
79 
80 }
81 
82 //
83 // -------------------------------------- get in the endLumi if needed --------------------------------------------
84 //
86 {
87  edm::LogInfo("TriggerRatesMonitorClient") << "TriggerRatesMonitorClient::endLumi " << std::endl;
88 }
89 
91 {
92 
94  desc.addUntracked<std::string>("dqmPath","HLT/Datasets");
95  descriptions.add("triggerRatesMonitorClient", desc);
96 
97 }
98 
99 // Define this as a plug-in
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
TH1F * getTH1F() const
MonitorElement * get(const std::string &path)
Definition: DQMStore.cc:302
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void dqmEndLuminosityBlock(DQMStore::IBooker &, DQMStore::IGetter &, edm::LuminosityBlock const &, edm::EventSetup const &) override
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
TH2F * getTH2F() const
TriggerRatesMonitorClient(const edm::ParameterSet &ps)
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:136
void add(std::string const &label, ParameterSetDescription const &psetDescription)
std::vector< TH2F * > m_hltXpd_counts
double rate(double x)
Definition: Constants.cc:3
std::vector< std::string > getSubdirs()
Definition: DQMStore.cc:320
dbl *** dir
Definition: mlp_gen.cc:35
void dqmEndJob(DQMStore::IBooker &, DQMStore::IGetter &) override
std::vector< std::string > getMEs()
Definition: DQMStore.cc:324