CMS 3D CMS Logo

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

#include <Histograms.h>

Inherits Histograms.

Public Member Functions

void Clear () override
 
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) override
 
void Fill (const reco::Particle::LorentzVector &recoP1, const int charge1, const reco::Particle::LorentzVector &recoP2, const int charge2, const double &recoMass, const double &genMass) override
 
 HMassResolutionVSPart (TFile *outputFile, const TString &name)
 
void Write () override
 
 ~HMassResolutionVSPart () override
 

Protected Attributes

std::map< TString, TH1 * > mapHisto_
 
std::unique_ptr< HDeltamuMinus
 
std::unique_ptr< HDeltamuPlus
 
TString nameSuffix_ [2]
 

Detailed Description

A set of histograms for resolution. The fill method requires the two selected muons, their charges and the reconstructed and generated masses. It evaluates the resolution on the mass vs: Pt of the pair DeltaEta of the pair DeltaPhi of the pair pt, eta and phi of the plus and minus muon separately

Definition at line 2265 of file Histograms.h.

Constructor & Destructor Documentation

HMassResolutionVSPart::HMassResolutionVSPart ( TFile *  outputFile,
const TString &  name 
)
inline

Definition at line 2267 of file Histograms.h.

References mps_fire::i, and Skims_PA_cff::name.

2267  : Histograms(outputFile, name) {
2268  // Kinematical variables
2269  nameSuffix_[0] = "Plus";
2270  nameSuffix_[1] = "Minus";
2271  TString titleSuffix[] = {" for mu+", " for mu-"};
2272 
2273  mapHisto_[name] = new TH1F(name, "#Delta M/M", 4000, -1, 1);
2274  mapHisto_[name + "VSPairPt"] =
2275  new HResolution(name + "VSPairPt", "resolution VS pt of the pair", 100, 0, 200, -1, 1, histoDir_);
2276  mapHisto_[name + "VSPairDeltaEta"] = new HResolution(
2277  name + "VSPairDeltaEta", "resolution VS #Delta#eta of the pair", 100, -0.1, 6.2, -1, 1, histoDir_);
2278  mapHisto_[name + "VSPairDeltaPhi"] = new HResolution(
2279  name + "VSPairDeltaPhi", "resolution VS #Delta#phi of the pair", 100, -0.1, 3.2, -1, 1, histoDir_);
2280 
2281  for (int i = 0; i < 2; ++i) {
2282  mapHisto_[name + "VSPt" + nameSuffix_[i]] = new HResolution(
2283  name + "VSPt" + nameSuffix_[i], "resolution VS pt" + titleSuffix[i], 100, 0, 200, -1, 1, histoDir_);
2284  mapHisto_[name + "VSEta" + nameSuffix_[i]] = new HResolution(
2285  name + "VSEta" + nameSuffix_[i], "resolution VS #eta" + titleSuffix[i], 100, -3, 3, -1, 1, histoDir_);
2286  mapHisto_[name + "VSPhi" + nameSuffix_[i]] = new HResolution(
2287  name + "VSPhi" + nameSuffix_[i], "resolution VS #phi" + titleSuffix[i], 100, -3.2, 3.2, -1, 1, histoDir_);
2288  }
2289 
2290  // single particles histograms
2291  muMinus.reset(new HDelta("muMinus"));
2292  muPlus.reset(new HDelta("muPlus"));
2293  }
std::unique_ptr< HDelta > muPlus
Definition: Histograms.h:2403
std::unique_ptr< HDelta > muMinus
Definition: Histograms.h:2402
TString nameSuffix_[2]
Definition: Histograms.h:2401
std::map< TString, TH1 * > mapHisto_
Definition: Histograms.h:2400
HMassResolutionVSPart::~HMassResolutionVSPart ( )
inlineoverride

Definition at line 2295 of file Histograms.h.

References timingPdfMaker::histo.

2295  {
2296  for (std::map<TString, TH1*>::const_iterator histo = mapHisto_.begin(); histo != mapHisto_.end(); histo++) {
2297  delete (*histo).second;
2298  }
2299  }
std::map< TString, TH1 * > mapHisto_
Definition: Histograms.h:2400

Member Function Documentation

void HMassResolutionVSPart::Clear ( )
inlineoverride

Definition at line 2377 of file Histograms.h.

References timingPdfMaker::histo.

2377  {
2378  for (std::map<TString, TH1*>::const_iterator histo = mapHisto_.begin(); histo != mapHisto_.end(); histo++) {
2379  (*histo).second->Clear();
2380  }
2381  muMinus->Clear();
2382  muPlus->Clear();
2383  }
std::unique_ptr< HDelta > muPlus
Definition: Histograms.h:2403
std::unique_ptr< HDelta > muMinus
Definition: Histograms.h:2402
std::map< TString, TH1 * > mapHisto_
Definition: Histograms.h:2400
void HMassResolutionVSPart::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 
)
inlineoverride

Definition at line 2301 of file Histograms.h.

References HcalObjRepresent::Fill().

