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 }
25 
27  dqmBaseFolder_(ps.getUntrackedParameter<std::string>("DQMBaseFolder"))
28 {}
29 
31 
33 {
34  if(!iGetter.dirExists(dqmBaseFolder_)) {
35  LogDebug("HLTTauDQMOffline") << "Folder " << dqmBaseFolder_ << " does not exist";
36  return;
37  }
39  for(const std::string& subfolder: iGetter.getSubdirs()) {
40  std::size_t pos = subfolder.rfind('/');
41  if(pos == std::string::npos)
42  continue;
43  ++pos; // position of first letter after /
44  if(subfolder.compare(pos, 4, "HLT_") == 0) { // start with HLT_
45  LogDebug("HLTTauDQMOffline") << "Processing path " << subfolder.substr(pos);
46  plotFilterEfficiencies(iBooker, iGetter, subfolder);
47  }
48  }
49 }
50 
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", "Efficiency to previous filter", eventsPerFilter->getNbinsX()-1,0,eventsPerFilter->getNbinsX()-1, 100,0,1);
62  efficiency->setAxisTitle("Efficiency", 2);
63  const TAxis *xaxis = eventsPerFilter->getTH1F()->GetXaxis();
64  for(int bin=1; bin < eventsPerFilter->getNbinsX(); ++bin) {
65  efficiency->setBinLabel(bin, xaxis->GetBinLabel(bin+1));
66  }
67 
68  // Fill efficiency TProfile
69  TProfile *prev = efficiency->getTProfile();
70  for(int i=2; i <= eventsPerFilter->getNbinsX(); ++i) {
71  if(eventsPerFilter->getBinContent(i-1) < eventsPerFilter->getBinContent(i)) {
72  LogDebug("HLTTauDQMOffline") << "HLTTauPostProcessor: Encountered denominator < numerator with efficiency plot EfficiencyRefPrevious in folder " << folder << ", bin " << i << " numerator " << eventsPerFilter->getBinContent(i) << " denominator " << eventsPerFilter->getBinContent(i-1);
73  continue;
74  }
75  const std::tuple<float, float> effErr = calcEfficiency(eventsPerFilter->getBinContent(i), eventsPerFilter->getBinContent(i-1));
76  const float efficiency = std::get<0>(effErr);
77  const float err = std::get<1>(effErr);
78 
79  prev->SetBinContent(i-1, efficiency);
80  prev->SetBinEntries(i-1, 1);
81  prev->SetBinError(i-1, std::sqrt(efficiency*efficiency + err*err));
82  }
83 }
84 
#define LogDebug(id)
TProfile * getTProfile() const
~HLTTauPostProcessor() override
MonitorElement * bookProfile(Args &&...args)
Definition: DQMStore.h:113
TH1F * getTH1F() const
void setBinLabel(int bin, const std::string &label, int axis=1)
set bin label for x, y or z axis (axis=1, 2, 3 respectively)
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:361
void dqmEndJob(DQMStore::IBooker &iBooker, DQMStore::IGetter &iGetter) override
HLTTauPostProcessor(const edm::ParameterSet &)
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
T sqrt(T t)
Definition: SSEVec.h:18
double f[11][100]
MonitorElement * get(std::string const &path)
Definition: DQMStore.cc:303
bin
set the eta bin as selection string.
void plotFilterEfficiencies(DQMStore::IBooker &iBooker, DQMStore::IGetter &iGetter, const std::string &folder) const
bool dirExists(std::string const &path)
Definition: DQMStore.cc:343
double getBinContent(int binx) const
get content of bin (1-D)
std::vector< std::string > getSubdirs()
Definition: DQMStore.cc:325
int getNbinsX() const
get # of bins in X-axis
void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
#define constexpr
const std::string dqmBaseFolder_