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 Member Functions | Protected Attributes
HFunctionResolution Class Reference

#include <Histograms.h>

Inheritance diagram for HFunctionResolution:
Histograms HFunctionResolutionVarianceCheck

Public Member Functions

virtual void Clear ()
 
virtual void Fill (const reco::Particle::LorentzVector &p4, const double &resValue, const int charge)
 
 HFunctionResolution (TFile *outputFile, const TString &name, const double &ptMax=100, const int totBinsY=300)
 
virtual void Write ()
 
 ~HFunctionResolution ()
 
- Public Member Functions inherited from Histograms
void declareHistograms ()
 
virtual void Fill (const reco::Particle::LorentzVector &p1, const reco::Particle::LorentzVector &p2)
 
virtual void Fill (const reco::Particle::LorentzVector &p1, const reco::Particle::LorentzVector &p2, const int charge, const double &weight=1.)
 
virtual void Fill (const CLHEP::HepLorentzVector &momentum1, const CLHEP::HepLorentzVector &momentum2)
 
virtual void Fill (const CLHEP::HepLorentzVector &momentum1, const CLHEP::HepLorentzVector &momentum2, const int charge, const double &weight=1.)
 
virtual void Fill (const CLHEP::HepLorentzVector &p1, const reco::Particle::LorentzVector &p2)
 
virtual void Fill (const reco::Particle::LorentzVector &p4, const double &weight=1.)
 
virtual void Fill (const reco::Particle::LorentzVector &p4, const int charge, const double &weight=1.)
 
virtual void Fill (const CLHEP::HepLorentzVector &momentum, const int charge, const double &weight=1.)
 
virtual void Fill (const reco::Particle::LorentzVector &p4, const double &genValue, const double recValue, const int charge)
 
virtual void Fill (const CLHEP::HepLorentzVector &p, const double &likeValue)
 
virtual void Fill (const unsigned int number)
 
virtual void Fill (const reco::Particle::LorentzVector &recoP1, const int charge1, const reco::Particle::LorentzVector &genP1, const reco::Particle::LorentzVector &recoP2, const int charge2, const reco::Particle::LorentzVector &genP2, const double &recoMass, const double &genMass)
 
virtual void Fill (const reco::Particle::LorentzVector &recoP1, const int charge1, const reco::Particle::LorentzVector &recoP2, const int charge2, const double &recoMass, const double &genMass)
 
virtual void Fill (const reco::Particle::LorentzVector &recoP1, const reco::Particle::LorentzVector &genP1, const reco::Particle::LorentzVector &recoP2, const reco::Particle::LorentzVector &genP2)
 
virtual void Fill (const double &x, const double &y)
 
virtual void Fill (const double &x, const double &y, const double &a, const double &b)
 
virtual void Fill (const reco::Particle::LorentzVector &p41, const reco::Particle::LorentzVector &p42, const reco::Particle::LorentzVector &p4Res, const double &weight=1.)
 
virtual void Fill (const CLHEP::HepLorentzVector &momentum1, const CLHEP::HepLorentzVector &momentum2, const CLHEP::HepLorentzVector &momentumRes, const double &weight=1.)
 
void fillEventInfo (int proc, int strk, int ntrkr)
 
void fillRecHistograms (const RecTrack_t &r)
 
void fillSimHistograms (const SimTrack_t &s)
 
void fillVzeroHistograms (const RecVzero_t &r, int part)
 
virtual double Get (const reco::Particle::LorentzVector &recoP1, const TString &covarianceName)
 
virtual TString GetName ()
 
 Histograms ()
 
 Histograms (const TString &name)
 
 Histograms (TFile *outputFile, const TString &name)
 
 Histograms (const edm::ParameterSet &pset)
 
virtual void SetWeight (double weight)
 
void writeHistograms ()
 
virtual ~Histograms ()
 
 ~Histograms ()
 

Protected Member Functions

int getXindex (const double &x) const
 
int getYindex (const double &y) const
 

Protected Attributes

double deltaX_
 
double deltaY_
 
TH1F * hReso_
 
TProfile * hResoVSEta_prof_
 
TProfile * hResoVSPhi_prof_
 
