CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_9_patch3/src/Validation/RecoEgamma/plugins/ElectronMcSignalPostValidator.cc

Go to the documentation of this file.
00001 
00002 #include "Validation/RecoEgamma/plugins/ElectronMcSignalPostValidator.h"
00003 #include "DQMServices/Core/interface/MonitorElement.h"
00004 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00005 
00006 ElectronMcSignalPostValidator::ElectronMcSignalPostValidator( const edm::ParameterSet & conf )
00007  : ElectronDqmAnalyzerBase(conf)
00008  {}
00009 
00010 ElectronMcSignalPostValidator::~ElectronMcSignalPostValidator()
00011  {}
00012 
00013 void ElectronMcSignalPostValidator::book()
00014  { setBookIndex(-1) ; }
00015 
00016 void ElectronMcSignalPostValidator::finalize()
00017  {
00018   setBookPrefix("h_ele") ;
00019 
00020   edm::LogInfo("ElectronMcSignalPostValidator::finalize") << "efficiency calculation" ;
00021   bookH1andDivide("etaEff","mc_Eta_matched","mc_Eta","#eta","Efficiency","");
00022   bookH1andDivide("zEff","mc_Z_matched","mc_Z","z (cm)","Efficiency","");
00023   bookH1andDivide("absetaEff","mc_AbsEta_matched","mc_AbsEta","|#eta|","Efficiency");
00024   bookH1andDivide("ptEff","mc_Pt_matched","mc_Pt","p_{T} (GeV/c)","Efficiency");
00025   bookH1andDivide("phiEff","mc_Phi_matched","mc_Phi","#phi (rad)","Efficiency");
00026   bookH2andDivide("ptEtaEff","mc_PtEta_matched","mc_PtEta","#eta","p_{T} (GeV/c)");
00027 
00028   edm::LogInfo("ElectronMcSignalPostValidator::finalize") << "q-misid calculation" ;
00029   bookH1andDivide("etaQmisid","mc_Eta_matched_qmisid","mc_Eta","#eta","q misId","");
00030   bookH1andDivide("zQmisid","mc_Z_matched_qmisid","mc_Z","z (cm)","q misId","");
00031   bookH1andDivide("absetaQmisid","mc_AbsEta_matched_qmisid","mc_AbsEta","|#eta|","q misId");
00032   bookH1andDivide("ptQmisid","mc_Pt_matched_qmisid","mc_Pt","p_{T} (GeV/c)","q misId");
00033 
00034   edm::LogInfo("ElectronMcSignalPostValidator::finalize") << "all reco electrons" ;
00035   bookH1andDivide("etaEff_all","vertexEta_all","h_mc_Eta","#eta","N_{rec}/N_{gen}","");
00036   bookH1andDivide("ptEff_all","vertexPt_all","h_mc_Pt","p_{T} (GeV/c)","N_{rec}/N_{gen}","");
00037 
00038   edm::LogInfo("ElectronMcSignalPostValidator::finalize") << "classes" ;
00039   bookH1andDivide("eta_goldenFrac","eta_golden","h_ele_eta","|#eta|","Fraction of electrons","fraction of golden electrons vs eta");
00040   bookH1andDivide("eta_bbremFrac" ,"eta_bbrem","h_ele_eta","|#eta|","Fraction of electrons","fraction of big brem electrons vs eta");
00041   bookH1andDivide("eta_showerFrac","eta_shower","h_ele_eta","|#eta|","Fraction of electrons","fraction of showering electrons vs eta");
00042 
00043   // fbrem
00044   MonitorElement * p1_ele_fbremVsEta_mean = get("fbremvsEtamean") ;
00045   TAxis * etaAxis = p1_ele_fbremVsEta_mean->getTProfile()->GetXaxis() ;
00046   MonitorElement * h1_ele_xOverX0VsEta = bookH1withSumw2("xOverx0VsEta","mean X/X_0 vs eta",etaAxis->GetNbins(),etaAxis->GetXmin(),etaAxis->GetXmax());
00047   for (int ibin=1;ibin<etaAxis->GetNbins()+1;ibin++) {
00048     double xOverX0 = 0.;
00049     if (p1_ele_fbremVsEta_mean->getBinContent(ibin)>0.)
00050      { xOverX0 = -log(p1_ele_fbremVsEta_mean->getBinContent(ibin)) ; }
00051     h1_ele_xOverX0VsEta->setBinContent(ibin,xOverX0) ;
00052   }
00053 
00054   // profiles from 2D histos
00055   profileX("PoPtrueVsEta","mean ele momentum / gen momentum vs eta","#eta","<P/P_{gen}>");
00056   profileX("PoPtrueVsPhi","mean ele momentum / gen momentum vs phi","#phi (rad)","<P/P_{gen}>");
00057   profileX("EoEtruePfVsEg","mean pflow sc energy / true energy vs e/g sc energy","E/E_{gen} (e/g)","<E/E_{gen}> (pflow)") ;
00058   profileY("EoEtruePfVsEg","mean e/g sc energy / true energy vs pflow sc energy","E/E_{gen} (pflow)","<E/E_{gen}> (eg)") ;
00059   profileX("EtaMnEtaTrueVsEta","mean ele eta - gen eta vs eta","#eta","<#eta_{rec} - #eta_{gen}>");
00060   profileX("EtaMnEtaTrueVsPhi","mean ele eta - gen eta vs phi","#phi (rad)","<#eta_{rec} - #eta_{gen}>");
00061   profileX("PhiMnPhiTrueVsEta","mean ele phi - gen phi vs eta","#eta","<#phi_{rec} - #phi_{gen}> (rad)");
00062   profileX("PhiMnPhiTrueVsPhi","mean ele phi - gen phi vs phi","#phi (rad)","");
00063   profileX("vertexPtVsEta","mean ele transverse momentum vs eta","#eta","<p_{T}> (GeV/c)");
00064   profileX("vertexPtVsPhi","mean ele transverse momentum vs phi","#phi (rad)","<p_{T}> (GeV/c)");
00065   profileX("EoPVsEta","mean ele E/p vs eta","#eta","<E/P_{vertex}>");
00066   profileX("EoPVsPhi","mean ele E/p vs phi","#phi (rad)","<E/P_{vertex}>");
00067   profileX("EoPoutVsEta","mean ele E/pout vs eta","#eta","<E_{seed}/P_{out}>");
00068   profileX("EoPoutVsPhi","mean ele E/pout vs phi","#phi (rad)","<E_{seed}/P_{out}>");
00069   profileX("EeleOPoutVsEta","mean ele Eele/pout vs eta","#eta","<E_{ele}/P_{out}>");
00070   profileX("EeleOPoutVsPhi","mean ele Eele/pout vs phi","#phi (rad)","<E_{ele}/P_{out}>");
00071   profileX("HoEVsEta","mean ele H/E vs eta","#eta","<H/E>");
00072   profileX("HoEVsPhi","mean ele H/E vs phi","#phi (rad)","<H/E>");
00073   profileX("chi2VsEta","mean ele track chi2 vs eta","#eta","<#Chi^{2}>");
00074   profileX("chi2VsPhi","mean ele track chi2 vs phi","#phi (rad)","<#Chi^{2}>");
00075   profileX("ambiguousTracksVsEta","mean ele # ambiguous tracks  vs eta","#eta","<N_{ambiguous}>");
00076   profileX("foundHitsVsEta","mean ele track # found hits vs eta","#eta","<N_{hits}>");
00077   profileX("foundHitsVsPhi","mean ele track # found hits vs phi","#phi (rad)","<N_{hits}>");
00078   profileX("lostHitsVsEta","mean ele track # lost hits vs eta","#eta","<N_{hits}>");
00079   profileX("lostHitsVsPhi","mean ele track # lost hits vs phi","#phi (rad)","<N_{hits}>");
00080   profileX("vertexTIPVsEta","mean tip (wrt gen vtx) vs eta","#eta","<TIP> (cm)");
00081   profileX("vertexTIPVsPhi","mean tip (wrt gen vtx) vs phi","#phi","<TIP> (cm)");
00082   profileX("vertexTIPVsPt","mean tip (wrt gen vtx) vs phi","p_{T} (GeV/c)","<TIP> (cm)");
00083   profileX("seedDphi2_VsEta","mean ele seed dphi 2nd layer vs eta","#eta","<#phi_{pred} - #phi_{hit}, 2nd layer> (rad)",-0.004,0.004);
00084   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);
00085   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);
00086   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);
00087   profileX("seedDphi2Pos_VsEta","mean ele seed dphi 2nd layer positron vs eta","#eta","<#phi_{pred} - #phi_{hit}, 2nd layer> (rad)",-0.004,0.004);
00088   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);
00089   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);
00090   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);
00091 
00092 //  // investigation
00093 //  TH2F * h2 = get("PoPtrueVsEta")->getTH2F() ;
00094 //  std::cout<<"H2   entries : "<<h2->GetEntries()<<std::endl ;
00095 //  std::cout<<"H2 effective entries : "<<h2->GetEffectiveEntries()<<std::endl ;
00096 //  Int_t ix, nx = h2->GetNbinsX(), iy, ny = h2->GetNbinsY(), is, nu = 0, no = 0, nb = 0 ;
00097 //  for ( iy = 0 ; iy<=(ny+1) ; ++iy )
00098 //    for ( ix = 0 ; ix<=(nx+1) ; ++ix )
00099 //     {
00100 //      is = iy*(nx+2) + ix ;
00101 //      if (h2->IsBinUnderflow(is)) ++nu ;
00102 //      if (h2->IsBinOverflow(is)) ++no ;
00103 //     }
00104 //  ix = 0 ;
00105 //  for ( iy = 0 ; iy<=(ny+1) ; ++iy )
00106 //   {
00107 //    is = iy*(nx+2) + ix ;
00108 //    nb += (*h2->GetSumw2())[is] ;
00109 //   }
00110 //  ix = nx+1 ;
00111 //  for ( iy = 0 ; iy<=(ny+1) ; ++iy )
00112 //   {
00113 //    is = iy*(nx+2) + ix ;
00114 //    nb += (*h2->GetSumw2())[is] ;
00115 //   }
00116 //  for ( ix = 1 ; ix<=nx ; ++ix )
00117 //   {
00118 //    iy = 0 ;
00119 //    is = iy*(nx+2) + ix ;
00120 //    nb += (*h2->GetSumw2())[is] ;
00121 //    iy = ny+1 ;
00122 //    is = iy*(nx+2) + ix ;
00123 //    nb += (*h2->GetSumw2())[is] ;
00124 //   }
00125 //  std::cout<<"H2   nx      : "<<nx<<std::endl ;
00126 //  std::cout<<"H2   ny      : "<<ny<<std::endl ;
00127 //  std::cout<<"H2   nsumw2  : "<<(*h2->GetSumw2()).fN<<std::endl ;
00128 //  std::cout<<"H2   nu      : "<<nu<<std::endl ;
00129 //  std::cout<<"H2   no      : "<<no<<std::endl ;
00130 //  std::cout<<"H2   outside : "<<nb<<std::endl ;
00131 //  std::cout<<"PFX  entries : "<<h2->ProfileX()->GetEntries()<<std::endl ;
00132 //  std::cout<<"PFX effective entries : "<<h2->ProfileX()->GetEffectiveEntries()<<std::endl ;
00133 }
00134 
00135