CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Protected Attributes | Private Attributes
HResolution Class Reference

#include <Histograms.h>

Inheritance diagram for HResolution:

Public Member Functions

void Fill (double p, double pt, double eta, double phi, double rp, double rpt, double reta, double rphi, double rcharge)
 
void Fill (double p, double pt, double eta, double phi, double rp, double rpt)
 
void Fill (double rp, double rpt, double reta, double rphi, double rcharge)
 
virtual Int_t Fill (Double_t x, Double_t y)
 
 HResolution (std::string dirName_, std::string name, std::string whereIs)
 
 HResolution (std::string name, TFile *file)
 
 HResolution (const TString &name, const TString &title, const int totBins, const double &xMin, const double &xMax, const double &yMin, const double &yMax, TDirectory *dir=0)
 
virtual Int_t Write (const char *name=0, Int_t option=0, Int_t bufsize=0)
 
 ~HResolution ()
 
 ~HResolution ()
 

Protected Attributes

TDirectory * diffDir_
 
TProfile * diffHisto_
 
TDirectory * dir2D_
 
TDirectory * dir_
 
TH2F * histo2D_
 
TH1F * resoHisto_
 

Private Attributes

DQMStoredbe_
 
MonitorElementh2Eta
 
MonitorElementh2EtaVsPhi
 
MonitorElementh2EtaVsPt
 
MonitorElementh2P
 
MonitorElementh2Phi
 
MonitorElementh2PhiVsEta
 
MonitorElementh2PhiVsPt
 
MonitorElementh2Pt
 
MonitorElementh2PtVsEta
 
MonitorElementh2PtVsPhi
 
MonitorElementhCharge
 
MonitorElementhEta
 
MonitorElementhP
 
MonitorElementhPhi
 
MonitorElementhPt
 
std::string theName
 
std::string where
 

Detailed Description

This histogram class can be used to evaluate the resolution of a variable. It has a TProfile, a TH2F and a TH1F. The TProfile is used to compute the rms of the distribution which is filled in the TH1F (the resolution histogram) in the Write method. If a TDirectory is passed to the constructor, the different histograms are placed in subdirectories.

Definition at line 1714 of file Histograms.h.

Constructor & Destructor Documentation

HResolution::HResolution ( const TString &  name,
const TString &  title,
const int  totBins,
const double &  xMin,
const double &  xMax,
const double &  yMin,
const double &  yMax,
TDirectory *  dir = 0 
)
inline

Definition at line 1716 of file Histograms.h.

References diffDir_, diffHisto_, dir, dir2D_, dir_, histo2D_, and resoHisto_.

1718  :
1719  dir_(dir),
1720  dir2D_(0),
1721  diffDir_(0)
1722  {
1723  if( dir_ != 0 ) {
1724  dir2D_ = (TDirectory*) dir_->Get("2D");
1725  if(dir2D_ == 0) dir2D_ = dir_->mkdir("2D");
1726  diffDir_ = (TDirectory*) dir_->Get("deltaXoverX");
1727  if(diffDir_ == 0) diffDir_ = dir->mkdir("deltaXoverX");
1728  }
1729  diffHisto_ = new TProfile(name+"_prof", title+" profile", totBins, xMin, xMax, yMin, yMax);
1730  histo2D_ = new TH2F(name+"2D", title, totBins, xMin, xMax, 4000, yMin, yMax);
1731  resoHisto_ = new TH1F(name, title, totBins, xMin, xMax);
1732  }
TDirectory * dir2D_
Definition: Histograms.h:1767
TDirectory * dir_
Definition: Histograms.h:1766
TDirectory * diffDir_
Definition: Histograms.h:1768
TH1F * resoHisto_
Definition: Histograms.h:1771
dbl *** dir
Definition: mlp_gen.cc:35
TH2F * histo2D_
Definition: Histograms.h:1770
TProfile * diffHisto_
Definition: Histograms.h:1769
HResolution::~HResolution ( )
inline

Definition at line 1733 of file Histograms.h.

References diffHisto_, histo2D_, and resoHisto_.

