CMS 3D CMS Logo

HLTTauPostProcessor.cc
Go to the documentation of this file.
4 
5 #include "TEfficiency.h"
6 #include "TProfile.h"
7 
8 #include <tuple>
9 
10 namespace {
11  std::tuple<float, float> calcEfficiency(float num, float denom) {
12  if (denom == 0.0f)
13  return std::make_tuple(0.0f, 0.0f);
14 
15  //float eff = num/denom;
16  constexpr double cl = 0.683f; // "1-sigma"
17  const float eff = num / denom;
18  const float errDown = TEfficiency::ClopperPearson(denom, num, cl, false);
19  const float errUp = TEfficiency::ClopperPearson(denom, num, cl, true);
20 
21  // Because of limitation of TProfile, just take max
22  return std::make_tuple(eff, std::max(eff - errDown, errUp - eff));
23  }
24 } // namespace
25 
27  : dqmBaseFolder_(ps.getUntrackedParameter<std::string>("DQMBaseFolder")) {}
28 
30 
32  if (!iGetter.dirExists(dqmBaseFolder_)) {
33  LogDebug("HLTTauDQMOffline") << "Folder " << dqmBaseFolder_ << " does not exist";
34  return;
35  }
37  for (const std::string& subfolder : iGetter.getSubdirs()) {
38  std::size_t pos = subfolder.rfind('/');
39  if (pos == std::string::npos)
40  continue;
41  ++pos; // position of first letter after /
42  if (subfolder.compare(pos, 4, "HLT_") == 0) { // start with HLT_
43  LogDebug("HLTTauDQMOffline") << "Processing path " << subfolder.substr(pos);
44  plotFilterEfficiencies(iBooker, iGetter, subfolder);
45  }
46  }
47 }
48 
50  DQMStore::IGetter& iGetter,
51  const std::string& folder) const {
52  // Get the source
53  const MonitorElement* eventsPerFilter = iGetter.get(folder + "/EventsPerFilter");
54  if (!eventsPerFilter) {
55  LogDebug("HLTTauDQMOffline") << "ME " << folder << "/EventsPerFilter not found";
56  return;
57  }
58 
59  // Book efficiency TProfile
60  iBooker.setCurrentFolder(folder);
61  MonitorElement* efficiency = iBooker.bookProfile("EfficiencyRefPrevious",
62  "Efficiency to previous filter",
63  eventsPerFilter->getNbinsX() - 1,
64  0,
65  eventsPerFilter->getNbinsX() - 1,
66  100,
67  0,
68  1);
69  efficiency->setAxisTitle("Efficiency", 2);
70  const TAxis* xaxis = eventsPerFilter->getTH1F()->GetXaxis();
71  for (int bin = 1; bin < eventsPerFilter->getNbinsX(); ++bin) {
72  efficiency->setBinLabel(bin, xaxis->GetBinLabel(bin + 1));
73  }
74 
75  // Fill efficiency TProfile
76  TProfile* prev = efficiency->getTProfile();
77  for (int i = 2; i <= eventsPerFilter->getNbinsX(); ++i) {
78  if (eventsPerFilter->getBinContent(i - 1) < eventsPerFilter->getBinContent(i)) {
79  LogDebug("HLTTauDQMOffline") << "HLTTauPostProcessor: Encountered denominator < numerator with efficiency plot "
80  "EfficiencyRefPrevious in folder "
81  << folder << ", bin " << i << " numerator " << eventsPerFilter->getBinContent(i)
82  << " denominator " << eventsPerFilter->getBinContent(i - 1);
83  continue;
84  }
85  const std::tuple<float, float> effErr =
86  calcEfficiency(eventsPerFilter->getBinContent(i), eventsPerFilter->getBinContent(i - 1));
87  const float efficiency = std::get<0>(effErr);
88  const float err = std::get<1>(effErr);
89 
90  prev->SetBinContent(i - 1, efficiency);
91  prev->SetBinEntries(i - 1, 1);
92  prev->SetBinError(i - 1, std::sqrt(efficiency * efficiency + err * err));
93  }
94 }
~HLTTauPostProcessor() override
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
virtual bool dirExists(std::string const &path) const
Definition: DQMStore.cc:747
void dqmEndJob(DQMStore::IBooker &iBooker, DQMStore::IGetter &iGetter) override
HLTTauPostProcessor(const edm::ParameterSet &)
MonitorElement * bookProfile(TString const &name, TString const &title, int nchX, double lowX, double highX, int, double lowY, double highY, char const *option="s", FUNC onbooking=NOOP())
Definition: DQMStore.h:399
T sqrt(T t)
Definition: SSEVec.h:19
double f[11][100]
void plotFilterEfficiencies(DQMStore::IBooker &iBooker, DQMStore::IGetter &iGetter, const std::string &folder) const
virtual TH1F * getTH1F() const
virtual MonitorElement * get(std::string const &fullpath) const
Definition: DQMStore.cc:690
virtual int getNbinsX() const
get # of bins in X-axis
#define LogDebug(id)
virtual double getBinContent(int binx) const
get content of bin (1-D)
const std::string dqmBaseFolder_
virtual DQM_DEPRECATED std::vector< std::string > getSubdirs() const
Definition: DQMStore.cc:717