CMS 3D CMS Logo

MtdTracksHarvester.cc
Go to the documentation of this file.
1 #include <string>
2 
7 
10 
12 
14 public:
15  explicit MtdTracksHarvester(const edm::ParameterSet& iConfig);
16  ~MtdTracksHarvester() override;
17 
18  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
19 
20 protected:
22 
23 private:
25 
27 
28  // --- Histograms
39 };
40 
41 // ------------ constructor and destructor --------------
43  : folder_(iConfig.getParameter<std::string>("folder")) {}
44 
46 
47 // auxiliary method to compute efficiency from the ratio of two 1D MonitorElement
49  for (int ibin = 1; ibin <= den->getNbinsX(); ibin++) {
50  double eff = num->getBinContent(ibin) / den->getBinContent(ibin);
51  double bin_err = sqrt((num->getBinContent(ibin) * (den->getBinContent(ibin) - num->getBinContent(ibin))) /
52  pow(den->getBinContent(ibin), 3));
53  if (den->getBinContent(ibin) == 0) {
54  eff = 0;
55  bin_err = 0;
56  }
57  result->setBinContent(ibin, eff);
58  result->setBinError(ibin, bin_err);
59  }
60 }
61 
62 // ------------ endjob tasks ----------------------------
64  // --- Get the monitoring histograms
65  MonitorElement* meBTLTrackEffEtaTot = igetter.get(folder_ + "TrackBTLEffEtaTot");
66  MonitorElement* meBTLTrackEffPhiTot = igetter.get(folder_ + "TrackBTLEffPhiTot");
67  MonitorElement* meBTLTrackEffPtTot = igetter.get(folder_ + "TrackBTLEffPtTot");
68  MonitorElement* meBTLTrackEffEtaMtd = igetter.get(folder_ + "TrackBTLEffEtaMtd");
69  MonitorElement* meBTLTrackEffPhiMtd = igetter.get(folder_ + "TrackBTLEffPhiMtd");
70  MonitorElement* meBTLTrackEffPtMtd = igetter.get(folder_ + "TrackBTLEffPtMtd");
71  MonitorElement* meETLTrackEffEtaTotZneg = igetter.get(folder_ + "TrackETLEffEtaTotZneg");
72  MonitorElement* meETLTrackEffPhiTotZneg = igetter.get(folder_ + "TrackETLEffPhiTotZneg");
73  MonitorElement* meETLTrackEffPtTotZneg = igetter.get(folder_ + "TrackETLEffPtTotZneg");
74  MonitorElement* meETLTrackEffEtaMtdZneg = igetter.get(folder_ + "TrackETLEffEtaMtdZneg");
75  MonitorElement* meETLTrackEffPhiMtdZneg = igetter.get(folder_ + "TrackETLEffPhiMtdZneg");
76  MonitorElement* meETLTrackEffPtMtdZneg = igetter.get(folder_ + "TrackETLEffPtMtdZneg");
77  MonitorElement* meETLTrackEffEtaTotZpos = igetter.get(folder_ + "TrackETLEffEtaTotZpos");
78  MonitorElement* meETLTrackEffPhiTotZpos = igetter.get(folder_ + "TrackETLEffPhiTotZpos");
79  MonitorElement* meETLTrackEffPtTotZpos = igetter.get(folder_ + "TrackETLEffPtTotZpos");
80  MonitorElement* meETLTrackEffEtaMtdZpos = igetter.get(folder_ + "TrackETLEffEtaMtdZpos");
81  MonitorElement* meETLTrackEffPhiMtdZpos = igetter.get(folder_ + "TrackETLEffPhiMtdZpos");
82  MonitorElement* meETLTrackEffPtMtdZpos = igetter.get(folder_ + "TrackETLEffPtMtdZpos");
83  MonitorElement* meMVATrackEffPtTot = igetter.get(folder_ + "MVAEffPtTot");
84  MonitorElement* meMVATrackMatchedEffPtTot = igetter.get(folder_ + "MVAMatchedEffPtTot");
85  MonitorElement* meMVATrackMatchedEffPtMtd = igetter.get(folder_ + "MVAMatchedEffPtMtd");
86  MonitorElement* meMVATrackEffEtaTot = igetter.get(folder_ + "MVAEffEtaTot");
87  MonitorElement* meMVATrackMatchedEffEtaTot = igetter.get(folder_ + "MVAMatchedEffEtaTot");
88  MonitorElement* meMVATrackMatchedEffEtaMtd = igetter.get(folder_ + "MVAMatchedEffEtaMtd");
89 
90  if (!meBTLTrackEffEtaTot || !meBTLTrackEffPhiTot || !meBTLTrackEffPtTot || !meBTLTrackEffEtaMtd ||
91  !meBTLTrackEffPhiMtd || !meBTLTrackEffPtMtd || !meETLTrackEffEtaTotZneg || !meETLTrackEffPhiTotZneg ||
92  !meETLTrackEffPtTotZneg || !meETLTrackEffEtaMtdZneg || !meETLTrackEffPhiMtdZneg || !meETLTrackEffPtMtdZneg ||
93  !meETLTrackEffEtaTotZpos || !meETLTrackEffPhiTotZpos || !meETLTrackEffPtTotZpos || !meETLTrackEffEtaMtdZpos ||
94  !meETLTrackEffPhiMtdZpos || !meETLTrackEffPtMtdZpos || !meMVATrackEffPtTot || !meMVATrackMatchedEffPtTot ||
95  !meMVATrackMatchedEffPtMtd || !meMVATrackEffEtaTot || !meMVATrackMatchedEffEtaTot ||
96  !meMVATrackMatchedEffEtaMtd) {
97  edm::LogError("MtdTracksHarvester") << "Monitoring histograms not found!" << std::endl;
98  return;
99  }
100 
101  // --- Book histograms
102  ibook.cd(folder_);
103  meBtlEtaEff_ = ibook.book1D("BtlEtaEff",
104  " Track Efficiency VS Eta;#eta;Efficiency",
105  meBTLTrackEffEtaTot->getNbinsX(),
106  meBTLTrackEffEtaTot->getTH1()->GetXaxis()->GetXmin(),
107  meBTLTrackEffEtaTot->getTH1()->GetXaxis()->GetXmax());
108  meBtlEtaEff_->getTH1()->SetMinimum(0.);
109  computeEfficiency1D(meBTLTrackEffEtaMtd, meBTLTrackEffEtaTot, meBtlEtaEff_);
110 
111  meBtlPhiEff_ = ibook.book1D("BtlPhiEff",
112  "Track Efficiency VS Phi;#phi [rad];Efficiency",
113  meBTLTrackEffPhiTot->getNbinsX(),
114  meBTLTrackEffPhiTot->getTH1()->GetXaxis()->GetXmin(),
115  meBTLTrackEffPhiTot->getTH1()->GetXaxis()->GetXmax());
116  meBtlPhiEff_->getTH1()->SetMinimum(0.);
117  computeEfficiency1D(meBTLTrackEffPhiMtd, meBTLTrackEffPhiTot, meBtlPhiEff_);
118 
119  meBtlPtEff_ = ibook.book1D("BtlPtEff",
120  "Track Efficiency VS Pt;Pt [GeV];Efficiency",
121  meBTLTrackEffPtTot->getNbinsX(),
122  meBTLTrackEffPtTot->getTH1()->GetXaxis()->GetXmin(),
123  meBTLTrackEffPtTot->getTH1()->GetXaxis()->GetXmax());
124  meBtlPtEff_->getTH1()->SetMinimum(0.);
125  computeEfficiency1D(meBTLTrackEffPtMtd, meBTLTrackEffPtTot, meBtlPtEff_);
126 
127  meEtlEtaEff_[0] = ibook.book1D("EtlEtaEffZneg",
128  " Track Efficiency VS Eta (-Z);#eta;Efficiency",
129  meETLTrackEffEtaTotZneg->getNbinsX(),
130  meETLTrackEffEtaTotZneg->getTH1()->GetXaxis()->GetXmin(),
131  meETLTrackEffEtaTotZneg->getTH1()->GetXaxis()->GetXmax());
132  meEtlEtaEff_[0]->getTH1()->SetMinimum(0.);
133  computeEfficiency1D(meETLTrackEffEtaMtdZneg, meETLTrackEffEtaTotZneg, meEtlEtaEff_[0]);
134 
135  meEtlPhiEff_[0] = ibook.book1D("EtlPhiEffZneg",
136  "Track Efficiency VS Phi (-Z);#phi [rad];Efficiency",
137  meETLTrackEffPhiTotZneg->getNbinsX(),
138  meETLTrackEffPhiTotZneg->getTH1()->GetXaxis()->GetXmin(),
139  meETLTrackEffPhiTotZneg->getTH1()->GetXaxis()->GetXmax());
140  meEtlPhiEff_[0]->getTH1()->SetMinimum(0.);
141  computeEfficiency1D(meETLTrackEffPhiMtdZneg, meETLTrackEffPhiTotZneg, meEtlPhiEff_[0]);
142 
143  meEtlPtEff_[0] = ibook.book1D("EtlPtEffZneg",
144  "Track Efficiency VS Pt (-Z);Pt [GeV];Efficiency",
145  meETLTrackEffPtTotZneg->getNbinsX(),
146  meETLTrackEffPtTotZneg->getTH1()->GetXaxis()->GetXmin(),
147  meETLTrackEffPtTotZneg->getTH1()->GetXaxis()->GetXmax());
148  meEtlPtEff_[0]->getTH1()->SetMinimum(0.);
149  computeEfficiency1D(meETLTrackEffPtMtdZneg, meETLTrackEffPtTotZneg, meEtlPtEff_[0]);
150 
151  meEtlEtaEff_[1] = ibook.book1D("EtlEtaEffZpos",
152  " Track Efficiency VS Eta (+Z);#eta;Efficiency",
153  meETLTrackEffEtaTotZpos->getNbinsX(),
154  meETLTrackEffEtaTotZpos->getTH1()->GetXaxis()->GetXmin(),
155  meETLTrackEffEtaTotZpos->getTH1()->GetXaxis()->GetXmax());
156  meEtlEtaEff_[1]->getTH1()->SetMinimum(0.);
157  computeEfficiency1D(meETLTrackEffEtaMtdZpos, meETLTrackEffEtaTotZpos, meEtlEtaEff_[1]);
158 
159  meEtlPhiEff_[1] = ibook.book1D("EtlPhiEffZpos",
160  "Track Efficiency VS Phi (+Z);#phi [rad];Efficiency",
161  meETLTrackEffPhiTotZpos->getNbinsX(),
162  meETLTrackEffPhiTotZpos->getTH1()->GetXaxis()->GetXmin(),
163  meETLTrackEffPhiTotZpos->getTH1()->GetXaxis()->GetXmax());
164  meEtlPhiEff_[1]->getTH1()->SetMinimum(0.);
165  computeEfficiency1D(meETLTrackEffPhiMtdZpos, meETLTrackEffPhiTotZpos, meEtlPhiEff_[1]);
166 
167  meEtlPtEff_[1] = ibook.book1D("EtlPtEffZpos",
168  "Track Efficiency VS Pt (+Z);Pt [GeV];Efficiency",
169  meETLTrackEffPtTotZpos->getNbinsX(),
170  meETLTrackEffPtTotZpos->getTH1()->GetXaxis()->GetXmin(),
171  meETLTrackEffPtTotZpos->getTH1()->GetXaxis()->GetXmax());
172  meEtlPtEff_[1]->getTH1()->SetMinimum(0.);
173  computeEfficiency1D(meETLTrackEffPtMtdZpos, meETLTrackEffPtTotZpos, meEtlPtEff_[1]);
174 
175  meMVAPtSelEff_ = ibook.book1D("MVAPtSelEff",
176  "Track selected efficiency VS Pt;Pt [GeV];Efficiency",
177  meMVATrackEffPtTot->getNbinsX(),
178  meMVATrackEffPtTot->getTH1()->GetXaxis()->GetXmin(),
179  meMVATrackEffPtTot->getTH1()->GetXaxis()->GetXmax());
180  meMVAPtSelEff_->getTH1()->SetMinimum(0.);
181  computeEfficiency1D(meMVATrackMatchedEffPtTot, meMVATrackEffPtTot, meMVAPtSelEff_);
182 
183  meMVAEtaSelEff_ = ibook.book1D("MVAEtaSelEff",
184  "Track selected efficiency VS Eta;Eta;Efficiency",
185  meMVATrackEffEtaTot->getNbinsX(),
186  meMVATrackEffEtaTot->getTH1()->GetXaxis()->GetXmin(),
187  meMVATrackEffEtaTot->getTH1()->GetXaxis()->GetXmax());
188  meMVAEtaSelEff_->getTH1()->SetMinimum(0.);
189  computeEfficiency1D(meMVATrackMatchedEffEtaTot, meMVATrackEffEtaTot, meMVAEtaSelEff_);
190 
191  meMVAPtMatchEff_ = ibook.book1D("MVAPtMatchEff",
192  "Track matched to GEN efficiency VS Pt;Pt [GeV];Efficiency",
193  meMVATrackMatchedEffPtTot->getNbinsX(),
194  meMVATrackMatchedEffPtTot->getTH1()->GetXaxis()->GetXmin(),
195  meMVATrackMatchedEffPtTot->getTH1()->GetXaxis()->GetXmax());
196  meMVAPtMatchEff_->getTH1()->SetMinimum(0.);
197  computeEfficiency1D(meMVATrackMatchedEffPtMtd, meMVATrackMatchedEffPtTot, meMVAPtMatchEff_);
198 
199  meMVAEtaMatchEff_ = ibook.book1D("MVAEtaMatchEff",
200  "Track matched to GEN efficiency VS Eta;Eta;Efficiency",
201  meMVATrackMatchedEffEtaTot->getNbinsX(),
202  meMVATrackMatchedEffEtaTot->getTH1()->GetXaxis()->GetXmin(),
203  meMVATrackMatchedEffEtaTot->getTH1()->GetXaxis()->GetXmax());
204  meMVAEtaMatchEff_->getTH1()->SetMinimum(0.);
205  computeEfficiency1D(meMVATrackMatchedEffEtaMtd, meMVATrackMatchedEffEtaTot, meMVAEtaMatchEff_);
206 
207  meBtlPhiEff_->getTH1()->SetMinimum(0.);
208  meBtlPtEff_->getTH1()->SetMinimum(0.);
209  for (int i = 0; i < 2; i++) {
210  meEtlEtaEff_[i]->getTH1()->SetMinimum(0.);
211  meEtlPhiEff_[i]->getTH1()->SetMinimum(0.);
212  meEtlPtEff_[i]->getTH1()->SetMinimum(0.);
213  }
214  meMVAPtSelEff_->getTH1()->SetMinimum(0.);
215  meMVAEtaSelEff_->getTH1()->SetMinimum(0.);
216  meMVAPtMatchEff_->getTH1()->SetMinimum(0.);
217  meMVAEtaMatchEff_->getTH1()->SetMinimum(0.);
218 }
219 
220 // ------------ method fills 'descriptions' with the allowed parameters for the module ----------
223 
224  desc.add<std::string>("folder", "MTD/Tracks/");
225 
226  descriptions.add("MtdTracksPostProcessor", desc);
227 }
228 
MtdTracksHarvester(const edm::ParameterSet &iConfig)
std::string folder_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
MonitorElement * meMVAEtaMatchEff_
Log< level::Error, false > LogError
void computeEfficiency1D(MonitorElement *num, MonitorElement *den, MonitorElement *result)
void dqmEndJob(DQMStore::IBooker &, DQMStore::IGetter &) override
MonitorElement * meMVAPtSelEff_
T sqrt(T t)
Definition: SSEVec.h:19
MonitorElement * meEtlPhiEff_[2]
MonitorElement * meMVAPtMatchEff_
MonitorElement * meEtlPtEff_[2]
const std::string folder_
MonitorElement * meBtlPtEff_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
MonitorElement * meBtlEtaEff_
virtual MonitorElement * get(std::string const &fullpath) const
Definition: DQMStore.cc:680
virtual TH1 * getTH1() const
~MtdTracksHarvester() override
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
MonitorElement * meBtlPhiEff_
virtual int getNbinsX() const
get # of bins in X-axis
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
MonitorElement * meMVAEtaSelEff_
MonitorElement * meEtlEtaEff_[2]
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
virtual double getBinContent(int binx) const
get content of bin (1-D)