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
HCovarianceVSParts Class Reference

#include <Histograms.h>

Inheritance diagram for HCovarianceVSParts:
Histograms

Public Member Functions

virtual void Clear ()
 
virtual void Fill (const reco::Particle::LorentzVector &recoP1, const reco::Particle::LorentzVector &genP1, const reco::Particle::LorentzVector &recoP2, const reco::Particle::LorentzVector &genP2)
 
virtual double Get (const reco::Particle::LorentzVector &recoP1, const TString &covarianceName)
 
 HCovarianceVSParts (TFile *outputFile, const TString &name, const double &ptMax)
 
 HCovarianceVSParts (const TString &inputFileName, const TString &name)
 Constructor reading the histograms from file. More...
 
virtual void Write ()
 
 ~HCovarianceVSParts ()
 
- Public Member Functions inherited from Histograms
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 &resValue, const int charge)
 
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 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.)
 
virtual TString GetName ()
 
 Histograms ()
 
 Histograms (const TString &name)
 
 Histograms (TFile *outputFile, const TString &name)
 
virtual void SetWeight (double weight)
 
virtual ~Histograms ()
 

Protected Attributes

std::map< TString,
HCovarianceVSxy * > 
mapHisto_
 
bool readMode_
 
- Protected Attributes inherited from Histograms
TDirectory * histoDir_
 
TString name_
 
TFile * outputFile_
 
double theWeight_
 

Detailed Description

This class uses the HCovariance histograms to compute the covariances between the two input muons kinematic quantities. The covariances are computed against pt and eta.

Definition at line 1983 of file Histograms.h.

Constructor & Destructor Documentation

HCovarianceVSParts::HCovarianceVSParts ( TFile *  outputFile,
const TString &  name,
const double &  ptMax 
)
inline

Definition at line 1986 of file Histograms.h.

References jptDQMConfig_cff::etaMax, jptDQMConfig_cff::etaMin, Histograms::histoDir_, mapHisto_, PtMinSelector_cfg::ptMin, and readMode_.

1986  : Histograms( outputFile, name ) {
1987  int totBinsX = 40;
1988  int totBinsY = 40;
1989  double etaMin = -3.;
1990  double etaMax = 3.;
1991  double ptMin = 0.;
1992 
1993  readMode_ = false;
1994 
1995  // Variances
1996  mapHisto_[name+"Pt"] = new HCovarianceVSxy(name+"Pt_", "Pt", totBinsX, ptMin, ptMax, totBinsY, etaMin, etaMax, histoDir_, true);
1997  mapHisto_[name+"CotgTheta"] = new HCovarianceVSxy(name+"CotgTheta_", "CotgTheta", totBinsX, ptMin, ptMax, totBinsY, etaMin, etaMax, histoDir_, true);
1998  mapHisto_[name+"Phi"] = new HCovarianceVSxy(name+"Phi_", "Phi", totBinsX, ptMin, ptMax, totBinsY, etaMin, etaMax, histoDir_, true);
1999  // Covariances
2000  mapHisto_[name+"Pt-CotgTheta"] = new HCovarianceVSxy(name+"Pt_CotgTheta_", "Pt-CotgTheta", totBinsX, ptMin, ptMax, totBinsY, etaMin, etaMax, histoDir_);
2001  mapHisto_[name+"Pt-Phi"] = new HCovarianceVSxy(name+"Pt_Phi_", "Pt-Phi", totBinsX, ptMin, ptMax, totBinsY, etaMin, etaMax, histoDir_);
2002  mapHisto_[name+"CotgTheta-Phi"] = new HCovarianceVSxy(name+"CotgTheta_Phi_", "CotgTheta-Phi", totBinsX, ptMin, ptMax, totBinsY, etaMin, etaMax, histoDir_);
2003  mapHisto_[name+"Pt1-Pt2"] = new HCovarianceVSxy(name+"Pt1_Pt2_", "Pt1-Pt2", totBinsX, ptMin, ptMax, totBinsY, etaMin, etaMax, histoDir_);
2004  mapHisto_[name+"CotgTheta1-CotgTheta2"] = new HCovarianceVSxy(name+"CotgTheta1_CotgTheta2_", "CotgTheta1-CotgTheta2", totBinsX, ptMin, ptMax, totBinsY, etaMin, etaMax, histoDir_);
2005  mapHisto_[name+"Phi1-Phi2"] = new HCovarianceVSxy(name+"Phi1_Phi2_", "Phi1-Phi2", totBinsX, ptMin, ptMax, totBinsY, etaMin, etaMax, histoDir_);
2006  mapHisto_[name+"Pt12-CotgTheta21"] = new HCovarianceVSxy(name+"Pt12_CotgTheta21_", "Pt12-CotgTheta21", totBinsX, ptMin, ptMax, totBinsY, etaMin, etaMax, histoDir_);
2007  mapHisto_[name+"Pt12-Phi21"] = new HCovarianceVSxy(name+"Pt12_Phi21_", "Pt12-Phi21", totBinsX, ptMin, ptMax, totBinsY, etaMin, etaMax, histoDir_);
2008  mapHisto_[name+"CotgTheta12-Phi21"] = new HCovarianceVSxy(name+"CotgTheta12_Phi21_", "CotgTheta12-Phi21", totBinsX, ptMin, ptMax, totBinsY, etaMin, etaMax, histoDir_);
2009  }
std::map< TString, HCovarianceVSxy * > mapHisto_
Definition: Histograms.h:2130
TDirectory * histoDir_
Definition: Histograms.h:123
HCovarianceVSParts::HCovarianceVSParts ( const TString &  inputFileName,
const TString &  name 
)
inline