TProfile * hResoVSPhiMinus_prof_
 
TProfile * hResoVSPhiPlus_prof_
 
TProfile * hResoVSPt_Bar_prof_
 
TProfile * hResoVSPt_Endc_17_prof_
 
TProfile * hResoVSPt_Endc_20_prof_
 
TProfile * hResoVSPt_Endc_24_prof_
 
TProfile * hResoVSPt_Ovlap_prof_
 
TProfile * hResoVSPt_prof_
 
TH2F * hResoVSPtEta_
 
int ** resoCount_
 
double ** resoVsPtEta_
 
int totBinsX_
 
int totBinsY_
 
double xMin_
 
double yMin_
 
- Protected Attributes inherited from Histograms
TDirectory * histoDir_
 
TString name_
 
TFile * outputFile_
 
double theWeight_
 

Detailed Description

This histogram class fills a TProfile with the resolution evaluated from the resolution functions for single muon quantities. The resolution functions are used by MuScleFit to evaluate the mass resolution, which is the value seen by minuit and through it, corrections are evaluated.
In the end we will compare the histograms filled by this class (from the resolution function, reflecting the parameters changes done by minuit) with those filled comparing recoMuons with genMuons (the real resolutions).

Definition at line 1472 of file Histograms.h.

Constructor & Destructor Documentation

HFunctionResolution::HFunctionResolution ( TFile *  outputFile,
const TString &  name,
const double &  ptMax = 100,
const int  totBinsY = 300 
)
inline

Definition at line 1475 of file Histograms.h.

References deltaX_, deltaY_, hReso_, hResoVSEta_prof_, hResoVSPhi_prof_, hResoVSPhiMinus_prof_, hResoVSPhiPlus_prof_, hResoVSPt_Bar_prof_, hResoVSPt_Endc_17_prof_, hResoVSPt_Endc_20_prof_, hResoVSPt_Endc_24_prof_, hResoVSPt_Ovlap_prof_, hResoVSPt_prof_, hResoVSPtEta_, i, j, mergeVDriftHistosByStation::name, Histograms::name_, jptDQMConfig_cff::ptMax, resoCount_, resoVsPtEta_, totBinsX_, totBinsY_, xMin_, and yMin_.

