CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
MuonME0SegHarvestor.cc
Go to the documentation of this file.
1 // system include files
2 #include <memory>
3 
4 // user include files
12 #include "TGraphAsymmErrors.h"
13 
17 
22 
24 
28 
30 
32  dbe_path_ = std::string("MuonME0RecHitsV/ME0SegmentsTask/");
33 }
34 
36 
37 TProfile *MuonME0SegHarvestor::ComputeEff(TH1F *num, TH1F *denum, std::string nameHist) {
38  std::string name = "eff_" + nameHist;
39  std::string title = "Segment Efficiency" + std::string(num->GetTitle());
40  TProfile *efficHist = new TProfile(name.c_str(),
41  title.c_str(),
42  denum->GetXaxis()->GetNbins(),
43  denum->GetXaxis()->GetXmin(),
44  denum->GetXaxis()->GetXmax());
45 
46  for (int i = 1; i <= denum->GetNbinsX(); i++) {
47  double nNum = num->GetBinContent(i);
48  double nDenum = denum->GetBinContent(i);
49  if (nDenum == 0 || nNum == 0) {
50  continue;
51  }
52  if (nNum > nDenum) {
53  double temp = nDenum;
54  nDenum = nNum;
55  nNum = temp;
56  edm::LogWarning("MuonME0SegHarvestor")
57  << "Alert! specific bin's num is bigger than denum " << i << " " << nNum << " " << nDenum;
58  }
59  const double effVal = nNum / nDenum;
60  efficHist->SetBinContent(i, effVal);
61  efficHist->SetBinEntries(i, 1);
62  efficHist->SetBinError(i, 0);
63  const double errLo = TEfficiency::ClopperPearson((int)nDenum, (int)nNum, 0.683, false);
64  const double errUp = TEfficiency::ClopperPearson((int)nDenum, (int)nNum, 0.683, true);
65  const double errVal = (effVal - errLo > errUp - effVal) ? effVal - errLo : errUp - effVal;
66  efficHist->SetBinError(i, sqrt(effVal * effVal + errVal * errVal));
67  }
68  return efficHist;
69 }
70 
72  DQMStore::IBooker &ibooker, DQMStore::IGetter &ig, std::string nameHist, TH1F *num, TH1F *den) {
73  if (num != nullptr && den != nullptr) {
74  TProfile *profile = ComputeEff(num, den, nameHist);
75 
76  TString x_axis_title = TString(num->GetXaxis()->GetTitle());
77  TString title = TString::Format("Segment Efficiency;%s;Eff.", x_axis_title.Data());
78 
79  profile->SetTitle(title.Data());
80  ibooker.bookProfile(profile->GetName(), profile);
81 
82  delete profile;
83 
84  } else {
85  edm::LogWarning("MuonME0SegHarvestor") << "Can not find histograms";
86  if (num == nullptr)
87  edm::LogWarning("MuonME0SegHarvestor") << "num not found";
88  if (den == nullptr)
89  edm::LogWarning("MuonME0SegHarvestor") << "den not found";
90  }
91  return;
92 }
93 
96 
97  TString eta_label_den = TString(dbe_path_) + "me0_simsegment_eta";
98  TString eta_label_num = TString(dbe_path_) + "me0_matchedsimsegment_eta";
99  TString pt_label_den = TString(dbe_path_) + "me0_simsegment_pt";
100  TString pt_label_num = TString(dbe_path_) + "me0_matchedsimsegment_pt";
101  TString phi_label_den = TString(dbe_path_) + "me0_simsegment_phi";
102  TString phi_label_num = TString(dbe_path_) + "me0_matchedsimsegment_phi";
103 
104  if (ig.get(eta_label_num.Data()) != nullptr && ig.get(eta_label_den.Data()) != nullptr) {
105  TH1F *num_vs_eta = (TH1F *)ig.get(eta_label_num.Data())->getTH1F()->Clone();
106  num_vs_eta->Sumw2();
107  TH1F *den_vs_eta = (TH1F *)ig.get(eta_label_den.Data())->getTH1F()->Clone();
108  den_vs_eta->Sumw2();
109 
110  ProcessBooking(ibooker, ig, "me0segment_eff_vs_eta", num_vs_eta, den_vs_eta);
111 
112  delete num_vs_eta;
113  delete den_vs_eta;
114 
115  } else
116  edm::LogWarning("MuonME0SegHarvestor") << "Can not find histograms: " << eta_label_num << " or " << eta_label_den;
117 
118  if (ig.get(pt_label_num.Data()) != nullptr && ig.get(pt_label_den.Data()) != nullptr) {
119  TH1F *num_vs_pt = (TH1F *)ig.get(pt_label_num.Data())->getTH1F()->Clone();
120  num_vs_pt->Sumw2();
121  TH1F *den_vs_pt = (TH1F *)ig.get(pt_label_den.Data())->getTH1F()->Clone();
122  den_vs_pt->Sumw2();
123 
124  ProcessBooking(ibooker, ig, "me0segment_eff_vs_pt", num_vs_pt, den_vs_pt);
125 
126  delete num_vs_pt;
127  delete den_vs_pt;
128 
129  } else
130  edm::LogWarning("MuonME0SegHarvestor") << "Can not find histograms: " << pt_label_num << " or " << pt_label_den;
131 
132  if (ig.get(phi_label_num.Data()) != nullptr && ig.get(phi_label_den.Data()) != nullptr) {
133  TH1F *num_vs_phi = (TH1F *)ig.get(phi_label_num.Data())->getTH1F()->Clone();
134  num_vs_phi->Sumw2();
135  TH1F *den_vs_phi = (TH1F *)ig.get(phi_label_den.Data())->getTH1F()->Clone();
136  den_vs_phi->Sumw2();
137 
138  ProcessBooking(ibooker, ig, "me0segment_eff_vs_phi", num_vs_phi, den_vs_phi);
139 
140  delete num_vs_phi;
141  delete den_vs_phi;
142 
143  } else
144  edm::LogWarning("MuonME0SegHarvestor") << "Can not find histograms: " << phi_label_num << " or " << phi_label_den;
145 }
146 
147 // define this as a plug-in
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
MuonME0SegHarvestor(const edm::ParameterSet &)
constructor
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:322
void ProcessBooking(DQMStore::IBooker &, DQMStore::IGetter &, std::string nameHist, TH1F *num, TH1F *den)
T sqrt(T t)
Definition: SSEVec.h:19
~MuonME0SegHarvestor() override
destructor
virtual MonitorElement * get(std::string const &fullpath) const
Definition: DQMStore.cc:673
TProfile * ComputeEff(TH1F *num, TH1F *denum, std::string nameHist)
Log< level::Warning, false > LogWarning
void dqmEndJob(DQMStore::IBooker &, DQMStore::IGetter &) override