Constructor reading the histograms from file.

Definition at line 2012 of file Histograms.h.

References analyzePatCleaning_cfg::inputFile, mapHisto_, mergeVDriftHistosByStation::name, Histograms::name_, and readMode_.

2013  {
2014  name_ = name;
2015  TFile * inputFile = new TFile(inputFileName, "READ");
2016  readMode_ = true;
2017 
2018  // Variances
2019  mapHisto_[name_+"Pt"] = new HCovarianceVSxy(inputFile, name_+"Pt_"+name_, name_);
2020  mapHisto_[name_+"CotgTheta"] = new HCovarianceVSxy(inputFile, name_+"CotgTheta_"+name_, name_);
2021  mapHisto_[name_+"Phi"] = new HCovarianceVSxy(inputFile, name_+"Phi_"+name_, name_);
2022  // Covariances
2023  mapHisto_[name_+"Pt-CotgTheta"] = new HCovarianceVSxy(inputFile, name_+"Pt_CotgTheta_"+name_, name_);
2024  mapHisto_[name_+"Pt-Phi"] = new HCovarianceVSxy(inputFile, name_+"Pt_Phi_"+name_, name_);
2025  mapHisto_[name_+"CotgTheta-Phi"] = new HCovarianceVSxy(inputFile, name_+"CotgTheta_Phi_"+name_, name_);
2026  mapHisto_[name_+"Pt1-Pt2"] = new HCovarianceVSxy(inputFile, name_+"Pt1_Pt2_"+name_, name_);
2027  mapHisto_[name_+"CotgTheta1-CotgTheta2"] = new HCovarianceVSxy(inputFile, name_+"CotgTheta1_CotgTheta2_"+name_, name_);
2028  mapHisto_[name_+"Phi1-Phi2"] = new HCovarianceVSxy(inputFile, name_+"Phi1_Phi2_"+name_, name_);
2029  mapHisto_[name_+"Pt12-CotgTheta21"] = new HCovarianceVSxy(inputFile, name_+"Pt12_CotgTheta21_"+name_, name_);
2030  mapHisto_[name_+"Pt12-Phi21"] = new HCovarianceVSxy(inputFile, name_+"Pt12_Phi21_"+name_, name_);
2031  mapHisto_[name_+"CotgTheta12-Phi21"] = new HCovarianceVSxy(inputFile, name_+"CotgTheta12_Phi21_"+name_, name_);
2032  }
TString name_
Definition: Histograms.h:121
std::map< TString, HCovarianceVSxy * > mapHisto_
Definition: Histograms.h:2130
HCovarianceVSParts::~HCovarianceVSParts ( )
inline

Definition at line 2034 of file Histograms.h.

References timingPdfMaker::histo, and mapHisto_.

2034  {
2035  for (std::map<TString, HCovarianceVSxy*>::const_iterator histo=mapHisto_.begin();
2036  histo!=mapHisto_.end(); histo++) {
2037  delete (*histo).second;
2038  }
2039  }
std::map< TString, HCovarianceVSxy * > mapHisto_
Definition: Histograms.h:2130

Member Function Documentation

