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 1716 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 1718 of file Histograms.h.

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

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

Definition at line 1735 of file Histograms.h.

References diffHisto_, histo2D_, and resoHisto_.

1735  {
1736  delete diffHisto_;
1737  delete histo2D_;
1738  delete resoHisto_;
1739  }
TH1F * resoHisto_
Definition: Histograms.h:1773
TH2F * histo2D_
Definition: Histograms.h:1772
TProfile * diffHisto_
Definition: Histograms.h:1771
HResolution::HResolution ( std::string  dirName_,
std::string  name,
std::string  whereIs 
)
inline

Definition at line 167 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(), DQMStore::setCurrentFolder(), AlCaHLTBitMon_QueryRunRegistry::string, theName, and where.

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

Definition at line 206 of file Histograms.h.

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

Definition at line 210 of file Histograms.h.

210 {}

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 213 of file Histograms.h.

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

214  {
215 
216  Fill(rp, rpt, reta, rphi, rcharge);
217 
218 
219  h2Eta->Fill(eta,reta);
220  h2Phi->Fill(phi,rphi);
221 
222  h2P->Fill(p,rp);
223  h2Pt->Fill(pt,rpt);
224 
225  h2PtVsEta->Fill(eta,rpt);
226  h2PtVsPhi->Fill(phi,rpt);
227 
228  h2EtaVsPt ->Fill(pt,reta);
229  h2EtaVsPhi->Fill(phi,reta);
230 
231  h2PhiVsPt ->Fill(pt,rphi);
232  h2PhiVsEta->Fill(eta,rphi);
233  }
MonitorElement * h2PhiVsPt
Definition: Histograms.h:291
MonitorElement * h2Eta
Definition: Histograms.h:279
MonitorElement * h2EtaVsPhi
Definition: Histograms.h:289
MonitorElement * h2P
Definition: Histograms.h:282
T eta() const
MonitorElement * h2PhiVsEta
Definition: Histograms.h:292
void Fill(long long x)
MonitorElement * h2Phi
Definition: Histograms.h:280
MonitorElement * h2Pt
Definition: Histograms.h:283
MonitorElement * h2EtaVsPt
Definition: Histograms.h:288
MonitorElement * h2PtVsPhi
Definition: Histograms.h:286
virtual Int_t Fill(Double_t x, Double_t y)
Definition: Histograms.h:1740
MonitorElement * h2PtVsEta
Definition: Histograms.h:285
Definition: DDAxes.h:10
void HResolution::Fill ( double  p,
double  pt,
double  eta,
double  phi,
double  rp,
double  rpt 
)
inline

Definition at line 236 of file Histograms.h.

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

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

Definition at line 251 of file Histograms.h.

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

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

Definition at line 1740 of file Histograms.h.

References diffHisto_, and histo2D_.

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

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

Definition at line 1745 of file Histograms.h.

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

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

Member Data Documentation

DQMStore* HResolution::dbe_
private

Definition at line 266 of file Histograms.h.

Referenced by HResolution().

TDirectory* HResolution::diffDir_
protected

Definition at line 1770 of file Histograms.h.

Referenced by HResolution(), and Write().

TProfile* HResolution::diffHisto_
protected

Definition at line 1771 of file Histograms.h.

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

TDirectory* HResolution::dir2D_
protected

Definition at line 1769 of file Histograms.h.

Referenced by HResolution(), and Write().

TDirectory* HResolution::dir_
protected

Definition at line 1768 of file Histograms.h.

Referenced by HResolution(), and Write().

MonitorElement* HResolution::h2Eta
private

Definition at line 279 of file Histograms.h.

Referenced by Fill(), and HResolution().

MonitorElement* HResolution::h2EtaVsPhi
private

Definition at line 289 of file Histograms.h.

Referenced by Fill(), and HResolution().

MonitorElement* HResolution::h2EtaVsPt
private

Definition at line 288 of file Histograms.h.

Referenced by Fill(), and HResolution().

MonitorElement* HResolution::h2P
private

Definition at line 282 of file Histograms.h.

Referenced by Fill(), and HResolution().

MonitorElement* HResolution::h2Phi
private

Definition at line 280 of file Histograms.h.

Referenced by Fill(), and HResolution().

MonitorElement* HResolution::h2PhiVsEta
private

Definition at line 292 of file Histograms.h.

Referenced by Fill(), and HResolution().

MonitorElement* HResolution::h2PhiVsPt
private

Definition at line 291 of file Histograms.h.

Referenced by Fill(), and HResolution().

MonitorElement* HResolution::h2Pt
private

Definition at line 283 of file Histograms.h.

Referenced by Fill(), and HResolution().

MonitorElement* HResolution::h2PtVsEta
private

Definition at line 285 of file Histograms.h.

Referenced by Fill(), and HResolution().

MonitorElement* HResolution::h2PtVsPhi
private

Definition at line 286 of file Histograms.h.

Referenced by Fill(), and HResolution().

MonitorElement* HResolution::hCharge
private

Definition at line 277 of file Histograms.h.

Referenced by Fill(), and HResolution().

MonitorElement* HResolution::hEta
private

Definition at line 271 of file Histograms.h.

Referenced by Fill(), and HResolution().

TH2F* HResolution::histo2D_
protected

Definition at line 1772 of file Histograms.h.

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

MonitorElement* HResolution::hP
private

Definition at line 274 of file Histograms.h.

Referenced by Fill(), and HResolution().

MonitorElement* HResolution::hPhi
private

Definition at line 272 of file Histograms.h.

Referenced by Fill(), and HResolution().

MonitorElement* HResolution::hPt
private

Definition at line 275 of file Histograms.h.

Referenced by Fill(), and HResolution().

TH1F* HResolution::resoHisto_
protected

Definition at line 1773 of file Histograms.h.

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

std::string HResolution::theName
private

Definition at line 268 of file Histograms.h.

Referenced by HResolution().

std::string HResolution::where
private

Definition at line 269 of file Histograms.h.

Referenced by HResolution().