CMS 3D CMS Logo

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_
 

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

Constructor & Destructor Documentation

◆ HFunctionResolution()

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

Definition at line 1550 of file Histograms.h.

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

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, Skims_PA_cff::name, AlignmentTrackSelector_cfi::ptMax, resoCount_, resoVsPtEta_, totBinsX_, totBinsY_, multiplicitycorr_cfi::xMax, xMin_, multiplicitycorr_cfi::yMax, and yMin_.

◆ ~HFunctionResolution()

HFunctionResolution::~HFunctionResolution ( )
inlineoverride

Definition at line 1592 of file Histograms.h.

1592  {
1593  delete hReso_;
1594  delete hResoVSPtEta_;
1595  // Free the resolution arrays
1596  for (int i = 0; i < totBinsX_; ++i) {
1597  delete[] resoVsPtEta_[i];
1598  delete[] resoCount_[i];
1599  }
1600  delete[] resoVsPtEta_;
1601  delete[] resoCount_;
1602  // Free the profile histograms
1603  delete hResoVSPt_prof_;
1604  delete hResoVSPt_Bar_prof_;
1605  delete hResoVSPt_Endc_17_prof_;
1606  delete hResoVSPt_Endc_20_prof_;
1607  delete hResoVSPt_Endc_24_prof_;
1608  delete hResoVSPt_Ovlap_prof_;
1609  delete hResoVSEta_prof_;
1610  //delete hResoVSTheta_prof_;
1611  delete hResoVSPhiPlus_prof_;
1612  delete hResoVSPhiMinus_prof_;
1613  delete hResoVSPhi_prof_;
1614  }

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_.

Member Function Documentation

◆ Clear()

void HFunctionResolution::Clear ( )
inlineoverride

Definition at line 1700 of file Histograms.h.

1700  {
1701  hReso_->Clear();
1702  hResoVSPtEta_->Clear();
1703  hResoVSPt_prof_->Clear();
1704  hResoVSPt_Bar_prof_->Clear();
1705  hResoVSPt_Endc_17_prof_->Clear();
1706  hResoVSPt_Endc_20_prof_->Clear();
1707  hResoVSPt_Endc_24_prof_->Clear();
1708  hResoVSPt_Ovlap_prof_->Clear();
1709  hResoVSEta_prof_->Clear();
1710  //hResoVSTheta_prof_->Clear();
1711  hResoVSPhiPlus_prof_->Clear();
1712  hResoVSPhiMinus_prof_->Clear();
1713  hResoVSPhi_prof_->Clear();
1714  }

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_.

◆ Fill()

void HFunctionResolution::Fill ( const reco::Particle::LorentzVector p4,
const double &  resValue,
const int  charge 
)
inlineoverride

Definition at line 1615 of file Histograms.h.

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

References ALCARECOTkAlJpsiMuMu_cff::charge, 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_, p4, resoCount_, resoVsPtEta_, totBinsX_, and totBinsY_.

Referenced by HFunctionResolutionVarianceCheck::Fill().

◆ getXindex()

int HFunctionResolution::getXindex ( const double &  x) const
inlineprotected

Definition at line 1717 of file Histograms.h.

1717 { return int((x - xMin_) / deltaX_ * totBinsX_); }

References deltaX_, createfilelist::int, totBinsX_, x, and xMin_.

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

◆ getYindex()

int HFunctionResolution::getYindex ( const double &  y) const
inlineprotected

Definition at line 1718 of file Histograms.h.

1718 { return int((y - yMin_) / deltaY_ * totBinsY_); }

References deltaY_, createfilelist::int, totBinsY_, y, and yMin_.

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

◆ Write()

void HFunctionResolution::Write ( )
inlineoverride

Definition at line 1657 of file Histograms.h.

