CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Histograms.h
Go to the documentation of this file.
1 #ifndef Validation_RecoMuon_Histograms_H
2 #define Validation_RecoMuon_Histograms_H
3 
10 #include "TH1F.h"
11 #include "TH2F.h"
12 #include "TString.h"
13 #include "TFile.h"
14 
18 
20 #include <iostream>
21 #include <vector>
22 #include <math.h>
23 
24 
26 public:
27 
28  HTrackVariables(DQMStore::IBooker &ibooker,std::string dirName_, std::string name,std::string whereIs =""):theName(name.c_str()),where(whereIs.c_str()){
29  ibooker.cd();
30  std::string dirName=dirName_;
31  //dirName+="/";
32  //dirName+=name.c_str();
33 
34  ibooker.setCurrentFolder(dirName.c_str());
35 
36  hEta = ibooker.book1D(theName+"_Eta_"+where,"#eta at "+where,120,-3.,3.);
37  hPhi = ibooker.book1D(theName+"_Phi_"+where,"#phi at "+where,100,-Geom::pi(),Geom::pi());
38  hP = ibooker.book1D(theName+"_P_"+where,"p_{t} at "+where,1000,0,2000);
39  hPt = ibooker.book1D(theName+"_Pt_"+where,"p_{t} at "+where,1000,0,2000);
40  hCharge = ibooker.book1D(theName+"_charge_"+where,"Charge at "+where,4,-2.,2.);
41 
42  hEtaVsGen = ibooker.book1D(theName+"_EtaVsGen_"+where,"#eta at "+where,120,-3.,3.);
43  hPhiVsGen = ibooker.book1D(theName+"_PhiVsGen_"+where,"#phi at "+where,100,-Geom::pi(),Geom::pi());
44  hPtVsGen = ibooker.book1D(theName+"_PtVsGen_"+where,"p_{t} at "+where,1000,0,2000);
45 
46  hDeltaR = ibooker.book1D(theName+"_DeltaR_"+where,"Delta R w.r.t. sim track for "+where,1000,0,20);
47 
48  theEntries = 0;
49  }
50 
51  /*
52  HTrackVariables(std::string name, TFile* file):theName(name.c_str()){
53 
54  hEta = dynamic_cast<TH1F*>( file->Get(theName+"_Eta_"+where));
55  hPhi = dynamic_cast<TH1F*>( file->Get(theName+"_Phi_"+where));
56  hP = dynamic_cast<TH1F*>( file->Get(theName+"_P_"+where));
57  hPt = dynamic_cast<TH1F*>( file->Get(theName+"_Pt_"+where));
58  hCharge = dynamic_cast<TH1F*>( file->Get(theName+"_charge_"+where));
59  }
60  */
61 
63 
64  MonitorElement *eta() {return hEta;}
65  MonitorElement *phi() {return hPhi;}
66  MonitorElement *p() {return hP;}
67  MonitorElement *pt() {return hPt;}
69  int entries() {return theEntries;}
70 
71  void Fill(double p, double pt, double eta, double phi, double charge){
72  hEta->Fill(eta);
73  hPhi->Fill(phi);
74  hP->Fill(p);
75  hPt->Fill(pt);
76  hCharge->Fill(charge);
77  ++theEntries;
78  }
79 
80  void Fill(double pt, double eta, double phi){
81  hEtaVsGen->Fill(eta);
82  hPhiVsGen->Fill(phi);
83  hPtVsGen->Fill(pt);
84  }
85 
86  void FillDeltaR(double deltaR){
87  hDeltaR->Fill(deltaR);
88  }
89 
91 
92  efficiencyHistos.push_back(computeEfficiency(hEtaVsGen,sim->eta(),ibooker));
93  efficiencyHistos.push_back(computeEfficiency(hPhiVsGen,sim->phi(),ibooker));
94  efficiencyHistos.push_back(computeEfficiency(hPtVsGen,sim->pt(),ibooker));
95 
96  double efficiency = 100*entries()/sim->entries();
97  return efficiency;
98  }
99 
101 
102  TH1F* hReco = reco->getTH1F();
103  TH1F* hSim = sim->getTH1F();
104 
105  std::string name = hReco->GetName();
106  std::string title = hReco->GetTitle();
107 
108  MonitorElement * me = ibooker.book1D(
109  "Eff_"+name,
110  "Efficiecny as a function of "+title,
111  hSim->GetNbinsX(),
112  hSim->GetXaxis()->GetXmin(),
113  hSim->GetXaxis()->GetXmax()
114  );
115 
116  me->getTH1F()->Divide(hReco,hSim,1.,1.,"b");
117 
118  // Set the error accordingly to binomial statistics
119  int nBinsEta = me->getTH1F()->GetNbinsX();
120  for(int bin = 1; bin <= nBinsEta; bin++) {
121  float nSimHit = hSim->GetBinContent(bin);
122  float eff = me->getTH1F()->GetBinContent(bin);
123  float error = 0;
124  if(nSimHit != 0 && eff <= 1) {
125  error = sqrt(eff*(1-eff)/nSimHit);
126  }
127  me->getTH1F()->SetBinError(bin, error);
128  }
129 
130  return me;
131  }
132 
133 
134  private:
135 
138 
140 
146 
150 
152 
153  std::vector<MonitorElement*> efficiencyHistos;
154 
155 
156 };
157 
158 
159 class HResolution {
160 public:
161 
162  HResolution(DQMStore::IBooker &ibooker, std::string dirName_,std::string name,std::string whereIs):theName(name.c_str()),where(whereIs.c_str()){
163 
164  ibooker.cd();
165  std::string dirName=dirName_;
166  //dirName+="/";
167  //dirName+=name.c_str();
168 
169  ibooker.setCurrentFolder(dirName.c_str());
170 
171  double eta = 15.; int neta = 800;
172  double phi = 12.; int nphi = 400;
173  double pt = 60.; int npt = 2000;
174 
175  hEta = ibooker.book1D(theName+"_Eta_"+where,"#eta "+theName,neta,-eta,eta); // 400
176  hPhi = ibooker.book1D(theName+"_Phi_"+where,"#phi "+theName,nphi,-phi,phi); // 100
177 
178  hP = ibooker.book1D(theName+"_P_"+where,"P "+theName,400,-4,4); // 200
179  hPt = ibooker.book1D(theName+"_Pt_"+where,"P_{t} "+theName,npt,-pt,pt); // 200
180 
181  hCharge = ibooker.book1D(theName+"_charge_"+where,"Charge "+theName,4,-2.,2.);
182 
183 
184  h2Eta = ibooker.book2D(theName+"_Eta_vs_Eta"+where,"#eta "+theName+" as a function of #eta",200,-2.5,2.5,neta,-eta,eta);
185  h2Phi = ibooker.book2D(theName+"_Phi_vs_Phi"+where,"#phi "+theName+" as a function of #phi",100,-Geom::pi(),Geom::pi(),nphi,-phi,phi);
186 
187  h2P = ibooker.book2D(theName+"_P_vs_P"+where,"P "+theName+" as a function of P",1000,0,2000,400,-4,4);
188  h2Pt = ibooker.book2D(theName+"_Pt_vs_Pt"+where,"P_{t} "+theName+" as a function of P_{t}",1000,0,2000,npt,-pt,pt);
189 
190  h2PtVsEta = ibooker.book2D(theName+"_Pt_vs_Eta"+where,"P_{t} "+theName+" as a function of #eta",200,-2.5,2.5,npt,-pt,pt);
191  h2PtVsPhi = ibooker.book2D(theName+"_Pt_vs_Phi"+where,"P_{t} "+theName+" as a function of #phi",100,-Geom::pi(),Geom::pi(),npt,-pt,pt);
192 
193  h2EtaVsPt = ibooker.book2D(theName+"_Eta_vs_Pt"+where,"#eta "+theName+" as a function of P_{t}",1000,0,2000,neta,-eta,eta);
194  h2EtaVsPhi = ibooker.book2D(theName+"_Eta_vs_Phi"+where,"#eta "+theName+" as a function of #phi",100,-Geom::pi(),Geom::pi(),neta,-eta,eta);
195 
196  h2PhiVsPt = ibooker.book2D(theName+"_Phi_vs_Pt"+where,"#phi "+theName+" as a function of P_{t}",1000,0,2000,nphi,-phi,phi);
197  h2PhiVsEta = ibooker.book2D(theName+"_Phi_vs_Eta"+where,"#phi "+theName+" as a function of #eta",200,-2.5,2.5,nphi,-phi,phi);
198  }
199 
200  HResolution(DQMStore::IBooker &ibooker, std::string name, TFile* file):theName(name.c_str()){
201  // dynamic_cast<TH1F*>( file->Get(theName+"") );
202  }
203 
205 
206 
207  void Fill(double p, double pt, double eta, double phi,
208  double rp, double rpt,double reta, double rphi, double rcharge){
209 
210  Fill(rp, rpt, reta, rphi, rcharge);
211 
212 
213  h2Eta->Fill(eta,reta);
214  h2Phi->Fill(phi,rphi);
215 
216  h2P->Fill(p,rp);
217  h2Pt->Fill(pt,rpt);
218 
219  h2PtVsEta->Fill(eta,rpt);
220  h2PtVsPhi->Fill(phi,rpt);
221 
222  h2EtaVsPt ->Fill(pt,reta);
223  h2EtaVsPhi->Fill(phi,reta);
224 
225  h2PhiVsPt ->Fill(pt,rphi);
226  h2PhiVsEta->Fill(eta,rphi);
227  }
228 
229 
230  void Fill(double p, double pt, double eta, double phi,
231  double rp, double rpt){
232 
233  hP->Fill(rp);
234  hPt->Fill(rpt);
235 
236  h2P->Fill(p,rp);
237  // h2PVsEta->Fill(eta,rp);
238  // h2PVsPhi->Fill(phi,rp);
239 
240  h2Pt->Fill(pt,rpt);
241  h2PtVsEta->Fill(eta,rpt);
242  h2PtVsPhi->Fill(phi,rpt);
243  }
244 
245  void Fill(double rp, double rpt,
246  double reta, double rphi, double rcharge){
247 
248  hEta->Fill(reta);
249  hPhi->Fill(rphi);
250 
251  hP->Fill(rp);
252  hPt->Fill(rpt);
253 
254  hCharge->Fill(rcharge);
255  }
256 
257 
258 
259 private:
260 
263 
266 
269 
271 
274 
277 
280 
283 
286 
287 };
288 
289 
291  public:
293 
294  // Position, sigma, residual, pull
295  hResX = ibooker.book1D (theName+"_X_Res", "X residual", 5000, -0.5,0.5);
296  hResY = ibooker.book1D (theName+"_Y_Res", "Y residual", 5000, -1.,1.);
297  hResZ = ibooker.book1D (theName+"_Z_Res", "Z residual", 5000, -1.5,1.5);
298 
299  hResXVsEta = ibooker.book2D(theName+"_XResVsEta", "X residual vs eta",
300  200, -2.5,2.5, 5000, -1.5,1.5);
301  hResYVsEta = ibooker.book2D(theName+"_YResVsEta", "Y residual vs eta",
302  200, -2.5,2.5, 5000, -1.,1.);
303  hResZVsEta = ibooker.book2D(theName+"_ZResVsEta", "Z residual vs eta",
304  200, -2.5,2.5, 5000, -1.5,1.5);
305 
306  hXResVsPhi = ibooker.book2D(theName+"_XResVsPhi", "X residual vs phi",
307  100,-Geom::pi(),Geom::pi(), 5000, -0.5,0.5);
308  hYResVsPhi = ibooker.book2D(theName+"_YResVsPhi", "Y residual vs phi",
309  100,-Geom::pi(),Geom::pi(), 5000, -1.,1.);
310  hZResVsPhi = ibooker.book2D(theName+"_ZResVsPhi", "Z residual vs phi",
311  100,-Geom::pi(),Geom::pi(), 5000, -1.5,1.5);
312 
313  hXResVsPos = ibooker.book2D(theName+"_XResVsPos", "X residual vs position",
314  10000, -750,750, 5000, -0.5,0.5);
315  hYResVsPos = ibooker.book2D(theName+"_YResVsPos", "Y residual vs position",
316  10000, -740,740, 5000, -1.,1.);
317  hZResVsPos = ibooker.book2D(theName+"_ZResVsPos", "Z residual vs position",
318  10000, -1100,1100, 5000, -1.5,1.5);
319 
320  hXPull = ibooker.book1D (theName+"_XPull", "X pull", 600, -2,2);
321  hYPull = ibooker.book1D (theName+"_YPull", "Y pull", 600, -2,2);
322  hZPull = ibooker.book1D (theName+"_ZPull", "Z pull", 600, -2,2);
323 
324  hXPullVsPos = ibooker.book2D (theName+"_XPullVsPos", "X pull vs position",10000, -750,750, 600, -2,2);
325  hYPullVsPos = ibooker.book2D (theName+"_YPullVsPos", "Y pull vs position",10000, -740,740, 600, -2,2);
326  hZPullVsPos = ibooker.book2D (theName+"_ZPullVsPos", "Z pull vs position",10000, -1100,1100, 600, -2,2);
327  }
328 
329  HResolution1DRecHit(DQMStore::IBooker &ibooker, const TString& name_, TFile* file){}
330 
332 
333  void Fill(double x, double y, double z,
334  double dx, double dy, double dz,
335  double errx, double erry, double errz,
336  double eta, double phi) {
337 
338  double rx = dx/x, ry = dy/y, rz = dz/z;
339 
340  hResX->Fill(rx);
341  hResY->Fill(ry);
342  hResZ->Fill(rz);
343 
344  hResXVsEta->Fill(eta,rx);
345  hResYVsEta->Fill(eta,ry);
346  hResZVsEta->Fill(eta,rz);
347 
348  hXResVsPhi->Fill(phi,rx);
349  hYResVsPhi->Fill(phi,ry);
350  hZResVsPhi->Fill(phi,rz);
351 
352  hXResVsPos->Fill(x,rx);
353  hYResVsPos->Fill(y,ry);
354  hZResVsPos->Fill(z,rz);
355 
356  if(errx < 1e-6)
357  std::cout << "NO proper error set for X: "<<errx<<std::endl;
358  else{
359  hXPull->Fill(dx/errx);
360  hXPullVsPos->Fill(x,dx/errx);
361  }
362 
363  if(erry < 1e-6)
364  std::cout << "NO proper error set for Y: "<<erry<<std::endl;
365  else{
366  hYPull->Fill(dy/erry);
367  hYPullVsPos->Fill(y,dy/erry);
368  }
369  if(errz < 1e-6)
370  std::cout << "NO proper error set for Z: "<<errz<<std::endl;
371  else{
372  hZPull->Fill(dz/errz);
373  hZPullVsPos->Fill(z,dz/errz);
374  }
375  }
376 
377 
378 
379  public:
381 
385 
389 
393 
397 
401 
405 
406 };
407 #endif
408 
MonitorElement * hResXVsEta
Definition: Histograms.h:386
MonitorElement * hZPull
Definition: Histograms.h:400
MonitorElement * h2PhiVsPt
Definition: Histograms.h:284
MonitorElement * h2Eta
Definition: Histograms.h:272
double computeEfficiency(HTrackVariables *sim, DQMStore::IBooker &ibooker)
Definition: Histograms.h:90
MonitorElement * h2EtaVsPhi
Definition: Histograms.h:282
void Fill(double p, double pt, double eta, double phi, double rp, double rpt, double reta, double rphi, double rcharge)
Definition: Histograms.h:207
MonitorElement * hCharge
Definition: Histograms.h:145
std::string where
Definition: Histograms.h:262
MonitorElement * hPt
Definition: Histograms.h:144
MonitorElement * hP
Definition: Histograms.h:267
void cd(void)
Definition: DQMStore.cc:266
HResolution(DQMStore::IBooker &ibooker, std::string dirName_, std::string name, std::string whereIs)
Definition: Histograms.h:162
MonitorElement * hPt
Definition: Histograms.h:268
MonitorElement * hPhi
Definition: Histograms.h:142
MonitorElement * hXPullVsPos
Definition: Histograms.h:402
MonitorElement * h2P
Definition: Histograms.h:275
HResolution1DRecHit(DQMStore::IBooker &ibooker, const TString &name_, TFile *file)
Definition: Histograms.h:329
MonitorElement * hCharge
Definition: Histograms.h:270
MonitorElement * p()
Definition: Histograms.h:66
HResolution(DQMStore::IBooker &ibooker, std::string name, TFile *file)
Definition: Histograms.h:200
HResolution1DRecHit(DQMStore::IBooker &ibooker, std::string name)
Definition: Histograms.h:292
MonitorElement * hPhi
Definition: Histograms.h:265
MonitorElement * h2PhiVsEta
Definition: Histograms.h:285
MonitorElement * eta()
Definition: Histograms.h:64
Definition: sim.h:19
void Fill(double pt, double eta, double phi)
Definition: Histograms.h:80
MonitorElement * hP
Definition: Histograms.h:143
MonitorElement * hXPull
Definition: Histograms.h:398
std::string theName
Definition: Histograms.h:261
void Fill(long long x)
MonitorElement * hYPullVsPos
Definition: Histograms.h:403
MonitorElement * hYResVsPhi
Definition: Histograms.h:391
std::vector< MonitorElement * > efficiencyHistos
Definition: Histograms.h:153
MonitorElement * hPhiVsGen
Definition: Histograms.h:148
std::string theName
Definition: Histograms.h:380
MonitorElement * hResYVsEta
Definition: Histograms.h:387
HTrackVariables(DQMStore::IBooker &ibooker, std::string dirName_, std::string name, std::string whereIs="")
Definition: Histograms.h:28
MonitorElement * phi()
Definition: Histograms.h:65
MonitorElement * hResZVsEta
Definition: Histograms.h:388
T sqrt(T t)
Definition: SSEVec.h:48
MonitorElement * h2Phi
Definition: Histograms.h:273
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
void Fill(double rp, double rpt, double reta, double rphi, double rcharge)
Definition: Histograms.h:245
MonitorElement * h2Pt
Definition: Histograms.h:276
MonitorElement * computeEfficiency(MonitorElement *reco, MonitorElement *sim, DQMStore::IBooker &ibooker)
Definition: Histograms.h:100
MonitorElement * hXResVsPos
Definition: Histograms.h:394
void Fill(double p, double pt, double eta, double phi, double rp, double rpt)
Definition: Histograms.h:230
MonitorElement * hZPullVsPos
Definition: Histograms.h:404
MonitorElement * hPtVsGen
Definition: Histograms.h:149
MonitorElement * charge()
Definition: Histograms.h:68
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
MonitorElement * pt()
Definition: Histograms.h:67
MonitorElement * hYPull
Definition: Histograms.h:399
void Fill(double p, double pt, double eta, double phi, double charge)
Definition: Histograms.h:71
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
MonitorElement * hXResVsPhi
Definition: Histograms.h:390
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:133
MonitorElement * hZResVsPhi
Definition: Histograms.h:392
TH1F * getTH1F(void) const
MonitorElement * hEtaVsGen
Definition: Histograms.h:147
MonitorElement * hYResVsPos
Definition: Histograms.h:395
double pi()
Definition: Pi.h:31
std::string theName
Definition: Histograms.h:136
MonitorElement * h2EtaVsPt
Definition: Histograms.h:281
MonitorElement * h2PtVsPhi
Definition: Histograms.h:279
tuple cout
Definition: gather_cfg.py:121
void Fill(double x, double y, double z, double dx, double dy, double dz, double errx, double erry, double errz, double eta, double phi)
Definition: Histograms.h:333
MonitorElement * hResY
Definition: Histograms.h:383
MonitorElement * hZResVsPos
Definition: Histograms.h:396
virtual Int_t Fill(Double_t x, Double_t y)
Definition: Histograms.h:1739
MonitorElement * hEta
Definition: Histograms.h:264
MonitorElement * hResZ
Definition: Histograms.h:384
MonitorElement * h2PtVsEta
Definition: Histograms.h:278
void FillDeltaR(double deltaR)
Definition: Histograms.h:86
MonitorElement * hDeltaR
Definition: Histograms.h:151
MonitorElement * hResX
Definition: Histograms.h:382
MonitorElement * hEta
Definition: Histograms.h:141
std::string where
Definition: Histograms.h:137