1733  {
1734  delete diffHisto_;
1735  delete histo2D_;
1736  delete resoHisto_;
1737  }
TH1F * resoHisto_
Definition: Histograms.h:1771
TH2F * histo2D_
Definition: Histograms.h:1770
TProfile * diffHisto_
Definition: Histograms.h:1769
HResolution::HResolution ( std::string  dirName_,
std::string  name,
std::string  whereIs 
)
inline

Definition at line 165 of file Histograms.h.

References DQMStore::book1D(), DQMStore::book2D(), DQMStore::cd(), dbe_, TrackerOfflineValidation_Dqm_cff::dirName, eta(), h2Eta, h2EtaVsPhi, h2EtaVsPt, h2P, h2Phi, h2PhiVsEta, h2PhiVsPt, h2Pt, h2PtVsEta, h2PtVsPhi, hCharge, hEta, hP, hPhi, hPt, cppFunctionSkipper::operator, phi, Geom::pi(), RecoTauCleanerPlugins::pt, DQMStore::setCurrentFolder(), AlCaHLTBitMon_QueryRunRegistry::string, theName, and where.

165  :theName(name.c_str()),where(whereIs.c_str()){
166 
168  dbe_->cd();
169  std::string dirName=dirName_;
170  //dirName+="/";
171  //dirName+=name.c_str();
172 
173  dbe_->setCurrentFolder(dirName.c_str());
174 
175  double eta = 15.; int neta = 800;
176  double phi = 12.; int nphi = 400;
177  double pt = 60.; int npt = 2000;
178 
179  hEta = dbe_->book1D(theName+"_Eta_"+where,"#eta "+theName,neta,-eta,eta); // 400
180  hPhi = dbe_->book1D(theName+"_Phi_"+where,"#phi "+theName,nphi,-phi,phi); // 100
181 
182  hP = dbe_->book1D(theName+"_P_"+where,"P "+theName,400,-4,4); // 200
183  hPt = dbe_->book1D(theName+"_Pt_"+where,"P_{t} "+theName,npt,-pt,pt); // 200
184 
185  hCharge = dbe_->book1D(theName+"_charge_"+where,"Charge "+theName,4,-2.,2.);
186 
187 
188  h2Eta = dbe_->book2D(theName+"_Eta_vs_Eta"+where,"#eta "+theName+" as a function of #eta",200,-2.5,2.5,neta,-eta,eta);
189  h2Phi = dbe_->book2D(theName+"_Phi_vs_Phi"+where,"#phi "+theName+" as a function of #phi",100,-Geom::pi(),Geom::pi(),nphi,-phi,phi);
190 
191  h2P = dbe_->book2D(theName+"_P_vs_P"+where,"P "+theName+" as a function of P",1000,0,2000,400,-4,4);
192  h2Pt = dbe_->book2D(theName+"_Pt_vs_Pt"+where,"P_{t} "+theName+" as a function of P_{t}",1000,0,2000,npt,-pt,pt);
193 
194  h2PtVsEta = dbe_->book2D(theName+"_Pt_vs_Eta"+where,"P_{t} "+theName+" as a function of #eta",200,-2.5,2.5,npt,-pt,pt);
195  h2PtVsPhi = dbe_->book2D(theName+"_Pt_vs_Phi"+where,"P_{t} "+theName+" as a function of #phi",100,-Geom::pi(),Geom::pi(),npt,-pt,pt);
196 
197  h2EtaVsPt = dbe_->book2D(theName+"_Eta_vs_Pt"+where,"#eta "+theName+" as a function of P_{t}",1000,0,2000,neta,-eta,eta);
198  h2EtaVsPhi = dbe_->book2D(theName+"_Eta_vs_Phi"+where,"#eta "+theName+" as a function of #phi",100,-Geom::pi(),Geom::pi(),neta,-eta,eta);
199 
200  h2PhiVsPt = dbe_->book2D(theName+"_Phi_vs_Pt"+where,"#phi "+theName+" as a function of P_{t}",1000,0,2000,nphi,-phi,phi);
201  h2PhiVsEta = dbe_->book2D(theName+"_Phi_vs_Eta"+where,"#phi "+theName+" as a function of #eta",200,-2.5,2.5,nphi,-phi,phi);
202  }
DQMStore * dbe_
Definition: Histograms.h:264
MonitorElement * h2PhiVsPt
Definition: Histograms.h:289
MonitorElement * h2Eta
Definition: Histograms.h:277
MonitorElement * h2EtaVsPhi
Definition: Histograms.h:287
std::string where
Definition: Histograms.h:267
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:872
MonitorElement * hP
Definition: Histograms.h:272
void cd(void)
go to top directory (ie. root)
Definition: DQMStore.cc:561
MonitorElement * hPt
Definition: Histograms.h:273
MonitorElement * h2P
Definition: Histograms.h:280
MonitorElement * hCharge
Definition: Histograms.h:275
MonitorElement * hPhi
Definition: Histograms.h:270
T eta() const
MonitorElement * h2PhiVsEta
Definition: Histograms.h:290
std::string theName
Definition: Histograms.h:266
MonitorElement * h2Phi
Definition: Histograms.h:278
MonitorElement * h2Pt
Definition: Histograms.h:281
double pi()
Definition: Pi.h:31
MonitorElement * h2EtaVsPt
Definition: Histograms.h:286
MonitorElement * h2PtVsPhi
Definition: Histograms.h:284
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:1000
MonitorElement * hEta
Definition: Histograms.h:269
MonitorElement * h2PtVsEta
Definition: Histograms.h:283
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:584
Definition: DDAxes.h:10
HResolution::HResolution ( std::string  name,
TFile *  file 
)
inline

