CMS 3D CMS Logo

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

#include <Histograms.h>

Inheritance diagram for HFunctionResolution:
HFunctionResolutionVarianceCheck

Public Member Functions

void Clear () override
 
void Fill (const reco::Particle::LorentzVector &p4, const double &resValue, const int charge) override
 
 HFunctionResolution (TFile *outputFile, const TString &name, const double &ptMax=100, const int totBinsY=300)
 
void Write () override
 
 ~HFunctionResolution () override
 

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_
 

Additional Inherited Members

- Public Attributes inherited from Histograms
dqm::reco::MonitorElementlumiVsLS
 
dqm::reco::MonitorElementnumberOfPixelClustersVsLS
 
dqm::reco::MonitorElementnumberOfPixelClustersVsLumi
 
dqm::reco::MonitorElementpixelLumiVsLS
 
dqm::reco::MonitorElementpixelLumiVsLumi
 
dqm::reco::MonitorElementpuVsLS
 

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 1550 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 1552 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_, mps_fire::i, dqmiolumiharvest::j, mergeVDriftHistosByStation::name, HLT_FULL_cff::ptMax, resoCount_, resoVsPtEta_, totBinsX_, totBinsY_, xMin_, HLT_FULL_cff::yMax, and yMin_.

1553  : Histograms(outputFile, name) {
1554  name_ = name;
1555  totBinsX_ = 300;
1556  totBinsY_ = totBinsY;
1557  xMin_ = 0.;
1558  yMin_ = -3.0;
1559  double xMax = ptMax;
1560  double yMax = 3.0;
1561  deltaX_ = xMax - xMin_;
1562  deltaY_ = yMax - yMin_;
1563  hReso_ = new TH1F(name + "_Reso", "resolution", 1000, 0, 1);
1564  hResoVSPtEta_ =
1565  new TH2F(name + "_ResoVSPtEta", "resolution vs pt and #eta", totBinsX_, xMin_, xMax, totBinsY_, yMin_, yMax);
1566  // Create and initialize the resolution arrays
1567  resoVsPtEta_ = new double*[totBinsX_];
1568  resoCount_ = new int*[totBinsX_];
1569  for (int i = 0; i < totBinsX_; ++i) {
1570  resoVsPtEta_[i] = new double[totBinsY_];
1571  resoCount_[i] = new int[totBinsY_];
1572  for (int j = 0; j < totBinsY_; ++j) {
1573  resoVsPtEta_[i][j] = 0;
1574  resoCount_[i][j] = 0;
1575  }
1576  }
1577  hResoVSPt_prof_ = new TProfile(name + "_ResoVSPt_prof", "resolution VS pt", totBinsX_, xMin_, xMax, yMin_, yMax);
1579  new TProfile(name + "_ResoVSPt_Bar_prof", "resolution VS pt Barrel", totBinsX_, xMin_, xMax, yMin_, yMax);
1580  hResoVSPt_Endc_17_prof_ = new TProfile(
1581  name + "_ResoVSPt_Endc_1.7_prof", "resolution VS pt Endcap (1.4<eta<1.7)", totBinsX_, xMin_, xMax, yMin_, yMax);
1582  hResoVSPt_Endc_20_prof_ = new TProfile(
1583  name + "_ResoVSPt_Endc_2.0_prof", "resolution VS pt Endcap (1.7<eta<2.0)", totBinsX_, xMin_, xMax, yMin_, yMax);
1584  hResoVSPt_Endc_24_prof_ = new TProfile(
1585  name + "_ResoVSPt_Endc_2.4_prof", "resolution VS pt Endcap (2.0<eta<2.4)", totBinsX_, xMin_, xMax, yMin_, yMax);
1587  new TProfile(name + "_ResoVSPt_Ovlap_prof", "resolution VS pt Overlap", totBinsX_, xMin_, xMax, yMin_, yMax);
1588  hResoVSEta_prof_ = new TProfile(name + "_ResoVSEta_prof", "resolution VS eta", totBinsY_, yMin_, yMax, 0, 1);
1589  //hResoVSTheta_prof_ = new TProfile( name+"_ResoVSTheta_prof", "resolution VS theta", 30, 0, TMath::Pi(), 0, 1);
1590  hResoVSPhiPlus_prof_ = new TProfile(name + "_ResoVSPhiPlus_prof", "resolution VS phi mu+", 16, -3.2, 3.2, 0, 1);
1591  hResoVSPhiMinus_prof_ = new TProfile(name + "_ResoVSPhiMinus_prof", "resolution VS phi mu-", 16, -3.2, 3.2, 0, 1);
1592  hResoVSPhi_prof_ = new TProfile(name + "_ResoVSPhi_prof", "resolution VS phi", 16, -3.2, 3.2, -1, 1);
1593  }
TProfile * hResoVSPt_Endc_20_prof_
Definition: Histograms.h:1728
TProfile * hResoVSPt_Endc_17_prof_
Definition: Histograms.h:1727
TProfile * hResoVSPt_Bar_prof_
Definition: Histograms.h:1726
TProfile * hResoVSPt_prof_
Definition: Histograms.h:1725
TProfile * hResoVSPt_Ovlap_prof_
Definition: Histograms.h:1730
TProfile * hResoVSPt_Endc_24_prof_
Definition: Histograms.h:1729
TProfile * hResoVSPhi_prof_
Definition: Histograms.h:1735
TProfile * hResoVSPhiMinus_prof_
Definition: Histograms.h:1733
TProfile * hResoVSEta_prof_
Definition: Histograms.h:1731
TProfile * hResoVSPhiPlus_prof_
Definition: Histograms.h:1734
double ** resoVsPtEta_
Definition: Histograms.h:1723
HFunctionResolution::~HFunctionResolution ( )
inlineoverride

