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

Constructor & Destructor Documentation

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

Definition at line 2147 of file Histograms.h.

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

2147  : Histograms( outputFile, name ) {
2148  // Kinematical variables
2149  nameSuffix_[0] = "Plus";
2150  nameSuffix_[1] = "Minus";
2151  TString titleSuffix[] = {" for mu+", " for mu-"};
2152 
2153  mapHisto_[name] = new TH1F (name, "#Delta M/M", 4000, -1, 1);
2154  mapHisto_[name+"VSPairPt"] = new HResolution (name+"VSPairPt", "resolution VS pt of the pair", 100, 0, 200, -1, 1, histoDir_);
2155  mapHisto_[name+"VSPairDeltaEta"] = new HResolution (name+"VSPairDeltaEta", "resolution VS #Delta#eta of the pair", 100, -0.1, 6.2, -1, 1, histoDir_);
2156  mapHisto_[name+"VSPairDeltaPhi"] = new HResolution (name+"VSPairDeltaPhi", "resolution VS #Delta#phi of the pair", 100, -0.1, 3.2, -1, 1, histoDir_);
2157 
2158  for( int i=0; i<2; ++i ) {
2159  mapHisto_[name+"VSPt"+nameSuffix_[i]] = new HResolution (name+"VSPt"+nameSuffix_[i], "resolution VS pt"+titleSuffix[i], 100, 0, 200, -1, 1, histoDir_);
2160  mapHisto_[name+"VSEta"+nameSuffix_[i]] = new HResolution (name+"VSEta"+nameSuffix_[i], "resolution VS #eta"+titleSuffix[i], 100, -3, 3, -1, 1, histoDir_);
2161  mapHisto_[name+"VSPhi"+nameSuffix_[i]] = new HResolution (name+"VSPhi"+nameSuffix_[i], "resolution VS #phi"+titleSuffix[i], 100, -3.2, 3.2, -1, 1, histoDir_);
2162  }
2163 
2164  // single particles histograms
2165  muMinus.reset( new HDelta("muMinus") );
2166  muPlus.reset( new HDelta("muPlus") );
2167  }
std::unique_ptr< HDelta > muPlus
Definition: Histograms.h:2274
std::unique_ptr< HDelta > muMinus
Definition: Histograms.h:2273
TString nameSuffix_[2]
Definition: Histograms.h:2272
std::map< TString, TH1 * > mapHisto_
Definition: Histograms.h:2271
HMassResolutionVSPart::~HMassResolutionVSPart ( )
inlineoverride

Definition at line 2169 of file Histograms.h.

References trackerHits::histo.

2169  {
2170  for (std::map<TString, TH1*>::const_iterator histo=mapHisto_.begin();
2171  histo!=mapHisto_.end(); histo++) {
2172  delete (*histo).second;
2173  }
2174  }
std::map< TString, TH1 * > mapHisto_
Definition: Histograms.h:2271

Member Function Documentation

void HMassResolutionVSPart::Clear ( )
inlineoverride

Definition at line 2247 of file Histograms.h.

References trackerHits::histo.

2247  {
2248  for (std::map<TString, TH1*>::const_iterator histo=mapHisto_.begin();
2249  histo!=mapHisto_.end(); histo++) {
2250  (*histo).second->Clear();
2251  }
2252  muMinus->Clear();
2253  muPlus->Clear();
2254  }
std::unique_ptr< HDelta > muPlus
Definition: Histograms.h:2274
std::unique_ptr< HDelta > muMinus
Definition: Histograms.h:2273
std::map< TString, TH1 * > mapHisto_
Definition: Histograms.h:2271
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 2176 of file Histograms.h.

References HcalObjRepresent::Fill().

2180  {
2181  muMinus->Fill(recoP1, genP1);
2182  muPlus->Fill(recoP2, genP2);
2183 
2184  Fill( recoP1, charge1, recoP2, charge2, recoMass, genMass );
2185  }
std::unique_ptr< HDelta > muPlus
Definition: Histograms.h:2274
std::unique_ptr< HDelta > muMinus
Definition: Histograms.h:2273
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:2176
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 2187 of file Histograms.h.

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

