CMS 3D CMS Logo

Public Member Functions | Protected Attributes

HMassResolutionVSPart Class Reference

#include <Histograms.h>

Inheritance diagram for HMassResolutionVSPart:
Histograms

List of all members.

Public Member Functions

virtual void Clear ()
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 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)
 HMassResolutionVSPart (TFile *outputFile, const TString &name)
virtual void Write ()
 ~HMassResolutionVSPart ()

Protected Attributes

std::map< TString, TH1 * > mapHisto_
std::auto_ptr< HDeltamuMinus
std::auto_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 1862 of file Histograms.h.


Constructor & Destructor Documentation

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

Definition at line 1865 of file Histograms.h.

References Histograms::histoDir_, i, mapHisto_, muMinus, muPlus, AlCaRecoCosmics_cfg::name, and nameSuffix_.

                                                                  : Histograms( outputFile, name ) {
    // Kinematical variables
    nameSuffix_[0] = "Plus";
    nameSuffix_[1] = "Minus";
    TString titleSuffix[] = {" for mu+", " for mu-"};

    mapHisto_[name]                  = new TH1F (name, "#Delta M/M", 4000, -1, 1);
    mapHisto_[name+"VSPairPt"]       = new HResolution (name+"VSPairPt", "resolution VS pt of the pair", 100, 0, 200, -1, 1, histoDir_);
    mapHisto_[name+"VSPairDeltaEta"] = new HResolution (name+"VSPairDeltaEta", "resolution VS #Delta#eta of the pair", 100, -0.1, 6.2, -1, 1, histoDir_);
    mapHisto_[name+"VSPairDeltaPhi"] = new HResolution (name+"VSPairDeltaPhi", "resolution VS #Delta#phi of the pair", 100, -0.1, 3.2, -1, 1, histoDir_);

    for( int i=0; i<2; ++i ) {
      mapHisto_[name+"VSPt"+nameSuffix_[i]]  = new HResolution (name+"VSPt"+nameSuffix_[i], "resolution VS pt"+titleSuffix[i], 100, 0, 200, -1, 1, histoDir_);
      mapHisto_[name+"VSEta"+nameSuffix_[i]] = new HResolution (name+"VSEta"+nameSuffix_[i], "resolution VS #eta"+titleSuffix[i], 100, -3, 3, -1, 1, histoDir_);
      mapHisto_[name+"VSPhi"+nameSuffix_[i]] = new HResolution (name+"VSPhi"+nameSuffix_[i], "resolution VS #phi"+titleSuffix[i], 100, -3.2, 3.2, -1, 1, histoDir_);
    }

    // single particles histograms
    muMinus.reset( new HDelta("muMinus") );
    muPlus.reset( new HDelta("muPlus") );
  }
HMassResolutionVSPart::~HMassResolutionVSPart ( ) [inline]

Definition at line 1887 of file Histograms.h.

References trackerHits::histo, and mapHisto_.

                          {
    for (std::map<TString, TH1*>::const_iterator histo=mapHisto_.begin(); 
         histo!=mapHisto_.end(); histo++) {
      delete (*histo).second;
    }
  }

Member Function Documentation

virtual void HMassResolutionVSPart::Clear ( ) [inline, virtual]

Implements Histograms.

Definition at line 1965 of file Histograms.h.

References trackerHits::histo, mapHisto_, muMinus, and muPlus.

                       {
    for (std::map<TString, TH1*>::const_iterator histo=mapHisto_.begin(); 
         histo!=mapHisto_.end(); histo++) {
      (*histo).second->Clear();
    }
    muMinus->Clear();
    muPlus->Clear();
  }
virtual 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 
) [inline, virtual]

Reimplemented from Histograms.

Definition at line 1905 of file Histograms.h.