Definition at line 204 of file Histograms.h.

204  :theName(name.c_str()){
205  // dynamic_cast<TH1F*>( file->Get(theName+"") );
206  }
std::string theName
Definition: Histograms.h:266
HResolution::~HResolution ( )
inline

Definition at line 208 of file Histograms.h.

208 {}

Member Function Documentation

void HResolution::Fill ( double  p,
double  pt,
double  eta,
double  phi,
double  rp,
double  rpt,
double  reta,
double  rphi,
double  rcharge 
)
inline

Definition at line 211 of file Histograms.h.

References MonitorElement::Fill(), Fill(), h2Eta, h2EtaVsPhi, h2EtaVsPt, h2P, h2Phi, h2PhiVsEta, h2PhiVsPt, h2Pt, h2PtVsEta, and h2PtVsPhi.

212  {
213 
214  Fill(rp, rpt, reta, rphi, rcharge);
215 
216 
217  h2Eta->Fill(eta,reta);
218  h2Phi->Fill(phi,rphi);
219 
220  h2P->Fill(p,rp);
221  h2Pt->Fill(pt,rpt);
222 
223  h2PtVsEta->Fill(eta,rpt);
224  h2PtVsPhi->Fill(phi,rpt);
225 
226  h2EtaVsPt ->Fill(pt,reta);
227  h2EtaVsPhi->Fill(phi,reta);
228 
229  h2PhiVsPt ->Fill(pt,rphi);
230  h2PhiVsEta->Fill(eta,rphi);
231  }
MonitorElement * h2PhiVsPt
Definition: Histograms.h:289
MonitorElement * h2Eta
Definition: Histograms.h:277
MonitorElement * h2EtaVsPhi
Definition: Histograms.h:287
MonitorElement * h2P
Definition: Histograms.h:280
T eta() const
MonitorElement * h2PhiVsEta
Definition: Histograms.h:290
void Fill(long long x)
MonitorElement * h2Phi
Definition: Histograms.h:278
MonitorElement * h2Pt
Definition: Histograms.h:281
MonitorElement * h2EtaVsPt
Definition: Histograms.h:286
MonitorElement * h2PtVsPhi
Definition: Histograms.h:284
virtual Int_t Fill(Double_t x, Double_t y)
Definition: Histograms.h:1738
MonitorElement * h2PtVsEta
Definition: Histograms.h:283
Definition: DDAxes.h:10
void HResolution::Fill ( double  p,
double  pt,
double  eta,
double  phi,
double  rp,
double  rpt 
)
inline

