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
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 &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.)
 
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 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 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 1986 of file Histograms.h.

Constructor & Destructor Documentation

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

Definition at line 1989 of file Histograms.h.

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

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

Constructor reading the histograms from file.

Definition at line 2015 of file Histograms.h.

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

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

Definition at line 2037 of file Histograms.h.

References interpolateCardsSimple::histo, and mapHisto_.

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

Member Function Documentation

virtual void HCovarianceVSParts::Clear ( )
inlinevirtual

Implements Histograms.

Definition at line 2126 of file Histograms.h.

References interpolateCardsSimple::histo, and mapHisto_.

2126  {
2127  for (std::map<TString, HCovarianceVSxy*>::const_iterator histo=mapHisto_.begin();
2128  histo!=mapHisto_.end(); histo++) {
2129  (*histo).second->Clear();
2130  }
2131  }
std::map< TString, HCovarianceVSxy * > mapHisto_
Definition: Histograms.h:2133
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 2048 of file Histograms.h.

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

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

Reimplemented from Histograms.

Definition at line 2044 of file Histograms.h.

References mapHisto_, and Histograms::name_.

2044  {
2045  return mapHisto_[name_+covarianceName]->Get(recoP1.pt(), recoP1.eta());
2046  }
TString name_
Definition: Histograms.h:124
std::map< TString, HCovarianceVSxy * > mapHisto_
Definition: Histograms.h:2133
virtual void HCovarianceVSParts::Write ( )
inlinevirtual

Implements Histograms.

Definition at line 2117 of file Histograms.h.

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

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

Member Data Documentation

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

Definition at line 2133 of file Histograms.h.

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

bool HCovarianceVSParts::readMode_
protected

Definition at line 2134 of file Histograms.h.

Referenced by HCovarianceVSParts(), and Write().