CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
PlotCombiner.cc
Go to the documentation of this file.
11 
12 #include "TH1F.h"
13 #include "TH2F.h"
14 
15 using namespace edm;
16 using namespace std;
17 
18 class PlotCombiner : public DQMEDHarvester {
19 public:
21  ~PlotCombiner() override;
22 
23 private:
25 
28 
29 protected:
30  void dqmEndJob(DQMStore::IBooker &, DQMStore::IGetter &) override; //performed in the endJob
31 };
32 
34  : myDQMrootFolder(pset.getUntrackedParameter<string>("MyDQMrootFolder")),
35  plots(pset.getUntrackedParameter<VParameterSet>("Plots")) {}
36 
38  for (VParameterSet::const_iterator pset = plots.begin(); pset != plots.end(); pset++) {
39  makePlot(*pset, ibooker_, igetter_);
40  }
41 }
42 
44  //get hold of MEs
45  vector<string> inputMEnames = pset.getUntrackedParameter<vector<string> >("InputMEnames");
46  vector<string> inputLabels = pset.getUntrackedParameter<vector<string> >("InputLabels");
47  if (inputMEnames.size() != inputLabels.size()) {
48  LogDebug("HLTriggerOfflineHeavyFlavor") << "Number of labels must match the histos[0]ber of InputMEnames" << endl;
49  return;
50  }
51  vector<TH1 *> histos;
52  vector<TString> labels;
53  for (size_t i = 0; i < inputMEnames.size(); i++) {
54  string MEname = myDQMrootFolder + "/" + inputMEnames[i];
55  MonitorElement *ME = igetter_.get(MEname);
56  if (ME == nullptr) {
57  LogDebug("HLTriggerOfflineHeavyFlavor") << "Could not find ME: " << MEname << endl;
58  continue;
59  }
60  histos.push_back(ME->getTH1());
61  labels.push_back(inputLabels[i]);
62  }
63  if (histos.empty()) {
64  return;
65  }
66  //figure out the output directory name
67  string outputMEname = pset.getUntrackedParameter<string>("OutputMEname");
68  ;
69  string outputDir = myDQMrootFolder;
70  string::size_type slashPos = outputMEname.rfind('/');
71  if (string::npos != slashPos) {
72  outputDir += "/" + outputMEname.substr(0, slashPos);
73  outputMEname.erase(0, slashPos + 1);
74  }
75  ibooker_.setCurrentFolder(outputDir);
76  //create output ME
77  TH2F *output;
78  if (histos[0]->GetXaxis()->GetXbins()->GetSize() == 0) {
79  output = new TH2F(outputMEname.c_str(),
80  outputMEname.c_str(),
81  histos[0]->GetXaxis()->GetNbins(),
82  histos[0]->GetXaxis()->GetXmin(),
83  histos[0]->GetXaxis()->GetXmax(),
84  histos.size(),
85  0,
86  histos.size());
87  } else {
88  output = new TH2F(outputMEname.c_str(),
89  outputMEname.c_str(),
90  histos[0]->GetXaxis()->GetNbins(),
91  histos[0]->GetXaxis()->GetXbins()->GetArray(),
92  histos.size(),
93  0,
94  histos.size());
95  }
96  output->SetTitle(outputMEname.c_str());
97  output->SetXTitle(histos[0]->GetXaxis()->GetTitle());
98  output->SetStats(kFALSE);
99  output->SetOption("colztexte");
100  for (size_t i = 0; i < histos.size(); i++) {
101  for (int j = 1; j <= histos[0]->GetNbinsX(); j++) {
102  output->SetBinContent(j, i + 1, histos[i]->GetBinContent(j));
103  output->SetBinError(j, i + 1, histos[i]->GetBinError(j));
104  }
105  output->GetYaxis()->SetBinLabel(i + 1, labels[i]);
106  }
107  ibooker_.book2D(outputMEname, output);
108  delete output;
109 }
110 
112 
113 //define this as a plug-in
T getUntrackedParameter(std::string const &, T const &) const
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
string myDQMrootFolder
Definition: PlotCombiner.cc:26
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::vector< ParameterSet > VParameterSet
Definition: ParameterSet.h:34
~PlotCombiner() override
uint16_t size_type
Definition: ME.h:11
virtual MonitorElement * get(std::string const &fullpath) const
Definition: DQMStore.cc:673
void dqmEndJob(DQMStore::IBooker &, DQMStore::IGetter &) override
Definition: PlotCombiner.cc:37
const VParameterSet plots
Definition: PlotCombiner.cc:27
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:212
void makePlot(const ParameterSet &pset, DQMStore::IBooker &, DQMStore::IGetter &)
Definition: PlotCombiner.cc:43
virtual TH1 * getTH1() const
PlotCombiner(const edm::ParameterSet &pset)
Definition: PlotCombiner.cc:33
#define LogDebug(id)