1657  {
1658  if (histoDir_ != nullptr)
1659  histoDir_->cd();
1660 
1661  hReso_->Write();
1662 
1663  for (int i = 0; i < totBinsX_; ++i) {
1664  for (int j = 0; j < totBinsY_; ++j) {
1665  int N = resoCount_[i][j];
1666  // Fill with the mean value
1667  if (N != 0)
1668  hResoVSPtEta_->SetBinContent(i + 1, j + 1, resoVsPtEta_[i][j] / N);
1669  else
1670  hResoVSPtEta_->SetBinContent(i + 1, j + 1, 0);
1671  }
1672  }
1673  hResoVSPtEta_->Write();
1674 
1675  hResoVSPt_prof_->Write();
1676  hResoVSPt_Bar_prof_->Write();
1677  hResoVSPt_Endc_17_prof_->Write();
1678  hResoVSPt_Endc_20_prof_->Write();
1679  hResoVSPt_Endc_24_prof_->Write();
1680  hResoVSPt_Ovlap_prof_->Write();
1681  hResoVSEta_prof_->Write();
1682  //hResoVSTheta_prof_->Write();
1683  hResoVSPhiMinus_prof_->Write();
1684  hResoVSPhiPlus_prof_->Write();
1685  hResoVSPhi_prof_->Write();
1686 
1687  TCanvas canvas(
1688  TString(hResoVSPtEta_->GetName()) + "_canvas", TString(hResoVSPtEta_->GetTitle()) + " canvas", 1000, 800);
1689  canvas.Divide(2);
1690  canvas.cd(1);
1691  hResoVSPtEta_->Draw("lego");
1692  canvas.cd(2);
1693  hResoVSPtEta_->Draw("surf5");
1694  canvas.Write();
1695  hResoVSPtEta_->Write();
1696 
1697  outputFile_->cd();
1698  }

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().

Member Data Documentation

◆ deltaX_

double HFunctionResolution::deltaX_
protected

Definition at line 1736 of file Histograms.h.

Referenced by getXindex(), and HFunctionResolution().

◆ deltaY_

double HFunctionResolution::deltaY_
protected

Definition at line 1736 of file Histograms.h.

Referenced by getYindex(), and HFunctionResolution().

◆ hReso_

TH1F* HFunctionResolution::hReso_
protected

Definition at line 1719 of file Histograms.h.

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

◆ hResoVSEta_prof_

TProfile* HFunctionResolution::hResoVSEta_prof_
protected

Definition at line 1729 of file Histograms.h.

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

◆ hResoVSPhi_prof_

TProfile* HFunctionResolution::hResoVSPhi_prof_
protected

Definition at line 1733 of file Histograms.h.

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

◆ hResoVSPhiMinus_prof_

TProfile* HFunctionResolution::hResoVSPhiMinus_prof_
protected

Definition at line 1731 of file Histograms.h.

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

◆ hResoVSPhiPlus_prof_

TProfile* HFunctionResolution::hResoVSPhiPlus_prof_
protected

Definition at line 1732 of file Histograms.h.

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

◆ hResoVSPt_Bar_prof_

TProfile* HFunctionResolution::hResoVSPt_Bar_prof_
protected

Definition at line 1724 of file Histograms.h.

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

◆ hResoVSPt_Endc_17_prof_

TProfile* HFunctionResolution::hResoVSPt_Endc_17_prof_
protected

Definition at line 1725 of file Histograms.h.

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

◆ hResoVSPt_Endc_20_prof_

TProfile* HFunctionResolution::hResoVSPt_Endc_20_prof_
protected

Definition at line 1726 of file Histograms.h.

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

◆ hResoVSPt_Endc_24_prof_

TProfile* HFunctionResolution::hResoVSPt_Endc_24_prof_
protected

Definition at line 1727 of file Histograms.h.

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

◆ hResoVSPt_Ovlap_prof_

TProfile* HFunctionResolution::hResoVSPt_Ovlap_prof_
protected

Definition at line 1728 of file Histograms.h.

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

◆ hResoVSPt_prof_

TProfile* HFunctionResolution::hResoVSPt_prof_
protected

Definition at line 1723 of file Histograms.h.

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

◆ hResoVSPtEta_

TH2F* HFunctionResolution::hResoVSPtEta_
protected

Definition at line 1720 of file Histograms.h.

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

◆ resoCount_

int** HFunctionResolution::resoCount_
protected

Definition at line 1722 of file Histograms.h.

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

◆ resoVsPtEta_

double** HFunctionResolution::resoVsPtEta_
protected

Definition at line 1721 of file Histograms.h.

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

◆ totBinsX_