Definition at line 234 of file Histograms.h.

References MonitorElement::Fill(), h2P, h2Pt, h2PtVsEta, h2PtVsPhi, hP, and hPt.

235  {
236 
237  hP->Fill(rp);
238  hPt->Fill(rpt);
239 
240  h2P->Fill(p,rp);
241  // h2PVsEta->Fill(eta,rp);
242  // h2PVsPhi->Fill(phi,rp);
243 
244  h2Pt->Fill(pt,rpt);
245  h2PtVsEta->Fill(eta,rpt);
246  h2PtVsPhi->Fill(phi,rpt);
247  }
MonitorElement * hP
Definition: Histograms.h:272
MonitorElement * hPt
Definition: Histograms.h:273
MonitorElement * h2P
Definition: Histograms.h:280
T eta() const
void Fill(long long x)
MonitorElement * h2Pt
Definition: Histograms.h:281
MonitorElement * h2PtVsPhi
Definition: Histograms.h:284
MonitorElement * h2PtVsEta
Definition: Histograms.h:283
Definition: DDAxes.h:10
void HResolution::Fill ( double  rp,
double  rpt,
double  reta,
double  rphi,
double  rcharge 
)
inline

Definition at line 249 of file Histograms.h.

References MonitorElement::Fill(), hCharge, hEta, hP, hPhi, and hPt.

250  {
251 
252  hEta->Fill(reta);
253  hPhi->Fill(rphi);
254 
255  hP->Fill(rp);
256  hPt->Fill(rpt);
257 
258  hCharge->Fill(rcharge);
259  }
MonitorElement * hP
Definition: Histograms.h:272
MonitorElement * hPt
Definition: Histograms.h:273
MonitorElement * hCharge
Definition: Histograms.h:275
MonitorElement * hPhi
Definition: Histograms.h:270
void Fill(long long x)
MonitorElement * hEta
Definition: Histograms.h:269
virtual Int_t HResolution::Fill ( Double_t  x,
Double_t  y 
)
inlinevirtual

Definition at line 1738 of file Histograms.h.

References diffHisto_, and histo2D_.

Referenced by HTrack::computePull(), HTrack::computeResolution(), HTrack::computeTDRResolution(), and Fill().

1738  {
1739  diffHisto_->Fill(x, y);
1740  histo2D_->Fill(x, y);
1741  return 0;
1742  }
Definition: DDAxes.h:10
TH2F * histo2D_
Definition: Histograms.h:1770
TProfile * diffHisto_
Definition: Histograms.h:1769
virtual Int_t HResolution::Write ( const char *  name = 0,
Int_t  option = 0,
Int_t  bufsize = 0 
)
inlinevirtual

Definition at line 1743 of file Histograms.h.

References diffDir_, diffHisto_, dir2D_, dir_, histo2D_, resoHisto_, and mathSSE::sqrt().

1743  {
1744  // Loop on all the bins and take the rms.
1745  // The TProfile bin error is by default the standard error on the mean, that is
1746  // rms/sqrt(N). If it is created with the "S" option (as we did NOT do), it would
1747  // already be the rms. Thus we take the error and multiply it by the sqrt of the
1748  // bin entries to get the rms.
1749  // bin 0 is the underflow, bin totBins+1 is the overflow.
1750  unsigned int totBins = diffHisto_->GetNbinsX();
1751  // std::cout << "totBins = " << totBins << std::endl;
1752  for( unsigned int iBin=1; iBin<=totBins; ++iBin ) {
1753  // std::cout << "iBin = " << iBin << ", " << diffHisto_->GetBinError(iBin)*sqrt(diffHisto_->GetBinEntries(iBin)) << std::endl;
1754  resoHisto_->SetBinContent( iBin, diffHisto_->GetBinError(iBin)*sqrt(diffHisto_->GetBinEntries(iBin)) );
1755  }
1756  if( dir_ != 0 ) dir_->cd();
1757  resoHisto_->Write();
1758  if( diffDir_ != 0 ) diffDir_->cd();
1759  diffHisto_->Write();
1760  if( dir2D_ != 0 ) dir2D_->cd();
1761  histo2D_->Write();
1762 
1763  return 0;
1764  }
TDirectory * dir2D_
Definition: Histograms.h:1767
T sqrt(T t)
Definition: SSEVec.h:48
TDirectory * dir_
Definition: Histograms.h:1766
TDirectory * diffDir_
Definition: Histograms.h:1768
TH1F * resoHisto_
Definition: Histograms.h:1771
TH2F * histo2D_
Definition: Histograms.h:1770
TProfile * diffHisto_
Definition: Histograms.h:1769