virtual void HCovarianceVSParts::Clear ( )
inlinevirtual

Implements Histograms.

Definition at line 2123 of file Histograms.h.

References timingPdfMaker::histo, and mapHisto_.

2123  {
2124  for (std::map<TString, HCovarianceVSxy*>::const_iterator histo=mapHisto_.begin();
2125  histo!=mapHisto_.end(); histo++) {
2126  (*histo).second->Clear();
2127  }
2128  }
std::map< TString, HCovarianceVSxy * > mapHisto_
Definition: Histograms.h:2130
virtual void HCovarianceVSParts::Fill ( const reco::Particle::LorentzVector recoP1,
const reco::Particle::LorentzVector genP1,
const reco::Particle::LorentzVector recoP2,
const reco::Particle::LorentzVector genP2 
)
inlinevirtual

Reimplemented from Histograms.

Definition at line 2045 of file Histograms.h.

References MuScleFitUtils::deltaPhiNoFabs(), mapHisto_, and Histograms::name_.

2048  {
2049 
2050  double pt1 = recoP1.pt();
2051  double eta1 = recoP1.eta();
2052  double pt2 = recoP2.pt();
2053  double eta2 = recoP2.eta();
2054 
2055  double diffPt1 = (pt1 - genP1.pt())/genP1.pt();
2056  double diffPt2 = (pt2 - genP2.pt())/genP2.pt();
2057 
2058  double genTheta1 = genP1.theta();
2059  double genTheta2 = genP2.theta();
2060  double recoTheta1 = recoP1.theta();
2061  double recoTheta2 = recoP2.theta();
2062 
2063  double genCotgTheta1 = TMath::Cos(genTheta1)/(TMath::Sin(genTheta1));
2064  double genCotgTheta2 = TMath::Cos(genTheta2)/(TMath::Sin(genTheta2));
2065  double recoCotgTheta1 = TMath::Cos(recoTheta1)/(TMath::Sin(recoTheta1));
2066  double recoCotgTheta2 = TMath::Cos(recoTheta2)/(TMath::Sin(recoTheta2));
2067 
2068  // double diffCotgTheta1 = (recoCotgTheta1 - genCotgTheta1)/genCotgTheta1;
2069  // double diffCotgTheta2 = (recoCotgTheta2 - genCotgTheta2)/genCotgTheta2;
2070  double diffCotgTheta1 = recoCotgTheta1 - genCotgTheta1;
2071  double diffCotgTheta2 = recoCotgTheta2 - genCotgTheta2;
2072 
2073  // double diffPhi1 = (recoP1.phi() - genP1.phi())/genP1.phi();
2074  // double diffPhi2 = (recoP2.phi() - genP2.phi())/genP2.phi();
2075  double diffPhi1 = MuScleFitUtils::deltaPhiNoFabs(recoP1.phi(), genP1.phi());
2076  double diffPhi2 = MuScleFitUtils::deltaPhiNoFabs(recoP2.phi(), genP2.phi());
2077 
2078  // Fill the variances
2079  mapHisto_[name_+"Pt"]->Fill(pt1, eta1, diffPt1, diffPt1);
2080  mapHisto_[name_+"Pt"]->Fill(pt2, eta2, diffPt2, diffPt2);
2081  mapHisto_[name_+"CotgTheta"]->Fill(pt1, eta1, diffCotgTheta1, diffCotgTheta1);
2082  mapHisto_[name_+"CotgTheta"]->Fill(pt2, eta2, diffCotgTheta2, diffCotgTheta2);
2083  mapHisto_[name_+"Phi"]->Fill(pt1, eta1, diffPhi1, diffPhi1);
2084  mapHisto_[name_+"Phi"]->Fill(pt2, eta2, diffPhi2, diffPhi2);
2085 
2086  // Fill these histograms with both muons
2087  mapHisto_[name_+"Pt-CotgTheta"]->Fill(pt1, eta1, diffPt1, diffCotgTheta1 );
2088  mapHisto_[name_+"Pt-CotgTheta"]->Fill(pt2, eta2, diffPt2, diffCotgTheta2 );
2089  mapHisto_[name_+"Pt-Phi"]->Fill(pt1, eta1, diffPt1, diffPhi1);
2090  mapHisto_[name_+"Pt-Phi"]->Fill(pt2, eta2, diffPt2, diffPhi2);
2091  mapHisto_[name_+"CotgTheta-Phi"]->Fill(pt1, eta1, diffCotgTheta1, diffPhi1);
2092  mapHisto_[name_+"CotgTheta-Phi"]->Fill(pt2, eta2, diffCotgTheta2, diffPhi2);
2093 
2094  // We fill two (pt, eta) bins for each pair of values. The bin of the
2095  // first and of the second muon. This should take account for the
2096  // assumed symmetry between the exchange of the first with the second muon.
2097  mapHisto_[name_+"Pt1-Pt2"]->Fill(pt1, eta1, diffPt1, diffPt2);
2098  mapHisto_[name_+"Pt1-Pt2"]->Fill(pt2, eta2, diffPt1, diffPt2);
2099  mapHisto_[name_+"CotgTheta1-CotgTheta2"]->Fill(pt1, eta1, diffCotgTheta1, diffCotgTheta2);
2100  mapHisto_[name_+"CotgTheta1-CotgTheta2"]->Fill(pt2, eta2, diffCotgTheta1, diffCotgTheta2);
2101  mapHisto_[name_+"Phi1-Phi2"]->Fill(pt1, eta1, diffPhi1, diffPhi2);
2102  mapHisto_[name_+"Phi1-Phi2"]->Fill(pt2, eta2, diffPhi1, diffPhi2);
2103 
2104  // Fill the following histograms again for each muon (pt, eta) bin. Same
2105  // reason as in the previous case. If the symmetry is true, it does not
2106  // make any difference the order by which we fill the pt and cotgTheta combinations.
2107  mapHisto_[name_+"Pt12-CotgTheta21"]->Fill(pt1, eta1, diffPt1, diffCotgTheta2);
2108  mapHisto_[name_+"Pt12-CotgTheta21"]->Fill(pt2, eta2, diffPt2, diffCotgTheta1);
2109  mapHisto_[name_+"Pt12-Phi21"]->Fill(pt1, eta1, diffPt1, diffPhi2);
2110  mapHisto_[name_+"Pt12-Phi21"]->Fill(pt2, eta2, diffPt2, diffPhi1);
2111  mapHisto_[name_+"CotgTheta12-Phi21"]->Fill(pt1, eta1, diffCotgTheta1, diffPhi2);
2112  mapHisto_[name_+"CotgTheta12-Phi21"]->Fill(pt2, eta2, diffCotgTheta2, diffPhi1);
2113  }
static double deltaPhiNoFabs(const double &phi1, const double &phi2)
Without fabs at the end, used to have a symmetric distribution for the resolution fits and variance c...
TString name_
Definition: Histograms.h:121
std::map< TString, HCovarianceVSxy * > mapHisto_
Definition: Histograms.h:2130
virtual double HCovarianceVSParts::Get ( const reco::Particle::LorentzVector recoP1,
const TString &  covarianceName 
)
inlinevirtual