References gather_cfg::cout, Geom::deltaPhi(), i, getHLTprescales::index, mapHisto_, Histograms::name_, and nameSuffix_.

                                                                       {

    if ( charge1 == charge2 ) std::cout << "Error: must get two opposite charge particles" << std::endl;

    double massRes = (recoMass - genMass)/genMass;

    reco::Particle::LorentzVector recoPair( recoP1 + recoP2 );
    double pairPt = recoPair.Pt();

    double recoPt[2]  = {recoP1.Pt(),  recoP2.Pt()};
    double recoEta[2] = {recoP1.Eta(), recoP2.Eta()};
    double recoPhi[2] = {recoP1.Phi(), recoP2.Phi()};

    // std::cout << "pairPt = " << pairPt << ", massRes = ("<<recoMass<<" - "<<genMass<<")/"<<genMass<<" = " << massRes
    //      << ", recoPt[0] = " << recoPt[0] << ", recoPt[1] = " << recoPt[1] << std::endl;

    // Index of the histogram. If the muons have charge1 = -1 and charge2 = 1, they already have the
    // correct histogram indeces. Otherwise, swap the indeces.
    // In any case the negative muon has index = 0  and the positive muon has index = 1.
    int id[2] = {0,1};
    if ( charge1 > 0 ) {
      id[0] = 1;
      id[1] = 0;
    }

    double pairDeltaEta = fabs(recoEta[0] - recoEta[1]);
    double pairDeltaPhi = MuScleFitUtils::deltaPhi(recoPhi[0], recoPhi[1]);

    mapHisto_[name_]->Fill(massRes);
    mapHisto_[name_+"VSPairPt"]->Fill(pairPt, massRes);
    mapHisto_[name_+"VSPairDeltaEta"]->Fill(pairDeltaEta, massRes);
    mapHisto_[name_+"VSPairDeltaPhi"]->Fill(pairDeltaPhi, massRes);

    // Resolution vs single muon quantities
    // ------------------------------------
    int index = 0;
    for( int i=0; i<2; ++i ) {
      index = id[i];
      mapHisto_[name_+"VSPt"+nameSuffix_[i]]->Fill(recoPt[i], massRes);
      mapHisto_[name_+"VSEta"+nameSuffix_[i]]->Fill(recoEta[i], massRes);
      mapHisto_[name_+"VSPhi"+nameSuffix_[i]]->Fill(recoPhi[i], massRes);
    }
  } 
virtual 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 
) [inline, virtual]

Reimplemented from Histograms.

Definition at line 1894 of file Histograms.h.

References muMinus, and muPlus.

                                                                       {
    muMinus->Fill(recoP1, genP1);
    muPlus->Fill(recoP2, genP2);

    Fill( recoP1, charge1, recoP2, charge2, recoMass, genMass );
  }
virtual void HMassResolutionVSPart::Write ( ) [inline, virtual]

Implements Histograms.

Definition at line 1953 of file Histograms.h.

References trackerHits::histo, Histograms::histoDir_, mapHisto_, muMinus, and muPlus.

                       {
    histoDir_->cd();
    for (std::map<TString, TH1*>::const_iterator histo=mapHisto_.begin(); 
         histo!=mapHisto_.end(); histo++) {
      (*histo).second->Write();
    }
    // Create the new dir and cd into it
    (histoDir_->mkdir("singleMuonsVSgen"))->cd();
    muMinus->Write();
    muPlus->Write();
  }

Member Data Documentation

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

Definition at line 1989 of file Histograms.h.

Referenced by Clear(), Fill(), HMassResolutionVSPart(), Write(), and ~HMassResolutionVSPart().

std::auto_ptr<HDelta> HMassResolutionVSPart::muMinus [protected]

Definition at line 1991 of file Histograms.h.

Referenced by Clear(), Fill(), HMassResolutionVSPart(), and Write().

std::auto_ptr<HDelta> HMassResolutionVSPart::muPlus [protected]

Definition at line 1992 of file Histograms.h.

Referenced by Clear(), Fill(), HMassResolutionVSPart(), and Write().

TString HMassResolutionVSPart::nameSuffix_[2] [protected]

Definition at line 1990 of file Histograms.h.

Referenced by Fill(), and HMassResolutionVSPart().