int HFunctionResolution::totBinsX_
protected

◆ totBinsY_

int HFunctionResolution::totBinsY_
protected

◆ xMin_

double HFunctionResolution::xMin_
protected

Definition at line 1735 of file Histograms.h.

Referenced by getXindex(), and HFunctionResolution().

◆ yMin_

double HFunctionResolution::yMin_
protected

Definition at line 1735 of file Histograms.h.

Referenced by getYindex(), and HFunctionResolution().

HFunctionResolution::hReso_
TH1F * hReso_
Definition: Histograms.h:1719
svgfig.canvas
def canvas(*sub, **attr)
Definition: svgfig.py:482
DDAxes::y
mps_fire.i
i
Definition: mps_fire.py:355
HFunctionResolution::resoCount_
int ** resoCount_
Definition: Histograms.h:1722
HFunctionResolution::deltaY_
double deltaY_
Definition: Histograms.h:1736
HFunctionResolution::hResoVSPt_Bar_prof_
TProfile * hResoVSPt_Bar_prof_
Definition: Histograms.h:1724
HFunctionResolution::hResoVSPt_Ovlap_prof_
TProfile * hResoVSPt_Ovlap_prof_
Definition: Histograms.h:1728
HFunctionResolution::hResoVSPt_Endc_24_prof_
TProfile * hResoVSPt_Endc_24_prof_
Definition: Histograms.h:1727
DDAxes::x
HFunctionResolution::hResoVSPt_prof_
TProfile * hResoVSPt_prof_
Definition: Histograms.h:1723
HFunctionResolution::hResoVSPhiPlus_prof_
TProfile * hResoVSPhiPlus_prof_
Definition: Histograms.h:1732
download_sqlite_cfg.outputFile
outputFile
Definition: download_sqlite_cfg.py:5
HFunctionResolution::hResoVSPt_Endc_20_prof_
TProfile * hResoVSPt_Endc_20_prof_
Definition: Histograms.h:1726
AlignmentTrackSelector_cfi.ptMax
ptMax
Definition: AlignmentTrackSelector_cfi.py:12
HFunctionResolution::xMin_
double xMin_
Definition: Histograms.h:1735
N
#define N
Definition: blowfish.cc:9
HFunctionResolution::hResoVSPhiMinus_prof_
TProfile * hResoVSPhiMinus_prof_
Definition: Histograms.h:1731
multiplicitycorr_cfi.yMax
yMax
Definition: multiplicitycorr_cfi.py:6
HFunctionResolution::getYindex
int getYindex(const double &y) const
Definition: Histograms.h:1718
ALCARECOTkAlJpsiMuMu_cff.charge
charge
Definition: ALCARECOTkAlJpsiMuMu_cff.py:47
HFunctionResolution::deltaX_
double deltaX_
Definition: Histograms.h:1736
createfilelist.int
int
Definition: createfilelist.py:10
HFunctionResolution::hResoVSPtEta_
TH2F * hResoVSPtEta_
Definition: Histograms.h:1720
p4
double p4[4]
Definition: TauolaWrapper.h:92
HFunctionResolution::hResoVSPhi_prof_
TProfile * hResoVSPhi_prof_
Definition: Histograms.h:1733
HFunctionResolution::totBinsX_
int totBinsX_
Definition: Histograms.h:1734
HFunctionResolution::resoVsPtEta_
double ** resoVsPtEta_
Definition: Histograms.h:1721
HFunctionResolution::getXindex
int getXindex(const double &x) const
Definition: Histograms.h:1717
HFunctionResolution::hResoVSEta_prof_
TProfile * hResoVSEta_prof_
Definition: Histograms.h:1729
Skims_PA_cff.name
name
Definition: Skims_PA_cff.py:17
HFunctionResolution::yMin_
double yMin_
Definition: Histograms.h:1735
multiplicitycorr_cfi.xMax
xMax
Definition: multiplicitycorr_cfi.py:5
HFunctionResolution::hResoVSPt_Endc_17_prof_
TProfile * hResoVSPt_Endc_17_prof_
Definition: Histograms.h:1725
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
HFunctionResolution::totBinsY_
int totBinsY_
Definition: Histograms.h:1734