Definition at line 1594 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_, mps_fire::i, resoCount_, resoVsPtEta_, and totBinsX_.

1594  {
1595  delete hReso_;
1596  delete hResoVSPtEta_;
1597  // Free the resolution arrays
1598  for (int i = 0; i < totBinsX_; ++i) {
1599  delete[] resoVsPtEta_[i];
1600  delete[] resoCount_[i];
1601  }
1602  delete[] resoVsPtEta_;
1603  delete[] resoCount_;
1604  // Free the profile histograms
1605  delete hResoVSPt_prof_;
1606  delete hResoVSPt_Bar_prof_;
1607  delete hResoVSPt_Endc_17_prof_;
1608  delete hResoVSPt_Endc_20_prof_;
1609  delete hResoVSPt_Endc_24_prof_;
1610  delete hResoVSPt_Ovlap_prof_;
1611  delete hResoVSEta_prof_;
1612  //delete hResoVSTheta_prof_;
1613  delete hResoVSPhiPlus_prof_;
1614  delete hResoVSPhiMinus_prof_;
1615  delete hResoVSPhi_prof_;
1616  }
TProfile * hResoVSPt_Endc_20_prof_
Definition: Histograms.h:1728
TProfile * hResoVSPt_Endc_17_prof_
Definition: Histograms.h:1727
TProfile * hResoVSPt_Bar_prof_
Definition: Histograms.h:1726
TProfile * hResoVSPt_prof_
Definition: Histograms.h:1725
TProfile * hResoVSPt_Ovlap_prof_
Definition: Histograms.h:1730
TProfile * hResoVSPt_Endc_24_prof_
Definition: Histograms.h:1729
TProfile * hResoVSPhi_prof_
Definition: Histograms.h:1735
TProfile * hResoVSPhiMinus_prof_
Definition: Histograms.h:1733
TProfile * hResoVSEta_prof_
Definition: Histograms.h:1731
TProfile * hResoVSPhiPlus_prof_
Definition: Histograms.h:1734
double ** resoVsPtEta_
Definition: Histograms.h:1723

