00001
00002
00003
00004
00005
00006
00007
00008
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
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
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
00107 tunePBin= parameters.getParameter<int>("tunePBin");
00108 tunePMax= parameters.getParameter<double>("tunePMax");
00109 tunePMin= parameters.getParameter<double>("tunePMin");
00110
00111 tunePResolution = dbe->book1D("Res_TuneP_pt", "Pt_{MuonBestTrack}-Pt_{tunePMuonBestTrack}/Pt_{MuonBestTrack}", tunePBin, tunePMin, tunePMax);
00112
00113
00114
00115
00116 phiBin = parameters.getParameter<int>("phiBin");
00117 phiMin = parameters.getParameter<double>("phiMin");
00118 phiMax = parameters.getParameter<double>("phiMax");
00119 phiGlbTrack.push_back(dbe->book1D(histname+"Glb_phi", "#phi_{GLB}", phiBin, phiMin, phiMax));
00120 phiGlbTrack[0]->setAxisTitle("rad");
00121 phiGlbTrack.push_back(dbe->book1D(histname+"Tk_phi", "#phi_{TKfromGLB}", phiBin, phiMin, phiMax));
00122 phiGlbTrack[1]->setAxisTitle("rad");
00123 phiGlbTrack.push_back(dbe->book1D(histname+"Sta_phi", "#phi_{STAfromGLB}", phiBin, phiMin, phiMax));
00124 phiGlbTrack[2]->setAxisTitle("rad");
00125 phiResolution.push_back(dbe->book1D("Res_TkGlb_phi", "#phi_{TKfromGLB} - #phi_{GLB}", phiBin*binFactor, phiMin/3000, phiMax/3000));
00126 phiResolution[0]->setAxisTitle("rad");
00127 phiResolution.push_back(dbe->book1D("Res_GlbSta_phi", "#phi_{GLB} - #phi_{STAfromGLB}", phiBin*binFactor, phiMin/100, phiMax/100));
00128 phiResolution[1]->setAxisTitle("rad");
00129 phiResolution.push_back(dbe->book1D("Res_TkSta_phi", "#phi_{TKfromGLB} - #phi_{STAfromGLB}", phiBin*binFactor, phiMin/100, phiMax/100));
00130 phiResolution[2]->setAxisTitle("rad");
00131 phiResolution.push_back(dbe->book2D("ResVsPhi_TkGlb_phi", "(#phi_{TKfromGLB} - #phi_{GLB}) vs #phi_GLB", phiBin, phiMin, phiMax, phiBin*binFactor, phiMin/3000, phiMax/3000));
00132 phiResolution[3]->setAxisTitle("rad",1);
00133 phiResolution[3]->setAxisTitle("rad",2);
00134 phiResolution.push_back(dbe->book2D("ResVsPhi_GlbSta_phi", "(#phi_{GLB} - #phi_{STAfromGLB}) vs #phi_{GLB}", phiBin, phiMin, phiMax, phiBin*binFactor, phiMin/100, phiMax/100));
00135 phiResolution[4]->setAxisTitle("rad",1);
00136 phiResolution[4]->setAxisTitle("rad",2);
00137 phiResolution.push_back(dbe->book2D("ResVsPhi_TkSta_phi", "(#phi_{TKfromGLB} - #phi_{STAfromGLB}) vs #phi_{TKfromGLB}", phiBin, phiMin, phiMax, phiBin*binFactor, phiMin/100, phiMax/100));
00138 phiResolution[5]->setAxisTitle("rad",1);
00139 phiResolution[5]->setAxisTitle("rad",2);
00140 phiPull = dbe->book1D("Pull_TkSta_phi", "#phi_{TKfromGLB} - #phi_{STAfromGLB} / error", 100,-10,10);
00141 phiTrack = dbe->book1D("TkMuon_phi", "#phi_{TK}", phiBin, phiMin, phiMax);
00142 phiTrack->setAxisTitle("rad");
00143 phiStaTrack = dbe->book1D("StaMuon_phi", "#phi_{STA}", phiBin, phiMin, phiMax);
00144 phiStaTrack->setAxisTitle("rad");
00145 phiEfficiency.push_back(dbe->book1D("StaPhi", "#phi_{STAfromGLB}", phiBin, phiMin, phiMax));
00146 phiEfficiency[0]->setAxisTitle("rad");
00147 phiEfficiency.push_back(dbe->book1D("StaPhi_ifCombinedAlso", "#phi_{STAfromGLB} if the isGlb=true", phiBin, phiMin, phiMax));
00148 phiEfficiency[1]->setAxisTitle("rad");
00149
00150
00151 chi2Bin = parameters.getParameter<int>("chi2Bin");
00152 chi2Min = parameters.getParameter<double>("chi2Min");
00153 chi2Max = parameters.getParameter<double>("chi2Max");
00154 chi2OvDFGlbTrack.push_back(dbe->book1D(histname+"Glb_chi2OverDf", "#chi_{2}OverDF_{GLB}", chi2Bin, chi2Min, chi2Max));
00155 chi2OvDFGlbTrack.push_back(dbe->book1D(histname+"Tk_chi2OverDf", "#chi_{2}OverDF_{TKfromGLB}", phiBin, chi2Min, chi2Max));
00156 chi2OvDFGlbTrack.push_back(dbe->book1D(histname+"Sta_chi2OverDf", "#chi_{2}OverDF_{STAfromGLB}", chi2Bin, chi2Min, chi2Max));
00157 chi2OvDFTrack = dbe->book1D("TkMuon_chi2OverDf", "#chi_{2}OverDF_{TK}", chi2Bin, chi2Min, chi2Max);
00158 chi2OvDFStaTrack = dbe->book1D("StaMuon_chi2OverDf", "#chi_{2}OverDF_{STA}", chi2Bin, chi2Min, chi2Max);
00159
00160 probchi2GlbTrack.push_back(dbe->book1D(histname+"Glb_probchi", "Prob #chi_{GLB}", 120, chi2Min, 1.20));
00161 probchi2GlbTrack.push_back(dbe->book1D(histname+"Tk_probchi", "Prob #chi_{TKfromGLB}", 120, chi2Min, 1.20));
00162 probchi2GlbTrack.push_back(dbe->book1D(histname+"Sta_probchi", "Prob #chi_{STAfromGLB}", 120, chi2Min, 1.20));
00163 probchi2Track=dbe->book1D("TkMuon_probchi", "Prob #chi_{TK}", 120, chi2Min, 1.20);
00164 probchi2StaTrack=dbe->book1D("StaMuon_probchi", "Prob #chi_{STA}", 120, chi2Min, 1.20);
00165
00166
00167 pBin = parameters.getParameter<int>("pBin");
00168 pMin = parameters.getParameter<double>("pMin");
00169 pMax = parameters.getParameter<double>("pMax");
00170 pGlbTrack.push_back(dbe->book1D(histname+"Glb_p", "p_{GLB}", pBin, pMin, pMax));
00171 pGlbTrack[0]->setAxisTitle("GeV");
00172 pGlbTrack.push_back(dbe->book1D(histname+"Tk_p", "p_{TKfromGLB}", pBin, pMin, pMax));
00173 pGlbTrack[1]->setAxisTitle("GeV");
00174 pGlbTrack.push_back(dbe->book1D(histname+"Sta_p", "p_{STAfromGLB}", pBin, pMin, pMax));
00175 pGlbTrack[2]->setAxisTitle("GeV");
00176 pTrack = dbe->book1D("TkMuon_p", "p_{TK}", pBin, pMin, pMax);
00177 pTrack->setAxisTitle("GeV");
00178 pStaTrack = dbe->book1D("StaMuon_p", "p_{STA}", pBin, pMin, pMax);
00179 pStaTrack->setAxisTitle("GeV");
00180
00181
00182 ptBin = parameters.getParameter<int>("ptBin");
00183 ptMin = parameters.getParameter<double>("ptMin");
00184 ptMax = parameters.getParameter<double>("ptMax");
00185 ptGlbTrack.push_back(dbe->book1D(histname+"Glb_pt", "pt_{GLB}", ptBin, ptMin, ptMax));
00186 ptGlbTrack[0]->setAxisTitle("GeV");
00187 ptGlbTrack.push_back(dbe->book1D(histname+"Tk_pt", "pt_{TKfromGLB}", ptBin, ptMin, ptMax));
00188 ptGlbTrack[1]->setAxisTitle("GeV");
00189 ptGlbTrack.push_back(dbe->book1D(histname+"Sta_pt", "pt_{STAfromGLB}", ptBin, ptMin, ptMax));
00190 ptGlbTrack[2]->setAxisTitle("GeV");
00191 ptTrack = dbe->book1D("TkMuon_pt", "pt_{TK}", ptBin, ptMin, ptMax);
00192 ptTrack->setAxisTitle("GeV");
00193 ptStaTrack = dbe->book1D("StaMuon_pt", "pt_{STA}", ptBin, ptMin, pMax);
00194 ptStaTrack->setAxisTitle("GeV");
00195
00196
00197 qGlbTrack.push_back(dbe->book1D(histname+"Glb_q", "q_{GLB}", 5, -2.5, 2.5));
00198 qGlbTrack.push_back(dbe->book1D(histname+"Tk_q", "q_{TKfromGLB}", 5, -2.5, 2.5));
00199 qGlbTrack.push_back(dbe->book1D(histname+"Sta_q", "q_{STAformGLB}", 5, -2.5, 2.5));
00200 qGlbTrack.push_back(dbe->book1D(histname+"qComparison", "comparison between q_{GLB} and q_{TKfromGLB}, q_{STAfromGLB}", 8, 0.5, 8.5));
00201 qGlbTrack[3]->setBinLabel(1,"qGlb=qSta");
00202 qGlbTrack[3]->setBinLabel(2,"qGlb!=qSta");
00203 qGlbTrack[3]->setBinLabel(3,"qGlb=qTk");
00204 qGlbTrack[3]->setBinLabel(4,"qGlb!=qTk");
00205 qGlbTrack[3]->setBinLabel(5,"qSta=qTk");
00206 qGlbTrack[3]->setBinLabel(6,"qSta!=qTk");
00207 qGlbTrack[3]->setBinLabel(7,"qGlb!=qSta,qGlb!=Tk");
00208 qGlbTrack[3]->setBinLabel(8,"qGlb=qSta,qGlb=Tk");
00209 qTrack = dbe->book1D("TkMuon_q", "q_{TK}", 5, -2.5, 2.5);
00210 qStaTrack = dbe->book1D("StaMuon_q", "q_{STA}", 5, -2.5, 2.5);
00211
00212
00213 pResBin = parameters.getParameter<int>("pResBin");
00214 pResMin = parameters.getParameter<double>("pResMin");
00215 pResMax = parameters.getParameter<double>("pResMax");
00216 qOverpResolution.push_back(dbe->book1D("Res_TkGlb_qOverp", "(q/p)_{TKfromGLB} - (q/p)_{GLB}", pResBin*binFactor*2, pResMin/10, pResMax/10));
00217 qOverpResolution[0]->setAxisTitle("GeV^{-1}");
00218 qOverpResolution.push_back(dbe->book1D("Res_GlbSta_qOverp", "(q/p)_{GLB} - (q/p)_{STAfromGLB}", pResBin*binFactor, pResMin, pResMax));
00219 qOverpResolution[1]->setAxisTitle("GeV^{-1}");
00220 qOverpResolution.push_back(dbe->book1D("Res_TkSta_qOverp", "(q/p)_{TKfromGLB} - (q/p)_{STAfromGLB}", pResBin*binFactor, pResMin, pResMax));
00221 qOverpResolution[2]->setAxisTitle("GeV^{-1}");
00222 qOverpPull = dbe->book1D("Pull_TkSta_qOverp", "(q/p)_{TKfromGLB} - (q/p)_{STAfromGLB} / error", 100,-10,10);
00223
00224 oneOverpResolution.push_back(dbe->book1D("Res_TkGlb_oneOverp", "(1/p)_{TKfromGLB} - (1/p)_{GLB}", pResBin*binFactor*2, pResMin/10, pResMax/10));
00225 oneOverpResolution[0]->setAxisTitle("GeV^{-1}");
00226 oneOverpResolution.push_back(dbe->book1D("Res_GlbSta_oneOverp", "(1/p)_{GLB} - (1/p)_{STAfromGLB}", pResBin*binFactor, pResMin, pResMax));
00227 oneOverpResolution[1]->setAxisTitle("GeV^{-1}");
00228 oneOverpResolution.push_back(dbe->book1D("Res_TkSta_oneOverp", "(q/p)_{TKfromGLB} - (q/p)_{STAfromGLB}", pResBin*binFactor, pResMin, pResMax));
00229 oneOverpResolution[2]->setAxisTitle("GeV^{-1}");
00230 oneOverpPull = dbe->book1D("Pull_TkSta_oneOverp", "(1/p)_{TKfromGLB} - (1/p)_{STAfromGLB} / error", 100,-10,10);
00231
00232
00233 qOverptResolution.push_back(dbe->book1D("Res_TkGlb_qOverpt", "(q/p_{t})_{TKfromGLB} - (q/p_{t})_{GLB}", pResBin*binFactor*2, pResMin/10, pResMax/10));
00234 qOverptResolution[0]->setAxisTitle("GeV^{-1}");
00235 qOverptResolution.push_back(dbe->book1D("Res_GlbSta_qOverpt", "(q/p_{t})_{GLB} - (q/p_{t})_{STAfromGLB}", pResBin*binFactor, pResMin, pResMax));
00236 qOverptResolution[1]->setAxisTitle("GeV^{-1}");
00237 qOverptResolution.push_back(dbe->book1D("Res_TkSta_qOverpt", "(q/p_{t})_{TKfromGLB} - (q/p_{t})_{STAfromGLB}", pResBin*binFactor, pResMin, pResMax));
00238 qOverptResolution[2]->setAxisTitle("GeV^{-1}");
00239 qOverptPull = dbe->book1D("Pull_TkSta_qOverpt", "(q/pt)_{TKfromGLB} - (q/pt)_{STAfromGLB} / error", 100,-10,10);
00240
00241 oneOverptResolution.push_back(dbe->book1D("Res_TkGlb_oneOverpt", "(1/p_{t})_{TKfromGLB} - (1/p_{t})_{GLB}", pResBin*binFactor*2, pResMin/10, pResMax/10));
00242 oneOverptResolution[0]->setAxisTitle("GeV^{-1}");
00243 oneOverptResolution.push_back(dbe->book1D("Res_GlbSta_oneOverpt", "(1/p_{t})_{GLB} - (1/p_{t})_{STAfromGLB}", pResBin*binFactor, pResMin, pResMax));
00244 oneOverptResolution[1]->setAxisTitle("GeV^{-1}");
00245 oneOverptResolution.push_back(dbe->book1D("Res_TkSta_oneOverpt", "(1/p_{t})_{TKfromGLB} - (1/p_{t})_{STAfromGLB}", pResBin*binFactor, pResMin, pResMax));
00246 oneOverptResolution[2]->setAxisTitle("GeV^{-1}");
00247 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));
00248 oneOverptResolution[3]->setAxisTitle("GeV^{-1}",2);
00249 oneOverptResolution.push_back(dbe->book2D("ResVsEta_GlbSta_oneOverpt", "(#eta_{GLB} - #eta_{STAfromGLB} vs (1/p_{t})_{GLB}", etaBin, etaMin, etaMax, pResBin*binFactor, pResMin, pResMax));
00250 oneOverptResolution[4]->setAxisTitle("GeV^{-1}",2);
00251 oneOverptResolution.push_back(dbe->book2D("ResVsEta_TkSta_oneOverpt", "(#eta_{TKfromGLB} - #eta_{STAfromGLB}) vs (1/p_{t})_{TKfromGLB}", etaBin, etaMin, etaMax, pResBin*binFactor, pResMin, pResMax));
00252 oneOverptResolution[5]->setAxisTitle("GeV^{-1}",2);
00253 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));
00254 oneOverptResolution[6]->setAxisTitle("rad",1);
00255 oneOverptResolution[6]->setAxisTitle("GeV^{-1}",2);
00256 oneOverptResolution.push_back(dbe->book2D("ResVsPhi_GlbSta_oneOverpt", "(#phi_{GLB} - #phi_{STAfromGLB} vs (1/p_{t})_{GLB}", phiBin, phiMin, phiMax, pResBin*binFactor, pResMin, pResMax));
00257 oneOverptResolution[7]->setAxisTitle("rad",1);
00258 oneOverptResolution[7]->setAxisTitle("GeV^{-1}",2);
00259 oneOverptResolution.push_back(dbe->book2D("ResVsPhi_TkSta_oneOverpt", "(#phi_{TKfromGLB} - #phi_{STAfromGLB}) vs (1/p_{t})_{TKfromGLB}", phiBin, phiMin, phiMax, pResBin*binFactor, pResMin, pResMax));
00260 oneOverptResolution[8]->setAxisTitle("rad",1);
00261 oneOverptResolution[8]->setAxisTitle("GeV^{-1}",2);
00262 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));
00263 oneOverptResolution[9]->setAxisTitle("GeV^{-1}",1);
00264 oneOverptResolution[9]->setAxisTitle("GeV^{-1}",2);
00265 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));
00266 oneOverptResolution[10]->setAxisTitle("GeV^{-1}",1);
00267 oneOverptResolution[10]->setAxisTitle("GeV^{-1}",2);
00268 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));
00269 oneOverptResolution[11]->setAxisTitle("GeV^{-1}",1);
00270 oneOverptResolution[11]->setAxisTitle("GeV^{-1}",2);
00271 oneOverptPull = dbe->book1D("Pull_TkSta_oneOverpt", "(1/pt)_{TKfromGLB} - (1/pt)_{STAfromGLB} / error", 100,-10,10);
00272
00273
00274
00275 rhBin=parameters.getParameter<int>("rhBin");
00276 rhMin=parameters.getParameter<double>("rhMin");
00277 rhMax=parameters.getParameter<double>("rhMax");
00278 rhAnalysis.push_back(dbe->book1D("StaRh_Frac_inGlb", "recHits_{STAinGLB} / recHits_{GLB}", rhBin, rhMin, rhMax));
00279 rhAnalysis.push_back(dbe->book1D("TkRh_Frac_inGlb", "recHits_{TKinGLB} / recHits_{GLB}", rhBin, rhMin, rhMax));
00280 rhAnalysis.push_back(dbe->book1D("StaRh_inGlb_Div_RhAssoSta", "recHits_{STAinGLB} / recHits_{STAfromGLB}", rhBin, rhMin, rhMax));
00281 rhAnalysis.push_back(dbe->book1D("TkRh_inGlb_Div_RhAssoTk", "recHits_{TKinGLB} / recHits_{TKfromGLB}", rhBin, rhMin, rhMax));
00282 rhAnalysis.push_back(dbe->book1D("GlbRh_Div_RhAssoStaTk", "recHits_{GLB} / (recHits_{TKfromGLB}+recHits_{STAfromGLB})", rhBin, rhMin, rhMax));
00283 rhAnalysis.push_back(dbe->book1D("invalidRh_Frac_inTk", "Invalid recHits / rechits_{GLB}", rhBin, rhMin, rhMax));
00284
00285
00286 muVStkSytemRotation.push_back(dbe->book2D("muVStkSytemRotation_posMu", "pT_{TK} / pT_{GLB} vs pT_{GLB} for #mu^{+}", 50,0,200,100,0.8,1.2));
00287 muVStkSytemRotation.push_back(dbe->book2D("muVStkSytemRotation_negMu", "pT_{TK} / pT_{GLB} vs pT_{GLB} for #mu^{-}", 50,0,200,100,0.8,1.2));
00288
00289 }
00290
00291
00292 void MuonRecoAnalyzer::GetRes( reco::TrackRef t1, reco::TrackRef t2, string par, float &res, float &pull){
00293
00294 float p1=0, p2=0, p1e=1, p2e=1;
00295
00296 if(par == "eta") {
00297 p1 = t1->eta(); p1e = t1->etaError();
00298 p2 = t2->eta(); p2e = t2->etaError();
00299 }
00300 else if(par == "theta") {
00301 p1 = t1->theta(); p1e = t1->thetaError();
00302 p2 = t2->theta(); p2e = t2->thetaError();
00303 }
00304 else if(par == "phi") {
00305 p1 = t1->phi(); p1e = t1->phiError();
00306 p2 = t2->phi(); p2e = t2->phiError();
00307 }
00308 else if(par == "qOverp") {
00309 p1 = t1->charge()/t1->p(); p1e = t1->qoverpError();
00310 p2 = t2->charge()/t2->p(); p2e = t2->qoverpError();
00311 }
00312 else if(par == "oneOverp") {
00313 p1 = 1./t1->p(); p1e = t1->qoverpError();
00314 p2 = 1./t2->p(); p2e = t2->qoverpError();
00315 }
00316 else if(par == "qOverpt") {
00317 p1 = t1->charge()/t1->pt(); p1e = t1->ptError()*p1*p1;
00318 p2 = t2->charge()/t2->pt(); p2e = t2->ptError()*p2*p2;
00319 }
00320 else if(par == "oneOverpt") {
00321 p1 = 1./t1->pt(); p1e = t1->ptError()*p1*p1;
00322 p2 = 1./t2->pt(); p2e = t2->ptError()*p2*p2;
00323 }
00324
00325 res = p1 - p2;
00326 if(p1e!=0 || p2e!=0) pull = res / sqrt(p1e*p1e + p2e*p2e);
00327 else pull = -99;
00328 return;
00329 }
00330
00331
00332
00333
00334
00335
00336 void MuonRecoAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup, const reco::Muon& recoMu) {
00337
00338 LogTrace(metname)<<"[MuonRecoAnalyzer] Analyze the mu";
00339
00340 float res=0, pull=0;
00341
00342 if(recoMu.isGlobalMuon()) {
00343
00344 LogTrace(metname)<<"[MuonRecoAnalyzer] The mu is global - filling the histos";
00345 if(recoMu.isTrackerMuon() && recoMu.isStandAloneMuon())
00346 muReco->Fill(1);
00347 if(!(recoMu.isTrackerMuon()) && recoMu.isStandAloneMuon())
00348 muReco->Fill(2);
00349 if(!recoMu.isStandAloneMuon())
00350 LogTrace(metname)<<"[MuonRecoAnalyzer] ERROR: the mu is global but not standalone!";
00351
00352
00353 reco::TrackRef recoCombinedGlbTrack = recoMu.combinedMuon();
00354
00355 reco::TrackRef recoTkGlbTrack = recoMu.track();
00356
00357 reco::TrackRef recoStaGlbTrack = recoMu.standAloneMuon();
00358
00359
00360 etaGlbTrack[0]->Fill(recoCombinedGlbTrack->eta());
00361 etaGlbTrack[1]->Fill(recoTkGlbTrack->eta());
00362 etaGlbTrack[2]->Fill(recoStaGlbTrack->eta());
00363
00364 GetRes(recoTkGlbTrack, recoCombinedGlbTrack, "eta", res, pull);
00365 etaResolution[0]->Fill(res);
00366 GetRes(recoCombinedGlbTrack, recoStaGlbTrack, "eta", res, pull);
00367 etaResolution[1]->Fill(res);
00368 GetRes(recoTkGlbTrack, recoStaGlbTrack, "eta", res, pull);
00369 etaResolution[2]->Fill(res);
00370 etaPull->Fill(pull);
00371
00372 etaResolution[3]->Fill(recoCombinedGlbTrack->eta(), recoTkGlbTrack->eta()-recoCombinedGlbTrack->eta());
00373 etaResolution[4]->Fill(recoCombinedGlbTrack->eta(), -recoStaGlbTrack->eta()+recoCombinedGlbTrack->eta());
00374 etaResolution[5]->Fill(recoCombinedGlbTrack->eta(), recoTkGlbTrack->eta()-recoStaGlbTrack->eta());
00375
00376
00377 thetaGlbTrack[0]->Fill(recoCombinedGlbTrack->theta());
00378 thetaGlbTrack[1]->Fill(recoTkGlbTrack->theta());
00379 thetaGlbTrack[2]->Fill(recoStaGlbTrack->theta());
00380 GetRes(recoTkGlbTrack, recoCombinedGlbTrack, "theta", res, pull);
00381 thetaResolution[0]->Fill(res);
00382
00383 GetRes(recoCombinedGlbTrack, recoStaGlbTrack, "theta", res, pull);
00384 thetaResolution[1]->Fill(res);
00385
00386 GetRes(recoTkGlbTrack, recoStaGlbTrack, "theta", res, pull);
00387 thetaResolution[2]->Fill(res);
00388 thetaPull->Fill(pull);
00389
00390 thetaResolution[3]->Fill(recoCombinedGlbTrack->theta(), recoTkGlbTrack->theta()-recoCombinedGlbTrack->theta());
00391 thetaResolution[4]->Fill(recoCombinedGlbTrack->theta(), -recoStaGlbTrack->theta()+recoCombinedGlbTrack->theta());
00392 thetaResolution[5]->Fill(recoCombinedGlbTrack->theta(), recoTkGlbTrack->theta()-recoStaGlbTrack->theta());
00393
00394
00395
00396 phiGlbTrack[0]->Fill(recoCombinedGlbTrack->phi());
00397 phiGlbTrack[1]->Fill(recoTkGlbTrack->phi());
00398 phiGlbTrack[2]->Fill(recoStaGlbTrack->phi());
00399 GetRes(recoTkGlbTrack, recoCombinedGlbTrack, "phi", res, pull);
00400 phiResolution[0]->Fill(res);
00401 GetRes(recoCombinedGlbTrack, recoStaGlbTrack, "phi", res, pull);
00402 phiResolution[1]->Fill(res);
00403 GetRes(recoTkGlbTrack, recoStaGlbTrack, "phi", res, pull);
00404 phiResolution[2]->Fill(res);
00405 phiPull->Fill(pull);
00406 phiResolution[3]->Fill(recoCombinedGlbTrack->phi(), recoTkGlbTrack->phi()-recoCombinedGlbTrack->phi());
00407 phiResolution[4]->Fill(recoCombinedGlbTrack->phi(), -recoStaGlbTrack->phi()+recoCombinedGlbTrack->phi());
00408 phiResolution[5]->Fill(recoCombinedGlbTrack->phi(), recoTkGlbTrack->phi()-recoStaGlbTrack->phi());
00409
00410 chi2OvDFGlbTrack[0]->Fill(recoCombinedGlbTrack->normalizedChi2());
00411 chi2OvDFGlbTrack[1]->Fill(recoTkGlbTrack->normalizedChi2());
00412 chi2OvDFGlbTrack[2]->Fill(recoStaGlbTrack->normalizedChi2());
00413
00414
00415
00416 probchi2GlbTrack[0]->Fill(TMath::Prob(recoCombinedGlbTrack->chi2(),recoCombinedGlbTrack->ndof()));
00417 probchi2GlbTrack[1]->Fill(TMath::Prob(recoTkGlbTrack->chi2(),recoTkGlbTrack->ndof()));
00418 probchi2GlbTrack[2]->Fill(TMath::Prob(recoStaGlbTrack->chi2(),recoStaGlbTrack->ndof()));
00419
00420
00421
00422 pGlbTrack[0]->Fill(recoCombinedGlbTrack->p());
00423 pGlbTrack[1]->Fill(recoTkGlbTrack->p());
00424 pGlbTrack[2]->Fill(recoStaGlbTrack->p());
00425
00426 ptGlbTrack[0]->Fill(recoCombinedGlbTrack->pt());
00427 ptGlbTrack[1]->Fill(recoTkGlbTrack->pt());
00428 ptGlbTrack[2]->Fill(recoStaGlbTrack->pt());
00429
00430 qGlbTrack[0]->Fill(recoCombinedGlbTrack->charge());
00431 qGlbTrack[1]->Fill(recoTkGlbTrack->charge());
00432 qGlbTrack[2]->Fill(recoStaGlbTrack->charge());
00433 if(recoCombinedGlbTrack->charge()==recoStaGlbTrack->charge()) qGlbTrack[3]->Fill(1);
00434 else qGlbTrack[3]->Fill(2);
00435 if(recoCombinedGlbTrack->charge()==recoTkGlbTrack->charge()) qGlbTrack[3]->Fill(3);
00436 else qGlbTrack[3]->Fill(4);
00437 if(recoStaGlbTrack->charge()==recoTkGlbTrack->charge()) qGlbTrack[3]->Fill(5);
00438 else qGlbTrack[3]->Fill(6);
00439 if(recoCombinedGlbTrack->charge()!=recoStaGlbTrack->charge() && recoCombinedGlbTrack->charge()!=recoTkGlbTrack->charge()) qGlbTrack[3]->Fill(7);
00440 if(recoCombinedGlbTrack->charge()==recoStaGlbTrack->charge() && recoCombinedGlbTrack->charge()==recoTkGlbTrack->charge()) qGlbTrack[3]->Fill(8);
00441
00442 GetRes(recoTkGlbTrack, recoCombinedGlbTrack, "qOverp", res, pull);
00443 qOverpResolution[0]->Fill(res);
00444 GetRes(recoCombinedGlbTrack, recoStaGlbTrack, "qOverp", res, pull);
00445 qOverpResolution[1]->Fill(res);
00446 GetRes(recoTkGlbTrack, recoStaGlbTrack, "qOverp", res, pull);
00447 qOverpResolution[2]->Fill(res);
00448 qOverpPull->Fill(pull);
00449
00450
00451 GetRes(recoTkGlbTrack, recoCombinedGlbTrack, "oneOverp", res, pull);
00452 oneOverpResolution[0]->Fill(res);
00453 GetRes(recoCombinedGlbTrack, recoStaGlbTrack, "oneOverp", res, pull);
00454 oneOverpResolution[1]->Fill(res);
00455 GetRes(recoTkGlbTrack, recoStaGlbTrack, "oneOverp", res, pull);
00456 oneOverpResolution[2]->Fill(res);
00457 oneOverpPull->Fill(pull);
00458
00459
00460 GetRes(recoTkGlbTrack, recoCombinedGlbTrack, "qOverpt", res, pull);
00461 qOverptResolution[0]->Fill(res);
00462 GetRes(recoCombinedGlbTrack, recoStaGlbTrack, "qOverpt", res, pull);
00463 qOverptResolution[1]->Fill(res);
00464 GetRes(recoTkGlbTrack, recoStaGlbTrack, "qOverpt", res, pull);
00465 qOverptResolution[2]->Fill(res);
00466 qOverptPull->Fill(pull);
00467
00468 GetRes(recoTkGlbTrack, recoCombinedGlbTrack, "oneOverpt", res, pull);
00469 oneOverptResolution[0]->Fill(res);
00470 GetRes(recoCombinedGlbTrack, recoStaGlbTrack, "oneOverpt", res, pull);
00471 oneOverptResolution[1]->Fill(res);
00472 GetRes(recoTkGlbTrack, recoStaGlbTrack, "oneOverpt", res, pull);
00473 oneOverptResolution[2]->Fill(res);
00474 oneOverptPull->Fill(pull);
00475
00476
00477
00478
00479 reco::TrackRef recoBestTrack = recoMu.muonBestTrack();
00480
00481 reco::TrackRef recoTunePBestTrack = recoMu.tunePMuonBestTrack();
00482
00483 double bestTrackPt = recoBestTrack->pt();
00484
00485 double tunePBestTrackPt = recoTunePBestTrack->pt();
00486
00487 double tunePBestTrackRes = (bestTrackPt - tunePBestTrackPt) / bestTrackPt;
00488
00489 tunePResolution->Fill(tunePBestTrackRes);
00490
00491
00492 oneOverptResolution[3]->Fill(recoCombinedGlbTrack->eta(),(1/recoTkGlbTrack->pt())-(1/recoCombinedGlbTrack->pt()));
00493 oneOverptResolution[4]->Fill(recoCombinedGlbTrack->eta(),-(1/recoStaGlbTrack->pt())+(1/recoCombinedGlbTrack->pt()));
00494 oneOverptResolution[5]->Fill(recoCombinedGlbTrack->eta(),(1/recoTkGlbTrack->pt())-(1/recoStaGlbTrack->pt()));
00495 oneOverptResolution[6]->Fill(recoCombinedGlbTrack->phi(),(1/recoTkGlbTrack->pt())-(1/recoCombinedGlbTrack->pt()));
00496 oneOverptResolution[7]->Fill(recoCombinedGlbTrack->phi(),-(1/recoStaGlbTrack->pt())+(1/recoCombinedGlbTrack->pt()));
00497 oneOverptResolution[8]->Fill(recoCombinedGlbTrack->phi(),(1/recoTkGlbTrack->pt())-(1/recoStaGlbTrack->pt()));
00498 oneOverptResolution[9]->Fill(recoCombinedGlbTrack->pt(),(1/recoTkGlbTrack->pt())-(1/recoCombinedGlbTrack->pt()));
00499 oneOverptResolution[10]->Fill(recoCombinedGlbTrack->pt(),-(1/recoStaGlbTrack->pt())+(1/recoCombinedGlbTrack->pt()));
00500 oneOverptResolution[11]->Fill(recoCombinedGlbTrack->pt(),(1/recoTkGlbTrack->pt())-(1/recoStaGlbTrack->pt()));
00501
00502
00503
00504
00505
00506
00507 double rhGlb = recoCombinedGlbTrack->found();
00508
00509 double rhGlb_StaProvenance=0;
00510
00511 double rhGlb_TkProvenance=0;
00512 for (trackingRecHit_iterator recHit = recoCombinedGlbTrack->recHitsBegin();
00513 recHit!=recoCombinedGlbTrack->recHitsEnd(); ++recHit){
00514 if((*recHit)->isValid()){
00515 DetId id = (*recHit)->geographicalId();
00516 if (id.det() == DetId::Muon)
00517 rhGlb_StaProvenance++;
00518 if (id.det() == DetId::Tracker)
00519 rhGlb_TkProvenance++;
00520 }
00521 }
00522
00523 double rhStaGlb = recoStaGlbTrack->recHitsSize();
00524
00525 double rhTkGlb = recoTkGlbTrack->found();
00526
00527 double rhTkGlb_notValid = recoTkGlbTrack->lost();
00528
00529
00530 rhAnalysis[0]->Fill(rhGlb_StaProvenance/rhGlb);
00531 rhAnalysis[1]->Fill(rhGlb_TkProvenance/rhGlb);
00532 rhAnalysis[2]->Fill(rhGlb_StaProvenance/rhStaGlb);
00533 rhAnalysis[3]->Fill(rhGlb_TkProvenance/rhTkGlb);
00534 rhAnalysis[4]->Fill(rhGlb/(rhStaGlb+rhTkGlb));
00535 rhAnalysis[5]->Fill(rhTkGlb_notValid/rhGlb);
00536
00537
00538 if(recoCombinedGlbTrack->charge()>0)
00539 muVStkSytemRotation[0]->Fill(recoCombinedGlbTrack->pt(),recoTkGlbTrack->pt()/recoCombinedGlbTrack->pt());
00540 else
00541 muVStkSytemRotation[1]->Fill(recoCombinedGlbTrack->pt(),recoTkGlbTrack->pt()/recoCombinedGlbTrack->pt());
00542
00543 }
00544
00545
00546 if(recoMu.isTrackerMuon() && !(recoMu.isGlobalMuon())) {
00547 LogTrace(metname)<<"[MuonRecoAnalyzer] The mu is tracker only - filling the histos";
00548 if(recoMu.isStandAloneMuon())
00549 muReco->Fill(3);
00550 if(!(recoMu.isStandAloneMuon()))
00551 muReco->Fill(4);
00552
00553
00554 reco::TrackRef recoTrack = recoMu.track();
00555
00556 etaTrack->Fill(recoTrack->eta());
00557 thetaTrack->Fill(recoTrack->theta());
00558 phiTrack->Fill(recoTrack->phi());
00559 chi2OvDFTrack->Fill(recoTrack->normalizedChi2());
00560 probchi2Track->Fill(TMath::Prob(recoTrack->chi2(),recoTrack->ndof()));
00561 pTrack->Fill(recoTrack->p());
00562 ptTrack->Fill(recoTrack->pt());
00563 qTrack->Fill(recoTrack->charge());
00564
00565 }
00566
00567 if(recoMu.isStandAloneMuon() && !(recoMu.isGlobalMuon())) {
00568 LogTrace(metname)<<"[MuonRecoAnalyzer] The mu is STA only - filling the histos";
00569 if(!(recoMu.isTrackerMuon()))
00570 muReco->Fill(5);
00571
00572
00573 reco::TrackRef recoStaTrack = recoMu.standAloneMuon();
00574
00575 etaStaTrack->Fill(recoStaTrack->eta());
00576 thetaStaTrack->Fill(recoStaTrack->theta());
00577 phiStaTrack->Fill(recoStaTrack->phi());
00578 chi2OvDFStaTrack->Fill(recoStaTrack->normalizedChi2());
00579 probchi2StaTrack->Fill(TMath::Prob(recoStaTrack->chi2(),recoStaTrack->ndof()));
00580 pStaTrack->Fill(recoStaTrack->p());
00581 ptStaTrack->Fill(recoStaTrack->pt());
00582 qStaTrack->Fill(recoStaTrack->charge());
00583
00584 }
00585
00586 if(recoMu.isCaloMuon() && !(recoMu.isGlobalMuon()) && !(recoMu.isTrackerMuon()) && !(recoMu.isStandAloneMuon()))
00587 muReco->Fill(6);
00588
00589
00590
00591
00592 reco::TrackRef recoStaGlbTrack = recoMu.standAloneMuon();
00593
00594 if(recoMu.isStandAloneMuon()){
00595 etaEfficiency[0]->Fill(recoStaGlbTrack->eta());
00596 phiEfficiency[0]->Fill(recoStaGlbTrack->phi());
00597 }
00598 if(recoMu.isStandAloneMuon() && recoMu.isGlobalMuon()){
00599 etaEfficiency[1]->Fill(recoStaGlbTrack->eta());
00600 phiEfficiency[1]->Fill(recoStaGlbTrack->phi());
00601 }
00602
00603
00604
00605
00606 }