Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:09:35

0001 #include "DQMOffline/L1Trigger/interface/L1TDiffHarvesting.h"
0002 
0003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0004 #include "FWCore/Framework/interface/MakerMacros.h"
0005 
0006 namespace dqmoffline {
0007   namespace l1t {
0008 
0009     L1TDiffHarvesting::L1TDiffPlotHandler::L1TDiffPlotHandler(const edm::ParameterSet &ps, std::string plotName)
0010         : dir1_(ps.getUntrackedParameter<std::string>("dir1")),
0011           dir2_(ps.getUntrackedParameter<std::string>("dir2")),
0012           outputDir_(ps.getUntrackedParameter<std::string>("outputDir", dir1_)),
0013           plotName_(plotName),
0014           h1_(),
0015           h2_(),
0016           h_diff_(),
0017           histType1_(),
0018           histType2_() {}
0019 
0020     L1TDiffHarvesting::L1TDiffPlotHandler::L1TDiffPlotHandler(const L1TDiffHarvesting::L1TDiffPlotHandler &handler)
0021         : dir1_(handler.dir1_),
0022           dir2_(handler.dir2_),
0023           outputDir_(handler.outputDir_),
0024           plotName_(handler.plotName_),
0025           h1_(),
0026           h2_(),
0027           h_diff_(),
0028           histType1_(),
0029           histType2_() {}
0030 
0031     void L1TDiffHarvesting::L1TDiffPlotHandler::computeDiff(DQMStore::IBooker &ibooker, DQMStore::IGetter &igetter) {
0032       loadHistograms(igetter);
0033       if (!isValid()) {
0034         return;
0035       }
0036       bookDiff(ibooker);
0037 
0038       TH1 *h_diff;
0039       TH1 *h1;
0040       TH1 *h2;
0041       bool is1D(histType1_ == MonitorElement::Kind::TH1F || histType1_ == MonitorElement::Kind::TH1D);
0042       bool is2D(histType1_ == MonitorElement::Kind::TH2F || histType1_ == MonitorElement::Kind::TH2D);
0043       bool isProfile(histType1_ == MonitorElement::Kind::TPROFILE);
0044 
0045       if (is1D) {
0046         h_diff = h_diff_->getTH1F();
0047         h1 = h1_->getTH1F();
0048         h2 = h2_->getTH1F();
0049       } else if (is2D) {
0050         h_diff = h_diff_->getTH2F();
0051         h1 = h1_->getTH2F();
0052         h2 = h2_->getTH2F();
0053       } else if (isProfile) {
0054         h_diff = h_diff_->getTProfile();
0055         h1 = h1_->getTProfile();
0056         h2 = h2_->getTProfile();
0057       } else {
0058         edm::LogWarning("L1TDiffHarvesting::L1TDiffPlotHandler::computeDiff")
0059             << "Unknown histogram type. Quitting booking" << std::endl;
0060 
0061         return;
0062       }
0063       h_diff->Add(h1);
0064       h_diff->Add(h2, -1);
0065       // if histograms are identical h_diff will have 0 entries -> not good to check if anything happened
0066       // let's fix it
0067       h_diff->SetEntries(h1->GetEntries() + h2->GetEntries());
0068     }
0069 
0070     void L1TDiffHarvesting::L1TDiffPlotHandler::loadHistograms(DQMStore::IGetter &igetter) {
0071       std::string h1Name = dir1_ + "/" + plotName_;
0072       std::string h2Name = dir2_ + "/" + plotName_;
0073       h1_ = igetter.get(h1Name);
0074       h2_ = igetter.get(h2Name);
0075 
0076       if (!h1_ || !h2_) {
0077         edm::LogWarning("L1TDiffHarvesting::L1TDiffPlotHandler::loadHistograms")
0078             << (!h1_ && !h2_ ? h1Name + " && " + h2Name
0079                 : !h1_       ? h1Name
0080                              : h2Name)
0081             << " not gettable. Quitting booking" << std::endl;
0082 
0083         return;
0084       }
0085 
0086       histType1_ = h1_->kind();
0087       histType2_ = h2_->kind();
0088     }
0089 
0090     bool L1TDiffHarvesting::L1TDiffPlotHandler::isValid() const {
0091       if (histType1_ == MonitorElement::Kind::INVALID) {
0092         edm::LogWarning("L1TDiffHarvesting::L1TDiffPlotHandler::isValid")
0093             << " Could not find a supported histogram type" << std::endl;
0094         return false;
0095       }
0096       if (histType1_ != histType2_) {
0097         edm::LogWarning("L1TDiffHarvesting::L1TDiffPlotHandler::isValid")
0098             << " Histogram 1 and 2 have different histogram types" << std::endl;
0099         return false;
0100       }
0101       return true;
0102     }
0103 
0104     void L1TDiffHarvesting::L1TDiffPlotHandler::bookDiff(DQMStore::IBooker &ibooker) {
0105       ibooker.setCurrentFolder(outputDir_);
0106 
0107       bool is1D(histType1_ == MonitorElement::Kind::TH1F || histType1_ == MonitorElement::Kind::TH1D);
0108       bool is2D(histType1_ == MonitorElement::Kind::TH2F || histType1_ == MonitorElement::Kind::TH2D);
0109       bool isProfile(histType1_ == MonitorElement::Kind::TPROFILE);
0110 
0111       if (is1D) {
0112         TH1F *h1 = h1_->getTH1F();
0113         double min = h1->GetXaxis()->GetXmin();
0114         double max = h1->GetXaxis()->GetXmax();
0115         int nBins = h1->GetNbinsX();
0116         h_diff_ = ibooker.book1D(plotName_, plotName_, nBins, min, max);
0117       } else if (is2D) {
0118         TH2F *h1 = h1_->getTH2F();
0119         double minX = h1->GetXaxis()->GetXmin();
0120         double maxX = h1->GetXaxis()->GetXmax();
0121         double minY = h1->GetYaxis()->GetXmin();
0122         double maxY = h1->GetYaxis()->GetXmax();
0123         int nBinsX = h1->GetNbinsX();
0124         int nBinsY = h1->GetNbinsY();
0125 
0126         h_diff_ = ibooker.book2D(plotName_, plotName_, nBinsX, minX, maxX, nBinsY, minY, maxY);
0127       } else if (isProfile) {
0128         TProfile *h1 = h1_->getTProfile();
0129         double minX = h1->GetXaxis()->GetXmin();
0130         double maxX = h1->GetXaxis()->GetXmax();
0131         double minY = h1->GetYaxis()->GetXmin();
0132         double maxY = h1->GetYaxis()->GetXmax();
0133         int nBins = h1->GetNbinsX();
0134         h_diff_ = ibooker.bookProfile(plotName_, plotName_, nBins, minX, maxX, minY, maxY);
0135       } else {
0136         edm::LogWarning("L1TDiffHarvesting::L1TDiffPlotHandler::bookDiff")
0137             << "Unknown histogram type. Quitting booking" << std::endl;
0138 
0139         return;
0140       }
0141     }
0142 
0143     L1TDiffHarvesting::L1TDiffHarvesting(const edm::ParameterSet &ps) : plotHandlers_() {
0144       using namespace std;
0145       for (const auto &plotConfig : ps.getUntrackedParameter<std::vector<edm::ParameterSet>>("plotCfgs")) {
0146         vector<string> plots = plotConfig.getUntrackedParameter<vector<string>>("plots");
0147         for (const auto &plot : plots) {
0148           plotHandlers_.push_back(L1TDiffPlotHandler(plotConfig, plot));
0149         }
0150       }
0151     }
0152 
0153     L1TDiffHarvesting::~L1TDiffHarvesting() {}
0154 
0155     void L1TDiffHarvesting::dqmEndJob(DQMStore::IBooker &ibooker, DQMStore::IGetter &igetter) {
0156       edm::LogInfo("L1TEfficiencyHarvesting") << "Called endRun." << std::endl;
0157 
0158       for (auto plotHandler : plotHandlers_) {
0159         plotHandler.computeDiff(ibooker, igetter);
0160       }
0161     }
0162 
0163     DEFINE_FWK_MODULE(L1TDiffHarvesting);
0164 
0165   }  // namespace l1t
0166 }  // namespace dqmoffline