1475  : Histograms(outputFile, name) {
1476  name_ = name;
1477  totBinsX_ = 300;
1478  totBinsY_ = totBinsY;
1479  xMin_ = 0.;
1480  yMin_ = -3.0;
1481  double xMax = ptMax;
1482  double yMax = 3.0;
1483  deltaX_ = xMax - xMin_;
1484  deltaY_ = yMax - yMin_;
1485  hReso_ = new TH1F( name+"_Reso", "resolution", 1000, 0, 1 );
1486  hResoVSPtEta_ = new TH2F( name+"_ResoVSPtEta", "resolution vs pt and #eta", totBinsX_, xMin_, xMax, totBinsY_, yMin_, yMax );
1487  // Create and initialize the resolution arrays
1488  resoVsPtEta_ = new double*[totBinsX_];
1489  resoCount_ = new int*[totBinsX_];
1490  for( int i=0; i<totBinsX_; ++i ) {
1491  resoVsPtEta_[i] = new double[totBinsY_];
1492  resoCount_[i] = new int[totBinsY_];
1493  for( int j=0; j<totBinsY_; ++j ) {
1494  resoVsPtEta_[i][j] = 0;
1495  resoCount_[i][j] = 0;
1496  }
1497  }
1498  hResoVSPt_prof_ = new TProfile( name+"_ResoVSPt_prof", "resolution VS pt", totBinsX_, xMin_, xMax, yMin_, yMax);
1499  hResoVSPt_Bar_prof_ = new TProfile( name+"_ResoVSPt_Bar_prof", "resolution VS pt Barrel", totBinsX_, xMin_, xMax, yMin_, yMax);
1500  hResoVSPt_Endc_17_prof_ = new TProfile( name+"_ResoVSPt_Endc_1.7_prof", "resolution VS pt Endcap (1.4<eta<1.7)", totBinsX_, xMin_, xMax, yMin_, yMax);
1501  hResoVSPt_Endc_20_prof_ = new TProfile( name+"_ResoVSPt_Endc_2.0_prof", "resolution VS pt Endcap (1.7<eta<2.0)", totBinsX_, xMin_, xMax, yMin_, yMax);
1502  hResoVSPt_Endc_24_prof_ = new TProfile( name+"_ResoVSPt_Endc_2.4_prof", "resolution VS pt Endcap (2.0<eta<2.4)", totBinsX_, xMin_, xMax, yMin_, yMax);
1503  hResoVSPt_Ovlap_prof_ = new TProfile( name+"_ResoVSPt_Ovlap_prof", "resolution VS pt Overlap", totBinsX_, xMin_, xMax, yMin_, yMax);
1504  hResoVSEta_prof_ = new TProfile( name+"_ResoVSEta_prof", "resolution VS eta", totBinsY_, yMin_, yMax, 0, 1);
1505  //hResoVSTheta_prof_ = new TProfile( name+"_ResoVSTheta_prof", "resolution VS theta", 30, 0, TMath::Pi(), 0, 1);
1506  hResoVSPhiPlus_prof_ = new TProfile( name+"_ResoVSPhiPlus_prof", "resolution VS phi mu+", 16, -3.2, 3.2, 0, 1);
1507  hResoVSPhiMinus_prof_ = new TProfile( name+"_ResoVSPhiMinus_prof", "resolution VS phi mu-", 16, -3.2, 3.2, 0, 1);
1508  hResoVSPhi_prof_ = new TProfile( name+"_ResoVSPhi_prof", "resolution VS phi", 16, -3.2, 3.2, -1, 1);
1509  }
int i
Definition: DBlmapReader.cc:9
TString name_
Definition: Histograms.h:123
TProfile * hResoVSPt_Endc_20_prof_
Definition: Histograms.h:1643
TProfile * hResoVSPt_Endc_17_prof_
Definition: Histograms.h:1642
TProfile * hResoVSPt_Bar_prof_
Definition: Histograms.h:1641
TProfile * hResoVSPt_prof_
Definition: Histograms.h:1640
int j
Definition: DBlmapReader.cc:9
TProfile * hResoVSPt_Ovlap_prof_
Definition: Histograms.h:1645
TProfile * hResoVSPt_Endc_24_prof_
Definition: Histograms.h:1644
TProfile * hResoVSPhi_prof_
Definition: Histograms.h:1650
TProfile * hResoVSPhiMinus_prof_
Definition: Histograms.h:1648
TProfile * hResoVSEta_prof_
Definition: Histograms.h:1646
TProfile * hResoVSPhiPlus_prof_
Definition: Histograms.h:1649
double ** resoVsPtEta_
Definition: Histograms.h:1638
HFunctionResolution::~HFunctionResolution ( )
inline

Definition at line 1510 of file Histograms.h.

References hReso_, hResoVSEta_prof_, hResoVSPhi_prof_, hResoVSPhiMinus_prof_, hResoVSPhiPlus_prof_, hResoVSPt_Bar_prof_, hResoVSPt_Endc_17_prof_, hResoVSPt_Endc_20_prof_, hResoVSPt_Endc_24_prof_, hResoVSPt_Ovlap_prof_, hResoVSPt_prof_, hResoVSPtEta_, i, resoCount_, resoVsPtEta_, and totBinsX_.

