CMS 3D CMS Logo

List of all members | Public Member Functions | Protected Attributes
HCovarianceVSParts Class Reference

#include <Histograms.h>

Inherits Histograms.

Public Member Functions

void Clear () override
 
void Fill (const reco::Particle::LorentzVector &recoP1, const reco::Particle::LorentzVector &genP1, const reco::Particle::LorentzVector &recoP2, const reco::Particle::LorentzVector &genP2) override
 
double Get (const reco::Particle::LorentzVector &recoP1, const TString &covarianceName) override
 
 HCovarianceVSParts (TFile *outputFile, const TString &name, const double &ptMax)
 
 HCovarianceVSParts (const TString &inputFileName, const TString &name)
 Constructor reading the histograms from file. More...
 
void Write () override
 
 ~HCovarianceVSParts () override
 

Protected Attributes

std::map< TString, HCovarianceVSxy * > mapHisto_
 
bool readMode_
 

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

Constructor & Destructor Documentation

◆ HCovarianceVSParts() [1/2]

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

Definition at line 2090 of file Histograms.h.

References ALCARECOTkAlBeamHalo_cff::etaMax, ALCARECOTkAlBeamHalo_cff::etaMin, mapHisto_, Skims_PA_cff::name, AlignmentTrackSelector_cfi::ptMax, ptMin, and readMode_.

2090  : Histograms(outputFile, name) {
2091  int totBinsX = 40;
2092  int totBinsY = 40;
2093  double etaMin = -3.;
2094  double etaMax = 3.;
2095  double ptMin = 0.;
2096 
2097  readMode_ = false;
2098 
2099  // Variances
2100  mapHisto_[name + "Pt"] =
2101  new HCovarianceVSxy(name + "Pt_", "Pt", totBinsX, ptMin, ptMax, totBinsY, etaMin, etaMax, histoDir_, true);
2102  mapHisto_[name + "CotgTheta"] = new HCovarianceVSxy(
2103  name + "CotgTheta_", "CotgTheta", totBinsX, ptMin, ptMax, totBinsY, etaMin, etaMax, histoDir_, true);
2104  mapHisto_[name + "Phi"] =
2105  new HCovarianceVSxy(name + "Phi_", "Phi", totBinsX, ptMin, ptMax, totBinsY, etaMin, etaMax, histoDir_, true);
2106  // Covariances
2107  mapHisto_[name + "Pt-CotgTheta"] = new HCovarianceVSxy(
2108  name + "Pt_CotgTheta_", "Pt-CotgTheta", totBinsX, ptMin, ptMax, totBinsY, etaMin, etaMax, histoDir_);
2109  mapHisto_[name + "Pt-Phi"] =
2110  new HCovarianceVSxy(name + "Pt_Phi_", "Pt-Phi", totBinsX, ptMin, ptMax, totBinsY, etaMin, etaMax, histoDir_);
2111  mapHisto_[name + "CotgTheta-Phi"] = new HCovarianceVSxy(
2112  name + "CotgTheta_Phi_", "CotgTheta-Phi", totBinsX, ptMin, ptMax, totBinsY, etaMin, etaMax, histoDir_);
2113  mapHisto_[name + "Pt1-Pt2"] =
2114  new HCovarianceVSxy(name + "Pt1_Pt2_", "Pt1-Pt2", totBinsX, ptMin, ptMax, totBinsY, etaMin, etaMax, histoDir_);
2115  mapHisto_[name + "CotgTheta1-CotgTheta2"] = new HCovarianceVSxy(name + "CotgTheta1_CotgTheta2_",
2116  "CotgTheta1-CotgTheta2",
2117  totBinsX,
2118  ptMin,
2119  ptMax,
2120  totBinsY,
2121  etaMin,
2122  etaMax,
2123  histoDir_);
2124  mapHisto_[name + "Phi1-Phi2"] = new HCovarianceVSxy(
2125  name + "Phi1_Phi2_", "Phi1-Phi2", totBinsX, ptMin, ptMax, totBinsY, etaMin, etaMax, histoDir_);
2126  mapHisto_[name + "Pt12-CotgTheta21"] = new HCovarianceVSxy(
2127  name + "Pt12_CotgTheta21_", "Pt12-CotgTheta21", totBinsX, ptMin, ptMax, totBinsY, etaMin, etaMax, histoDir_);
2128  mapHisto_[name + "Pt12-Phi21"] = new HCovarianceVSxy(
2129  name + "Pt12_Phi21_", "Pt12-Phi21", totBinsX, ptMin, ptMax, totBinsY, etaMin, etaMax, histoDir_);
2130  mapHisto_[name + "CotgTheta12-Phi21"] = new HCovarianceVSxy(
2131  name + "CotgTheta12_Phi21_", "CotgTheta12-Phi21", totBinsX, ptMin, ptMax, totBinsY, etaMin, etaMax, histoDir_);
2132  }
constexpr float ptMin
std::map< TString, HCovarianceVSxy * > mapHisto_
Definition: Histograms.h:2254