Reimplemented from Histograms.

Definition at line 2041 of file Histograms.h.

References mapHisto_, and Histograms::name_.

2041  {
2042  return mapHisto_[name_+covarianceName]->Get(recoP1.pt(), recoP1.eta());
2043  }
TString name_
Definition: Histograms.h:121
std::map< TString, HCovarianceVSxy * > mapHisto_
Definition: Histograms.h:2130
virtual void HCovarianceVSParts::Write ( )
inlinevirtual

Implements Histograms.

Definition at line 2114 of file Histograms.h.

References timingPdfMaker::histo, Histograms::histoDir_, mapHisto_, and readMode_.

2114  {
2115  if( !readMode_ ) {
2116  histoDir_->cd();
2117  for (std::map<TString, HCovarianceVSxy*>::const_iterator histo=mapHisto_.begin();
2118  histo!=mapHisto_.end(); histo++) {
2119  (*histo).second->Write();
2120  }
2121  }
2122  }
std::map< TString, HCovarianceVSxy * > mapHisto_
Definition: Histograms.h:2130
TDirectory * histoDir_
Definition: Histograms.h:123

Member Data Documentation

std::map<TString, HCovarianceVSxy*> HCovarianceVSParts::mapHisto_
protected

Definition at line 2130 of file Histograms.h.

Referenced by Clear(), Fill(), Get(), HCovarianceVSParts(), Write(), and ~HCovarianceVSParts().

bool HCovarianceVSParts::readMode_
protected

Definition at line 2131 of file Histograms.h.

Referenced by HCovarianceVSParts(), and Write().