1510  {
1511  delete hReso_;
1512  delete hResoVSPtEta_;
1513  // Free the resolution arrays
1514  for( int i=0; i<totBinsX_; ++i ) {
1515  delete[] resoVsPtEta_[i];
1516  delete[] resoCount_[i];
1517  }
1518  delete[] resoVsPtEta_;
1519  delete[] resoCount_;
1520  // Free the profile histograms
1521  delete hResoVSPt_prof_;
1522  delete hResoVSPt_Bar_prof_;
1523  delete hResoVSPt_Endc_17_prof_;
1524  delete hResoVSPt_Endc_20_prof_;
1525  delete hResoVSPt_Endc_24_prof_;
1526  delete hResoVSPt_Ovlap_prof_;
1527  delete hResoVSEta_prof_;
1528  //delete hResoVSTheta_prof_;
1529  delete hResoVSPhiPlus_prof_;
1530  delete hResoVSPhiMinus_prof_;
1531  delete hResoVSPhi_prof_;
1532  }
int i
Definition: DBlmapReader.cc:9
TProfile * hResoVSPt_Endc_20_prof_
Definition: Histograms.h:1643
TProfile * hResoVSPt_Endc_17_prof_
Definition: Histograms.h:1642
TProfile * hResoVSPt_Bar_prof_
Definition: Histograms.h:1641
TProfile * hResoVSPt_prof_
Definition: Histograms.h:1640
TProfile * hResoVSPt_Ovlap_prof_
Definition: Histograms.h:1645
TProfile * hResoVSPt_Endc_24_prof_
Definition: Histograms.h:1644
TProfile * hResoVSPhi_prof_
Definition: Histograms.h:1650
TProfile * hResoVSPhiMinus_prof_
Definition: Histograms.h:1648
TProfile * hResoVSEta_prof_
Definition: Histograms.h:1646
TProfile * hResoVSPhiPlus_prof_
Definition: Histograms.h:1649
double ** resoVsPtEta_
Definition: Histograms.h:1638

Member Function Documentation

virtual void HFunctionResolution::Clear ( )
inlinevirtual

Implements Histograms.

Definition at line 1613 of file Histograms.h.

References hReso_, hResoVSEta_prof_, hResoVSPhi_prof_, hResoVSPhiMinus_prof_, hResoVSPhiPlus_prof_, hResoVSPt_Bar_prof_, hResoVSPt_Endc_17_prof_, hResoVSPt_Endc_20_prof_, hResoVSPt_Endc_24_prof_, hResoVSPt_Ovlap_prof_, hResoVSPt_prof_, and hResoVSPtEta_.

1613  {
1614  hReso_->Clear();
1615  hResoVSPtEta_->Clear();
1616  hResoVSPt_prof_->Clear();
1617  hResoVSPt_Bar_prof_->Clear();
1618  hResoVSPt_Endc_17_prof_->Clear();
1619  hResoVSPt_Endc_20_prof_->Clear();
1620  hResoVSPt_Endc_24_prof_->Clear();
1621  hResoVSPt_Ovlap_prof_->Clear();
1622  hResoVSEta_prof_->Clear();
1623  //hResoVSTheta_prof_->Clear();
1624  hResoVSPhiPlus_prof_->Clear();
1625  hResoVSPhiMinus_prof_->Clear();
1626  hResoVSPhi_prof_->Clear();
1627  }
TProfile * hResoVSPt_Endc_20_prof_
Definition: Histograms.h:1643
TProfile * hResoVSPt_Endc_17_prof_
Definition: Histograms.h:1642
TProfile * hResoVSPt_Bar_prof_
Definition: Histograms.h:1641
TProfile * hResoVSPt_prof_
Definition: Histograms.h:1640
TProfile * hResoVSPt_Ovlap_prof_
Definition: Histograms.h:1645
TProfile * hResoVSPt_Endc_24_prof_
Definition: Histograms.h:1644
TProfile * hResoVSPhi_prof_
Definition: Histograms.h:1650
TProfile * hResoVSPhiMinus_prof_
Definition: Histograms.h:1648
TProfile * hResoVSEta_prof_
Definition: Histograms.h:1646
TProfile * hResoVSPhiPlus_prof_
Definition: Histograms.h:1649
virtual void HFunctionResolution::Fill ( const reco::Particle::LorentzVector p4,
const double &  resValue,
const int  charge 
)
inlinevirtual

Reimplemented from Histograms.

Reimplemented in HFunctionResolutionVarianceCheck.

Definition at line 1533 of file Histograms.h.

References getXindex(), getYindex(), hReso_, hResoVSEta_prof_, hResoVSPhi_prof_, hResoVSPhiMinus_prof_, hResoVSPhiPlus_prof_, hResoVSPt_Bar_prof_, hResoVSPt_Endc_17_prof_, hResoVSPt_Endc_20_prof_, hResoVSPt_Endc_24_prof_, hResoVSPt_Ovlap_prof_, hResoVSPt_prof_, resoCount_, resoVsPtEta_, totBinsX_, and totBinsY_.