◆ HCovarianceVSParts() [2/2]

HCovarianceVSParts::HCovarianceVSParts ( const TString &  inputFileName,
const TString &  name 
)
inline

Constructor reading the histograms from file.

Definition at line 2135 of file Histograms.h.

References dtResolutionTest_cfi::inputFile, InefficientDoubleROC::inputFileName, mapHisto_, Skims_PA_cff::name, and readMode_.

2135  {
2136  name_ = name;
2137  TFile* inputFile = new TFile(inputFileName, "READ");
2138  readMode_ = true;
2139 
2140  // Variances
2141  mapHisto_[name_ + "Pt"] = new HCovarianceVSxy(inputFile, name_ + "Pt_" + name_, name_);
2142  mapHisto_[name_ + "CotgTheta"] = new HCovarianceVSxy(inputFile, name_ + "CotgTheta_" + name_, name_);
2143  mapHisto_[name_ + "Phi"] = new HCovarianceVSxy(inputFile, name_ + "Phi_" + name_, name_);
2144  // Covariances
2145  mapHisto_[name_ + "Pt-CotgTheta"] = new HCovarianceVSxy(inputFile, name_ + "Pt_CotgTheta_" + name_, name_);
2146  mapHisto_[name_ + "Pt-Phi"] = new HCovarianceVSxy(inputFile, name_ + "Pt_Phi_" + name_, name_);
2147  mapHisto_[name_ + "CotgTheta-Phi"] = new HCovarianceVSxy(inputFile, name_ + "CotgTheta_Phi_" + name_, name_);
2148  mapHisto_[name_ + "Pt1-Pt2"] = new HCovarianceVSxy(inputFile, name_ + "Pt1_Pt2_" + name_, name_);
2149  mapHisto_[name_ + "CotgTheta1-CotgTheta2"] =
2150  new HCovarianceVSxy(inputFile, name_ + "CotgTheta1_CotgTheta2_" + name_, name_);
2151  mapHisto_[name_ + "Phi1-Phi2"] = new HCovarianceVSxy(inputFile, name_ + "Phi1_Phi2_" + name_, name_);
2152  mapHisto_[name_ + "Pt12-CotgTheta21"] = new HCovarianceVSxy(inputFile, name_ + "Pt12_CotgTheta21_" + name_, name_);
2153  mapHisto_[name_ + "Pt12-Phi21"] = new HCovarianceVSxy(inputFile, name_ + "Pt12_Phi21_" + name_, name_);
2154  mapHisto_[name_ + "CotgTheta12-Phi21"] =
2155  new HCovarianceVSxy(inputFile, name_ + "CotgTheta12_Phi21_" + name_, name_);
2156  }
std::map< TString, HCovarianceVSxy * > mapHisto_
Definition: Histograms.h:2254

◆ ~HCovarianceVSParts()

HCovarianceVSParts::~HCovarianceVSParts ( )
inlineoverride

Definition at line 2158 of file Histograms.h.

References timingPdfMaker::histo, and mapHisto_.

2158  {
2159  for (std::map<TString, HCovarianceVSxy*>::const_iterator histo = mapHisto_.begin(); histo != mapHisto_.end();
2160  histo++) {
2161  delete (*histo).second;
2162  }
2163  }
std::map< TString, HCovarianceVSxy * > mapHisto_
Definition: Histograms.h:2254

Member Function Documentation

