Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include "DQMOffline/RecoB/plugins/BTagPerformanceHarvester.h"
0002 #include "FWCore/Utilities/interface/InputTag.h"
0003 #include "DQMOffline/RecoB/interface/Tools.h"
0004 #include "DQMOffline/RecoB/interface/TagInfoPlotterFactory.h"
0005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0006 #include "FWCore/Framework/interface/MakerMacros.h"
0007 
0008 using namespace edm;
0009 using namespace std;
0010 using namespace RecoBTag;
0011 
0012 BTagPerformanceHarvester::BTagPerformanceHarvester(const edm::ParameterSet& pSet)
0013     : etaRanges(pSet.getParameter<vector<double>>("etaRanges")),
0014       ptRanges(pSet.getParameter<vector<double>>("ptRanges")),
0015       produceEps(pSet.getParameter<bool>("produceEps")),
0016       producePs(pSet.getParameter<bool>("producePs")),
0017       moduleConfig(pSet.getParameter<vector<edm::ParameterSet>>("tagConfig")),
0018       flavPlots_(pSet.getParameter<std::string>("flavPlots")),
0019       makeDiffPlots_(pSet.getParameter<bool>("differentialPlots")) {
0020   //mcPlots_ : 1=b+c+l+ni; 2=all+1; 3=1+d+u+s+g; 4=3+all . Default is 2. Don't use 0.
0021   if (flavPlots_.find("dusg") < 15) {
0022     if (flavPlots_.find("all") < 15)
0023       mcPlots_ = 4;
0024     else
0025       mcPlots_ = 3;
0026   } else if (flavPlots_.find("bcl") < 15) {
0027     if (flavPlots_.find("all") < 15)
0028       mcPlots_ = 2;
0029     else
0030       mcPlots_ = 1;
0031   } else
0032     mcPlots_ = 0;
0033 
0034   if (etaRanges.size() <= 1)
0035     etaRanges = {pSet.getParameter<double>("etaMin"), pSet.getParameter<double>("etaMax")};
0036   if (ptRanges.size() <= 1)
0037     ptRanges = {pSet.getParameter<double>("ptRecJetMin"), pSet.getParameter<double>("ptRecJetMax")};
0038 
0039   for (vector<edm::ParameterSet>::const_iterator iModule = moduleConfig.begin(); iModule != moduleConfig.end();
0040        ++iModule) {
0041     const string& dataFormatType = iModule->exists("type") ? iModule->getParameter<string>("type") : "JetTag";
0042     if (dataFormatType == "JetTag") {
0043       const InputTag& moduleLabel = iModule->getParameter<InputTag>("label");
0044       jetTagInputTags.push_back(moduleLabel);
0045       binJetTagPlotters.push_back(vector<std::shared_ptr<JetTagPlotter>>());
0046       if (mcPlots_ && makeDiffPlots_) {
0047         differentialPlots.push_back(vector<std::unique_ptr<BTagDifferentialPlot>>());
0048       }
0049     } else if (dataFormatType == "TagCorrelation") {
0050       const InputTag& label1 = iModule->getParameter<InputTag>("label1");
0051       const InputTag& label2 = iModule->getParameter<InputTag>("label2");
0052       tagCorrelationInputTags.push_back(std::pair<edm::InputTag, edm::InputTag>(label1, label2));
0053       binTagCorrelationPlotters.push_back(vector<std::unique_ptr<TagCorrelationPlotter>>());
0054     } else {
0055       tagInfoInputTags.push_back(vector<edm::InputTag>());
0056       tiDataFormatType.push_back(dataFormatType);
0057       binTagInfoPlotters.push_back(vector<std::unique_ptr<BaseTagInfoPlotter>>());
0058     }
0059   }
0060 }
0061 
0062 EtaPtBin BTagPerformanceHarvester::getEtaPtBin(const int& iEta, const int& iPt) {
0063   // DEFINE BTagBin:
0064   bool etaActive_, ptActive_;
0065   double etaMin_, etaMax_, ptMin_, ptMax_;
0066 
0067   if (iEta != -1) {
0068     etaActive_ = true;
0069     etaMin_ = etaRanges[iEta];
0070     etaMax_ = etaRanges[iEta + 1];
0071   } else {
0072     etaActive_ = false;
0073     etaMin_ = etaRanges[0];
0074     etaMax_ = etaRanges[etaRanges.size() - 1];
0075   }
0076 
0077   if (iPt != -1) {
0078     ptActive_ = true;
0079     ptMin_ = ptRanges[iPt];
0080     ptMax_ = ptRanges[iPt + 1];
0081   } else {
0082     ptActive_ = false;
0083     ptMin_ = ptRanges[0];
0084     ptMax_ = ptRanges[ptRanges.size() - 1];
0085   }
0086   return EtaPtBin(etaActive_, etaMin_, etaMax_, ptActive_, ptMin_, ptMax_);
0087 }
0088 
0089 BTagPerformanceHarvester::~BTagPerformanceHarvester() {}
0090 
0091 void BTagPerformanceHarvester::dqmEndJob(DQMStore::IBooker& ibook, DQMStore::IGetter& iget) {
0092   // Book all histograms.
0093 
0094   // iterate over ranges:
0095   const int iEtaStart = -1;  // this will be the inactive one
0096   const int iEtaEnd = etaRanges.size() > 2 ? etaRanges.size() - 1
0097                                            : 0;  // if there is only one bin defined, leave it as the inactive one
0098   const int iPtStart = -1;                       // this will be the inactive one
0099   const int iPtEnd =
0100       ptRanges.size() > 2 ? ptRanges.size() - 1 : 0;  // if there is only one bin defined, leave it as the inactive one
0101   setTDRStyle();
0102 
0103   TagInfoPlotterFactory theFactory;
0104   int iTag = -1;
0105   int iTagCorr = -1;
0106   int iInfoTag = -1;
0107   for (vector<edm::ParameterSet>::const_iterator iModule = moduleConfig.begin(); iModule != moduleConfig.end();
0108        ++iModule) {
0109     const string& dataFormatType = iModule->exists("type") ? iModule->getParameter<string>("type") : "JetTag";
0110     const bool& doCTagPlots = iModule->exists("doCTagPlots") ? iModule->getParameter<bool>("doCTagPlots") : false;
0111 
0112     if (dataFormatType == "JetTag") {
0113       iTag++;
0114       const string& folderName = iModule->getParameter<string>("folder");
0115 
0116       // Contains plots for each bin of rapidity and pt.
0117       auto differentialPlotsConstantEta = std::make_unique<std::vector<std::unique_ptr<BTagDifferentialPlot>>>();
0118       auto differentialPlotsConstantPt = std::make_unique<std::vector<std::unique_ptr<BTagDifferentialPlot>>>();
0119       if (mcPlots_ && makeDiffPlots_) {
0120         // the constant b-efficiency for the differential plots versus pt and eta
0121         const double& effBConst =
0122             iModule->getParameter<edm::ParameterSet>("parameters").getParameter<double>("effBConst");
0123 
0124         // the objects for the differential plots vs. eta,pt for
0125         for (int iEta = iEtaStart; iEta < iEtaEnd; iEta++) {
0126           std::unique_ptr<BTagDifferentialPlot> etaConstDifferentialPlot =
0127               std::make_unique<BTagDifferentialPlot>(effBConst, BTagDifferentialPlot::constETA, folderName, mcPlots_);
0128           differentialPlotsConstantEta->push_back(std::move(etaConstDifferentialPlot));
0129         }
0130 
0131         for (int iPt = iPtStart; iPt < iPtEnd; iPt++) {
0132           // differentialPlots for this pt bin
0133           std::unique_ptr<BTagDifferentialPlot> ptConstDifferentialPlot =
0134               std::make_unique<BTagDifferentialPlot>(effBConst, BTagDifferentialPlot::constPT, folderName, mcPlots_);
0135           differentialPlotsConstantPt->push_back(std::move(ptConstDifferentialPlot));
0136         }
0137       }
0138       // eta loop
0139       for (int iEta = iEtaStart; iEta < iEtaEnd; iEta++) {
0140         // pt loop
0141         for (int iPt = iPtStart; iPt < iPtEnd; iPt++) {
0142           const EtaPtBin& etaPtBin = getEtaPtBin(iEta, iPt);
0143 
0144           // Instantiate the generic b tag plotter
0145           bool doDifferentialPlots =
0146               iModule->exists("differentialPlots") && iModule->getParameter<bool>("differentialPlots") == true;
0147           std::shared_ptr<JetTagPlotter> jetTagPlotter =
0148               std::make_shared<JetTagPlotter>(folderName,
0149                                               etaPtBin,
0150                                               iModule->getParameter<edm::ParameterSet>("parameters"),
0151                                               mcPlots_,
0152                                               true,
0153                                               ibook,
0154                                               doCTagPlots,
0155                                               doDifferentialPlots);
0156           binJetTagPlotters.at(iTag).push_back(jetTagPlotter);
0157 
0158           // Add to the corresponding differential plotters
0159           if (mcPlots_ && makeDiffPlots_) {
0160             (*differentialPlotsConstantEta)[iEta + 1]->addBinPlotter(jetTagPlotter);
0161             (*differentialPlotsConstantPt)[iPt + 1]->addBinPlotter(jetTagPlotter);
0162           }
0163         }
0164       }
0165       // the objects for the differential plots vs. eta, pt: collect all from constant eta and constant pt
0166       if (mcPlots_ && makeDiffPlots_) {
0167         differentialPlots.at(iTag).reserve(differentialPlotsConstantEta->size() + differentialPlotsConstantPt->size());
0168         differentialPlots.at(iTag).insert(differentialPlots.at(iTag).end(),
0169                                           std::make_move_iterator(differentialPlotsConstantEta->begin()),
0170                                           std::make_move_iterator(differentialPlotsConstantEta->end()));
0171         differentialPlots.at(iTag).insert(differentialPlots.at(iTag).end(),
0172                                           std::make_move_iterator(differentialPlotsConstantPt->begin()),
0173                                           std::make_move_iterator(differentialPlotsConstantPt->end()));
0174 
0175         edm::LogInfo("Info") << "====>>>> ## sizeof differentialPlots = " << differentialPlots.size();
0176       }
0177     } else if (dataFormatType == "TagCorrelation") {
0178       iTagCorr++;
0179       const InputTag& label1 = iModule->getParameter<InputTag>("label1");
0180       const InputTag& label2 = iModule->getParameter<InputTag>("label2");
0181 
0182       // eta loop
0183       for (int iEta = iEtaStart; iEta != iEtaEnd; ++iEta) {
0184         // pt loop
0185         for (int iPt = iPtStart; iPt != iPtEnd; ++iPt) {
0186           const EtaPtBin& etaPtBin = getEtaPtBin(iEta, iPt);
0187           // Instantiate the generic b tag correlation plotter
0188           std::unique_ptr<TagCorrelationPlotter> tagCorrelationPlotter =
0189               std::make_unique<TagCorrelationPlotter>(label1.label(),
0190                                                       label2.label(),
0191                                                       etaPtBin,
0192                                                       iModule->getParameter<edm::ParameterSet>("parameters"),
0193                                                       mcPlots_,
0194                                                       doCTagPlots,
0195                                                       true,
0196                                                       ibook);
0197           binTagCorrelationPlotters.at(iTagCorr).push_back(std::move(tagCorrelationPlotter));
0198         }
0199       }
0200     } else {
0201       iInfoTag++;
0202       // tag info retrievel is deferred(needs availability of EventSetup)
0203       const InputTag& moduleLabel = iModule->getParameter<InputTag>("label");
0204       const string& folderName = iModule->getParameter<string>("folder");
0205       // eta loop
0206       for (int iEta = iEtaStart; iEta < iEtaEnd; iEta++) {
0207         // pt loop
0208         for (int iPt = iPtStart; iPt < iPtEnd; iPt++) {
0209           const EtaPtBin& etaPtBin = getEtaPtBin(iEta, iPt);
0210 
0211           // Instantiate the tagInfo plotter
0212 
0213           std::unique_ptr<BaseTagInfoPlotter> jetTagPlotter =
0214               theFactory.buildPlotter(dataFormatType,
0215                                       moduleLabel.label(),
0216                                       etaPtBin,
0217                                       iModule->getParameter<edm::ParameterSet>("parameters"),
0218                                       folderName,
0219                                       mcPlots_,
0220                                       true,
0221                                       ibook);
0222           binTagInfoPlotters.at(iInfoTag).push_back(std::move(jetTagPlotter));
0223         }
0224       }
0225       edm::LogInfo("Info") << "====>>>> ## sizeof differentialPlots = " << differentialPlots.size();
0226     }
0227   }
0228 
0229   setTDRStyle();
0230   for (unsigned int iJetLabel = 0; iJetLabel != binJetTagPlotters.size(); ++iJetLabel) {
0231     int plotterSize = binJetTagPlotters[iJetLabel].size();
0232     for (int iPlotter = 0; iPlotter != plotterSize; ++iPlotter) {
0233       binJetTagPlotters[iJetLabel][iPlotter]->finalize(ibook, iget);
0234       if (producePs)
0235         (*binJetTagPlotters[iJetLabel][iPlotter]).psPlot(psBaseName);
0236       if (produceEps)
0237         (*binJetTagPlotters[iJetLabel][iPlotter]).epsPlot(epsBaseName);
0238     }
0239 
0240     if (makeDiffPlots_) {
0241       for (auto& iPlotter : differentialPlots[iJetLabel]) {
0242         iPlotter->process(ibook);
0243         if (producePs)
0244           iPlotter->psPlot(psBaseName);
0245         if (produceEps)
0246           iPlotter->epsPlot(epsBaseName);
0247       }
0248     }
0249   }
0250 
0251   for (auto& iJetLabel : binTagInfoPlotters) {
0252     for (auto& iPlotter : iJetLabel) {
0253       iPlotter->finalize(ibook, iget);
0254       if (producePs)
0255         iPlotter->psPlot(psBaseName);
0256       if (produceEps)
0257         iPlotter->epsPlot(epsBaseName);
0258     }
0259   }
0260   for (unsigned int iJetLabel = 0; iJetLabel != binTagCorrelationPlotters.size(); ++iJetLabel) {
0261     int plotterSize = binTagCorrelationPlotters[iJetLabel].size();
0262     for (int iPlotter = 0; iPlotter != plotterSize; ++iPlotter) {
0263       binTagCorrelationPlotters[iJetLabel][iPlotter]->finalize(ibook, iget);
0264       if (producePs)
0265         (*binTagCorrelationPlotters[iJetLabel][iPlotter]).psPlot(psBaseName);
0266       if (produceEps)
0267         (*binTagCorrelationPlotters[iJetLabel][iPlotter]).epsPlot(epsBaseName);
0268     }
0269   }
0270 }
0271 
0272 //define this as a plug-in
0273 DEFINE_FWK_MODULE(BTagPerformanceHarvester);