CMS 3D CMS Logo

ElectronMcSignalPostValidator.cc
Go to the documentation of this file.
1 
5 
8  // histos bining and limits
9 
10  edm::ParameterSet histosSet = conf.getParameter<edm::ParameterSet>("histosCfg");
11 
12  set_EfficiencyFlag = histosSet.getParameter<bool>("EfficiencyFlag");
13  set_StatOverflowFlag = histosSet.getParameter<bool>("StatOverflowFlag");
14 }
15 
17 
19  setBookIndex(-1);
20  setBookPrefix("h_ele");
23 
24  edm::LogInfo("ElectronMcSignalPostValidator::finalize") << "efficiency calculation";
26  iBooker, iGetter, "etaEff", "mc_Eta_matched", "mc_Eta", "#eta", "Efficiency", "Efficiency vs gen eta");
27  bookH1andDivide(iBooker,
28  iGetter,
29  "etaEff_Extended",
30  "mc_Eta_Extended_matched",
31  "mc_Eta_Extended",
32  "#eta",
33  "Efficiency",
34  "Efficiency vs gen eta"); //Efficiency vs gen eta --- Eta of matched electrons
35  bookH1andDivide(iBooker, iGetter, "zEff", "mc_Z_matched", "mc_Z", "z (cm)", "Efficiency", "");
36  bookH1andDivide(iBooker, iGetter, "absetaEff", "mc_AbsEta_matched", "mc_AbsEta", "|#eta|", "Efficiency", "");
37  bookH1andDivide(iBooker,
38  iGetter,
39  "absetaEff_Extended",
40  "mc_AbsEta_Extended_matched",
41  "mc_AbsEta_Extended",
42  "|#eta|",
43  "Efficiency",
44  "");
45  bookH1andDivide(iBooker, iGetter, "ptEff", "mc_Pt_matched", "mc_Pt", "p_{T} (GeV/c)", "Efficiency", "");
46  bookH1andDivide(iBooker, iGetter, "phiEff", "mc_Phi_matched", "mc_Phi", "#phi (rad)", "Efficiency", "");
47  //bookH2andDivide(iBooker, iGetter, "ptEtaEff", "mc_PtEta_matched", "mc_PtEta", "#eta", "p_{T} (GeV/c)", "");
48 
49  edm::LogInfo("ElectronMcSignalPostValidator::finalize") << "q-misid calculation";
50  bookH1andDivide(iBooker, iGetter, "etaQmisid", "mc_Eta_matched_qmisid", "mc_Eta", "#eta", "q misId", "");
51  bookH1andDivide(iBooker, iGetter, "zQmisid", "mc_Z_matched_qmisid", "mc_Z", "z (cm)", "q misId", "");
52  bookH1andDivide(iBooker, iGetter, "absetaQmisid", "mc_AbsEta_matched_qmisid", "mc_AbsEta", "|#eta|", "q misId", "");
53  bookH1andDivide(iBooker, iGetter, "ptQmisid", "mc_Pt_matched_qmisid", "mc_Pt", "p_{T} (GeV/c)", "q misId", "");
54 
55  edm::LogInfo("ElectronMcSignalPostValidator::finalize") << "all reco electrons";
56  bookH1andDivide(iBooker, iGetter, "etaEff_all", "vertexEta_all", "h_mc_Eta", "#eta", "N_{rec}/N_{gen}", "");
57  bookH1andDivide(iBooker, iGetter, "ptEff_all", "vertexPt_all", "h_mc_Pt", "p_{T} (GeV/c)", "N_{rec}/N_{gen}", "");
58 
59  edm::LogInfo("ElectronMcSignalPostValidator::finalize") << "classes";
60  bookH1andDivide(iBooker,
61  iGetter,
62  "eta_goldenFrac",
63  "eta_golden",
64  "h_ele_eta",
65  "|#eta|",
66  "Fraction of electrons",
67  "fraction of golden electrons vs eta");
68  bookH1andDivide(iBooker,
69  iGetter,
70  "eta_bbremFrac",
71  "eta_bbrem",
72  "h_ele_eta",
73  "|#eta|",
74  "Fraction of electrons",
75  "fraction of big brem electrons vs eta");
76  bookH1andDivide(iBooker,
77  iGetter,
78  "eta_showerFrac",
79  "eta_shower",
80  "h_ele_eta",
81  "|#eta|",
82  "Fraction of electrons",
83  "fraction of showering electrons vs eta");
84 
85  // fbrem
86  MonitorElement* p1_ele_fbremVsEta_mean = get(iGetter, "fbremvsEtamean");
87  TAxis* etaAxis = p1_ele_fbremVsEta_mean->getTProfile()->GetXaxis();
89  iBooker, "xOverx0VsEta", "mean X/X_0 vs eta", etaAxis->GetNbins(), etaAxis->GetXmin(), etaAxis->GetXmax());
90  for (int ibin = 1; ibin < etaAxis->GetNbins() + 1; ibin++) {
91  double xOverX0 = 0.;
92  if (p1_ele_fbremVsEta_mean->getBinContent(ibin) > 0.) {
93  xOverX0 = -log(p1_ele_fbremVsEta_mean->getBinContent(ibin));
94  }
95  h1_ele_xOverX0VsEta->setBinContent(ibin, xOverX0);
96  }
97 
98  MonitorElement* h1_ele_provenance = get(iGetter, "provenance");
99  if (h1_ele_provenance->getBinContent(3) > 0) {
100  h1_ele_provenance->getTH1F()->Scale(1. / h1_ele_provenance->getBinContent(3));
101  }
102  MonitorElement* h1_ele_provenance_barrel = get(iGetter, "provenance_barrel");
103  if (h1_ele_provenance_barrel->getBinContent(3) > 0) {
104  h1_ele_provenance_barrel->getTH1F()->Scale(1. / h1_ele_provenance_barrel->getBinContent(3));
105  }
106  MonitorElement* h1_ele_provenance_endcaps = get(iGetter, "provenance_endcaps");
107  if (h1_ele_provenance_endcaps->getBinContent(3) > 0) {
108  h1_ele_provenance_endcaps->getTH1F()->Scale(1. / h1_ele_provenance_endcaps->getBinContent(3));
109  }
110 
111  MonitorElement* h1_ele_provenance_Extended = get(iGetter, "provenance_Extended");
112  if (h1_ele_provenance_Extended->getBinContent(3) > 0) {
113  h1_ele_provenance_Extended->getTH1F()->Scale(1. / h1_ele_provenance_Extended->getBinContent(3));
114  }
115 
116  // profiles from 2D histos
117  profileX(iBooker,
118  iGetter,
119  "scl_EoEtrueVsrecOfflineVertices",
120  "E/Etrue vs number of primary vertices",
121  "N_{primary vertices}",
122  "E/E_{true}",
123  0.8);
124  profileX(iBooker,
125  iGetter,
126  "scl_EoEtrueVsrecOfflineVertices_Extended",
127  "E/Etrue vs number of primary vertices, 2.5<|eta|<3",
128  "N_{primary vertices}",
129  "E/E_{true}",
130  0.8);
131  profileX(iBooker,
132  iGetter,
133  "scl_EoEtrueVsrecOfflineVertices_barrel",
134  "E/Etrue vs number of primary vertices , barrel",
135  "N_{primary vertices}",
136  "E/E_{true}",
137  0.8);
138  profileX(iBooker,
139  iGetter,
140  "scl_EoEtrueVsrecOfflineVertices_endcaps",
141  "E/Etrue vs number of primary vertices , endcaps",
142  "N_{primary vertices}",
143  "E/E_{true}",
144  0.8);
145 
146  profileX(iBooker, iGetter, "PoPtrueVsEta", "mean ele momentum / gen momentum vs eta", "#eta", "<P/P_{gen}>");
147  profileX(iBooker, iGetter, "PoPtrueVsEta_Extended", "mean ele momentum / gen momentum vs eta", "#eta", "<P/P_{gen}>");
148  profileX(iBooker, iGetter, "PoPtrueVsPhi", "mean ele momentum / gen momentum vs phi", "#phi (rad)", "<P/P_{gen}>");
149  profileX(iBooker, iGetter, "sigmaIetaIetaVsPt", "SigmaIetaIeta vs pt", "p_{T} (GeV/c)", "SigmaIetaIeta");
150  profileX(iBooker,
151  iGetter,
152  "EoEtruePfVsEg",
153  "mean mustache SC/true energy vs final SC/true energy",
154  "E_{final SC}/E_{gen}",
155  "E_{mustache}/E_{gen}");
156  profileY(iBooker,
157  iGetter,
158  "EoEtruePfVsEg",
159  "mean mustache SC/true energy vs final SC/true energy",
160  "E_{final SC}/E_{gen}",
161  "E_{mustache}/E_{gen}");
162  profileX(iBooker, iGetter, "EtaMnEtaTrueVsEta", "mean ele eta - gen eta vs eta", "#eta", "<#eta_{rec} - #eta_{gen}>");
163  profileX(
164  iBooker, iGetter, "EtaMnEtaTrueVsPhi", "mean ele eta - gen eta vs phi", "#phi (rad)", "<#eta_{rec} - #eta_{gen}>");
165  profileX(
166  iBooker, iGetter, "PhiMnPhiTrueVsEta", "mean ele phi - gen phi vs eta", "#eta", "<#phi_{rec} - #phi_{gen}> (rad)");
167  profileX(iBooker, iGetter, "PhiMnPhiTrueVsPhi", "mean ele phi - gen phi vs phi", "#phi (rad)", "");
168  profileX(iBooker, iGetter, "vertexPtVsEta", "mean ele transverse momentum vs eta", "#eta", "<p_{T}> (GeV/c)");
169  profileX(iBooker, iGetter, "vertexPtVsPhi", "mean ele transverse momentum vs phi", "#phi (rad)", "<p_{T}> (GeV/c)");
170  profileX(iBooker, iGetter, "EoPVsEta", "mean ele E/p vs eta", "#eta", "<E/P_{vertex}>");
171  profileX(iBooker, iGetter, "EoPVsEta_Extended", "mean ele E/p vs eta", "#eta", "<E/P_{vertex}>");
172  profileX(iBooker, iGetter, "EoPVsPhi", "mean ele E/p vs phi", "#phi (rad)", "<E/P_{vertex}>");
173  profileX(iBooker, iGetter, "EoPoutVsEta", "mean ele E/pout vs eta", "#eta", "<E_{seed}/P_{out}>");
174  profileX(iBooker, iGetter, "EoPoutVsPhi", "mean ele E/pout vs phi", "#phi (rad)", "<E_{seed}/P_{out}>");
175  profileX(iBooker, iGetter, "EeleOPoutVsEta", "mean ele Eele/pout vs eta", "#eta", "<E_{ele}/P_{out}>");
176  profileX(iBooker, iGetter, "EeleOPoutVsPhi", "mean ele Eele/pout vs phi", "#phi (rad)", "<E_{ele}/P_{out}>");
177  profileX(iBooker, iGetter, "HoEVsEta", "mean ele H/E vs eta", "#eta", "<H/E>");
178  profileX(iBooker, iGetter, "HoEVsPhi", "mean ele H/E vs phi", "#phi (rad)", "<H/E>");
179  profileX(iBooker, iGetter, "chi2VsEta", "mean ele track chi2 vs eta", "#eta", "<#Chi^{2}>");
180  profileX(iBooker, iGetter, "chi2VsPhi", "mean ele track chi2 vs phi", "#phi (rad)", "<#Chi^{2}>");
181  profileX(iBooker, iGetter, "ambiguousTracksVsEta", "mean ele # ambiguous tracks vs eta", "#eta", "<N_{ambiguous}>");
182  profileX(iBooker, iGetter, "foundHitsVsEta", "mean ele track # found hits vs eta", "#eta", "<N_{hits}>");
183  profileX(iBooker, iGetter, "foundHitsVsEta_Extended", "mean ele track # found hits vs eta", "#eta", "<N_{hits}>");
184  profileX(iBooker, iGetter, "foundHitsVsEta_mAOD", "mean ele track # found hits vs eta", "#eta", "<N_{hits}>");
185  profileX(iBooker, iGetter, "foundHitsVsPhi", "mean ele track # found hits vs phi", "#phi (rad)", "<N_{hits}>");
186  profileX(iBooker, iGetter, "lostHitsVsEta", "mean ele track # lost hits vs eta", "#eta", "<N_{hits}>");
187  profileX(iBooker, iGetter, "lostHitsVsPhi", "mean ele track # lost hits vs phi", "#phi (rad)", "<N_{hits}>");
188  profileX(iBooker, iGetter, "vertexTIPVsEta", "mean tip (wrt gen vtx) vs eta", "#eta", "<TIP> (cm)");
189  profileX(iBooker, iGetter, "vertexTIPVsPhi", "mean tip (wrt gen vtx) vs phi", "#phi", "<TIP> (cm)");
190  profileX(iBooker, iGetter, "vertexTIPVsPt", "mean tip (wrt gen vtx) vs phi", "p_{T} (GeV/c)", "<TIP> (cm)");
191  profileX(iBooker,
192  iGetter,
193  "seedDphi2_VsEta",
194  "mean ele seed dphi 2nd layer vs eta",
195  "#eta",
196  "<#phi_{pred} - #phi_{hit}, 2nd layer> (rad)",
197  -0.004,
198  0.004);
199  profileX(iBooker,
200  iGetter,
201  "seedDphi2_VsPt",
202  "mean ele seed dphi 2nd layer vs pt",
203  "p_{T} (GeV/c)",
204  "<#phi_{pred} - #phi_{hit}, 2nd layer> (rad)",
205  -0.004,
206  0.004);
207  profileX(iBooker,
208  iGetter,
209  "seedDrz2_VsEta",
210  "mean ele seed dr(dz) 2nd layer vs eta",
211  "#eta",
212  "<r(z)_{pred} - r(z)_{hit}, 2nd layer> (cm)",
213  -0.15,
214  0.15);
215  profileX(iBooker,
216  iGetter,
217  "seedDrz2_VsPt",
218  "mean ele seed dr(dz) 2nd layer vs pt",
219  "p_{T} (GeV/c)",
220  "<r(z)_{pred} - r(z)_{hit}, 2nd layer> (cm)",
221  -0.15,
222  0.15);
223  profileX(iBooker,
224  iGetter,
225  "seedDphi2Pos_VsEta",
226  "mean ele seed dphi 2nd layer positron vs eta",
227  "#eta",
228  "<#phi_{pred} - #phi_{hit}, 2nd layer> (rad)",
229  -0.004,
230  0.004);
231  profileX(iBooker,
232  iGetter,
233  "seedDphi2Pos_VsPt",
234  "mean ele seed dphi 2nd layer positron vs pt",
235  "p_{T} (GeV/c)",
236  "<#phi_{pred} - #phi_{hit}, 2nd layer> (rad)",
237  -0.004,
238  0.004);
239  profileX(iBooker,
240  iGetter,
241  "seedDrz2Pos_VsEta",
242  "mean ele seed dr(dz) 2nd layer positron vs eta",
243  "#eta",
244  "<r(z)_{pred} - r(z)_{hit}, 2nd layer> (cm)",
245  -0.15,
246  0.15);
247  profileX(iBooker,
248  iGetter,
249  "seedDrz2Pos_VsPt",
250  "mean ele seed dr(dz) 2nd layer positron vs pt",
251  "p_{T} (GeV/c)",
252  "<r(z)_{pred} - r(z)_{hit}, 2nd layer> (cm)",
253  -0.15,
254  0.15);
255 
256 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
MonitorElement * profileX(DQMStore::IBooker &iBooker, DQMStore::IGetter &iGetter, MonitorElement *me2d, const std::string &title="", const std::string &titleX="", const std::string &titleY="", Double_t minimum=-1111, Double_t maximum=-1111)
virtual TProfile * getTProfile() const
MonitorElement * bookH1withSumw2(DQMStore::IBooker &, const std::string &name, const std::string &title, int nchX, double lowX, double highX, const std::string &titleX="", const std::string &titleY="Events", Option_t *option="E1 P")
void setBookPrefix(const std::string &)
ElectronMcSignalPostValidator(const edm::ParameterSet &conf)
MonitorElement * bookH1andDivide(DQMStore::IBooker &iBooker, DQMStore::IGetter &, const std::string &name, MonitorElement *num, MonitorElement *denom, const std::string &titleX, const std::string &titleY, const std::string &title="")
void finalize(DQMStore::IBooker &iBooker, DQMStore::IGetter &iGetter) override
MonitorElement * profileY(DQMStore::IBooker &iBooker, DQMStore::IGetter &iGetter, MonitorElement *me2d, const std::string &title="", const std::string &titleX="", const std::string &titleY="", Double_t minimum=-1111, Double_t maximum=-1111)
void setBookStatOverflowFlag(const bool &)
void setBookEfficiencyFlag(const bool &)
Log< level::Info, false > LogInfo
virtual void setBinContent(int binx, double content)
set content of bin (1-D)
virtual TH1F * getTH1F() const
virtual double getBinContent(int binx) const
get content of bin (1-D)