Referenced by HFunctionResolutionVarianceCheck::Fill().

1533  {
1534  if( resValue != resValue ) return;
1535  hReso_->Fill(resValue);
1536 
1537  // Fill the arrays with the resolution value and count
1538  int xIndex = getXindex(p4.Pt());
1539  int yIndex = getYindex(p4.Eta());
1540  if ( 0 <= xIndex && xIndex < totBinsX_ && 0 <= yIndex && yIndex < totBinsY_ ) {
1541  resoVsPtEta_[xIndex][yIndex] += resValue;
1542  // ATTENTION: we count only for positive muons because we are considering di-muon resonances
1543  // and we use this counter to compute the mean in the end. The resoVsPtEta value is filled with each muon,
1544  // but they must be considered independently (analogous to a factor 2) so in the end we would have
1545  // to divide by N/2, that is why we only increase the counter for half the muons.
1546  // if( charge > 0 )
1547  // No more. Changing it here influences also other uses of this class. The macro FunctionTerms.cc
1548  // multiplies the terms by the 2 factor.
1549 
1550  resoCount_[xIndex][yIndex] += 1;
1551 
1552  // hResoVSPtEta->Fill(p4.Pt(), p4.Eta(), resValue);
1553  hResoVSPt_prof_->Fill(p4.Pt(),resValue);
1554  if(fabs(p4.Eta())<=0.9)
1555  hResoVSPt_Bar_prof_->Fill(p4.Pt(),resValue);
1556  else if(fabs(p4.Eta())>0.9 && fabs(p4.Eta())<=1.4 )
1557  hResoVSPt_Ovlap_prof_->Fill(p4.Pt(),resValue);
1558  else if(fabs(p4.Eta())>1.4 && fabs(p4.Eta())<=1.7 )
1559  hResoVSPt_Endc_17_prof_->Fill(p4.Pt(),resValue);
1560  else if(fabs(p4.Eta())>1.7 && fabs(p4.Eta())<=2.0 )
1561  hResoVSPt_Endc_20_prof_->Fill(p4.Pt(),resValue);
1562  else
1563  hResoVSPt_Endc_24_prof_->Fill(p4.Pt(),resValue);
1564  hResoVSEta_prof_->Fill(p4.Eta(),resValue);
1565  //hResoVSTheta_prof_->Fill(p4.Theta(),resValue);
1566  if(charge>0)
1567  hResoVSPhiPlus_prof_->Fill(p4.Phi(),resValue);
1568  else if(charge<0)
1569  hResoVSPhiMinus_prof_->Fill(p4.Phi(),resValue);
1570  hResoVSPhi_prof_->Fill(p4.Phi(),resValue);
1571  }
1572  }
double charge(const std::vector< uint8_t > &Ampls)
int getYindex(const double &y) const
Definition: Histograms.h:1633
TProfile * hResoVSPt_Endc_20_prof_
Definition: Histograms.h:1643
TProfile * hResoVSPt_Endc_17_prof_
Definition: Histograms.h:1642
TProfile * hResoVSPt_Bar_prof_
Definition: Histograms.h:1641
double p4[4]
Definition: TauolaWrapper.h:92
TProfile * hResoVSPt_prof_
Definition: Histograms.h:1640
TProfile * hResoVSPt_Ovlap_prof_
Definition: Histograms.h:1645
TProfile * hResoVSPt_Endc_24_prof_
Definition: Histograms.h:1644
TProfile * hResoVSPhi_prof_
Definition: Histograms.h:1650
TProfile * hResoVSPhiMinus_prof_
Definition: Histograms.h:1648
int getXindex(const double &x) const
Definition: Histograms.h:1630
TProfile * hResoVSEta_prof_
Definition: Histograms.h:1646
TProfile * hResoVSPhiPlus_prof_
Definition: Histograms.h:1649
double ** resoVsPtEta_
Definition: Histograms.h:1638
int HFunctionResolution::getXindex ( const double &  x) const
inlineprotected

Definition at line 1630 of file Histograms.h.

