CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MuonRecoAnalyzer.cc
Go to the documentation of this file.
1 
2 /*
3  * See header file for a description of this class.
4  *
5  * $Date: 2010/01/15 17:58:23 $
6  * $Revision: 1.22 $
7  * \author G. Mila - INFN Torino
8  * updated: G. Hesketh, CERN
9  */
10 
12 
17 
20 
22 
23 #include <string>
24 #include "TMath.h"
25 using namespace std;
26 using namespace edm;
27 
28 
29 
31  parameters = pSet;
32 }
33 
34 
36 
37 
39 
40  metname = "muRecoAnalyzer";
41 
42  LogTrace(metname)<<"[MuonRecoAnalyzer] Parameters initialization";
43  dbe->setCurrentFolder("Muons/MuonRecoAnalyzer");
44 
45  muReco = dbe->book1D("muReco", "muon reconstructed tracks", 6, 1, 7);
46  muReco->setBinLabel(1,"glb+tk+sta");
47  muReco->setBinLabel(2,"glb+sta");
48  muReco->setBinLabel(3,"tk+sta");
49  muReco->setBinLabel(4,"tk");
50  muReco->setBinLabel(5,"sta");
51  muReco->setBinLabel(6,"calo");
52 
53  int binFactor = 4;
54 
55  // monitoring of eta parameter
56  etaBin = parameters.getParameter<int>("etaBin");
57  etaMin = parameters.getParameter<double>("etaMin");
58  etaMax = parameters.getParameter<double>("etaMax");
59  std::string histname = "GlbMuon_";
60  etaGlbTrack.push_back(dbe->book1D(histname+"Glb_eta", "#eta_{GLB}", etaBin, etaMin, etaMax));
61  etaGlbTrack.push_back(dbe->book1D(histname+"Tk_eta", "#eta_{TKfromGLB}", etaBin, etaMin, etaMax));
62  etaGlbTrack.push_back(dbe->book1D(histname+"Sta_eta", "#eta_{STAfromGLB}", etaBin, etaMin, etaMax));
63  etaResolution.push_back(dbe->book1D("Res_TkGlb_eta", "#eta_{TKfromGLB} - #eta_{GLB}", etaBin*binFactor, etaMin/3000, etaMax/3000));
64  etaResolution.push_back(dbe->book1D("Res_GlbSta_eta", "#eta_{GLB} - #eta_{STAfromGLB}", etaBin*binFactor, etaMin/100, etaMax/100));
65  etaResolution.push_back(dbe->book1D("Res_TkSta_eta", "#eta_{TKfromGLB} - #eta_{STAfromGLB}", etaBin*binFactor, etaMin/100, etaMax/100));
66  etaResolution.push_back(dbe->book2D("ResVsEta_TkGlb_eta", "(#eta_{TKfromGLB} - #eta_{GLB}) vs #eta_{GLB}", etaBin, etaMin, etaMax, etaBin*binFactor, etaMin/3000, etaMax/3000));
67  etaResolution.push_back(dbe->book2D("ResVsEta_GlbSta_eta", "(#eta_{GLB} - #eta_{STAfromGLB}) vs #eta_{GLB}", etaBin, etaMin, etaMax, etaBin*binFactor, etaMin/100, etaMax/100));
68  etaResolution.push_back(dbe->book2D("ResVsEta_TkSta_eta", "(#eta_{TKfromGLB} - #eta_{STAfromGLB}) vs #eta_{TKfromGLB}", etaBin, etaMin, etaMax, etaBin*binFactor, etaMin/100, etaMax/100));
69  etaPull = dbe->book1D("Pull_TkSta_eta", "#eta_{TKfromGLB} - #eta_{GLB} / error", 100,-10,10);
70  etaTrack = dbe->book1D("TkMuon_eta", "#eta_{TK}", etaBin, etaMin, etaMax);
71  etaStaTrack = dbe->book1D("StaMuon_eta", "#eta_{STA}", etaBin, etaMin, etaMax);
72  etaEfficiency.push_back(dbe->book1D("StaEta", "#eta_{STAfromGLB}", etaBin, etaMin, etaMax));
73  etaEfficiency.push_back(dbe->book1D("StaEta_ifCombinedAlso", "#eta_{STAfromGLB} if isGlb=true", etaBin, etaMin, etaMax));
74 
75  // monitoring of theta parameter
76  thetaBin = parameters.getParameter<int>("thetaBin");
77  thetaMin = parameters.getParameter<double>("thetaMin");
78  thetaMax = parameters.getParameter<double>("thetaMax");
79  thetaGlbTrack.push_back(dbe->book1D(histname+"Glb_theta", "#theta_{GLB}", thetaBin, thetaMin, thetaMax));
80  thetaGlbTrack[0]->setAxisTitle("rad");
81  thetaGlbTrack.push_back(dbe->book1D(histname+"Tk_theta", "#theta_{TKfromGLB}", thetaBin, thetaMin, thetaMax));
82  thetaGlbTrack[1]->setAxisTitle("rad");
83  thetaGlbTrack.push_back(dbe->book1D(histname+"Sta_theta", "#theta_{STAfromGLB}", thetaBin, thetaMin, thetaMax));
84  thetaGlbTrack[2]->setAxisTitle("rad");
85  thetaResolution.push_back(dbe->book1D("Res_TkGlb_theta", "#theta_{TKfromGLB} - #theta_{GLB}", thetaBin*binFactor, -(thetaMax/3000), thetaMax/3000));
86  thetaResolution[0]->setAxisTitle("rad");
87  thetaResolution.push_back(dbe->book1D("Res_GlbSta_theta", "#theta_{GLB} - #theta_{STAfromGLB}", thetaBin*binFactor,-(thetaMax/100), thetaMax/100));
88  thetaResolution[1]->setAxisTitle("rad");
89  thetaResolution.push_back(dbe->book1D("Res_TkSta_theta", "#theta_{TKfromGLB} - #theta_{STAfromGLB}", thetaBin*binFactor, -(thetaMax/100), thetaMax/100));
90  thetaResolution[2]->setAxisTitle("rad");
91  thetaResolution.push_back(dbe->book2D("ResVsTheta_TkGlb_theta", "(#theta_{TKfromGLB} - #theta_{GLB}) vs #theta_{GLB}", thetaBin, thetaMin, thetaMax, thetaBin*binFactor, -(thetaMax/3000), thetaMax/3000));
92  thetaResolution[3]->setAxisTitle("rad",1);
93  thetaResolution[3]->setAxisTitle("rad",2);
94  thetaResolution.push_back(dbe->book2D("ResVsTheta_GlbSta_theta", "(#theta_{GLB} - #theta_{STAfromGLB}) vs #theta_{GLB}", thetaBin, thetaMin, thetaMax, thetaBin*binFactor, -(thetaMax/100), thetaMax/100));
95  thetaResolution[4]->setAxisTitle("rad",1);
96  thetaResolution[4]->setAxisTitle("rad",2);
97  thetaResolution.push_back(dbe->book2D("ResVsTheta_TkSta_theta", "(#theta_{TKfromGLB} - #theta_{STAfromGLB}) vs #theta_{TKfromGLB}", thetaBin, thetaMin, thetaMax, thetaBin*binFactor, -(thetaMax/100), thetaMax/100));
98  thetaResolution[5]->setAxisTitle("rad",1);
99  thetaResolution[5]->setAxisTitle("rad",2);
100  thetaPull = dbe->book1D("Pull_TkSta_theta", "#theta_{TKfromGLB} - #theta_{STAfromGLB} / error", 100,-10,10);
101  thetaTrack = dbe->book1D("TkMuon_theta", "#theta_{TK}", thetaBin, thetaMin, thetaMax);
102  thetaTrack->setAxisTitle("rad");
103  thetaStaTrack = dbe->book1D("StaMuon_theta", "#theta_{STA}", thetaBin, thetaMin, thetaMax);
104  thetaStaTrack->setAxisTitle("rad");
105 
106  // monitoring of phi paramater
107  phiBin = parameters.getParameter<int>("phiBin");
108  phiMin = parameters.getParameter<double>("phiMin");
109  phiMax = parameters.getParameter<double>("phiMax");
110  phiGlbTrack.push_back(dbe->book1D(histname+"Glb_phi", "#phi_{GLB}", phiBin, phiMin, phiMax));
111  phiGlbTrack[0]->setAxisTitle("rad");
112  phiGlbTrack.push_back(dbe->book1D(histname+"Tk_phi", "#phi_{TKfromGLB}", phiBin, phiMin, phiMax));
113  phiGlbTrack[1]->setAxisTitle("rad");
114  phiGlbTrack.push_back(dbe->book1D(histname+"Sta_phi", "#phi_{STAfromGLB}", phiBin, phiMin, phiMax));
115  phiGlbTrack[2]->setAxisTitle("rad");
116  phiResolution.push_back(dbe->book1D("Res_TkGlb_phi", "#phi_{TKfromGLB} - #phi_{GLB}", phiBin*binFactor, phiMin/3000, phiMax/3000));
117  phiResolution[0]->setAxisTitle("rad");
118  phiResolution.push_back(dbe->book1D("Res_GlbSta_phi", "#phi_{GLB} - #phi_{STAfromGLB}", phiBin*binFactor, phiMin/100, phiMax/100));
119  phiResolution[1]->setAxisTitle("rad");
120  phiResolution.push_back(dbe->book1D("Res_TkSta_phi", "#phi_{TKfromGLB} - #phi_{STAfromGLB}", phiBin*binFactor, phiMin/100, phiMax/100));
121  phiResolution[2]->setAxisTitle("rad");
122  phiResolution.push_back(dbe->book2D("ResVsPhi_TkGlb_phi", "(#phi_{TKfromGLB} - #phi_{GLB}) vs #phi_GLB", phiBin, phiMin, phiMax, phiBin*binFactor, phiMin/3000, phiMax/3000));
123  phiResolution[3]->setAxisTitle("rad",1);
124  phiResolution[3]->setAxisTitle("rad",2);
125  phiResolution.push_back(dbe->book2D("ResVsPhi_GlbSta_phi", "(#phi_{GLB} - #phi_{STAfromGLB}) vs #phi_{GLB}", phiBin, phiMin, phiMax, phiBin*binFactor, phiMin/100, phiMax/100));
126  phiResolution[4]->setAxisTitle("rad",1);
127  phiResolution[4]->setAxisTitle("rad",2);
128  phiResolution.push_back(dbe->book2D("ResVsPhi_TkSta_phi", "(#phi_{TKfromGLB} - #phi_{STAfromGLB}) vs #phi_{TKfromGLB}", phiBin, phiMin, phiMax, phiBin*binFactor, phiMin/100, phiMax/100));
129  phiResolution[5]->setAxisTitle("rad",1);
130  phiResolution[5]->setAxisTitle("rad",2);
131  phiPull = dbe->book1D("Pull_TkSta_phi", "#phi_{TKfromGLB} - #phi_{STAfromGLB} / error", 100,-10,10);
132  phiTrack = dbe->book1D("TkMuon_phi", "#phi_{TK}", phiBin, phiMin, phiMax);
133  phiTrack->setAxisTitle("rad");
134  phiStaTrack = dbe->book1D("StaMuon_phi", "#phi_{STA}", phiBin, phiMin, phiMax);
135  phiStaTrack->setAxisTitle("rad");
136  phiEfficiency.push_back(dbe->book1D("StaPhi", "#phi_{STAfromGLB}", phiBin, phiMin, phiMax));
137  phiEfficiency[0]->setAxisTitle("rad");
138  phiEfficiency.push_back(dbe->book1D("StaPhi_ifCombinedAlso", "#phi_{STAfromGLB} if the isGlb=true", phiBin, phiMin, phiMax));
139  phiEfficiency[1]->setAxisTitle("rad");
140 
141  // monitoring of the chi2 parameter
142  chi2Bin = parameters.getParameter<int>("chi2Bin");
143  chi2Min = parameters.getParameter<double>("chi2Min");
144  chi2Max = parameters.getParameter<double>("chi2Max");
145  chi2OvDFGlbTrack.push_back(dbe->book1D(histname+"Glb_chi2OverDf", "#chi_{2}OverDF_{GLB}", chi2Bin, chi2Min, chi2Max));
146  chi2OvDFGlbTrack.push_back(dbe->book1D(histname+"Tk_chi2OverDf", "#chi_{2}OverDF_{TKfromGLB}", phiBin, chi2Min, chi2Max));
147  chi2OvDFGlbTrack.push_back(dbe->book1D(histname+"Sta_chi2OverDf", "#chi_{2}OverDF_{STAfromGLB}", chi2Bin, chi2Min, chi2Max));
148  chi2OvDFTrack = dbe->book1D("TkMuon_chi2OverDf", "#chi_{2}OverDF_{TK}", chi2Bin, chi2Min, chi2Max);
149  chi2OvDFStaTrack = dbe->book1D("StaMuon_chi2OverDf", "#chi_{2}OverDF_{STA}", chi2Bin, chi2Min, chi2Max);
150 //--------------------------
151  probchi2GlbTrack.push_back(dbe->book1D(histname+"Glb_probchi", "Prob #chi_{GLB}", 120, chi2Min, 1.20));
152  probchi2GlbTrack.push_back(dbe->book1D(histname+"Tk_probchi", "Prob #chi_{TKfromGLB}", 120, chi2Min, 1.20));
153  probchi2GlbTrack.push_back(dbe->book1D(histname+"Sta_probchi", "Prob #chi_{STAfromGLB}", 120, chi2Min, 1.20));
154  probchi2Track=dbe->book1D("TkMuon_probchi", "Prob #chi_{TK}", 120, chi2Min, 1.20);
155  probchi2StaTrack=dbe->book1D("StaMuon_probchi", "Prob #chi_{STA}", 120, chi2Min, 1.20);
156 
157  // monitoring of the momentum
158  pBin = parameters.getParameter<int>("pBin");
159  pMin = parameters.getParameter<double>("pMin");
160  pMax = parameters.getParameter<double>("pMax");
161  pGlbTrack.push_back(dbe->book1D(histname+"Glb_p", "p_{GLB}", pBin, pMin, pMax));
162  pGlbTrack[0]->setAxisTitle("GeV");
163  pGlbTrack.push_back(dbe->book1D(histname+"Tk_p", "p_{TKfromGLB}", pBin, pMin, pMax));
164  pGlbTrack[1]->setAxisTitle("GeV");
165  pGlbTrack.push_back(dbe->book1D(histname+"Sta_p", "p_{STAfromGLB}", pBin, pMin, pMax));
166  pGlbTrack[2]->setAxisTitle("GeV");
167  pTrack = dbe->book1D("TkMuon_p", "p_{TK}", pBin, pMin, pMax);
168  pTrack->setAxisTitle("GeV");
169  pStaTrack = dbe->book1D("StaMuon_p", "p_{STA}", pBin, pMin, pMax);
170  pStaTrack->setAxisTitle("GeV");
171 
172  // monitoring of the transverse momentum
173  ptBin = parameters.getParameter<int>("ptBin");
174  ptMin = parameters.getParameter<double>("ptMin");
175  ptMax = parameters.getParameter<double>("ptMax");
176  ptGlbTrack.push_back(dbe->book1D(histname+"Glb_pt", "pt_{GLB}", ptBin, ptMin, ptMax));
177  ptGlbTrack[0]->setAxisTitle("GeV");
178  ptGlbTrack.push_back(dbe->book1D(histname+"Tk_pt", "pt_{TKfromGLB}", ptBin, ptMin, ptMax));
179  ptGlbTrack[1]->setAxisTitle("GeV");
180  ptGlbTrack.push_back(dbe->book1D(histname+"Sta_pt", "pt_{STAfromGLB}", ptBin, ptMin, ptMax));
181  ptGlbTrack[2]->setAxisTitle("GeV");
182  ptTrack = dbe->book1D("TkMuon_pt", "pt_{TK}", ptBin, ptMin, ptMax);
183  ptTrack->setAxisTitle("GeV");
184  ptStaTrack = dbe->book1D("StaMuon_pt", "pt_{STA}", ptBin, ptMin, pMax);
185  ptStaTrack->setAxisTitle("GeV");
186 
187  // monitoring of the muon charge
188  qGlbTrack.push_back(dbe->book1D(histname+"Glb_q", "q_{GLB}", 5, -2.5, 2.5));
189  qGlbTrack.push_back(dbe->book1D(histname+"Tk_q", "q_{TKfromGLB}", 5, -2.5, 2.5));
190  qGlbTrack.push_back(dbe->book1D(histname+"Sta_q", "q_{STAformGLB}", 5, -2.5, 2.5));
191  qGlbTrack.push_back(dbe->book1D(histname+"qComparison", "comparison between q_{GLB} and q_{TKfromGLB}, q_{STAfromGLB}", 8, 0.5, 8.5));
192  qGlbTrack[3]->setBinLabel(1,"qGlb=qSta");
193  qGlbTrack[3]->setBinLabel(2,"qGlb!=qSta");
194  qGlbTrack[3]->setBinLabel(3,"qGlb=qTk");
195  qGlbTrack[3]->setBinLabel(4,"qGlb!=qTk");
196  qGlbTrack[3]->setBinLabel(5,"qSta=qTk");
197  qGlbTrack[3]->setBinLabel(6,"qSta!=qTk");
198  qGlbTrack[3]->setBinLabel(7,"qGlb!=qSta,qGlb!=Tk");
199  qGlbTrack[3]->setBinLabel(8,"qGlb=qSta,qGlb=Tk");
200  qTrack = dbe->book1D("TkMuon_q", "q_{TK}", 5, -2.5, 2.5);
201  qStaTrack = dbe->book1D("StaMuon_q", "q_{STA}", 5, -2.5, 2.5);
202 
203  // monitoring of the momentum resolution
204  pResBin = parameters.getParameter<int>("pResBin");
205  pResMin = parameters.getParameter<double>("pResMin");
206  pResMax = parameters.getParameter<double>("pResMax");
207  qOverpResolution.push_back(dbe->book1D("Res_TkGlb_qOverp", "(q/p)_{TKfromGLB} - (q/p)_{GLB}", pResBin*binFactor*2, pResMin/10, pResMax/10));
208  qOverpResolution[0]->setAxisTitle("GeV^{-1}");
209  qOverpResolution.push_back(dbe->book1D("Res_GlbSta_qOverp", "(q/p)_{GLB} - (q/p)_{STAfromGLB}", pResBin*binFactor, pResMin, pResMax));
210  qOverpResolution[1]->setAxisTitle("GeV^{-1}");
211  qOverpResolution.push_back(dbe->book1D("Res_TkSta_qOverp", "(q/p)_{TKfromGLB} - (q/p)_{STAfromGLB}", pResBin*binFactor, pResMin, pResMax));
212  qOverpResolution[2]->setAxisTitle("GeV^{-1}");
213  qOverpPull = dbe->book1D("Pull_TkSta_qOverp", "(q/p)_{TKfromGLB} - (q/p)_{STAfromGLB} / error", 100,-10,10);
214 
215  oneOverpResolution.push_back(dbe->book1D("Res_TkGlb_oneOverp", "(1/p)_{TKfromGLB} - (1/p)_{GLB}", pResBin*binFactor*2, pResMin/10, pResMax/10));
216  oneOverpResolution[0]->setAxisTitle("GeV^{-1}");
217  oneOverpResolution.push_back(dbe->book1D("Res_GlbSta_oneOverp", "(1/p)_{GLB} - (1/p)_{STAfromGLB}", pResBin*binFactor, pResMin, pResMax));
218  oneOverpResolution[1]->setAxisTitle("GeV^{-1}");
219  oneOverpResolution.push_back(dbe->book1D("Res_TkSta_oneOverp", "(q/p)_{TKfromGLB} - (q/p)_{STAfromGLB}", pResBin*binFactor, pResMin, pResMax));
220  oneOverpResolution[2]->setAxisTitle("GeV^{-1}");
221  oneOverpPull = dbe->book1D("Pull_TkSta_oneOverp", "(1/p)_{TKfromGLB} - (1/p)_{STAfromGLB} / error", 100,-10,10);
222 
223 
224  qOverptResolution.push_back(dbe->book1D("Res_TkGlb_qOverpt", "(q/p_{t})_{TKfromGLB} - (q/p_{t})_{GLB}", pResBin*binFactor*2, pResMin/10, pResMax/10));
225  qOverptResolution[0]->setAxisTitle("GeV^{-1}");
226  qOverptResolution.push_back(dbe->book1D("Res_GlbSta_qOverpt", "(q/p_{t})_{GLB} - (q/p_{t})_{STAfromGLB}", pResBin*binFactor, pResMin, pResMax));
227  qOverptResolution[1]->setAxisTitle("GeV^{-1}");
228  qOverptResolution.push_back(dbe->book1D("Res_TkSta_qOverpt", "(q/p_{t})_{TKfromGLB} - (q/p_{t})_{STAfromGLB}", pResBin*binFactor, pResMin, pResMax));
229  qOverptResolution[2]->setAxisTitle("GeV^{-1}");
230  qOverptPull = dbe->book1D("Pull_TkSta_qOverpt", "(q/pt)_{TKfromGLB} - (q/pt)_{STAfromGLB} / error", 100,-10,10);
231 
232  oneOverptResolution.push_back(dbe->book1D("Res_TkGlb_oneOverpt", "(1/p_{t})_{TKfromGLB} - (1/p_{t})_{GLB}", pResBin*binFactor*2, pResMin/10, pResMax/10));
233  oneOverptResolution[0]->setAxisTitle("GeV^{-1}");
234  oneOverptResolution.push_back(dbe->book1D("Res_GlbSta_oneOverpt", "(1/p_{t})_{GLB} - (1/p_{t})_{STAfromGLB}", pResBin*binFactor, pResMin, pResMax));
235  oneOverptResolution[1]->setAxisTitle("GeV^{-1}");
236  oneOverptResolution.push_back(dbe->book1D("Res_TkSta_oneOverpt", "(1/p_{t})_{TKfromGLB} - (1/p_{t})_{STAfromGLB}", pResBin*binFactor, pResMin, pResMax));
237  oneOverptResolution[2]->setAxisTitle("GeV^{-1}");
238  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));
239  oneOverptResolution[3]->setAxisTitle("GeV^{-1}",2);
240  oneOverptResolution.push_back(dbe->book2D("ResVsEta_GlbSta_oneOverpt", "(#eta_{GLB} - #eta_{STAfromGLB} vs (1/p_{t})_{GLB}", etaBin, etaMin, etaMax, pResBin*binFactor, pResMin, pResMax));
241  oneOverptResolution[4]->setAxisTitle("GeV^{-1}",2);
242  oneOverptResolution.push_back(dbe->book2D("ResVsEta_TkSta_oneOverpt", "(#eta_{TKfromGLB} - #eta_{STAfromGLB}) vs (1/p_{t})_{TKfromGLB}", etaBin, etaMin, etaMax, pResBin*binFactor, pResMin, pResMax));
243  oneOverptResolution[5]->setAxisTitle("GeV^{-1}",2);
244  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));
245  oneOverptResolution[6]->setAxisTitle("rad",1);
246  oneOverptResolution[6]->setAxisTitle("GeV^{-1}",2);
247  oneOverptResolution.push_back(dbe->book2D("ResVsPhi_GlbSta_oneOverpt", "(#phi_{GLB} - #phi_{STAfromGLB} vs (1/p_{t})_{GLB}", phiBin, phiMin, phiMax, pResBin*binFactor, pResMin, pResMax));
248  oneOverptResolution[7]->setAxisTitle("rad",1);
249  oneOverptResolution[7]->setAxisTitle("GeV^{-1}",2);
250  oneOverptResolution.push_back(dbe->book2D("ResVsPhi_TkSta_oneOverpt", "(#phi_{TKfromGLB} - #phi_{STAfromGLB}) vs (1/p_{t})_{TKfromGLB}", phiBin, phiMin, phiMax, pResBin*binFactor, pResMin, pResMax));
251  oneOverptResolution[8]->setAxisTitle("rad",1);
252  oneOverptResolution[8]->setAxisTitle("GeV^{-1}",2);
253  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));
254  oneOverptResolution[9]->setAxisTitle("GeV^{-1}",1);
255  oneOverptResolution[9]->setAxisTitle("GeV^{-1}",2);
256  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));
257  oneOverptResolution[10]->setAxisTitle("GeV^{-1}",1);
258  oneOverptResolution[10]->setAxisTitle("GeV^{-1}",2);
259  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));
260  oneOverptResolution[11]->setAxisTitle("GeV^{-1}",1);
261  oneOverptResolution[11]->setAxisTitle("GeV^{-1}",2);
262  oneOverptPull = dbe->book1D("Pull_TkSta_oneOverpt", "(1/pt)_{TKfromGLB} - (1/pt)_{STAfromGLB} / error", 100,-10,10);
263 
264 
265  // monitoring of the recHits provenance
266  rhBin=parameters.getParameter<int>("rhBin");
267  rhMin=parameters.getParameter<double>("rhMin");
268  rhMax=parameters.getParameter<double>("rhMax");
269  rhAnalysis.push_back(dbe->book1D("StaRh_Frac_inGlb", "recHits_{STAinGLB} / recHits_{GLB}", rhBin, rhMin, rhMax));
270  rhAnalysis.push_back(dbe->book1D("TkRh_Frac_inGlb", "recHits_{TKinGLB} / recHits_{GLB}", rhBin, rhMin, rhMax));
271  rhAnalysis.push_back(dbe->book1D("StaRh_inGlb_Div_RhAssoSta", "recHits_{STAinGLB} / recHits_{STAfromGLB}", rhBin, rhMin, rhMax));
272  rhAnalysis.push_back(dbe->book1D("TkRh_inGlb_Div_RhAssoTk", "recHits_{TKinGLB} / recHits_{TKfromGLB}", rhBin, rhMin, rhMax));
273  rhAnalysis.push_back(dbe->book1D("GlbRh_Div_RhAssoStaTk", "recHits_{GLB} / (recHits_{TKfromGLB}+recHits_{STAfromGLB})", rhBin, rhMin, rhMax));
274  rhAnalysis.push_back(dbe->book1D("invalidRh_Frac_inTk", "Invalid recHits / rechits_{GLB}", rhBin, rhMin, rhMax));
275 
276  // monitoring of the muon system rotation w.r.t. tracker
277  muVStkSytemRotation.push_back(dbe->book2D("muVStkSytemRotation_posMu", "pT_{TK} / pT_{GLB} vs pT_{GLB} for #mu^{+}", 50,0,200,100,0.8,1.2));
278  muVStkSytemRotation.push_back(dbe->book2D("muVStkSytemRotation_negMu", "pT_{TK} / pT_{GLB} vs pT_{GLB} for #mu^{-}", 50,0,200,100,0.8,1.2));
279 
280 }
281 
282 
283 void MuonRecoAnalyzer::GetRes( reco::TrackRef t1, reco::TrackRef t2, string par, float &res, float &pull){
284 
285  float p1=0, p2=0, p1e=1, p2e=1;
286 
287  if(par == "eta") {
288  p1 = t1->eta(); p1e = t1->etaError();
289  p2 = t2->eta(); p2e = t2->etaError();
290  }
291  else if(par == "theta") {
292  p1 = t1->theta(); p1e = t1->thetaError();
293  p2 = t2->theta(); p2e = t2->thetaError();
294  }
295  else if(par == "phi") {
296  p1 = t1->phi(); p1e = t1->phiError();
297  p2 = t2->phi(); p2e = t2->phiError();
298  }
299  else if(par == "qOverp") {
300  p1 = t1->charge()/t1->p(); p1e = t1->qoverpError();
301  p2 = t2->charge()/t2->p(); p2e = t2->qoverpError();
302  }
303  else if(par == "oneOverp") {
304  p1 = 1./t1->p(); p1e = t1->qoverpError();
305  p2 = 1./t2->p(); p2e = t2->qoverpError();
306  }
307  else if(par == "qOverpt") {
308  p1 = t1->charge()/t1->pt(); p1e = t1->ptError()*p1*p1;
309  p2 = t2->charge()/t2->pt(); p2e = t2->ptError()*p2*p2;
310  }
311  else if(par == "oneOverpt") {
312  p1 = 1./t1->pt(); p1e = t1->ptError()*p1*p1;
313  p2 = 1./t2->pt(); p2e = t2->ptError()*p2*p2;
314  }
315 
316  res = p1 - p2;
317  if(p1e!=0 || p2e!=0) pull = res / sqrt(p1e*p1e + p2e*p2e);
318  else pull = -99;
319  return;
320 }
321 
322 
323 
324 
325 
326 void MuonRecoAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup, const reco::Muon& recoMu) {
327 
328  LogTrace(metname)<<"[MuonRecoAnalyzer] Analyze the mu";
329 
330  float res=0, pull=0;
331 
332  if(recoMu.isGlobalMuon()) {
333 
334  LogTrace(metname)<<"[MuonRecoAnalyzer] The mu is global - filling the histos";
335  if(recoMu.isTrackerMuon() && recoMu.isStandAloneMuon())
336  muReco->Fill(1);
337  if(!(recoMu.isTrackerMuon()) && recoMu.isStandAloneMuon())
338  muReco->Fill(2);
339  if(!recoMu.isStandAloneMuon())
340  LogTrace(metname)<<"[MuonRecoAnalyzer] ERROR: the mu is global but not standalone!";
341 
342  // get the track combinig the information from both the Tracker and the Spectrometer
343  reco::TrackRef recoCombinedGlbTrack = recoMu.combinedMuon();
344  // get the track using only the tracker data
345  reco::TrackRef recoTkGlbTrack = recoMu.track();
346  // get the track using only the mu spectrometer data
347  reco::TrackRef recoStaGlbTrack = recoMu.standAloneMuon();
348 
349 
350  etaGlbTrack[0]->Fill(recoCombinedGlbTrack->eta());
351  etaGlbTrack[1]->Fill(recoTkGlbTrack->eta());
352  etaGlbTrack[2]->Fill(recoStaGlbTrack->eta());
353 
354  GetRes(recoTkGlbTrack, recoCombinedGlbTrack, "eta", res, pull);
355  etaResolution[0]->Fill(res);
356  GetRes(recoCombinedGlbTrack, recoStaGlbTrack, "eta", res, pull);
357  etaResolution[1]->Fill(res);
358  GetRes(recoTkGlbTrack, recoStaGlbTrack, "eta", res, pull);
359  etaResolution[2]->Fill(res);
360  etaPull->Fill(pull);
361 
362  etaResolution[3]->Fill(recoCombinedGlbTrack->eta(), recoTkGlbTrack->eta()-recoCombinedGlbTrack->eta());
363  etaResolution[4]->Fill(recoCombinedGlbTrack->eta(), -recoStaGlbTrack->eta()+recoCombinedGlbTrack->eta());
364  etaResolution[5]->Fill(recoCombinedGlbTrack->eta(), recoTkGlbTrack->eta()-recoStaGlbTrack->eta());
365 
366 
367 
368  thetaGlbTrack[0]->Fill(recoCombinedGlbTrack->theta());
369  thetaGlbTrack[1]->Fill(recoTkGlbTrack->theta());
370  thetaGlbTrack[2]->Fill(recoStaGlbTrack->theta());
371  GetRes(recoTkGlbTrack, recoCombinedGlbTrack, "theta", res, pull);
372  thetaResolution[0]->Fill(res);
373 
374  GetRes(recoCombinedGlbTrack, recoStaGlbTrack, "theta", res, pull);
375  thetaResolution[1]->Fill(res);
376 
377  GetRes(recoTkGlbTrack, recoStaGlbTrack, "theta", res, pull);
378  thetaResolution[2]->Fill(res);
379  thetaPull->Fill(pull);
380 
381  thetaResolution[3]->Fill(recoCombinedGlbTrack->theta(), recoTkGlbTrack->theta()-recoCombinedGlbTrack->theta());
382  thetaResolution[4]->Fill(recoCombinedGlbTrack->theta(), -recoStaGlbTrack->theta()+recoCombinedGlbTrack->theta());
383  thetaResolution[5]->Fill(recoCombinedGlbTrack->theta(), recoTkGlbTrack->theta()-recoStaGlbTrack->theta());
384 
385 
386 
387  phiGlbTrack[0]->Fill(recoCombinedGlbTrack->phi());
388  phiGlbTrack[1]->Fill(recoTkGlbTrack->phi());
389  phiGlbTrack[2]->Fill(recoStaGlbTrack->phi());
390  GetRes(recoTkGlbTrack, recoCombinedGlbTrack, "phi", res, pull);
391  phiResolution[0]->Fill(res);
392  GetRes(recoCombinedGlbTrack, recoStaGlbTrack, "phi", res, pull);
393  phiResolution[1]->Fill(res);
394  GetRes(recoTkGlbTrack, recoStaGlbTrack, "phi", res, pull);
395  phiResolution[2]->Fill(res);
396  phiPull->Fill(pull);
397  phiResolution[3]->Fill(recoCombinedGlbTrack->phi(), recoTkGlbTrack->phi()-recoCombinedGlbTrack->phi());
398  phiResolution[4]->Fill(recoCombinedGlbTrack->phi(), -recoStaGlbTrack->phi()+recoCombinedGlbTrack->phi());
399  phiResolution[5]->Fill(recoCombinedGlbTrack->phi(), recoTkGlbTrack->phi()-recoStaGlbTrack->phi());
400 
401  chi2OvDFGlbTrack[0]->Fill(recoCombinedGlbTrack->normalizedChi2());
402  chi2OvDFGlbTrack[1]->Fill(recoTkGlbTrack->normalizedChi2());
403  chi2OvDFGlbTrack[2]->Fill(recoStaGlbTrack->normalizedChi2());
404 //-------------------------
405 // double probchi = TMath::Prob(recoCombinedGlbTrack->normalizedChi2(),recoCombinedGlbTrack->ndof());
406 // cout << "rellenando histos."<<endl;
407  probchi2GlbTrack[0]->Fill(TMath::Prob(recoCombinedGlbTrack->chi2(),recoCombinedGlbTrack->ndof()));
408  probchi2GlbTrack[1]->Fill(TMath::Prob(recoTkGlbTrack->chi2(),recoTkGlbTrack->ndof()));
409  probchi2GlbTrack[2]->Fill(TMath::Prob(recoStaGlbTrack->chi2(),recoStaGlbTrack->ndof()));
410  // cout << "rellenados histos."<<endl;
411 //-------------------------
412 
413  pGlbTrack[0]->Fill(recoCombinedGlbTrack->p());
414  pGlbTrack[1]->Fill(recoTkGlbTrack->p());
415  pGlbTrack[2]->Fill(recoStaGlbTrack->p());
416 
417  ptGlbTrack[0]->Fill(recoCombinedGlbTrack->pt());
418  ptGlbTrack[1]->Fill(recoTkGlbTrack->pt());
419  ptGlbTrack[2]->Fill(recoStaGlbTrack->pt());
420 
421  qGlbTrack[0]->Fill(recoCombinedGlbTrack->charge());
422  qGlbTrack[1]->Fill(recoTkGlbTrack->charge());
423  qGlbTrack[2]->Fill(recoStaGlbTrack->charge());
424  if(recoCombinedGlbTrack->charge()==recoStaGlbTrack->charge()) qGlbTrack[3]->Fill(1);
425  else qGlbTrack[3]->Fill(2);
426  if(recoCombinedGlbTrack->charge()==recoTkGlbTrack->charge()) qGlbTrack[3]->Fill(3);
427  else qGlbTrack[3]->Fill(4);
428  if(recoStaGlbTrack->charge()==recoTkGlbTrack->charge()) qGlbTrack[3]->Fill(5);
429  else qGlbTrack[3]->Fill(6);
430  if(recoCombinedGlbTrack->charge()!=recoStaGlbTrack->charge() && recoCombinedGlbTrack->charge()!=recoTkGlbTrack->charge()) qGlbTrack[3]->Fill(7);
431  if(recoCombinedGlbTrack->charge()==recoStaGlbTrack->charge() && recoCombinedGlbTrack->charge()==recoTkGlbTrack->charge()) qGlbTrack[3]->Fill(8);
432 
433  GetRes(recoTkGlbTrack, recoCombinedGlbTrack, "qOverp", res, pull);
434  qOverpResolution[0]->Fill(res);
435  GetRes(recoCombinedGlbTrack, recoStaGlbTrack, "qOverp", res, pull);
436  qOverpResolution[1]->Fill(res);
437  GetRes(recoTkGlbTrack, recoStaGlbTrack, "qOverp", res, pull);
438  qOverpResolution[2]->Fill(res);
439  qOverpPull->Fill(pull);
440 
441 
442  GetRes(recoTkGlbTrack, recoCombinedGlbTrack, "oneOverp", res, pull);
443  oneOverpResolution[0]->Fill(res);
444  GetRes(recoCombinedGlbTrack, recoStaGlbTrack, "oneOverp", res, pull);
445  oneOverpResolution[1]->Fill(res);
446  GetRes(recoTkGlbTrack, recoStaGlbTrack, "oneOverp", res, pull);
447  oneOverpResolution[2]->Fill(res);
448  oneOverpPull->Fill(pull);
449 
450 
451  GetRes(recoTkGlbTrack, recoCombinedGlbTrack, "qOverpt", res, pull);
452  qOverptResolution[0]->Fill(res);
453  GetRes(recoCombinedGlbTrack, recoStaGlbTrack, "qOverpt", res, pull);
454  qOverptResolution[1]->Fill(res);
455  GetRes(recoTkGlbTrack, recoStaGlbTrack, "qOverpt", res, pull);
456  qOverptResolution[2]->Fill(res);
457  qOverptPull->Fill(pull);
458 
459  GetRes(recoTkGlbTrack, recoCombinedGlbTrack, "oneOverpt", res, pull);
460  oneOverptResolution[0]->Fill(res);
461  GetRes(recoCombinedGlbTrack, recoStaGlbTrack, "oneOverpt", res, pull);
462  oneOverptResolution[1]->Fill(res);
463  GetRes(recoTkGlbTrack, recoStaGlbTrack, "oneOverpt", res, pull);
464  oneOverptResolution[2]->Fill(res);
465  oneOverptPull->Fill(pull);
466 
467  oneOverptResolution[3]->Fill(recoCombinedGlbTrack->eta(),(1/recoTkGlbTrack->pt())-(1/recoCombinedGlbTrack->pt()));
468  oneOverptResolution[4]->Fill(recoCombinedGlbTrack->eta(),-(1/recoStaGlbTrack->pt())+(1/recoCombinedGlbTrack->pt()));
469  oneOverptResolution[5]->Fill(recoCombinedGlbTrack->eta(),(1/recoTkGlbTrack->pt())-(1/recoStaGlbTrack->pt()));
470  oneOverptResolution[6]->Fill(recoCombinedGlbTrack->phi(),(1/recoTkGlbTrack->pt())-(1/recoCombinedGlbTrack->pt()));
471  oneOverptResolution[7]->Fill(recoCombinedGlbTrack->phi(),-(1/recoStaGlbTrack->pt())+(1/recoCombinedGlbTrack->pt()));
472  oneOverptResolution[8]->Fill(recoCombinedGlbTrack->phi(),(1/recoTkGlbTrack->pt())-(1/recoStaGlbTrack->pt()));
473  oneOverptResolution[9]->Fill(recoCombinedGlbTrack->pt(),(1/recoTkGlbTrack->pt())-(1/recoCombinedGlbTrack->pt()));
474  oneOverptResolution[10]->Fill(recoCombinedGlbTrack->pt(),-(1/recoStaGlbTrack->pt())+(1/recoCombinedGlbTrack->pt()));
475  oneOverptResolution[11]->Fill(recoCombinedGlbTrack->pt(),(1/recoTkGlbTrack->pt())-(1/recoStaGlbTrack->pt()));
476 
477  // valid hits Glb track
478  double rhGlb = recoCombinedGlbTrack->found();
479  // valid hits Glb track from Tracker
480  double rhGlb_StaProvenance=0;
481  // valid hits Glb track from Sta system
482  double rhGlb_TkProvenance=0;
483  for (trackingRecHit_iterator recHit = recoCombinedGlbTrack->recHitsBegin();
484  recHit!=recoCombinedGlbTrack->recHitsEnd(); ++recHit){
485  if((*recHit)->isValid()){
486  DetId id = (*recHit)->geographicalId();
487  if (id.det() == DetId::Muon)
488  rhGlb_StaProvenance++;
489  if (id.det() == DetId::Tracker)
490  rhGlb_TkProvenance++;
491  }
492  }
493  // valid hits Sta track associated to Glb track
494  double rhStaGlb = recoStaGlbTrack->recHitsSize();
495  // valid hits Traker track associated to Glb track
496  double rhTkGlb = recoTkGlbTrack->found();
497  // invalid hits Traker track associated to Glb track
498  double rhTkGlb_notValid = recoTkGlbTrack->lost();
499 
500  // fill the histos
501  rhAnalysis[0]->Fill(rhGlb_StaProvenance/rhGlb);
502  rhAnalysis[1]->Fill(rhGlb_TkProvenance/rhGlb);
503  rhAnalysis[2]->Fill(rhGlb_StaProvenance/rhStaGlb);
504  rhAnalysis[3]->Fill(rhGlb_TkProvenance/rhTkGlb);
505  rhAnalysis[4]->Fill(rhGlb/(rhStaGlb+rhTkGlb));
506  rhAnalysis[5]->Fill(rhTkGlb_notValid/rhGlb);
507 
508  // aligment plots (mu system w.r.t. tracker rotation)
509  if(recoCombinedGlbTrack->charge()>0)
510  muVStkSytemRotation[0]->Fill(recoCombinedGlbTrack->pt(),recoTkGlbTrack->pt()/recoCombinedGlbTrack->pt());
511  else
512  muVStkSytemRotation[1]->Fill(recoCombinedGlbTrack->pt(),recoTkGlbTrack->pt()/recoCombinedGlbTrack->pt());
513 
514  }
515 
516 
517  if(recoMu.isTrackerMuon() && !(recoMu.isGlobalMuon())) {
518  LogTrace(metname)<<"[MuonRecoAnalyzer] The mu is tracker only - filling the histos";
519  if(recoMu.isStandAloneMuon())
520  muReco->Fill(3);
521  if(!(recoMu.isStandAloneMuon()))
522  muReco->Fill(4);
523 
524  // get the track using only the tracker data
525  reco::TrackRef recoTrack = recoMu.track();
526 
527  etaTrack->Fill(recoTrack->eta());
528  thetaTrack->Fill(recoTrack->theta());
529  phiTrack->Fill(recoTrack->phi());
530  chi2OvDFTrack->Fill(recoTrack->normalizedChi2());
531  probchi2Track->Fill(TMath::Prob(recoTrack->chi2(),recoTrack->ndof()));
532  pTrack->Fill(recoTrack->p());
533  ptTrack->Fill(recoTrack->pt());
534  qTrack->Fill(recoTrack->charge());
535 
536  }
537 
538  if(recoMu.isStandAloneMuon() && !(recoMu.isGlobalMuon())) {
539  LogTrace(metname)<<"[MuonRecoAnalyzer] The mu is STA only - filling the histos";
540  if(!(recoMu.isTrackerMuon()))
541  muReco->Fill(5);
542 
543  // get the track using only the mu spectrometer data
544  reco::TrackRef recoStaTrack = recoMu.standAloneMuon();
545 
546  etaStaTrack->Fill(recoStaTrack->eta());
547  thetaStaTrack->Fill(recoStaTrack->theta());
548  phiStaTrack->Fill(recoStaTrack->phi());
549  chi2OvDFStaTrack->Fill(recoStaTrack->normalizedChi2());
550  probchi2StaTrack->Fill(TMath::Prob(recoStaTrack->chi2(),recoStaTrack->ndof()));
551  pStaTrack->Fill(recoStaTrack->p());
552  ptStaTrack->Fill(recoStaTrack->pt());
553  qStaTrack->Fill(recoStaTrack->charge());
554 
555  }
556 
557  if(recoMu.isCaloMuon() && !(recoMu.isGlobalMuon()) && !(recoMu.isTrackerMuon()) && !(recoMu.isStandAloneMuon()))
558  muReco->Fill(6);
559 
560  //efficiency plots
561 
562  // get the track using only the mu spectrometer data
563  reco::TrackRef recoStaGlbTrack = recoMu.standAloneMuon();
564 
565  if(recoMu.isStandAloneMuon()){
566  etaEfficiency[0]->Fill(recoStaGlbTrack->eta());
567  phiEfficiency[0]->Fill(recoStaGlbTrack->phi());
568  }
569  if(recoMu.isStandAloneMuon() && recoMu.isGlobalMuon()){
570  etaEfficiency[1]->Fill(recoStaGlbTrack->eta());
571  phiEfficiency[1]->Fill(recoStaGlbTrack->phi());
572  }
573 
574 
575 
576 
577 }
T getParameter(std::string const &) const
virtual ~MuonRecoAnalyzer()
Destructor.
MonitorElement * muReco
std::vector< MonitorElement * > oneOverpResolution
MonitorElement * qOverptPull
std::vector< MonitorElement * > ptGlbTrack
std::vector< MonitorElement * > chi2OvDFGlbTrack
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:519
void analyze(const edm::Event &, const edm::EventSetup &, const reco::Muon &recoMu)
Get the analysis.
bool isTrackerMuon() const
Definition: Muon.h:149
MonitorElement * etaStaTrack
void setBinLabel(int bin, const std::string &label, int axis=1)
set bin label for x, y or z axis (axis=1, 2, 3 respectively)
virtual TrackRef track() const
reference to a Track
Definition: Muon.h:39
bool isGlobalMuon() const
Definition: Muon.h:148
std::vector< MonitorElement * > phiGlbTrack
void beginJob(DQMStore *dbe)
Inizialize parameters for histo binning.
bool isStandAloneMuon() const
Definition: Muon.h:150
MonitorElement * chi2OvDFStaTrack
MonitorElement * phiStaTrack
MonitorElement * qOverpPull
std::vector< MonitorElement * > phiEfficiency
bool isCaloMuon() const
Definition: Muon.h:151
void Fill(long long x)
MonitorElement * etaTrack
int iEvent
Definition: GenABIO.cc:243
MonitorElement * chi2OvDFTrack
std::vector< MonitorElement * > qOverptResolution
std::vector< MonitorElement * > etaGlbTrack
std::vector< MonitorElement * > pGlbTrack
std::vector< MonitorElement * > qOverpResolution
MonitorElement * oneOverpPull
T sqrt(T t)
Definition: SSEVec.h:28
std::string metname
std::vector< MonitorElement * > qGlbTrack
MonitorElement * pTrack
std::vector< MonitorElement * > probchi2GlbTrack
std::vector< MonitorElement * > etaResolution
MonitorElement * oneOverptPull
double p2[4]
Definition: TauolaWrapper.h:90
std::vector< MonitorElement * > etaEfficiency
#define LogTrace(id)
void GetRes(reco::TrackRef t1, reco::TrackRef t2, std::string par, float &res, float &pull)
MonitorElement * thetaStaTrack
Definition: DetId.h:20
std::vector< MonitorElement * > phiResolution
virtual TrackRef combinedMuon() const
reference to a stand-alone muon Track
Definition: Muon.h:45
edm::ParameterSet parameters
std::vector< MonitorElement * > thetaGlbTrack
MonitorElement * phiPull
MonitorElement * phiTrack
std::vector< MonitorElement * > oneOverptResolution
MonitorElement * probchi2Track
std::vector< MonitorElement * > muVStkSytemRotation
MonitorElement * ptTrack
MuonRecoAnalyzer(const edm::ParameterSet &, MuonServiceProxy *theService)
Constructor.
double p1[4]
Definition: TauolaWrapper.h:89
MonitorElement * probchi2StaTrack
MonitorElement * ptStaTrack
std::vector< MonitorElement * > rhAnalysis
MonitorElement * etaPull
MonitorElement * pStaTrack
MonitorElement * book2D(const char *name, const char *title, int nchX, double lowX, double highX, int nchY, double lowY, double highY)
Book 2D histogram.
Definition: DQMStore.cc:647
MonitorElement * thetaTrack
std::vector< MonitorElement * > thetaResolution
void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
MonitorElement * qTrack
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:237
MonitorElement * thetaPull
MonitorElement * qStaTrack
const double par[8 *NPar][4]
virtual TrackRef standAloneMuon() const
reference to a stand-alone muon Track
Definition: Muon.h:42