Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 11:10:02

0001 #include "FWCore/Framework/interface/Frameworkfwd.h"
0002 #include "FWCore/Framework/interface/Event.h"
0003 #include "FWCore/Framework/interface/MakerMacros.h"
0004 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0005 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0006 #include "FWCore/Utilities/interface/InputTag.h"
0007 #include "DQMServices/Core/interface/DQMEDHarvester.h"
0008 #include "DQMServices/Core/interface/DQMStore.h"
0009 
0010 #include "TH1F.h"
0011 #include "TH2F.h"
0012 
0013 using namespace edm;
0014 using namespace std;
0015 
0016 class PlotCombiner : public DQMEDHarvester {
0017 public:
0018   PlotCombiner(const edm::ParameterSet &pset);
0019   ~PlotCombiner() override;
0020 
0021 private:
0022   void makePlot(const ParameterSet &pset, DQMStore::IBooker &, DQMStore::IGetter &);
0023 
0024   string myDQMrootFolder;
0025   const VParameterSet plots;
0026 
0027 protected:
0028   void dqmEndJob(DQMStore::IBooker &, DQMStore::IGetter &) override;  //performed in the endJob
0029 };
0030 
0031 PlotCombiner::PlotCombiner(const edm::ParameterSet &pset)
0032     : myDQMrootFolder(pset.getUntrackedParameter<string>("MyDQMrootFolder")),
0033       plots(pset.getUntrackedParameter<VParameterSet>("Plots")) {}
0034 
0035 void PlotCombiner::dqmEndJob(DQMStore::IBooker &ibooker_, DQMStore::IGetter &igetter_) {
0036   for (VParameterSet::const_iterator pset = plots.begin(); pset != plots.end(); pset++) {
0037     makePlot(*pset, ibooker_, igetter_);
0038   }
0039 }
0040 
0041 void PlotCombiner::makePlot(const ParameterSet &pset, DQMStore::IBooker &ibooker_, DQMStore::IGetter &igetter_) {
0042   //get hold of MEs
0043   vector<string> inputMEnames = pset.getUntrackedParameter<vector<string> >("InputMEnames");
0044   vector<string> inputLabels = pset.getUntrackedParameter<vector<string> >("InputLabels");
0045   if (inputMEnames.size() != inputLabels.size()) {
0046     LogDebug("HLTriggerOfflineHeavyFlavor") << "Number of labels must match the histos[0]ber of InputMEnames" << endl;
0047     return;
0048   }
0049   vector<TH1 *> histos;
0050   vector<TString> labels;
0051   for (size_t i = 0; i < inputMEnames.size(); i++) {
0052     string MEname = myDQMrootFolder + "/" + inputMEnames[i];
0053     MonitorElement *ME = igetter_.get(MEname);
0054     if (ME == nullptr) {
0055       LogDebug("HLTriggerOfflineHeavyFlavor") << "Could not find ME: " << MEname << endl;
0056       continue;
0057     }
0058     histos.push_back(ME->getTH1());
0059     labels.push_back(inputLabels[i]);
0060   }
0061   if (histos.empty()) {
0062     return;
0063   }
0064   //figure out the output directory name
0065   string outputMEname = pset.getUntrackedParameter<string>("OutputMEname");
0066   ;
0067   string outputDir = myDQMrootFolder;
0068   string::size_type slashPos = outputMEname.rfind('/');
0069   if (string::npos != slashPos) {
0070     outputDir += "/" + outputMEname.substr(0, slashPos);
0071     outputMEname.erase(0, slashPos + 1);
0072   }
0073   ibooker_.setCurrentFolder(outputDir);
0074   //create output ME
0075   TH2F *output;
0076   if (histos[0]->GetXaxis()->GetXbins()->GetSize() == 0) {
0077     output = new TH2F(outputMEname.c_str(),
0078                       outputMEname.c_str(),
0079                       histos[0]->GetXaxis()->GetNbins(),
0080                       histos[0]->GetXaxis()->GetXmin(),
0081                       histos[0]->GetXaxis()->GetXmax(),
0082                       histos.size(),
0083                       0,
0084                       histos.size());
0085   } else {
0086     output = new TH2F(outputMEname.c_str(),
0087                       outputMEname.c_str(),
0088                       histos[0]->GetXaxis()->GetNbins(),
0089                       histos[0]->GetXaxis()->GetXbins()->GetArray(),
0090                       histos.size(),
0091                       0,
0092                       histos.size());
0093   }
0094   output->SetTitle(outputMEname.c_str());
0095   output->SetXTitle(histos[0]->GetXaxis()->GetTitle());
0096   output->SetStats(kFALSE);
0097   output->SetOption("colztexte");
0098   for (size_t i = 0; i < histos.size(); i++) {
0099     for (int j = 1; j <= histos[0]->GetNbinsX(); j++) {
0100       output->SetBinContent(j, i + 1, histos[i]->GetBinContent(j));
0101       output->SetBinError(j, i + 1, histos[i]->GetBinError(j));
0102     }
0103     output->GetYaxis()->SetBinLabel(i + 1, labels[i]);
0104   }
0105   ibooker_.book2D(outputMEname, output);
0106   delete output;
0107 }
0108 
0109 PlotCombiner::~PlotCombiner() {}
0110 
0111 //define this as a plug-in
0112 DEFINE_FWK_MODULE(PlotCombiner);