References deltaX_, totBinsX_, and xMin_.

Referenced by Fill(), and HFunctionResolutionVarianceCheck::Fill().

1630  {
1631  return int((x-xMin_)/deltaX_*totBinsX_);
1632  }
Definition: DDAxes.h:10
int HFunctionResolution::getYindex ( const double &  y) const
inlineprotected

Definition at line 1633 of file Histograms.h.

References deltaY_, totBinsY_, and yMin_.

Referenced by Fill(), and HFunctionResolutionVarianceCheck::Fill().

1633  {
1634  return int((y-yMin_)/deltaY_*totBinsY_);
1635  }
virtual void HFunctionResolution::Write ( )
inlinevirtual

Implements Histograms.

Reimplemented in HFunctionResolutionVarianceCheck.

Definition at line 1574 of file Histograms.h.

References svgfig::canvas(), Histograms::histoDir_, hReso_, hResoVSEta_prof_, hResoVSPhi_prof_, hResoVSPhiMinus_prof_, hResoVSPhiPlus_prof_, hResoVSPt_Bar_prof_, hResoVSPt_Endc_17_prof_, hResoVSPt_Endc_20_prof_, hResoVSPt_Endc_24_prof_, hResoVSPt_Ovlap_prof_, hResoVSPt_prof_, hResoVSPtEta_, i, j, N, Histograms::outputFile_, resoCount_, resoVsPtEta_, totBinsX_, and totBinsY_.

Referenced by HFunctionResolutionVarianceCheck::Write().

1574  {
1575  if(histoDir_ != 0) histoDir_->cd();
1576 
1577  hReso_->Write();
1578 
1579  for( int i=0; i<totBinsX_; ++i ) {
1580  for( int j=0; j<totBinsY_; ++j ) {
1581  int N = resoCount_[i][j];
1582  // Fill with the mean value
1583  if( N != 0 ) hResoVSPtEta_->SetBinContent( i+1, j+1, resoVsPtEta_[i][j]/N );
1584  else hResoVSPtEta_->SetBinContent( i+1, j+1, 0 );
1585  }
1586  }
1587  hResoVSPtEta_->Write();
1588 
1589  hResoVSPt_prof_->Write();
1590  hResoVSPt_Bar_prof_->Write();
1591  hResoVSPt_Endc_17_prof_->Write();
1592  hResoVSPt_Endc_20_prof_->Write();
1593  hResoVSPt_Endc_24_prof_->Write();
1594  hResoVSPt_Ovlap_prof_->Write();
1595  hResoVSEta_prof_->Write();
1596  //hResoVSTheta_prof_->Write();
1597  hResoVSPhiMinus_prof_->Write();
1598  hResoVSPhiPlus_prof_->Write();
1599  hResoVSPhi_prof_->Write();
1600 
1601  TCanvas canvas(TString(hResoVSPtEta_->GetName())+"_canvas", TString(hResoVSPtEta_->GetTitle())+" canvas", 1000, 800);
1602  canvas.Divide(2);
1603  canvas.cd(1);
1604  hResoVSPtEta_->Draw("lego");
1605  canvas.cd(2);
1606  hResoVSPtEta_->Draw("surf5");
1607  canvas.Write();
1608  hResoVSPtEta_->Write();
1609 
1610  outputFile_->cd();
1611  }
int i
Definition: DBlmapReader.cc:9
TFile * outputFile_
Definition: Histograms.h:124
def canvas
Definition: svgfig.py:481
TProfile * hResoVSPt_Endc_20_prof_
Definition: Histograms.h:1643
TProfile * hResoVSPt_Endc_17_prof_
Definition: Histograms.h:1642
TProfile * hResoVSPt_Bar_prof_
Definition: Histograms.h:1641
TProfile * hResoVSPt_prof_
Definition: Histograms.h:1640
TDirectory * histoDir_
Definition: Histograms.h:125
int j
Definition: DBlmapReader.cc:9
TProfile * hResoVSPt_Ovlap_prof_
Definition: Histograms.h:1645
TProfile * hResoVSPt_Endc_24_prof_
Definition: Histograms.h:1644
TProfile * hResoVSPhi_prof_
Definition: Histograms.h:1650
#define N
Definition: blowfish.cc:9
TProfile * hResoVSPhiMinus_prof_
Definition: Histograms.h:1648
TProfile * hResoVSEta_prof_
Definition: Histograms.h:1646
TProfile * hResoVSPhiPlus_prof_
Definition: Histograms.h:1649
double ** resoVsPtEta_
Definition: Histograms.h:1638

