CMS 3D CMS Logo

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