◆ Clear()

void HCovarianceVSParts::Clear ( )
inlineoverride

Definition at line 2246 of file Histograms.h.

References timingPdfMaker::histo, and mapHisto_.

2246  {
2247  for (std::map<TString, HCovarianceVSxy*>::const_iterator histo = mapHisto_.begin(); histo != mapHisto_.end();
2248  histo++) {
2249  (*histo).second->Clear();
2250  }
2251  }
std::map< TString, HCovarianceVSxy * > mapHisto_
Definition: Histograms.h:2254

◆ Fill()

void HCovarianceVSParts::Fill ( const reco::Particle::LorentzVector recoP1,
const reco::Particle::LorentzVector genP1,
const reco::Particle::LorentzVector recoP2,
const reco::Particle::LorentzVector genP2 
)
inlineoverride

Definition at line 2169 of file Histograms.h.

References MuScleFitUtils::deltaPhiNoFabs(), HLT_2022v12_cff::eta1, HLT_2022v12_cff::eta2, mapHisto_, HLT_2022v12_cff::pt1, and HLT_2022v12_cff::pt2.

2172  {
2173  double pt1 = recoP1.pt();
2174  double eta1 = recoP1.eta();
2175  double pt2 = recoP2.pt();
2176  double eta2 = recoP2.eta();
2177 
2178  double diffPt1 = (pt1 - genP1.pt()) / genP1.pt();
2179  double diffPt2 = (pt2 - genP2.pt()) / genP2.pt();
2180 
2181  double genTheta1 = genP1.theta();
2182  double genTheta2 = genP2.theta();
2183  double recoTheta1 = recoP1.theta();
2184  double recoTheta2 = recoP2.theta();
2185 
2186  double genCotgTheta1 = TMath::Cos(genTheta1) / (TMath::Sin(genTheta1));
2187  double genCotgTheta2 = TMath::Cos(genTheta2) / (TMath::Sin(genTheta2));
2188  double recoCotgTheta1 = TMath::Cos(recoTheta1) / (TMath::Sin(recoTheta1));
2189  double recoCotgTheta2 = TMath::Cos(recoTheta2) / (TMath::Sin(recoTheta2));
2190 
2191  // double diffCotgTheta1 = (recoCotgTheta1 - genCotgTheta1)/genCotgTheta1;
2192  // double diffCotgTheta2 = (recoCotgTheta2 - genCotgTheta2)/genCotgTheta2;
2193  double diffCotgTheta1 = recoCotgTheta1 - genCotgTheta1;
2194  double diffCotgTheta2 = recoCotgTheta2 - genCotgTheta2;
2195 
2196  // double diffPhi1 = (recoP1.phi() - genP1.phi())/genP1.phi();
2197  // double diffPhi2 = (recoP2.phi() - genP2.phi())/genP2.phi();
2198  double diffPhi1 = MuScleFitUtils::deltaPhiNoFabs(recoP1.phi(), genP1.phi());
2199  double diffPhi2 = MuScleFitUtils::deltaPhiNoFabs(recoP2.phi(), genP2.phi());
2200 
2201  // Fill the variances
2202  mapHisto_[name_ + "Pt"]->Fill(pt1, eta1, diffPt1, diffPt1);
2203  mapHisto_[name_ + "Pt"]->Fill(pt2, eta2, diffPt2, diffPt2);
2204  mapHisto_[name_ + "CotgTheta"]->Fill(pt1, eta1, diffCotgTheta1, diffCotgTheta1);
2205  mapHisto_[name_ + "CotgTheta"]->Fill(pt2, eta2, diffCotgTheta2, diffCotgTheta2);
2206  mapHisto_[name_ + "Phi"]->Fill(pt1, eta1, diffPhi1, diffPhi1);
2207  mapHisto_[name_ + "Phi"]->Fill(pt2, eta2, diffPhi2, diffPhi2);
2208 
2209  // Fill these histograms with both muons
2210  mapHisto_[name_ + "Pt-CotgTheta"]->Fill(pt1, eta1, diffPt1, diffCotgTheta1);
2211  mapHisto_[name_ + "Pt-CotgTheta"]->Fill(pt2, eta2, diffPt2, diffCotgTheta2);
2212  mapHisto_[name_ + "Pt-Phi"]->Fill(pt1, eta1, diffPt1, diffPhi1);
2213  mapHisto_[name_ + "Pt-Phi"]->Fill(pt2, eta2, diffPt2, diffPhi2);
2214  mapHisto_[name_ + "CotgTheta-Phi"]->Fill(pt1, eta1, diffCotgTheta1, diffPhi1);
2215  mapHisto_[name_ + "CotgTheta-Phi"]->Fill(pt2, eta2, diffCotgTheta2, diffPhi2);
2216 
2217  // We fill two (pt, eta) bins for each pair of values. The bin of the
2218  // first and of the second muon. This should take account for the
2219  // assumed symmetry between the exchange of the first with the second muon.
2220  mapHisto_[name_ + "Pt1-Pt2"]->Fill(pt1, eta1, diffPt1, diffPt2);
2221  mapHisto_[name_ + "Pt1-Pt2"]->Fill(pt2, eta2, diffPt1, diffPt2);
2222  mapHisto_[name_ + "CotgTheta1-CotgTheta2"]->Fill(pt1, eta1, diffCotgTheta1, diffCotgTheta2);
2223  mapHisto_[name_ + "CotgTheta1-CotgTheta2"]->Fill(pt2, eta2, diffCotgTheta1, diffCotgTheta2);
2224  mapHisto_[name_ + "Phi1-Phi2"]->Fill(pt1, eta1, diffPhi1, diffPhi2);
2225  mapHisto_[name_ + "Phi1-Phi2"]->Fill(pt2, eta2, diffPhi1, diffPhi2);
2226 
2227  // Fill the following histograms again for each muon (pt, eta) bin. Same
2228  // reason as in the previous case. If the symmetry is true, it does not
2229  // make any difference the order by which we fill the pt and cotgTheta combinations.
2230  mapHisto_[name_ + "Pt12-CotgTheta21"]->Fill(pt1, eta1, diffPt1, diffCotgTheta2);
2231  mapHisto_[name_ + "Pt12-CotgTheta21"]->Fill(pt2, eta2, diffPt2, diffCotgTheta1);
2232  mapHisto_[name_ + "Pt12-Phi21"]->Fill(pt1, eta1, diffPt1, diffPhi2);
2233  mapHisto_[name_ + "Pt12-Phi21"]->Fill(pt2, eta2, diffPt2, diffPhi1);
2234  mapHisto_[name_ + "CotgTheta12-Phi21"]->Fill(pt1, eta1, diffCotgTheta1, diffPhi2);
2235  mapHisto_[name_ + "CotgTheta12-Phi21"]->Fill(pt2, eta2, diffCotgTheta2, diffPhi1);
2236  }
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...
std::map< TString, HCovarianceVSxy * > mapHisto_
Definition: Histograms.h:2254