2191  {
2192 
2193  if ( charge1 == charge2 ) std::cout << "Error: must get two opposite charge particles" << std::endl;
2194 
2195  double massRes = (recoMass - genMass)/genMass;
2196 
2197  reco::Particle::LorentzVector recoPair( recoP1 + recoP2 );
2198  double pairPt = recoPair.Pt();
2199 
2200  double recoPt[2] = {recoP1.Pt(), recoP2.Pt()};
2201  double recoEta[2] = {recoP1.Eta(), recoP2.Eta()};
2202  double recoPhi[2] = {recoP1.Phi(), recoP2.Phi()};
2203 
2204  // std::cout << "pairPt = " << pairPt << ", massRes = ("<<recoMass<<" - "<<genMass<<")/"<<genMass<<" = " << massRes
2205  // << ", recoPt[0] = " << recoPt[0] << ", recoPt[1] = " << recoPt[1] << std::endl;
2206 
2207  // Index of the histogram. If the muons have charge1 = -1 and charge2 = 1, they already have the
2208  // correct histogram indeces. Otherwise, swap the indeces.
2209  // In any case the negative muon has index = 0 and the positive muon has index = 1.
2210  int id[2] = {0,1};
2211  if ( charge1 > 0 ) {
2212  id[0] = 1;
2213  id[1] = 0;
2214  }
2215 
2216  double pairDeltaEta = fabs(recoEta[0] - recoEta[1]);
2217  double pairDeltaPhi = MuScleFitUtils::deltaPhi(recoPhi[0], recoPhi[1]);
2218 
2219  mapHisto_[name_]->Fill(massRes);
2220  mapHisto_[name_+"VSPairPt"]->Fill(pairPt, massRes);
2221  mapHisto_[name_+"VSPairDeltaEta"]->Fill(pairDeltaEta, massRes);
2222  mapHisto_[name_+"VSPairDeltaPhi"]->Fill(pairDeltaPhi, massRes);
2223 
2224  // Resolution vs single muon quantities
2225  // ------------------------------------
2226  int index = 0;
2227  for( int i=0; i<2; ++i ) {
2228  index = id[i];
2229  mapHisto_[name_+"VSPt"+nameSuffix_[i]]->Fill(recoPt[index], massRes); // EM [index] or [i]???
2230  mapHisto_[name_+"VSEta"+nameSuffix_[i]]->Fill(recoEta[index], massRes);
2231  mapHisto_[name_+"VSPhi"+nameSuffix_[i]]->Fill(recoPhi[index], massRes);
2232  }
2233  }
TString nameSuffix_[2]
Definition: Histograms.h:2272
std::map< TString, TH1 * > mapHisto_
Definition: Histograms.h:2271
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 2235 of file Histograms.h.

References trackerHits::histo.

2235  {
2236  histoDir_->cd();
2237  for (std::map<TString, TH1*>::const_iterator histo=mapHisto_.begin();
2238  histo!=mapHisto_.end(); histo++) {
2239  (*histo).second->Write();
2240  }
2241  // Create the new dir and cd into it
2242  (histoDir_->mkdir("singleMuonsVSgen"))->cd();
2243  muMinus->Write();
2244  muPlus->Write();
2245  }
std::unique_ptr< HDelta > muPlus
Definition: Histograms.h:2274
std::unique_ptr< HDelta > muMinus
Definition: Histograms.h:2273
std::map< TString, TH1 * > mapHisto_
Definition: Histograms.h:2271

Member Data Documentation

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

Definition at line 2271 of file Histograms.h.

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

Definition at line 2273 of file Histograms.h.

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

Definition at line 2274 of file Histograms.h.

TString HMassResolutionVSPart::nameSuffix_[2]
protected

Definition at line 2272 of file Histograms.h.