CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ElectronMcFakePostValidator.cc
Go to the documentation of this file.
1 
5 
8  {}
9 
11  {}
12 
14  { setBookIndex(-1) ; }
15 
17  {
18  setBookPrefix("h_ele") ;
19 
20  edm::LogInfo("ElectronMcFakePostValidator::finalize") << "efficiency calculation " ;
21  bookH1andDivide("etaEff","matchingObjectEta_matched","matchingObject_eta","#eta","Efficiency");
22  bookH1andDivide("zEff","matchingObjectZ_matched","matchingObject_z","z (cm)","Efficiency");
23  bookH1andDivide("absetaEff","matchingObjectAbsEta_matched","matchingObject_abseta","|#eta|","Efficiency");
24  bookH1andDivide("ptEff","matchingObjectPt_matched","matchingObject_Pt","p_{T} (GeV/c)","Efficiency");
25  bookH1andDivide("phiEff","matchingObjectPhi_matched","matchingObject_phi","#phi (rad)","Efficiency");
26 // bookH2andDivide("ptEtaEff","matchingObjectPtEta_matched","matchingObjectPtEta","#eta","p_{T} (GeV/c)");
27 //
28 // std::cout << "[ElectronMcFakePostValidator] q-misid calculation " << std::endl;
29 // bookH1andDivide("etaQmisid","matchingObjectEta_matched_qmisid","h_simEta","#eta","q misId","",true);
30 // bookH1andDivide("zQmisid","matchingObjectZ_matched_qmisid","h_simZ","z (cm)","q misId","",true);
31 // bookH1andDivide("absetaQmisid","matchingObjectAbsEta_matched_qmisid","h_simAbsEta","|#eta|","q misId");
32 // bookH1andDivide("ptQmisid","matchingObjectPt_matched_qmisid","h_simPt","p_{T} (GeV/c)","q misId");
33 
34  edm::LogInfo("ElectronMcFakePostValidator::finalize") << "all reco electrons " ;
35  bookH1andDivide("etaEff_all","vertexEta_all","matchingObject_eta","#eta","N_{rec}/N_{gen}");
36  bookH1andDivide("ptEff_all", "vertexPt_all","matchingObject_Pt","p_{T} (GeV/c)","N_{rec}/N_{gen}");
37 
38  edm::LogInfo("ElectronMcFakePostValidator::finalize") << "classes" ;
39  bookH1andDivide("eta_goldenFrac","eta_golden","h_ele_eta","|#eta|","Fraction of electrons","fraction of golden electrons vs eta");
40  bookH1andDivide("eta_bbremFrac" ,"eta_bbrem", "h_ele_eta","|#eta|","Fraction of electrons","fraction of big brem electrons vs eta");
41 // bookH1andDivide("eta_narrowFrac","eta_narrow","h_ele_eta","|#eta|","Fraction of electrons","fraction of narrow electrons vs eta");
42  bookH1andDivide("eta_showerFrac","eta_shower","h_ele_eta","|#eta|","Fraction of electrons","fraction of showering electrons vs eta");
43 
44  // fbrem
45  MonitorElement * p1_ele_fbremVsEta_mean = get("fbremvsEtamean") ;
46  TAxis * etaAxis = p1_ele_fbremVsEta_mean->getTProfile()->GetXaxis() ;
47  MonitorElement * h1_ele_xOverX0VsEta = bookH1withSumw2("xOverx0VsEta","mean X/X_0 vs eta",etaAxis->GetNbins(),etaAxis->GetXmin(),etaAxis->GetXmax());
48  for (int ibin=1;ibin<etaAxis->GetNbins()+1;ibin++) {
49  double xOverX0 = 0.;
50  if (p1_ele_fbremVsEta_mean->getBinContent(ibin)>0.)
51  { xOverX0 = -log(p1_ele_fbremVsEta_mean->getBinContent(ibin)) ; }
52  h1_ele_xOverX0VsEta->setBinContent(ibin,xOverX0) ;
53  }
54 
55  // profiles from 2D histos
56  profileX("PoPmatchingObjectVsEta","","#eta","<P/P_{gen}>");
57  profileX("PoPmatchingObjectVsPhi","","#phi (rad)","<P/P_{gen}>");
58  profileX("EtaMnEtamatchingObjectVsEta","","#eta","<#eta_{rec} - #eta_{gen}>");
59  profileX("EtaMnEtamatchingObjectVsPhi","","#phi (rad)","<#eta_{rec} - #eta_{gen}>");
60  profileX("PhiMnPhimatchingObjectVsEta","","#eta","<#phi_{rec} - #phi_{gen}> (rad)");
61  profileX("PhiMnPhimatchingObjectVsPhi","","#phi (rad)","");
62  profileX("vertexPtVsEta","mean ele transverse momentum vs eta","#eta","<p_{T}> (GeV/c)");
63  profileX("vertexPtVsPhi","mean ele transverse momentum vs phi","#phi (rad)","<p_{T}> (GeV/c)");
64  profileX("EoPVsEta","mean ele E/p vs eta","#eta","<E/P_{vertex}>");
65  profileX("EoPVsPhi","mean ele E/p vs phi","#phi (rad)","<E/P_{vertex}>");
66  profileX("EoPoutVsEta","mean ele E/pout vs eta","#eta","<E_{seed}/P_{out}>");
67  profileX("EoPoutVsPhi","mean ele E/pout vs phi","#phi (rad)","<E_{seed}/P_{out}>");
68  profileX("EeleOPoutVsEta","mean ele Eele/pout vs eta","#eta","<E_{ele}/P_{out}>");
69  profileX("EeleOPoutVsPhi","mean ele Eele/pout vs phi","#phi (rad)","<E_{ele}/P_{out}>");
70  profileX("HoEVsEta","mean ele H/E vs eta","#eta","<H/E>");
71  profileX("HoEVsPhi","mean ele H/E vs phi","#phi (rad)","<H/E>");
72  profileX("chi2VsEta","mean ele track chi2 vs eta","#eta","<#Chi^{2}>");
73  profileX("chi2VsPhi","mean ele track chi2 vs phi","#phi (rad)","<#Chi^{2}>");
74  profileX("foundHitsVsEta","mean ele track # found hits vs eta","#eta","<N_{hits}>");
75  profileX("foundHitsVsPhi","mean ele track # found hits vs phi","#phi (rad)","<N_{hits}>");
76  profileX("lostHitsVsEta","mean ele track # lost hits vs eta","#eta","<N_{hits}>");
77  profileX("lostHitsVsPhi","mean ele track # lost hits vs phi","#phi (rad)","<N_{hits}>");
78  profileX("vertexTIPVsEta","mean tip (wrt gen vtx) vs eta","#eta","<TIP> (cm)");
79  profileX("vertexTIPVsPhi","mean tip (wrt gen vtx) vs phi","#phi","<TIP> (cm)");
80  profileX("vertexTIPVsPt","mean tip (wrt gen vtx) vs phi","p_{T} (GeV/c)","<TIP> (cm)");
81  profileX("seedDphi2_VsEta","mean ele seed dphi 2nd layer vs eta","#eta","<#phi_{pred} - #phi_{hit}, 2nd layer> (rad)",-0.004,0.004);
82  profileX("seedDphi2_VsPt","mean ele seed dphi 2nd layer vs pt","p_{T} (GeV/c)","<#phi_{pred} - #phi_{hit}, 2nd layer> (rad)",-0.004,0.004);
83  profileX("seedDrz2_VsEta","mean ele seed dr(dz) 2nd layer vs eta","#eta","<r(z)_{pred} - r(z)_{hit}, 2nd layer> (cm)",-0.15,0.15);
84  profileX("seedDrz2_VsPt","mean ele seed dr(dz) 2nd layer vs pt","p_{T} (GeV/c)","<r(z)_{pred} - r(z)_{hit}, 2nd layer> (cm)",-0.15,0.15);
85  profileX("seedDphi2Pos_VsEta","mean ele seed dphi 2nd layer positron vs eta","#eta","<#phi_{pred} - #phi_{hit}, 2nd layer> (rad)",-0.004,0.004);
86  profileX("seedDphi2Pos_VsPt","mean ele seed dphi 2nd layer positron vs pt","p_{T} (GeV/c)","<#phi_{pred} - #phi_{hit}, 2nd layer> (rad)",-0.004,0.004);
87  profileX("seedDrz2Pos_VsEta","mean ele seed dr(dz) 2nd layer positron vs eta","#eta","<r(z)_{pred} - r(z)_{hit}, 2nd layer> (cm)",-0.15,0.15);
88  profileX("seedDrz2Pos_VsPt","mean ele seed dr(dz) 2nd layer positron vs pt","p_{T} (GeV/c)","<r(z)_{pred} - r(z)_{hit}, 2nd layer> (cm)",-0.15,0.15);
89  }
90 
91 
void setBinContent(int binx, double content)
set content of bin (1-D)
ElectronMcFakePostValidator(const edm::ParameterSet &conf)
MonitorElement * bookH1withSumw2(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")
MonitorElement * profileX(MonitorElement *me2d, const std::string &title="", const std::string &titleX="", const std::string &titleY="", Double_t minimum=-1111, Double_t maximum=-1111)
void setBookPrefix(const std::string &)
tuple conf
Definition: dbtoconf.py:185
double getBinContent(int binx) const
get content of bin (1-D)
TProfile * getTProfile(void) const
MonitorElement * bookH1andDivide(const std::string &name, MonitorElement *num, MonitorElement *denom, const std::string &titleX, const std::string &titleY, const std::string &title="")