CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_13_patch3/src/DQMOffline/Muon/src/MuonRecoAnalyzer.cc

Go to the documentation of this file.
00001 
00002 /*
00003  *  See header file for a description of this class.
00004  *
00005  *  $Date: 2010/02/15 11:24:36 $
00006  *  $Revision: 1.23 $
00007  *  \author G. Mila - INFN Torino
00008  *  updated: G. Hesketh, CERN
00009  */
00010 
00011 #include "DQMOffline/Muon/src/MuonRecoAnalyzer.h"
00012 
00013 #include "DataFormats/Common/interface/Handle.h"
00014 #include "DataFormats/MuonReco/interface/Muon.h"
00015 #include "DataFormats/MuonReco/interface/MuonFwd.h" 
00016 #include "DataFormats/MuonReco/interface/MuonEnergy.h"
00017 
00018 #include "DataFormats/TrackReco/interface/Track.h"
00019 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00020 
00021 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00022 
00023 #include <string>
00024 #include "TMath.h"
00025 using namespace std;
00026 using namespace edm;
00027 
00028 
00029 
00030 MuonRecoAnalyzer::MuonRecoAnalyzer(const edm::ParameterSet& pSet, MuonServiceProxy *theService):MuonAnalyzerBase(theService) {
00031   parameters = pSet;
00032 }
00033 
00034 
00035 MuonRecoAnalyzer::~MuonRecoAnalyzer() { }
00036 
00037 
00038 void MuonRecoAnalyzer::beginJob(DQMStore * dbe) {
00039 
00040   metname = "muRecoAnalyzer";
00041 
00042   LogTrace(metname)<<"[MuonRecoAnalyzer] Parameters initialization";
00043   dbe->setCurrentFolder("Muons/MuonRecoAnalyzer");
00044 
00045   muReco = dbe->book1D("muReco", "muon reconstructed tracks", 6, 1, 7);
00046   muReco->setBinLabel(1,"glb+tk+sta");
00047   muReco->setBinLabel(2,"glb+sta");
00048   muReco->setBinLabel(3,"tk+sta");
00049   muReco->setBinLabel(4,"tk");
00050   muReco->setBinLabel(5,"sta");
00051   muReco->setBinLabel(6,"calo");
00052 
00053   int binFactor = 4;
00054 
00055   // monitoring of eta parameter
00056   etaBin = parameters.getParameter<int>("etaBin");
00057   etaMin = parameters.getParameter<double>("etaMin");
00058   etaMax = parameters.getParameter<double>("etaMax");
00059   std::string histname = "GlbMuon_";
00060   etaGlbTrack.push_back(dbe->book1D(histname+"Glb_eta", "#eta_{GLB}", etaBin, etaMin, etaMax));
00061   etaGlbTrack.push_back(dbe->book1D(histname+"Tk_eta", "#eta_{TKfromGLB}", etaBin, etaMin, etaMax));
00062   etaGlbTrack.push_back(dbe->book1D(histname+"Sta_eta", "#eta_{STAfromGLB}", etaBin, etaMin, etaMax));
00063   etaResolution.push_back(dbe->book1D("Res_TkGlb_eta", "#eta_{TKfromGLB} - #eta_{GLB}", etaBin*binFactor, etaMin/3000, etaMax/3000));
00064   etaResolution.push_back(dbe->book1D("Res_GlbSta_eta", "#eta_{GLB} - #eta_{STAfromGLB}", etaBin*binFactor, etaMin/100, etaMax/100));
00065   etaResolution.push_back(dbe->book1D("Res_TkSta_eta", "#eta_{TKfromGLB} - #eta_{STAfromGLB}", etaBin*binFactor, etaMin/100, etaMax/100));
00066   etaResolution.push_back(dbe->book2D("ResVsEta_TkGlb_eta", "(#eta_{TKfromGLB} - #eta_{GLB}) vs #eta_{GLB}", etaBin, etaMin, etaMax, etaBin*binFactor, etaMin/3000, etaMax/3000));
00067   etaResolution.push_back(dbe->book2D("ResVsEta_GlbSta_eta", "(#eta_{GLB} - #eta_{STAfromGLB}) vs #eta_{GLB}", etaBin, etaMin, etaMax, etaBin*binFactor, etaMin/100, etaMax/100));
00068   etaResolution.push_back(dbe->book2D("ResVsEta_TkSta_eta", "(#eta_{TKfromGLB} - #eta_{STAfromGLB}) vs #eta_{TKfromGLB}", etaBin, etaMin, etaMax, etaBin*binFactor, etaMin/100, etaMax/100));
00069   etaPull = dbe->book1D("Pull_TkSta_eta", "#eta_{TKfromGLB} - #eta_{GLB} / error", 100,-10,10);
00070   etaTrack = dbe->book1D("TkMuon_eta", "#eta_{TK}", etaBin, etaMin, etaMax);
00071   etaStaTrack = dbe->book1D("StaMuon_eta", "#eta_{STA}", etaBin, etaMin, etaMax);
00072   etaEfficiency.push_back(dbe->book1D("StaEta", "#eta_{STAfromGLB}", etaBin, etaMin, etaMax));
00073   etaEfficiency.push_back(dbe->book1D("StaEta_ifCombinedAlso", "#eta_{STAfromGLB} if isGlb=true", etaBin, etaMin, etaMax));
00074 
00075   // monitoring of theta parameter
00076   thetaBin = parameters.getParameter<int>("thetaBin");
00077   thetaMin = parameters.getParameter<double>("thetaMin");
00078   thetaMax = parameters.getParameter<double>("thetaMax");
00079   thetaGlbTrack.push_back(dbe->book1D(histname+"Glb_theta", "#theta_{GLB}", thetaBin, thetaMin, thetaMax));
00080   thetaGlbTrack[0]->setAxisTitle("rad");
00081   thetaGlbTrack.push_back(dbe->book1D(histname+"Tk_theta", "#theta_{TKfromGLB}", thetaBin, thetaMin, thetaMax));
00082   thetaGlbTrack[1]->setAxisTitle("rad");
00083   thetaGlbTrack.push_back(dbe->book1D(histname+"Sta_theta", "#theta_{STAfromGLB}", thetaBin, thetaMin, thetaMax));
00084   thetaGlbTrack[2]->setAxisTitle("rad");
00085   thetaResolution.push_back(dbe->book1D("Res_TkGlb_theta", "#theta_{TKfromGLB} - #theta_{GLB}", thetaBin*binFactor, -(thetaMax/3000), thetaMax/3000));
00086   thetaResolution[0]->setAxisTitle("rad");
00087   thetaResolution.push_back(dbe->book1D("Res_GlbSta_theta", "#theta_{GLB} - #theta_{STAfromGLB}", thetaBin*binFactor,-(thetaMax/100), thetaMax/100));
00088   thetaResolution[1]->setAxisTitle("rad");
00089   thetaResolution.push_back(dbe->book1D("Res_TkSta_theta", "#theta_{TKfromGLB} - #theta_{STAfromGLB}", thetaBin*binFactor, -(thetaMax/100), thetaMax/100));
00090   thetaResolution[2]->setAxisTitle("rad");
00091   thetaResolution.push_back(dbe->book2D("ResVsTheta_TkGlb_theta", "(#theta_{TKfromGLB} - #theta_{GLB}) vs #theta_{GLB}", thetaBin, thetaMin, thetaMax, thetaBin*binFactor, -(thetaMax/3000), thetaMax/3000));
00092   thetaResolution[3]->setAxisTitle("rad",1);
00093   thetaResolution[3]->setAxisTitle("rad",2);
00094   thetaResolution.push_back(dbe->book2D("ResVsTheta_GlbSta_theta", "(#theta_{GLB} - #theta_{STAfromGLB}) vs #theta_{GLB}", thetaBin, thetaMin, thetaMax, thetaBin*binFactor, -(thetaMax/100), thetaMax/100));
00095   thetaResolution[4]->setAxisTitle("rad",1);
00096   thetaResolution[4]->setAxisTitle("rad",2);
00097   thetaResolution.push_back(dbe->book2D("ResVsTheta_TkSta_theta", "(#theta_{TKfromGLB} - #theta_{STAfromGLB}) vs #theta_{TKfromGLB}", thetaBin, thetaMin, thetaMax, thetaBin*binFactor, -(thetaMax/100), thetaMax/100));
00098   thetaResolution[5]->setAxisTitle("rad",1);
00099   thetaResolution[5]->setAxisTitle("rad",2);
00100   thetaPull = dbe->book1D("Pull_TkSta_theta", "#theta_{TKfromGLB} - #theta_{STAfromGLB} / error", 100,-10,10);
00101   thetaTrack = dbe->book1D("TkMuon_theta", "#theta_{TK}", thetaBin, thetaMin, thetaMax);
00102   thetaTrack->setAxisTitle("rad");
00103   thetaStaTrack = dbe->book1D("StaMuon_theta", "#theta_{STA}", thetaBin, thetaMin, thetaMax);
00104   thetaStaTrack->setAxisTitle("rad");
00105 
00106   // monitoring of phi paramater
00107   phiBin = parameters.getParameter<int>("phiBin");
00108   phiMin = parameters.getParameter<double>("phiMin");
00109   phiMax = parameters.getParameter<double>("phiMax");
00110   phiGlbTrack.push_back(dbe->book1D(histname+"Glb_phi", "#phi_{GLB}", phiBin, phiMin, phiMax));
00111   phiGlbTrack[0]->setAxisTitle("rad");
00112   phiGlbTrack.push_back(dbe->book1D(histname+"Tk_phi", "#phi_{TKfromGLB}", phiBin, phiMin, phiMax));
00113   phiGlbTrack[1]->setAxisTitle("rad");
00114   phiGlbTrack.push_back(dbe->book1D(histname+"Sta_phi", "#phi_{STAfromGLB}", phiBin, phiMin, phiMax));
00115   phiGlbTrack[2]->setAxisTitle("rad");
00116   phiResolution.push_back(dbe->book1D("Res_TkGlb_phi", "#phi_{TKfromGLB} - #phi_{GLB}", phiBin*binFactor, phiMin/3000, phiMax/3000));
00117   phiResolution[0]->setAxisTitle("rad");
00118   phiResolution.push_back(dbe->book1D("Res_GlbSta_phi", "#phi_{GLB} - #phi_{STAfromGLB}", phiBin*binFactor, phiMin/100, phiMax/100));
00119   phiResolution[1]->setAxisTitle("rad");
00120   phiResolution.push_back(dbe->book1D("Res_TkSta_phi", "#phi_{TKfromGLB} - #phi_{STAfromGLB}", phiBin*binFactor, phiMin/100, phiMax/100));
00121   phiResolution[2]->setAxisTitle("rad");
00122   phiResolution.push_back(dbe->book2D("ResVsPhi_TkGlb_phi", "(#phi_{TKfromGLB} - #phi_{GLB}) vs #phi_GLB", phiBin, phiMin, phiMax, phiBin*binFactor, phiMin/3000, phiMax/3000));
00123   phiResolution[3]->setAxisTitle("rad",1);
00124   phiResolution[3]->setAxisTitle("rad",2);
00125   phiResolution.push_back(dbe->book2D("ResVsPhi_GlbSta_phi", "(#phi_{GLB} - #phi_{STAfromGLB}) vs #phi_{GLB}", phiBin, phiMin, phiMax, phiBin*binFactor, phiMin/100, phiMax/100));
00126   phiResolution[4]->setAxisTitle("rad",1);
00127   phiResolution[4]->setAxisTitle("rad",2);
00128   phiResolution.push_back(dbe->book2D("ResVsPhi_TkSta_phi", "(#phi_{TKfromGLB} - #phi_{STAfromGLB}) vs #phi_{TKfromGLB}", phiBin, phiMin, phiMax, phiBin*binFactor, phiMin/100, phiMax/100));
00129   phiResolution[5]->setAxisTitle("rad",1);
00130   phiResolution[5]->setAxisTitle("rad",2);
00131   phiPull = dbe->book1D("Pull_TkSta_phi", "#phi_{TKfromGLB} - #phi_{STAfromGLB} / error", 100,-10,10);
00132   phiTrack = dbe->book1D("TkMuon_phi", "#phi_{TK}", phiBin, phiMin, phiMax);
00133   phiTrack->setAxisTitle("rad"); 
00134   phiStaTrack = dbe->book1D("StaMuon_phi", "#phi_{STA}", phiBin, phiMin, phiMax);
00135   phiStaTrack->setAxisTitle("rad"); 
00136   phiEfficiency.push_back(dbe->book1D("StaPhi", "#phi_{STAfromGLB}", phiBin, phiMin, phiMax));
00137   phiEfficiency[0]->setAxisTitle("rad"); 
00138   phiEfficiency.push_back(dbe->book1D("StaPhi_ifCombinedAlso", "#phi_{STAfromGLB} if the isGlb=true", phiBin, phiMin, phiMax));
00139   phiEfficiency[1]->setAxisTitle("rad"); 
00140 
00141   // monitoring of the chi2 parameter
00142   chi2Bin = parameters.getParameter<int>("chi2Bin");
00143   chi2Min = parameters.getParameter<double>("chi2Min");
00144   chi2Max = parameters.getParameter<double>("chi2Max");
00145   chi2OvDFGlbTrack.push_back(dbe->book1D(histname+"Glb_chi2OverDf", "#chi_{2}OverDF_{GLB}", chi2Bin, chi2Min, chi2Max));
00146   chi2OvDFGlbTrack.push_back(dbe->book1D(histname+"Tk_chi2OverDf", "#chi_{2}OverDF_{TKfromGLB}", phiBin, chi2Min, chi2Max));
00147   chi2OvDFGlbTrack.push_back(dbe->book1D(histname+"Sta_chi2OverDf", "#chi_{2}OverDF_{STAfromGLB}", chi2Bin, chi2Min, chi2Max));
00148   chi2OvDFTrack = dbe->book1D("TkMuon_chi2OverDf", "#chi_{2}OverDF_{TK}", chi2Bin, chi2Min, chi2Max);
00149   chi2OvDFStaTrack = dbe->book1D("StaMuon_chi2OverDf", "#chi_{2}OverDF_{STA}", chi2Bin, chi2Min, chi2Max);
00150 //--------------------------
00151   probchi2GlbTrack.push_back(dbe->book1D(histname+"Glb_probchi", "Prob #chi_{GLB}", 120, chi2Min, 1.20));
00152   probchi2GlbTrack.push_back(dbe->book1D(histname+"Tk_probchi", "Prob #chi_{TKfromGLB}", 120, chi2Min, 1.20));
00153   probchi2GlbTrack.push_back(dbe->book1D(histname+"Sta_probchi", "Prob #chi_{STAfromGLB}", 120, chi2Min, 1.20));
00154   probchi2Track=dbe->book1D("TkMuon_probchi", "Prob #chi_{TK}", 120, chi2Min, 1.20);
00155   probchi2StaTrack=dbe->book1D("StaMuon_probchi", "Prob #chi_{STA}", 120, chi2Min, 1.20);
00156 
00157   // monitoring of the momentum
00158   pBin = parameters.getParameter<int>("pBin");
00159   pMin = parameters.getParameter<double>("pMin");
00160   pMax = parameters.getParameter<double>("pMax");
00161   pGlbTrack.push_back(dbe->book1D(histname+"Glb_p", "p_{GLB}", pBin, pMin, pMax));
00162   pGlbTrack[0]->setAxisTitle("GeV"); 
00163   pGlbTrack.push_back(dbe->book1D(histname+"Tk_p", "p_{TKfromGLB}", pBin, pMin, pMax));
00164   pGlbTrack[1]->setAxisTitle("GeV");
00165   pGlbTrack.push_back(dbe->book1D(histname+"Sta_p", "p_{STAfromGLB}", pBin, pMin, pMax));
00166   pGlbTrack[2]->setAxisTitle("GeV");
00167   pTrack = dbe->book1D("TkMuon_p", "p_{TK}", pBin, pMin, pMax);
00168   pTrack->setAxisTitle("GeV"); 
00169   pStaTrack = dbe->book1D("StaMuon_p", "p_{STA}", pBin, pMin, pMax);
00170   pStaTrack->setAxisTitle("GeV"); 
00171 
00172   // monitoring of the transverse momentum
00173   ptBin = parameters.getParameter<int>("ptBin");
00174   ptMin = parameters.getParameter<double>("ptMin");
00175   ptMax = parameters.getParameter<double>("ptMax");
00176   ptGlbTrack.push_back(dbe->book1D(histname+"Glb_pt", "pt_{GLB}", ptBin, ptMin, ptMax));
00177   ptGlbTrack[0]->setAxisTitle("GeV"); 
00178   ptGlbTrack.push_back(dbe->book1D(histname+"Tk_pt", "pt_{TKfromGLB}", ptBin, ptMin, ptMax));
00179   ptGlbTrack[1]->setAxisTitle("GeV"); 
00180   ptGlbTrack.push_back(dbe->book1D(histname+"Sta_pt", "pt_{STAfromGLB}", ptBin, ptMin, ptMax));
00181   ptGlbTrack[2]->setAxisTitle("GeV"); 
00182   ptTrack = dbe->book1D("TkMuon_pt", "pt_{TK}", ptBin, ptMin, ptMax);
00183   ptTrack->setAxisTitle("GeV"); 
00184   ptStaTrack = dbe->book1D("StaMuon_pt", "pt_{STA}", ptBin, ptMin, pMax);
00185   ptStaTrack->setAxisTitle("GeV"); 
00186 
00187   // monitoring of the muon charge
00188   qGlbTrack.push_back(dbe->book1D(histname+"Glb_q", "q_{GLB}", 5, -2.5, 2.5));
00189   qGlbTrack.push_back(dbe->book1D(histname+"Tk_q", "q_{TKfromGLB}", 5, -2.5, 2.5));
00190   qGlbTrack.push_back(dbe->book1D(histname+"Sta_q", "q_{STAformGLB}", 5, -2.5, 2.5));
00191   qGlbTrack.push_back(dbe->book1D(histname+"qComparison", "comparison between q_{GLB} and q_{TKfromGLB}, q_{STAfromGLB}", 8, 0.5, 8.5));
00192   qGlbTrack[3]->setBinLabel(1,"qGlb=qSta");
00193   qGlbTrack[3]->setBinLabel(2,"qGlb!=qSta");
00194   qGlbTrack[3]->setBinLabel(3,"qGlb=qTk");
00195   qGlbTrack[3]->setBinLabel(4,"qGlb!=qTk");
00196   qGlbTrack[3]->setBinLabel(5,"qSta=qTk");
00197   qGlbTrack[3]->setBinLabel(6,"qSta!=qTk");
00198   qGlbTrack[3]->setBinLabel(7,"qGlb!=qSta,qGlb!=Tk");
00199   qGlbTrack[3]->setBinLabel(8,"qGlb=qSta,qGlb=Tk");
00200   qTrack = dbe->book1D("TkMuon_q", "q_{TK}", 5, -2.5, 2.5);
00201   qStaTrack = dbe->book1D("StaMuon_q", "q_{STA}", 5, -2.5, 2.5);
00202 
00203   // monitoring of the momentum resolution
00204   pResBin = parameters.getParameter<int>("pResBin");
00205   pResMin = parameters.getParameter<double>("pResMin");
00206   pResMax = parameters.getParameter<double>("pResMax");
00207   qOverpResolution.push_back(dbe->book1D("Res_TkGlb_qOverp", "(q/p)_{TKfromGLB} - (q/p)_{GLB}", pResBin*binFactor*2, pResMin/10, pResMax/10));
00208   qOverpResolution[0]->setAxisTitle("GeV^{-1}"); 
00209   qOverpResolution.push_back(dbe->book1D("Res_GlbSta_qOverp", "(q/p)_{GLB} - (q/p)_{STAfromGLB}", pResBin*binFactor, pResMin, pResMax));
00210   qOverpResolution[1]->setAxisTitle("GeV^{-1}"); 
00211   qOverpResolution.push_back(dbe->book1D("Res_TkSta_qOverp", "(q/p)_{TKfromGLB} - (q/p)_{STAfromGLB}", pResBin*binFactor, pResMin, pResMax));
00212   qOverpResolution[2]->setAxisTitle("GeV^{-1}"); 
00213   qOverpPull = dbe->book1D("Pull_TkSta_qOverp", "(q/p)_{TKfromGLB} - (q/p)_{STAfromGLB} / error", 100,-10,10);
00214  
00215   oneOverpResolution.push_back(dbe->book1D("Res_TkGlb_oneOverp", "(1/p)_{TKfromGLB} - (1/p)_{GLB}", pResBin*binFactor*2, pResMin/10, pResMax/10));
00216   oneOverpResolution[0]->setAxisTitle("GeV^{-1}"); 
00217   oneOverpResolution.push_back(dbe->book1D("Res_GlbSta_oneOverp", "(1/p)_{GLB} - (1/p)_{STAfromGLB}", pResBin*binFactor, pResMin, pResMax));
00218   oneOverpResolution[1]->setAxisTitle("GeV^{-1}");
00219   oneOverpResolution.push_back(dbe->book1D("Res_TkSta_oneOverp", "(q/p)_{TKfromGLB} - (q/p)_{STAfromGLB}", pResBin*binFactor, pResMin, pResMax));
00220   oneOverpResolution[2]->setAxisTitle("GeV^{-1}");
00221   oneOverpPull = dbe->book1D("Pull_TkSta_oneOverp", "(1/p)_{TKfromGLB} - (1/p)_{STAfromGLB} / error", 100,-10,10);
00222 
00223 
00224   qOverptResolution.push_back(dbe->book1D("Res_TkGlb_qOverpt", "(q/p_{t})_{TKfromGLB} - (q/p_{t})_{GLB}", pResBin*binFactor*2, pResMin/10, pResMax/10));
00225   qOverptResolution[0]->setAxisTitle("GeV^{-1}");  
00226   qOverptResolution.push_back(dbe->book1D("Res_GlbSta_qOverpt", "(q/p_{t})_{GLB} - (q/p_{t})_{STAfromGLB}", pResBin*binFactor, pResMin, pResMax));
00227   qOverptResolution[1]->setAxisTitle("GeV^{-1}"); 
00228   qOverptResolution.push_back(dbe->book1D("Res_TkSta_qOverpt", "(q/p_{t})_{TKfromGLB} - (q/p_{t})_{STAfromGLB}", pResBin*binFactor, pResMin, pResMax));
00229   qOverptResolution[2]->setAxisTitle("GeV^{-1}"); 
00230   qOverptPull = dbe->book1D("Pull_TkSta_qOverpt", "(q/pt)_{TKfromGLB} - (q/pt)_{STAfromGLB} / error", 100,-10,10);
00231 
00232   oneOverptResolution.push_back(dbe->book1D("Res_TkGlb_oneOverpt", "(1/p_{t})_{TKfromGLB} - (1/p_{t})_{GLB}", pResBin*binFactor*2, pResMin/10, pResMax/10));
00233   oneOverptResolution[0]->setAxisTitle("GeV^{-1}");  
00234   oneOverptResolution.push_back(dbe->book1D("Res_GlbSta_oneOverpt", "(1/p_{t})_{GLB} - (1/p_{t})_{STAfromGLB}", pResBin*binFactor, pResMin, pResMax));
00235   oneOverptResolution[1]->setAxisTitle("GeV^{-1}");
00236   oneOverptResolution.push_back(dbe->book1D("Res_TkSta_oneOverpt", "(1/p_{t})_{TKfromGLB} - (1/p_{t})_{STAfromGLB}", pResBin*binFactor, pResMin, pResMax));
00237   oneOverptResolution[2]->setAxisTitle("GeV^{-1}");
00238   oneOverptResolution.push_back(dbe->book2D("ResVsEta_TkGlb_oneOverpt", "(#eta_{TKfromGLB} - #eta_{GLB}) vs (1/p_{t})_{GLB}", etaBin, etaMin, etaMax, pResBin*binFactor*2, pResMin/10, pResMax/10));
00239   oneOverptResolution[3]->setAxisTitle("GeV^{-1}",2);
00240   oneOverptResolution.push_back(dbe->book2D("ResVsEta_GlbSta_oneOverpt", "(#eta_{GLB} - #eta_{STAfromGLB} vs (1/p_{t})_{GLB}", etaBin, etaMin, etaMax, pResBin*binFactor, pResMin, pResMax));
00241   oneOverptResolution[4]->setAxisTitle("GeV^{-1}",2);
00242   oneOverptResolution.push_back(dbe->book2D("ResVsEta_TkSta_oneOverpt", "(#eta_{TKfromGLB} - #eta_{STAfromGLB}) vs (1/p_{t})_{TKfromGLB}", etaBin, etaMin, etaMax, pResBin*binFactor, pResMin, pResMax));
00243   oneOverptResolution[5]->setAxisTitle("GeV^{-1}",2);
00244   oneOverptResolution.push_back(dbe->book2D("ResVsPhi_TkGlb_oneOverpt", "(#phi_{TKfromGLB} - #phi_{GLB}) vs (1/p_{t})_{GLB}", phiBin, phiMin, phiMax, pResBin*binFactor*2, pResMin/10, pResMax/10));
00245   oneOverptResolution[6]->setAxisTitle("rad",1);
00246   oneOverptResolution[6]->setAxisTitle("GeV^{-1}",2);
00247   oneOverptResolution.push_back(dbe->book2D("ResVsPhi_GlbSta_oneOverpt", "(#phi_{GLB} - #phi_{STAfromGLB} vs (1/p_{t})_{GLB}", phiBin, phiMin, phiMax, pResBin*binFactor, pResMin, pResMax));
00248   oneOverptResolution[7]->setAxisTitle("rad",1);
00249   oneOverptResolution[7]->setAxisTitle("GeV^{-1}",2);
00250   oneOverptResolution.push_back(dbe->book2D("ResVsPhi_TkSta_oneOverpt", "(#phi_{TKfromGLB} - #phi_{STAfromGLB}) vs (1/p_{t})_{TKfromGLB}", phiBin, phiMin, phiMax, pResBin*binFactor, pResMin, pResMax));
00251   oneOverptResolution[8]->setAxisTitle("rad",1);
00252   oneOverptResolution[8]->setAxisTitle("GeV^{-1}",2);
00253   oneOverptResolution.push_back(dbe->book2D("ResVsPt_TkGlb_oneOverpt", "((1/p_{t})_{TKfromGLB} - (1/p_{t})_{GLB}) vs (1/p_{t})_{GLB}", ptBin/5, ptMin, ptMax/100, pResBin*binFactor*2, pResMin/10, pResMax/10));
00254   oneOverptResolution[9]->setAxisTitle("GeV^{-1}",1);
00255   oneOverptResolution[9]->setAxisTitle("GeV^{-1}",2);
00256   oneOverptResolution.push_back(dbe->book2D("ResVsPt_GlbSta_oneOverpt", "((1/p_{t})_{GLB} - (1/p_{t})_{STAfromGLB} vs (1/p_{t})_{GLB}", ptBin/5, ptMin, ptMax/100, pResBin*binFactor, pResMin, pResMax));
00257   oneOverptResolution[10]->setAxisTitle("GeV^{-1}",1);
00258   oneOverptResolution[10]->setAxisTitle("GeV^{-1}",2);
00259   oneOverptResolution.push_back(dbe->book2D("ResVsPt_TkSta_oneOverpt", "((1/p_{t})_{TKfromGLB} - (1/p_{t})_{STAfromGLB}) vs (1/p_{t})_{TKfromGLB}", ptBin/5, ptMin, ptMax/100, pResBin*binFactor, pResMin, pResMax));
00260   oneOverptResolution[11]->setAxisTitle("GeV^{-1}",1);
00261   oneOverptResolution[11]->setAxisTitle("GeV^{-1}",2);
00262   oneOverptPull = dbe->book1D("Pull_TkSta_oneOverpt", "(1/pt)_{TKfromGLB} - (1/pt)_{STAfromGLB} / error", 100,-10,10);
00263 
00264 
00265   // monitoring of the recHits provenance
00266   rhBin=parameters.getParameter<int>("rhBin");
00267   rhMin=parameters.getParameter<double>("rhMin");
00268   rhMax=parameters.getParameter<double>("rhMax");
00269   rhAnalysis.push_back(dbe->book1D("StaRh_Frac_inGlb", "recHits_{STAinGLB} / recHits_{GLB}", rhBin, rhMin, rhMax));
00270   rhAnalysis.push_back(dbe->book1D("TkRh_Frac_inGlb", "recHits_{TKinGLB} / recHits_{GLB}", rhBin, rhMin, rhMax));
00271   rhAnalysis.push_back(dbe->book1D("StaRh_inGlb_Div_RhAssoSta", "recHits_{STAinGLB} / recHits_{STAfromGLB}", rhBin, rhMin, rhMax));
00272   rhAnalysis.push_back(dbe->book1D("TkRh_inGlb_Div_RhAssoTk", "recHits_{TKinGLB} / recHits_{TKfromGLB}", rhBin, rhMin, rhMax));
00273   rhAnalysis.push_back(dbe->book1D("GlbRh_Div_RhAssoStaTk", "recHits_{GLB} / (recHits_{TKfromGLB}+recHits_{STAfromGLB})", rhBin, rhMin, rhMax));
00274   rhAnalysis.push_back(dbe->book1D("invalidRh_Frac_inTk", "Invalid recHits / rechits_{GLB}", rhBin, rhMin, rhMax));
00275 
00276   // monitoring of the muon system rotation w.r.t. tracker
00277   muVStkSytemRotation.push_back(dbe->book2D("muVStkSytemRotation_posMu", "pT_{TK} / pT_{GLB} vs pT_{GLB} for #mu^{+}", 50,0,200,100,0.8,1.2));
00278   muVStkSytemRotation.push_back(dbe->book2D("muVStkSytemRotation_negMu", "pT_{TK} / pT_{GLB} vs pT_{GLB} for #mu^{-}", 50,0,200,100,0.8,1.2));
00279 
00280 }
00281 
00282 
00283 void MuonRecoAnalyzer::GetRes( reco::TrackRef t1, reco::TrackRef t2, string par, float &res, float &pull){
00284   
00285   float p1=0, p2=0, p1e=1, p2e=1;
00286 
00287   if(par == "eta") {
00288     p1 = t1->eta(); p1e = t1->etaError();
00289     p2 = t2->eta(); p2e = t2->etaError();
00290   }
00291   else if(par == "theta") {
00292     p1 = t1->theta(); p1e = t1->thetaError();
00293     p2 = t2->theta(); p2e = t2->thetaError();
00294   }
00295   else if(par == "phi") {
00296     p1 = t1->phi(); p1e = t1->phiError();
00297     p2 = t2->phi(); p2e = t2->phiError();
00298   }
00299   else if(par == "qOverp") {
00300     p1 = t1->charge()/t1->p(); p1e = t1->qoverpError();
00301     p2 = t2->charge()/t2->p(); p2e = t2->qoverpError();
00302   }
00303   else if(par == "oneOverp") {
00304     p1 = 1./t1->p(); p1e = t1->qoverpError();
00305     p2 = 1./t2->p(); p2e = t2->qoverpError();
00306   }
00307   else if(par == "qOverpt") {
00308     p1 = t1->charge()/t1->pt(); p1e = t1->ptError()*p1*p1;
00309     p2 = t2->charge()/t2->pt(); p2e = t2->ptError()*p2*p2;
00310   }
00311   else if(par == "oneOverpt") {
00312     p1 = 1./t1->pt(); p1e = t1->ptError()*p1*p1;
00313     p2 = 1./t2->pt(); p2e = t2->ptError()*p2*p2;
00314   }
00315 
00316   res = p1 - p2;
00317   if(p1e!=0 || p2e!=0) pull = res / sqrt(p1e*p1e + p2e*p2e);
00318   else pull = -99;
00319   return;
00320 }
00321 
00322 
00323 
00324 
00325 
00326 void MuonRecoAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup, const reco::Muon& recoMu) {
00327 
00328   LogTrace(metname)<<"[MuonRecoAnalyzer] Analyze the mu";
00329 
00330   float res=0, pull=0;
00331 
00332   if(recoMu.isGlobalMuon()) {
00333 
00334     LogTrace(metname)<<"[MuonRecoAnalyzer] The mu is global - filling the histos";
00335     if(recoMu.isTrackerMuon() && recoMu.isStandAloneMuon())
00336       muReco->Fill(1);
00337     if(!(recoMu.isTrackerMuon()) && recoMu.isStandAloneMuon())
00338       muReco->Fill(2);
00339     if(!recoMu.isStandAloneMuon())
00340       LogTrace(metname)<<"[MuonRecoAnalyzer] ERROR: the mu is global but not standalone!";
00341 
00342     // get the track combinig the information from both the Tracker and the Spectrometer
00343     reco::TrackRef recoCombinedGlbTrack = recoMu.combinedMuon();
00344     // get the track using only the tracker data
00345     reco::TrackRef recoTkGlbTrack = recoMu.track();
00346     // get the track using only the mu spectrometer data
00347     reco::TrackRef recoStaGlbTrack = recoMu.standAloneMuon();
00348 
00349   
00350     etaGlbTrack[0]->Fill(recoCombinedGlbTrack->eta());
00351     etaGlbTrack[1]->Fill(recoTkGlbTrack->eta());
00352     etaGlbTrack[2]->Fill(recoStaGlbTrack->eta());
00353 
00354     GetRes(recoTkGlbTrack, recoCombinedGlbTrack, "eta", res, pull);
00355     etaResolution[0]->Fill(res);
00356     GetRes(recoCombinedGlbTrack, recoStaGlbTrack, "eta", res, pull);
00357     etaResolution[1]->Fill(res);
00358     GetRes(recoTkGlbTrack, recoStaGlbTrack, "eta", res, pull);
00359     etaResolution[2]->Fill(res);
00360     etaPull->Fill(pull);
00361 
00362     etaResolution[3]->Fill(recoCombinedGlbTrack->eta(), recoTkGlbTrack->eta()-recoCombinedGlbTrack->eta());
00363     etaResolution[4]->Fill(recoCombinedGlbTrack->eta(), -recoStaGlbTrack->eta()+recoCombinedGlbTrack->eta());
00364     etaResolution[5]->Fill(recoCombinedGlbTrack->eta(), recoTkGlbTrack->eta()-recoStaGlbTrack->eta());
00365  
00366 
00367 
00368     thetaGlbTrack[0]->Fill(recoCombinedGlbTrack->theta());
00369     thetaGlbTrack[1]->Fill(recoTkGlbTrack->theta());
00370     thetaGlbTrack[2]->Fill(recoStaGlbTrack->theta());
00371     GetRes(recoTkGlbTrack, recoCombinedGlbTrack, "theta", res, pull);
00372     thetaResolution[0]->Fill(res);
00373 
00374     GetRes(recoCombinedGlbTrack, recoStaGlbTrack, "theta", res, pull);
00375     thetaResolution[1]->Fill(res);
00376 
00377     GetRes(recoTkGlbTrack, recoStaGlbTrack, "theta", res, pull);
00378     thetaResolution[2]->Fill(res);
00379     thetaPull->Fill(pull);
00380 
00381     thetaResolution[3]->Fill(recoCombinedGlbTrack->theta(), recoTkGlbTrack->theta()-recoCombinedGlbTrack->theta());
00382     thetaResolution[4]->Fill(recoCombinedGlbTrack->theta(), -recoStaGlbTrack->theta()+recoCombinedGlbTrack->theta());
00383     thetaResolution[5]->Fill(recoCombinedGlbTrack->theta(), recoTkGlbTrack->theta()-recoStaGlbTrack->theta());
00384 
00385 
00386      
00387     phiGlbTrack[0]->Fill(recoCombinedGlbTrack->phi());
00388     phiGlbTrack[1]->Fill(recoTkGlbTrack->phi());
00389     phiGlbTrack[2]->Fill(recoStaGlbTrack->phi());
00390     GetRes(recoTkGlbTrack, recoCombinedGlbTrack, "phi", res, pull);
00391     phiResolution[0]->Fill(res);
00392     GetRes(recoCombinedGlbTrack, recoStaGlbTrack, "phi", res, pull);
00393     phiResolution[1]->Fill(res);
00394     GetRes(recoTkGlbTrack, recoStaGlbTrack, "phi", res, pull);
00395     phiResolution[2]->Fill(res);
00396     phiPull->Fill(pull);
00397     phiResolution[3]->Fill(recoCombinedGlbTrack->phi(), recoTkGlbTrack->phi()-recoCombinedGlbTrack->phi());
00398     phiResolution[4]->Fill(recoCombinedGlbTrack->phi(), -recoStaGlbTrack->phi()+recoCombinedGlbTrack->phi());
00399     phiResolution[5]->Fill(recoCombinedGlbTrack->phi(), recoTkGlbTrack->phi()-recoStaGlbTrack->phi());
00400     
00401     chi2OvDFGlbTrack[0]->Fill(recoCombinedGlbTrack->normalizedChi2());
00402     chi2OvDFGlbTrack[1]->Fill(recoTkGlbTrack->normalizedChi2());
00403     chi2OvDFGlbTrack[2]->Fill(recoStaGlbTrack->normalizedChi2());
00404 //-------------------------
00405 //    double probchi = TMath::Prob(recoCombinedGlbTrack->normalizedChi2(),recoCombinedGlbTrack->ndof());
00406 //    cout << "rellenando histos."<<endl;
00407     probchi2GlbTrack[0]->Fill(TMath::Prob(recoCombinedGlbTrack->chi2(),recoCombinedGlbTrack->ndof()));
00408     probchi2GlbTrack[1]->Fill(TMath::Prob(recoTkGlbTrack->chi2(),recoTkGlbTrack->ndof()));
00409     probchi2GlbTrack[2]->Fill(TMath::Prob(recoStaGlbTrack->chi2(),recoStaGlbTrack->ndof()));
00410     //    cout << "rellenados histos."<<endl;
00411 //-------------------------
00412 
00413     pGlbTrack[0]->Fill(recoCombinedGlbTrack->p());
00414     pGlbTrack[1]->Fill(recoTkGlbTrack->p());
00415     pGlbTrack[2]->Fill(recoStaGlbTrack->p());
00416 
00417     ptGlbTrack[0]->Fill(recoCombinedGlbTrack->pt());
00418     ptGlbTrack[1]->Fill(recoTkGlbTrack->pt());
00419     ptGlbTrack[2]->Fill(recoStaGlbTrack->pt());
00420 
00421     qGlbTrack[0]->Fill(recoCombinedGlbTrack->charge());
00422     qGlbTrack[1]->Fill(recoTkGlbTrack->charge());
00423     qGlbTrack[2]->Fill(recoStaGlbTrack->charge());
00424     if(recoCombinedGlbTrack->charge()==recoStaGlbTrack->charge()) qGlbTrack[3]->Fill(1);
00425     else qGlbTrack[3]->Fill(2);
00426     if(recoCombinedGlbTrack->charge()==recoTkGlbTrack->charge()) qGlbTrack[3]->Fill(3);
00427     else qGlbTrack[3]->Fill(4);
00428     if(recoStaGlbTrack->charge()==recoTkGlbTrack->charge()) qGlbTrack[3]->Fill(5);
00429     else qGlbTrack[3]->Fill(6);
00430     if(recoCombinedGlbTrack->charge()!=recoStaGlbTrack->charge() && recoCombinedGlbTrack->charge()!=recoTkGlbTrack->charge()) qGlbTrack[3]->Fill(7);
00431     if(recoCombinedGlbTrack->charge()==recoStaGlbTrack->charge() && recoCombinedGlbTrack->charge()==recoTkGlbTrack->charge()) qGlbTrack[3]->Fill(8);
00432     
00433     GetRes(recoTkGlbTrack, recoCombinedGlbTrack, "qOverp", res, pull);
00434     qOverpResolution[0]->Fill(res);
00435     GetRes(recoCombinedGlbTrack, recoStaGlbTrack, "qOverp", res, pull);
00436     qOverpResolution[1]->Fill(res);
00437     GetRes(recoTkGlbTrack, recoStaGlbTrack, "qOverp", res, pull);
00438     qOverpResolution[2]->Fill(res);
00439     qOverpPull->Fill(pull);
00440 
00441 
00442     GetRes(recoTkGlbTrack, recoCombinedGlbTrack, "oneOverp", res, pull);
00443     oneOverpResolution[0]->Fill(res);
00444     GetRes(recoCombinedGlbTrack, recoStaGlbTrack, "oneOverp", res, pull);
00445     oneOverpResolution[1]->Fill(res);
00446     GetRes(recoTkGlbTrack, recoStaGlbTrack, "oneOverp", res, pull);
00447     oneOverpResolution[2]->Fill(res);
00448     oneOverpPull->Fill(pull);
00449     
00450 
00451     GetRes(recoTkGlbTrack, recoCombinedGlbTrack, "qOverpt", res, pull);
00452     qOverptResolution[0]->Fill(res);
00453     GetRes(recoCombinedGlbTrack, recoStaGlbTrack, "qOverpt", res, pull);
00454     qOverptResolution[1]->Fill(res);
00455     GetRes(recoTkGlbTrack, recoStaGlbTrack, "qOverpt", res, pull);
00456     qOverptResolution[2]->Fill(res);
00457     qOverptPull->Fill(pull);
00458 
00459     GetRes(recoTkGlbTrack, recoCombinedGlbTrack, "oneOverpt", res, pull);
00460     oneOverptResolution[0]->Fill(res);
00461     GetRes(recoCombinedGlbTrack, recoStaGlbTrack, "oneOverpt", res, pull);
00462     oneOverptResolution[1]->Fill(res);
00463     GetRes(recoTkGlbTrack, recoStaGlbTrack, "oneOverpt", res, pull);
00464     oneOverptResolution[2]->Fill(res);
00465     oneOverptPull->Fill(pull);
00466 
00467     oneOverptResolution[3]->Fill(recoCombinedGlbTrack->eta(),(1/recoTkGlbTrack->pt())-(1/recoCombinedGlbTrack->pt()));
00468     oneOverptResolution[4]->Fill(recoCombinedGlbTrack->eta(),-(1/recoStaGlbTrack->pt())+(1/recoCombinedGlbTrack->pt()));
00469     oneOverptResolution[5]->Fill(recoCombinedGlbTrack->eta(),(1/recoTkGlbTrack->pt())-(1/recoStaGlbTrack->pt()));
00470     oneOverptResolution[6]->Fill(recoCombinedGlbTrack->phi(),(1/recoTkGlbTrack->pt())-(1/recoCombinedGlbTrack->pt()));
00471     oneOverptResolution[7]->Fill(recoCombinedGlbTrack->phi(),-(1/recoStaGlbTrack->pt())+(1/recoCombinedGlbTrack->pt()));
00472     oneOverptResolution[8]->Fill(recoCombinedGlbTrack->phi(),(1/recoTkGlbTrack->pt())-(1/recoStaGlbTrack->pt()));
00473     oneOverptResolution[9]->Fill(recoCombinedGlbTrack->pt(),(1/recoTkGlbTrack->pt())-(1/recoCombinedGlbTrack->pt()));
00474     oneOverptResolution[10]->Fill(recoCombinedGlbTrack->pt(),-(1/recoStaGlbTrack->pt())+(1/recoCombinedGlbTrack->pt()));
00475     oneOverptResolution[11]->Fill(recoCombinedGlbTrack->pt(),(1/recoTkGlbTrack->pt())-(1/recoStaGlbTrack->pt()));
00476     
00477     // valid hits Glb track
00478     double rhGlb = recoCombinedGlbTrack->found();
00479     // valid hits Glb track from Tracker
00480     double rhGlb_StaProvenance=0;
00481     // valid hits Glb track from Sta system
00482     double rhGlb_TkProvenance=0;
00483     for (trackingRecHit_iterator recHit = recoCombinedGlbTrack->recHitsBegin();
00484          recHit!=recoCombinedGlbTrack->recHitsEnd(); ++recHit){
00485       if((*recHit)->isValid()){
00486         DetId id = (*recHit)->geographicalId();
00487         if (id.det() == DetId::Muon)
00488           rhGlb_StaProvenance++;
00489         if (id.det() == DetId::Tracker)
00490           rhGlb_TkProvenance++;
00491       }
00492     }
00493     // valid hits Sta track associated to Glb track
00494     double rhStaGlb = recoStaGlbTrack->recHitsSize();
00495     // valid hits Traker track associated to Glb track
00496     double rhTkGlb = recoTkGlbTrack->found();
00497     // invalid hits Traker track associated to Glb track
00498     double rhTkGlb_notValid = recoTkGlbTrack->lost();
00499    
00500     // fill the histos
00501     rhAnalysis[0]->Fill(rhGlb_StaProvenance/rhGlb);
00502     rhAnalysis[1]->Fill(rhGlb_TkProvenance/rhGlb);
00503     rhAnalysis[2]->Fill(rhGlb_StaProvenance/rhStaGlb);
00504     rhAnalysis[3]->Fill(rhGlb_TkProvenance/rhTkGlb);
00505     rhAnalysis[4]->Fill(rhGlb/(rhStaGlb+rhTkGlb));
00506     rhAnalysis[5]->Fill(rhTkGlb_notValid/rhGlb);
00507 
00508     // aligment plots (mu system w.r.t. tracker rotation)
00509     if(recoCombinedGlbTrack->charge()>0)
00510       muVStkSytemRotation[0]->Fill(recoCombinedGlbTrack->pt(),recoTkGlbTrack->pt()/recoCombinedGlbTrack->pt());
00511     else
00512       muVStkSytemRotation[1]->Fill(recoCombinedGlbTrack->pt(),recoTkGlbTrack->pt()/recoCombinedGlbTrack->pt());
00513     
00514   }
00515 
00516 
00517   if(recoMu.isTrackerMuon() && !(recoMu.isGlobalMuon())) {
00518     LogTrace(metname)<<"[MuonRecoAnalyzer] The mu is tracker only - filling the histos";
00519      if(recoMu.isStandAloneMuon())
00520       muReco->Fill(3);
00521     if(!(recoMu.isStandAloneMuon()))
00522       muReco->Fill(4);
00523 
00524     // get the track using only the tracker data
00525     reco::TrackRef recoTrack = recoMu.track();
00526 
00527     etaTrack->Fill(recoTrack->eta());
00528     thetaTrack->Fill(recoTrack->theta());
00529     phiTrack->Fill(recoTrack->phi());
00530     chi2OvDFTrack->Fill(recoTrack->normalizedChi2());
00531     probchi2Track->Fill(TMath::Prob(recoTrack->chi2(),recoTrack->ndof()));
00532     pTrack->Fill(recoTrack->p());
00533     ptTrack->Fill(recoTrack->pt());
00534     qTrack->Fill(recoTrack->charge());
00535     
00536   }
00537 
00538   if(recoMu.isStandAloneMuon() && !(recoMu.isGlobalMuon())) {
00539     LogTrace(metname)<<"[MuonRecoAnalyzer] The mu is STA only - filling the histos";
00540     if(!(recoMu.isTrackerMuon()))
00541       muReco->Fill(5);
00542      
00543     // get the track using only the mu spectrometer data
00544     reco::TrackRef recoStaTrack = recoMu.standAloneMuon();
00545 
00546     etaStaTrack->Fill(recoStaTrack->eta());
00547     thetaStaTrack->Fill(recoStaTrack->theta());
00548     phiStaTrack->Fill(recoStaTrack->phi());
00549     chi2OvDFStaTrack->Fill(recoStaTrack->normalizedChi2());
00550     probchi2StaTrack->Fill(TMath::Prob(recoStaTrack->chi2(),recoStaTrack->ndof()));
00551     pStaTrack->Fill(recoStaTrack->p());
00552     ptStaTrack->Fill(recoStaTrack->pt());
00553     qStaTrack->Fill(recoStaTrack->charge());
00554 
00555   }
00556     
00557   if(recoMu.isCaloMuon() && !(recoMu.isGlobalMuon()) && !(recoMu.isTrackerMuon()) && !(recoMu.isStandAloneMuon()))
00558     muReco->Fill(6);
00559   
00560   //efficiency plots
00561   
00562   // get the track using only the mu spectrometer data
00563   reco::TrackRef recoStaGlbTrack = recoMu.standAloneMuon();
00564   
00565   if(recoMu.isStandAloneMuon()){
00566     etaEfficiency[0]->Fill(recoStaGlbTrack->eta());
00567     phiEfficiency[0]->Fill(recoStaGlbTrack->phi());
00568   }
00569   if(recoMu.isStandAloneMuon() && recoMu.isGlobalMuon()){
00570     etaEfficiency[1]->Fill(recoStaGlbTrack->eta());
00571     phiEfficiency[1]->Fill(recoStaGlbTrack->phi());
00572   }
00573 
00574 
00575 
00576 
00577 }