Member Function Documentation

void HFunctionResolution::Clear ( )
inlineoverride

Definition at line 1702 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_.

1702  {
1703  hReso_->Clear();
1704  hResoVSPtEta_->Clear();
1705  hResoVSPt_prof_->Clear();
1706  hResoVSPt_Bar_prof_->Clear();
1707  hResoVSPt_Endc_17_prof_->Clear();
1708  hResoVSPt_Endc_20_prof_->Clear();
1709  hResoVSPt_Endc_24_prof_->Clear();
1710  hResoVSPt_Ovlap_prof_->Clear();
1711  hResoVSEta_prof_->Clear();
1712  //hResoVSTheta_prof_->Clear();
1713  hResoVSPhiPlus_prof_->Clear();
1714  hResoVSPhiMinus_prof_->Clear();
1715  hResoVSPhi_prof_->Clear();
1716  }
TProfile * hResoVSPt_Endc_20_prof_
Definition: Histograms.h:1728
TProfile * hResoVSPt_Endc_17_prof_
Definition: Histograms.h:1727
TProfile * hResoVSPt_Bar_prof_
Definition: Histograms.h:1726
TProfile * hResoVSPt_prof_
Definition: Histograms.h:1725
TProfile * hResoVSPt_Ovlap_prof_
Definition: Histograms.h:1730
TProfile * hResoVSPt_Endc_24_prof_
Definition: Histograms.h:1729
TProfile * hResoVSPhi_prof_
Definition: Histograms.h:1735
TProfile * hResoVSPhiMinus_prof_
Definition: Histograms.h:1733
TProfile * hResoVSEta_prof_
Definition: Histograms.h:1731
TProfile * hResoVSPhiPlus_prof_
Definition: Histograms.h:1734
void HFunctionResolution::Fill ( const reco::Particle::LorentzVector p4,
const double &  resValue,
const int  charge 
)
inlineoverride

Definition at line 1617 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().

1617  {
1618  if (resValue != resValue)
1619  return;
1620  hReso_->Fill(resValue);
1621 
1622  // Fill the arrays with the resolution value and count
1623  int xIndex = getXindex(p4.Pt());
1624  int yIndex = getYindex(p4.Eta());
1625  if (0 <= xIndex && xIndex < totBinsX_ && 0 <= yIndex && yIndex < totBinsY_) {
1626  resoVsPtEta_[xIndex][yIndex] += resValue;
1627  // ATTENTION: we count only for positive muons because we are considering di-muon resonances
1628  // and we use this counter to compute the mean in the end. The resoVsPtEta value is filled with each muon,
1629  // but they must be considered independently (analogous to a factor 2) so in the end we would have
1630  // to divide by N/2, that is why we only increase the counter for half the muons.
1631  // if( charge > 0 )
1632  // No more. Changing it here influences also other uses of this class. The macro FunctionTerms.cc
1633  // multiplies the terms by the 2 factor.
1634 
1635  resoCount_[xIndex][yIndex] += 1;
1636 
1637  // hResoVSPtEta->Fill(p4.Pt(), p4.Eta(), resValue);
1638  hResoVSPt_prof_->Fill(p4.Pt(), resValue);
1639  if (fabs(p4.Eta()) <= 0.9)
1640  hResoVSPt_Bar_prof_->Fill(p4.Pt(), resValue);
1641  else if (fabs(p4.Eta()) > 0.9 && fabs(p4.Eta()) <= 1.4)
1642  hResoVSPt_Ovlap_prof_->Fill(p4.Pt(), resValue);
1643  else if (fabs(p4.Eta()) > 1.4 && fabs(p4.Eta()) <= 1.7)
1644  hResoVSPt_Endc_17_prof_->Fill(p4.Pt(), resValue);
1645  else if (fabs(p4.Eta()) > 1.7 && fabs(p4.Eta()) <= 2.0)
1646  hResoVSPt_Endc_20_prof_->Fill(p4.Pt(), resValue);
1647  else
1648  hResoVSPt_Endc_24_prof_->Fill(p4.Pt(), resValue);
1649  hResoVSEta_prof_->Fill(p4.Eta(), resValue);
1650  //hResoVSTheta_prof_->Fill(p4.Theta(),resValue);
1651  if (charge > 0)
1652  hResoVSPhiPlus_prof_->Fill(p4.Phi(), resValue);
1653  else if (charge < 0)
1654  hResoVSPhiMinus_prof_->Fill(p4.Phi(), resValue);
1655  hResoVSPhi_prof_->Fill(p4.Phi(), resValue);
1656  }
1657  }
int getYindex(const double &y) const
Definition: Histograms.h:1720
TProfile * hResoVSPt_Endc_20_prof_
Definition: Histograms.h:1728
TProfile * hResoVSPt_Endc_17_prof_
Definition: Histograms.h:1727
TProfile * hResoVSPt_Bar_prof_
Definition: Histograms.h:1726
TProfile * hResoVSPt_prof_
Definition: Histograms.h:1725
TProfile * hResoVSPt_Ovlap_prof_
Definition: Histograms.h:1730
TProfile * hResoVSPt_Endc_24_prof_
Definition: Histograms.h:1729
TProfile * hResoVSPhi_prof_
Definition: Histograms.h:1735
TProfile * hResoVSPhiMinus_prof_
Definition: Histograms.h:1733
int getXindex(const double &x) const
Definition: Histograms.h:1719
TProfile * hResoVSEta_prof_
Definition: Histograms.h:1731
TProfile * hResoVSPhiPlus_prof_
Definition: Histograms.h:1734
double ** resoVsPtEta_
Definition: Histograms.h:1723
int HFunctionResolution::getXindex ( const double &  x) const
inlineprotected