Member Data Documentation

double HFunctionResolution::deltaX_
protected

Definition at line 1653 of file Histograms.h.

Referenced by getXindex(), and HFunctionResolution().

double HFunctionResolution::deltaY_
protected

Definition at line 1653 of file Histograms.h.

Referenced by getYindex(), and HFunctionResolution().

TH1F* HFunctionResolution::hReso_
protected

Definition at line 1636 of file Histograms.h.

Referenced by Clear(), Fill(), HFunctionResolution(), Write(), and ~HFunctionResolution().

TProfile* HFunctionResolution::hResoVSEta_prof_
protected

Definition at line 1646 of file Histograms.h.

Referenced by Clear(), Fill(), HFunctionResolution(), Write(), and ~HFunctionResolution().

TProfile* HFunctionResolution::hResoVSPhi_prof_
protected

Definition at line 1650 of file Histograms.h.

Referenced by Clear(), Fill(), HFunctionResolution(), Write(), and ~HFunctionResolution().

TProfile* HFunctionResolution::hResoVSPhiMinus_prof_
protected

Definition at line 1648 of file Histograms.h.

Referenced by Clear(), Fill(), HFunctionResolution(), Write(), and ~HFunctionResolution().

TProfile* HFunctionResolution::hResoVSPhiPlus_prof_
protected

Definition at line 1649 of file Histograms.h.

Referenced by Clear(), Fill(), HFunctionResolution(), Write(), and ~HFunctionResolution().

TProfile* HFunctionResolution::hResoVSPt_Bar_prof_
protected

Definition at line 1641 of file Histograms.h.

Referenced by Clear(), Fill(), HFunctionResolution(), Write(), and ~HFunctionResolution().

TProfile* HFunctionResolution::hResoVSPt_Endc_17_prof_
protected

Definition at line 1642 of file Histograms.h.

Referenced by Clear(), Fill(), HFunctionResolution(), Write(), and ~HFunctionResolution().

TProfile* HFunctionResolution::hResoVSPt_Endc_20_prof_
protected

Definition at line 1643 of file Histograms.h.

Referenced by Clear(), Fill(), HFunctionResolution(), Write(), and ~HFunctionResolution().

TProfile* HFunctionResolution::hResoVSPt_Endc_24_prof_
protected

Definition at line 1644 of file Histograms.h.

Referenced by Clear(), Fill(), HFunctionResolution(), Write(), and ~HFunctionResolution().

TProfile* HFunctionResolution::hResoVSPt_Ovlap_prof_
protected

Definition at line 1645 of file Histograms.h.

Referenced by Clear(), Fill(), HFunctionResolution(), Write(), and ~HFunctionResolution().

TProfile* HFunctionResolution::hResoVSPt_prof_
protected

Definition at line 1640 of file Histograms.h.

Referenced by Clear(), Fill(), HFunctionResolution(), Write(), and ~HFunctionResolution().

TH2F* HFunctionResolution::hResoVSPtEta_
protected

Definition at line 1637 of file Histograms.h.

Referenced by Clear(), HFunctionResolution(), Write(), and ~HFunctionResolution().

int** HFunctionResolution::resoCount_
protected

Definition at line 1639 of file Histograms.h.

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

double** HFunctionResolution::resoVsPtEta_
protected

Definition at line 1638 of file Histograms.h.

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

int HFunctionResolution::totBinsX_
protected
int HFunctionResolution::totBinsY_
protected
double HFunctionResolution::xMin_
protected

Definition at line 1652 of file Histograms.h.

Referenced by getXindex(), and HFunctionResolution().

double HFunctionResolution::yMin_
protected

Definition at line 1652 of file Histograms.h.

Referenced by getYindex(), and HFunctionResolution().