Member Data Documentation

DQMStore* HResolution::dbe_
private

Definition at line 264 of file Histograms.h.

Referenced by HResolution().

TDirectory* HResolution::diffDir_
protected

Definition at line 1768 of file Histograms.h.

Referenced by HResolution(), and Write().

TProfile* HResolution::diffHisto_
protected

Definition at line 1769 of file Histograms.h.

Referenced by Fill(), HResolution(), Write(), and ~HResolution().

TDirectory* HResolution::dir2D_
protected

Definition at line 1767 of file Histograms.h.

Referenced by HResolution(), and Write().

TDirectory* HResolution::dir_
protected

Definition at line 1766 of file Histograms.h.

Referenced by HResolution(), and Write().

MonitorElement* HResolution::h2Eta
private

Definition at line 277 of file Histograms.h.

Referenced by Fill(), and HResolution().

MonitorElement* HResolution::h2EtaVsPhi
private

Definition at line 287 of file Histograms.h.

Referenced by Fill(), and HResolution().

MonitorElement* HResolution::h2EtaVsPt
private

Definition at line 286 of file Histograms.h.

Referenced by Fill(), and HResolution().

MonitorElement* HResolution::h2P
private

Definition at line 280 of file Histograms.h.

Referenced by Fill(), and HResolution().

MonitorElement* HResolution::h2Phi
private

Definition at line 278 of file Histograms.h.

Referenced by Fill(), and HResolution().

MonitorElement* HResolution::h2PhiVsEta
private

Definition at line 290 of file Histograms.h.

Referenced by Fill(), and HResolution().

MonitorElement* HResolution::h2PhiVsPt
private

Definition at line 289 of file Histograms.h.

Referenced by Fill(), and HResolution().

MonitorElement* HResolution::h2Pt
private

Definition at line 281 of file Histograms.h.

Referenced by Fill(), and HResolution().

MonitorElement* HResolution::h2PtVsEta
private

Definition at line 283 of file Histograms.h.

Referenced by Fill(), and HResolution().

MonitorElement* HResolution::h2PtVsPhi
private

Definition at line 284 of file Histograms.h.

Referenced by Fill(), and HResolution().

MonitorElement* HResolution::hCharge
private

Definition at line 275 of file Histograms.h.

Referenced by Fill(), and HResolution().

MonitorElement* HResolution::hEta
private

Definition at line 269 of file Histograms.h.

Referenced by Fill(), and HResolution().

TH2F* HResolution::histo2D_
protected

Definition at line 1770 of file Histograms.h.

Referenced by Fill(), HResolution(), Write(), and ~HResolution().

MonitorElement* HResolution::hP
private

Definition at line 272 of file Histograms.h.

Referenced by Fill(), and HResolution().

MonitorElement* HResolution::hPhi
private

Definition at line 270 of file Histograms.h.

Referenced by Fill(), and HResolution().

MonitorElement* HResolution::hPt
private

Definition at line 273 of file Histograms.h.

Referenced by Fill(), and HResolution().

TH1F* HResolution::resoHisto_
protected

Definition at line 1771 of file Histograms.h.

Referenced by HResolution(), Write(), and ~HResolution().

std::string HResolution::theName
private

Definition at line 266 of file Histograms.h.

Referenced by HResolution().

std::string HResolution::where
private

Definition at line 267 of file Histograms.h.

Referenced by HResolution().