Definition at line 1719 of file Histograms.h.

References deltaX_, totBinsX_, and xMin_.

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

1719 { return int((x - xMin_) / deltaX_ * totBinsX_); }
int HFunctionResolution::getYindex ( const double &  y) const
inlineprotected

Definition at line 1720 of file Histograms.h.

References deltaY_, totBinsY_, and yMin_.

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

1720 { return int((y - yMin_) / deltaY_ * totBinsY_); }
void HFunctionResolution::Write ( )
inlineoverride

Definition at line 1659 of file Histograms.h.

References svgfig::canvas(), 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_, mps_fire::i, dqmiolumiharvest::j, N, resoCount_, resoVsPtEta_, totBinsX_, and totBinsY_.

Referenced by HFunctionResolutionVarianceCheck::Write().

1659  {
1660  if (histoDir_ != nullptr)
1661  histoDir_->cd();
1662 
1663  hReso_->Write();
1664 
1665  for (int i = 0; i < totBinsX_; ++i) {
1666  for (int j = 0; j < totBinsY_; ++j) {
1667  int N = resoCount_[i][j];
1668  // Fill with the mean value
1669  if (N != 0)
1670  hResoVSPtEta_->SetBinContent(i + 1, j + 1, resoVsPtEta_[i][j] / N);
1671  else
1672  hResoVSPtEta_->SetBinContent(i + 1, j + 1, 0);
1673  }
1674  }
1675  hResoVSPtEta_->Write();
1676 
1677  hResoVSPt_prof_->Write();
1678  hResoVSPt_Bar_prof_->Write();
1679  hResoVSPt_Endc_17_prof_->Write();
1680  hResoVSPt_Endc_20_prof_->Write();
1681  hResoVSPt_Endc_24_prof_->Write();
1682  hResoVSPt_Ovlap_prof_->Write();
1683  hResoVSEta_prof_->Write();
1684  //hResoVSTheta_prof_->Write();
1685  hResoVSPhiMinus_prof_->Write();
1686  hResoVSPhiPlus_prof_->Write();
1687  hResoVSPhi_prof_->Write();
1688 
1689  TCanvas canvas(
1690  TString(hResoVSPtEta_->GetName()) + "_canvas", TString(hResoVSPtEta_->GetTitle()) + " canvas", 1000, 800);
1691  canvas.Divide(2);
1692  canvas.cd(1);
1693  hResoVSPtEta_->Draw("lego");
1694  canvas.cd(2);
1695  hResoVSPtEta_->Draw("surf5");
1696  canvas.Write();
1697  hResoVSPtEta_->Write();
1698 
1699  outputFile_->cd();
1700  }
def canvas
Definition: svgfig.py:482
TProfile * hResoVSPt_Endc_20_prof_
Definition: Histograms.h:1728
TProfile * hResoVSPt_Endc_17_prof_
Definition: Histograms.h:1727
TProfile * hResoVSPt_Bar_prof_
Definition: Histograms.h:1726
TProfile * hResoVSPt_prof_
Definition: Histograms.h:1725
TProfile * hResoVSPt_Ovlap_prof_
Definition: Histograms.h:1730
TProfile * hResoVSPt_Endc_24_prof_
Definition: Histograms.h:1729
TProfile * hResoVSPhi_prof_
Definition: Histograms.h:1735
#define N
Definition: blowfish.cc:9
TProfile * hResoVSPhiMinus_prof_
Definition: Histograms.h:1733
TProfile * hResoVSEta_prof_
Definition: Histograms.h:1731
TProfile * hResoVSPhiPlus_prof_
Definition: Histograms.h:1734
double ** resoVsPtEta_
Definition: Histograms.h:1723