2308  {
2309  muMinus->Fill(recoP1, genP1);
2310  muPlus->Fill(recoP2, genP2);
2311 
2312  Fill(recoP1, charge1, recoP2, charge2, recoMass, genMass);
2313  }
std::unique_ptr< HDelta > muPlus
Definition: Histograms.h:2403
std::unique_ptr< HDelta > muMinus
Definition: Histograms.h:2402
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) override
Definition: Histograms.h:2301
void HMassResolutionVSPart::Fill ( const reco::Particle::LorentzVector recoP1,
const int  charge1,
const reco::Particle::LorentzVector recoP2,
const int  charge2,
const double &  recoMass,
const double &  genMass 
)
inlineoverride

Definition at line 2315 of file Histograms.h.

References gather_cfg::cout, MuScleFitUtils::deltaPhi(), and mps_fire::i.

2322  {
2323  if (charge1 == charge2)
2324  std::cout << "Error: must get two opposite charge particles" << std::endl;
2325 
2326  double massRes = (recoMass - genMass) / genMass;
2327 
2328  reco::Particle::LorentzVector recoPair(recoP1 + recoP2);
2329  double pairPt = recoPair.Pt();
2330 
2331  double recoPt[2] = {recoP1.Pt(), recoP2.Pt()};
2332  double recoEta[2] = {recoP1.Eta(), recoP2.Eta()};
2333  double recoPhi[2] = {recoP1.Phi(), recoP2.Phi()};
2334 
2335  // std::cout << "pairPt = " << pairPt << ", massRes = ("<<recoMass<<" - "<<genMass<<")/"<<genMass<<" = " << massRes
2336  // << ", recoPt[0] = " << recoPt[0] << ", recoPt[1] = " << recoPt[1] << std::endl;
2337 
2338  // Index of the histogram. If the muons have charge1 = -1 and charge2 = 1, they already have the
2339  // correct histogram indeces. Otherwise, swap the indeces.
2340  // In any case the negative muon has index = 0 and the positive muon has index = 1.
2341  int id[2] = {0, 1};
2342  if (charge1 > 0) {
2343  id[0] = 1;
2344  id[1] = 0;
2345  }
2346 
2347  double pairDeltaEta = fabs(recoEta[0] - recoEta[1]);
2348  double pairDeltaPhi = MuScleFitUtils::deltaPhi(recoPhi[0], recoPhi[1]);
2349 
2350  mapHisto_[name_]->Fill(massRes);
2351  mapHisto_[name_ + "VSPairPt"]->Fill(pairPt, massRes);
2352  mapHisto_[name_ + "VSPairDeltaEta"]->Fill(pairDeltaEta, massRes);
2353  mapHisto_[name_ + "VSPairDeltaPhi"]->Fill(pairDeltaPhi, massRes);
2354 
2355  // Resolution vs single muon quantities
2356  // ------------------------------------
2357  int index = 0;
2358  for (int i = 0; i < 2; ++i) {
2359  index = id[i];
2360  mapHisto_[name_ + "VSPt" + nameSuffix_[i]]->Fill(recoPt[index], massRes); // EM [index] or [i]???
2361  mapHisto_[name_ + "VSEta" + nameSuffix_[i]]->Fill(recoEta[index], massRes);
2362  mapHisto_[name_ + "VSPhi" + nameSuffix_[i]]->Fill(recoPhi[index], massRes);
2363  }
2364  }
TString nameSuffix_[2]
Definition: Histograms.h:2401
std::map< TString, TH1 * > mapHisto_
Definition: Histograms.h:2400
static double deltaPhi(const double &phi1, const double &phi2)
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:21
void HMassResolutionVSPart::Write ( )
inlineoverride

Definition at line 2366 of file Histograms.h.

References hippyaddtobaddatafiles::cd(), and timingPdfMaker::histo.

2366  {
2367  histoDir_->cd();
2368  for (std::map<TString, TH1*>::const_iterator histo = mapHisto_.begin(); histo != mapHisto_.end(); histo++) {
2369  (*histo).second->Write();
2370  }
2371  // Create the new dir and cd into it
2372  (histoDir_->mkdir("singleMuonsVSgen"))->cd();
2373  muMinus->Write();
2374  muPlus->Write();
2375  }
std::unique_ptr< HDelta > muPlus
Definition: Histograms.h:2403
std::unique_ptr< HDelta > muMinus
Definition: Histograms.h:2402
std::map< TString, TH1 * > mapHisto_
Definition: Histograms.h:2400

Member Data Documentation

std::map<TString, TH1*> HMassResolutionVSPart::mapHisto_
protected

Definition at line 2400 of file Histograms.h.

std::unique_ptr<HDelta> HMassResolutionVSPart::muMinus
protected

Definition at line 2402 of file Histograms.h.

std::unique_ptr<HDelta> HMassResolutionVSPart::muPlus
protected

Definition at line 2403 of file Histograms.h.

TString HMassResolutionVSPart::nameSuffix_[2]
protected

Definition at line 2401 of file Histograms.h.