CMS 3D CMS Logo

BTagPerformanceHarvester.cc
Go to the documentation of this file.
7 
8 using namespace edm;
9 using namespace std;
10 using namespace RecoBTag;
11 
13  etaRanges(pSet.getParameter< vector<double> >("etaRanges")),
14  ptRanges(pSet.getParameter< vector<double> >("ptRanges")),
15  produceEps(pSet.getParameter< bool >("produceEps")),
16  producePs(pSet.getParameter< bool >("producePs")),
17  moduleConfig(pSet.getParameter< vector<edm::ParameterSet> >("tagConfig")),
18  flavPlots_(pSet.getParameter< std::string >("flavPlots")),
19  makeDiffPlots_(pSet.getParameter< bool >("differentialPlots"))
20 {
21  //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.
22  if (flavPlots_.find("dusg")<15) {
23  if (flavPlots_.find("all")<15)
24  mcPlots_ = 4;
25  else
26  mcPlots_ = 3;
27  } else if (flavPlots_.find("bcl")<15) {
28  if (flavPlots_.find("all")<15)
29  mcPlots_ = 2;
30  else
31  mcPlots_ = 1;
32  } else
33  mcPlots_ = 0;
34 
35  if (etaRanges.size() <= 1)
36  etaRanges = { pSet.getParameter<double>("etaMin"), pSet.getParameter<double>("etaMax") };
37  if (ptRanges.size() <= 1)
38  ptRanges = { pSet.getParameter<double>("ptRecJetMin"), pSet.getParameter<double>("ptRecJetMax") };
39 
40  for (vector<edm::ParameterSet>::const_iterator iModule = moduleConfig.begin();
41  iModule != moduleConfig.end(); ++iModule) {
42  const string& dataFormatType = iModule->exists("type") ?
43  iModule->getParameter<string>("type") :
44  "JetTag";
45  if (dataFormatType == "JetTag") {
46  const InputTag& moduleLabel = iModule->getParameter<InputTag>("label");
47  jetTagInputTags.push_back(moduleLabel);
48  binJetTagPlotters.push_back(vector<std::shared_ptr<JetTagPlotter>>());
49  if (mcPlots_ && makeDiffPlots_) {
50  differentialPlots.push_back(vector<std::unique_ptr<BTagDifferentialPlot>>());
51  }
52  } else if (dataFormatType == "TagCorrelation") {
53  const InputTag& label1 = iModule->getParameter<InputTag>("label1");
54  const InputTag& label2 = iModule->getParameter<InputTag>("label2");
55  tagCorrelationInputTags.push_back(std::pair<edm::InputTag, edm::InputTag>(label1, label2));
56  binTagCorrelationPlotters.push_back(vector<std::unique_ptr<TagCorrelationPlotter>>());
57  } else {
58  tagInfoInputTags.push_back(vector<edm::InputTag>());
59  tiDataFormatType.push_back(dataFormatType);
60  binTagInfoPlotters.push_back(vector<std::unique_ptr<BaseTagInfoPlotter>>());
61  }
62  }
63 
64 }
65 
66 EtaPtBin BTagPerformanceHarvester::getEtaPtBin(const int& iEta, const int& iPt)
67 {
68  // DEFINE BTagBin:
69  bool etaActive_ , ptActive_;
70  double etaMin_, etaMax_, ptMin_, ptMax_;
71 
72  if (iEta != -1) {
73  etaActive_ = true;
74  etaMin_ = etaRanges[iEta];
75  etaMax_ = etaRanges[iEta+1];
76  }
77  else {
78  etaActive_ = false;
79  etaMin_ = etaRanges[0];
80  etaMax_ = etaRanges[etaRanges.size() - 1];
81  }
82 
83  if (iPt != -1) {
84  ptActive_ = true;
85  ptMin_ = ptRanges[iPt];
86  ptMax_ = ptRanges[iPt+1];
87  }
88  else {
89  ptActive_ = false;
90  ptMin_ = ptRanges[0];
91  ptMax_ = ptRanges[ptRanges.size() - 1];
92  }
93  return EtaPtBin(etaActive_, etaMin_, etaMax_, ptActive_, ptMin_, ptMax_);
94 }
95 
97 
99 {
100  // Book all histograms.
101 
102  // iterate over ranges:
103  const int iEtaStart = -1 ; // this will be the inactive one
104  const int iEtaEnd = etaRanges.size() > 2 ? etaRanges.size() - 1 : 0; // if there is only one bin defined, leave it as the inactive one
105  const int iPtStart = -1 ; // this will be the inactive one
106  const int iPtEnd = ptRanges.size() > 2 ? ptRanges.size() - 1 : 0; // if there is only one bin defined, leave it as the inactive one
107  setTDRStyle();
108 
109  TagInfoPlotterFactory theFactory;
110  int iTag = -1; int iTagCorr = -1; int iInfoTag = -1;
111  for (vector<edm::ParameterSet>::const_iterator iModule = moduleConfig.begin();
112  iModule != moduleConfig.end(); ++iModule) {
113 
114  const string& dataFormatType = iModule->exists("type") ?
115  iModule->getParameter<string>("type") :
116  "JetTag";
117  const bool& doCTagPlots = iModule->exists("doCTagPlots") ?
118  iModule->getParameter<bool>("doCTagPlots") :
119  false;
120 
121  if (dataFormatType == "JetTag") {
122  iTag++;
123  const string& folderName = iModule->getParameter<string>("folder");
124 
125  // Contains plots for each bin of rapidity and pt.
126  auto differentialPlotsConstantEta = std::make_unique< std::vector<std::unique_ptr<BTagDifferentialPlot>> >();
127  auto differentialPlotsConstantPt = std::make_unique< std::vector<std::unique_ptr<BTagDifferentialPlot>> >();
128  if (mcPlots_ && makeDiffPlots_) {
129  // the constant b-efficiency for the differential plots versus pt and eta
130  const double& effBConst =
131  iModule->getParameter<edm::ParameterSet>("parameters").getParameter<double>("effBConst");
132 
133  // the objects for the differential plots vs. eta,pt for
134  for (int iEta = iEtaStart; iEta < iEtaEnd; iEta++) {
135  std::unique_ptr<BTagDifferentialPlot> etaConstDifferentialPlot = std::make_unique<BTagDifferentialPlot>(effBConst, BTagDifferentialPlot::constETA, folderName, mcPlots_);
136  differentialPlotsConstantEta->push_back(std::move(etaConstDifferentialPlot));
137  }
138 
139  for (int iPt = iPtStart; iPt < iPtEnd; iPt++) {
140  // differentialPlots for this pt bin
141  std::unique_ptr<BTagDifferentialPlot> ptConstDifferentialPlot = std::make_unique<BTagDifferentialPlot>(effBConst, BTagDifferentialPlot::constPT, folderName, mcPlots_);
142  differentialPlotsConstantPt->push_back(std::move(ptConstDifferentialPlot));
143  }
144  }
145  // eta loop
146  for (int iEta = iEtaStart; iEta < iEtaEnd; iEta++) {
147  // pt loop
148  for (int iPt = iPtStart; iPt < iPtEnd; iPt++) {
149 
150  const EtaPtBin& etaPtBin = getEtaPtBin(iEta, iPt);
151 
152  // Instantiate the generic b tag plotter
153  bool doDifferentialPlots = iModule->exists("differentialPlots") && iModule->getParameter<bool>("differentialPlots") == true;
154  std::shared_ptr<JetTagPlotter> jetTagPlotter = std::make_shared<JetTagPlotter>(folderName, etaPtBin,
155  iModule->getParameter<edm::ParameterSet>("parameters"),mcPlots_,true, ibook, doCTagPlots, doDifferentialPlots);
156  binJetTagPlotters.at(iTag).push_back(jetTagPlotter);
157 
158  // Add to the corresponding differential plotters
159  if (mcPlots_ && makeDiffPlots_) {
160  (*differentialPlotsConstantEta)[iEta+1]->addBinPlotter(jetTagPlotter);
161  (*differentialPlotsConstantPt)[iPt+1] ->addBinPlotter(jetTagPlotter);
162  }
163  }
164  }
165  // the objects for the differential plots vs. eta, pt: collect all from constant eta and constant pt
166  if (mcPlots_ && makeDiffPlots_) {
167  differentialPlots.at(iTag).reserve(differentialPlotsConstantEta->size() + differentialPlotsConstantPt->size());
168  differentialPlots.at(iTag).insert(differentialPlots.at(iTag).end(), std::make_move_iterator(differentialPlotsConstantEta->begin()), std::make_move_iterator(differentialPlotsConstantEta->end()));
169  differentialPlots.at(iTag).insert(differentialPlots.at(iTag).end(), std::make_move_iterator(differentialPlotsConstantPt->begin()), std::make_move_iterator(differentialPlotsConstantPt->end()));
170 
171  edm::LogInfo("Info")
172  << "====>>>> ## sizeof differentialPlots = " << differentialPlots.size();
173  }
174  } else if (dataFormatType == "TagCorrelation") {
175  iTagCorr++;
176  const InputTag& label1 = iModule->getParameter<InputTag>("label1");
177  const InputTag& label2 = iModule->getParameter<InputTag>("label2");
178 
179  // eta loop
180  for (int iEta = iEtaStart; iEta != iEtaEnd; ++iEta) {
181  // pt loop
182  for (int iPt = iPtStart; iPt != iPtEnd; ++iPt) {
183  const EtaPtBin& etaPtBin = getEtaPtBin(iEta, iPt);
184  // Instantiate the generic b tag correlation plotter
185  std::unique_ptr<TagCorrelationPlotter> tagCorrelationPlotter = std::make_unique<TagCorrelationPlotter>(label1.label(), label2.label(), etaPtBin,
186  iModule->getParameter<edm::ParameterSet>("parameters"),
187  mcPlots_, doCTagPlots, true, ibook);
188  binTagCorrelationPlotters.at(iTagCorr).push_back(std::move(tagCorrelationPlotter));
189  }
190  }
191  } else {
192  iInfoTag++;
193  // tag info retrievel is deferred(needs availability of EventSetup)
194  const InputTag& moduleLabel = iModule->getParameter<InputTag>("label");
195  const string& folderName = iModule->getParameter<string>("folder");
196  // eta loop
197  for (int iEta = iEtaStart; iEta < iEtaEnd; iEta++) {
198  // pt loop
199  for (int iPt = iPtStart; iPt < iPtEnd; iPt++) {
200  const EtaPtBin& etaPtBin = getEtaPtBin(iEta, iPt);
201 
202  // Instantiate the tagInfo plotter
203 
204  std::unique_ptr<BaseTagInfoPlotter> jetTagPlotter = theFactory.buildPlotter(dataFormatType, moduleLabel.label(),
205  etaPtBin, iModule->getParameter<edm::ParameterSet>("parameters"), folderName,
206  mcPlots_, true, ibook);
207  binTagInfoPlotters.at(iInfoTag).push_back(std::move(jetTagPlotter));
208  }
209  }
210  edm::LogInfo("Info")
211  << "====>>>> ## sizeof differentialPlots = " << differentialPlots.size();
212  }
213  }
214 
215  setTDRStyle();
216  for (unsigned int iJetLabel = 0; iJetLabel != binJetTagPlotters.size(); ++iJetLabel) {
217  int plotterSize = binJetTagPlotters[iJetLabel].size();
218  for (int iPlotter = 0; iPlotter != plotterSize; ++iPlotter) {
219  binJetTagPlotters[iJetLabel][iPlotter]->finalize(ibook, iget);
220  if (producePs) (*binJetTagPlotters[iJetLabel][iPlotter]).psPlot(psBaseName);
221  if (produceEps)(*binJetTagPlotters[iJetLabel][iPlotter]).epsPlot(epsBaseName);
222  }
223 
224  if (makeDiffPlots_) {
225  for (auto& iPlotter: differentialPlots[iJetLabel]) {
226  iPlotter->process(ibook);
227  if (producePs) iPlotter->psPlot(psBaseName);
228  if (produceEps) iPlotter->epsPlot(epsBaseName);
229  }
230  }
231  }
232 
233  for (auto& iJetLabel: binTagInfoPlotters) {
234  for (auto& iPlotter: iJetLabel) {
235  iPlotter->finalize(ibook, iget);
236  if (producePs) iPlotter->psPlot(psBaseName);
237  if (produceEps) iPlotter->epsPlot(epsBaseName);
238  }
239  }
240  for (unsigned int iJetLabel = 0; iJetLabel != binTagCorrelationPlotters.size(); ++iJetLabel) {
241  int plotterSize = binTagCorrelationPlotters[iJetLabel].size();
242  for (int iPlotter = 0; iPlotter != plotterSize; ++iPlotter) {
243  binTagCorrelationPlotters[iJetLabel][iPlotter]->finalize(ibook, iget);
244  if (producePs) (*binTagCorrelationPlotters[iJetLabel][iPlotter]).psPlot(psBaseName);
245  if (produceEps)(*binTagCorrelationPlotters[iJetLabel][iPlotter]).epsPlot(epsBaseName);
246  }
247  }
248 
249 }
250 
251 
252 //define this as a plug-in
T getParameter(std::string const &) const
std::vector< std::pair< edm::InputTag, edm::InputTag > > tagCorrelationInputTags
BTagPerformanceHarvester(const edm::ParameterSet &pSet)
std::vector< std::vector< edm::InputTag > > tagInfoInputTags
std::vector< std::vector< std::shared_ptr< JetTagPlotter > > > binJetTagPlotters
void dqmEndJob(DQMStore::IBooker &, DQMStore::IGetter &) override
EtaPtBin getEtaPtBin(const int &iEta, const int &iPt)
std::vector< edm::InputTag > jetTagInputTags
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
def setTDRStyle()
Definition: plotscripts.py:89
Definition: Tools.h:23
std::string const & label() const
Definition: InputTag.h:36
std::vector< std::vector< std::unique_ptr< BaseTagInfoPlotter > > > binTagInfoPlotters
HLT enums.
std::vector< std::string > tiDataFormatType
std::vector< edm::ParameterSet > moduleConfig
std::unique_ptr< BaseTagInfoPlotter > buildPlotter(const std::string &dataFormatType, const std::string &tagName, const EtaPtBin &etaPtBin, const edm::ParameterSet &pSet, const std::string &folderName, unsigned int mc, bool wf, DQMStore::IBooker &ibook)
std::vector< std::vector< std::unique_ptr< BTagDifferentialPlot > > > differentialPlots
std::vector< double > etaRanges
std::vector< std::vector< std::unique_ptr< TagCorrelationPlotter > > > binTagCorrelationPlotters
def move(src, dest)
Definition: eostools.py:511
std::vector< double > ptRanges