Member Data Documentation

double HFunctionResolution::deltaX_
protected

Definition at line 1738 of file Histograms.h.

Referenced by getXindex(), and HFunctionResolution().

double HFunctionResolution::deltaY_
protected

Definition at line 1738 of file Histograms.h.

Referenced by getYindex(), and HFunctionResolution().

TH1F* HFunctionResolution::hReso_
protected

Definition at line 1721 of file Histograms.h.

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

TProfile* HFunctionResolution::hResoVSEta_prof_
protected

Definition at line 1731 of file Histograms.h.

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

TProfile* HFunctionResolution::hResoVSPhi_prof_
protected

Definition at line 1735 of file Histograms.h.

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

TProfile* HFunctionResolution::hResoVSPhiMinus_prof_
protected

Definition at line 1733 of file Histograms.h.

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

TProfile* HFunctionResolution::hResoVSPhiPlus_prof_
protected

Definition at line 1734 of file Histograms.h.

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

TProfile* HFunctionResolution::hResoVSPt_Bar_prof_
protected

Definition at line 1726 of file Histograms.h.

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

TProfile* HFunctionResolution::hResoVSPt_Endc_17_prof_
protected

Definition at line 1727 of file Histograms.h.

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

TProfile* HFunctionResolution::hResoVSPt_Endc_20_prof_
protected

Definition at line 1728 of file Histograms.h.

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

TProfile* HFunctionResolution::hResoVSPt_Endc_24_prof_
protected

Definition at line 1729 of file Histograms.h.

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

TProfile* HFunctionResolution::hResoVSPt_Ovlap_prof_
protected

Definition at line 1730 of file Histograms.h.

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

TProfile* HFunctionResolution::hResoVSPt_prof_
protected

Definition at line 1725 of file Histograms.h.

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

TH2F* HFunctionResolution::hResoVSPtEta_
protected

Definition at line 1722 of file Histograms.h.

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

int** HFunctionResolution::resoCount_
protected

Definition at line 1724 of file Histograms.h.

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

double** HFunctionResolution::resoVsPtEta_
protected

Definition at line 1723 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 1737 of file Histograms.h.

Referenced by getXindex(), and HFunctionResolution().

double HFunctionResolution::yMin_
protected

Definition at line 1737 of file Histograms.h.

Referenced by getYindex(), and HFunctionResolution().