◆ Get()

double HCovarianceVSParts::Get ( const reco::Particle::LorentzVector recoP1,
const TString &  covarianceName 
)
inlineoverride

Definition at line 2165 of file Histograms.h.

References mapHisto_.

2165  {
2166  return mapHisto_[name_ + covarianceName]->Get(recoP1.pt(), recoP1.eta());
2167  }
std::map< TString, HCovarianceVSxy * > mapHisto_
Definition: Histograms.h:2254

◆ Write()

void HCovarianceVSParts::Write ( )
inlineoverride

Definition at line 2237 of file Histograms.h.

References timingPdfMaker::histo, mapHisto_, and readMode_.

2237  {
2238  if (!readMode_) {
2239  histoDir_->cd();
2240  for (std::map<TString, HCovarianceVSxy*>::const_iterator histo = mapHisto_.begin(); histo != mapHisto_.end();
2241  histo++) {
2242  (*histo).second->Write();
2243  }
2244  }
2245  }
std::map< TString, HCovarianceVSxy * > mapHisto_
Definition: Histograms.h:2254

Member Data Documentation

◆ mapHisto_

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

Definition at line 2254 of file Histograms.h.

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

◆ readMode_

bool HCovarianceVSParts::readMode_
protected

Definition at line 2255 of file Histograms.h.

Referenced by HCovarianceVSParts(), and Write().