CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

Calculate Class Reference

Helper class for the calculation of a top and a W boson mass estime. More...

#include <DQM/Physics/interface/TopDQMHelpers.h>

List of all members.

Public Member Functions

 Calculate (int maxNJets, double wMass)
 default constructor
double massTopQuark (const std::vector< reco::Jet > &jets)
 calculate W boson mass estimate
double massWBoson (const std::vector< reco::Jet > &jets)
 calculate W boson mass estimate
 ~Calculate ()
 default destructor

Private Member Functions

void operator() (const std::vector< reco::Jet > &jets)

Private Attributes

bool failed_
 indicate failed associations
double massTopQuark_
 cache of top quark mass estimate
double massWBoson_
 cache of w boson mass estimate
int maxNJets_
 max. number of jets to be considered
double wMass_
 paramater of the w boson mass

Detailed Description

Helper class for the calculation of a top and a W boson mass estime.

Helper class for the calculation of a top and a W boson mass estimate. The core implementation originates from the plugin TtSemiLepHypMaxSumPtWMass in TopQuarkAnalysis/TopJetCombination package. It may be extended to include b tag information.

Definition at line 64 of file TopDQMHelpers.h.


Constructor & Destructor Documentation

Calculate::Calculate ( int  maxNJets,
double  wMass 
)

default constructor

Definition at line 3 of file TopDQMHelpers.cc.

                                              : 
  failed_(false), maxNJets_(maxNJets), wMass_(wMass), massWBoson_(-1.), massTopQuark_(-1.)
{
}
Calculate::~Calculate ( ) [inline]

default destructor

Definition at line 69 of file TopDQMHelpers.h.

{};

Member Function Documentation

double Calculate::massTopQuark ( const std::vector< reco::Jet > &  jets)

calculate W boson mass estimate

Definition at line 15 of file TopDQMHelpers.cc.

References failed_, massTopQuark_, and operator()().

Referenced by TopSingleLepton::MonitorEnsemble::fill().

double Calculate::massWBoson ( const std::vector< reco::Jet > &  jets)

calculate W boson mass estimate

Definition at line 9 of file TopDQMHelpers.cc.

References failed_, massWBoson_, and operator()().

Referenced by TopSingleLepton::MonitorEnsemble::fill().

void Calculate::operator() ( const std::vector< reco::Jet > &  jets) [private]

do the calculation; this is called only once per event by the first function call to return a mass estimate. The once calculated values are cached afterwards

Definition at line 21 of file TopDQMHelpers.cc.

References failed_, massTopQuark_, massWBoson_, maxNJets_, and wMass_.

Referenced by massTopQuark(), and massWBoson().

{
  if(maxNJets_<0) maxNJets_=jets.size();
  failed_= jets.size()<(unsigned int) maxNJets_;
  if( failed_){ return; }

  // associate those jets with maximum pt of the vectorial 
  // sum to the hadronic decay chain
  double maxPt=-1.;
  std::vector<int> maxPtIndices;
  maxPtIndices.push_back(-1);
  maxPtIndices.push_back(-1);
  maxPtIndices.push_back(-1);
  for(int idx=0; idx<maxNJets_; ++idx){
    for(int jdx=0; jdx<maxNJets_; ++jdx){ if(jdx<=idx) continue;
      for(int kdx=0; kdx<maxNJets_; ++kdx){ if(kdx==idx || kdx==jdx) continue;
        reco::Particle::LorentzVector sum = jets[idx].p4()+jets[jdx].p4()+jets[kdx].p4();
        if( maxPt<0. || maxPt<sum.pt() ){
          maxPt=sum.pt();
          maxPtIndices.clear();
          maxPtIndices.push_back(idx);
          maxPtIndices.push_back(jdx);
          maxPtIndices.push_back(kdx);
        }
      }
    }
  }
  massTopQuark_= (jets[maxPtIndices[0]].p4()+
                  jets[maxPtIndices[1]].p4()+
                  jets[maxPtIndices[2]].p4()).mass();

  // associate those jets that get closest to the W mass
  // with their invariant mass to the W boson
  double wDist =-1.;
  std::vector<int> wMassIndices;
  wMassIndices.push_back(-1);
  wMassIndices.push_back(-1);
  for(unsigned idx=0; idx<maxPtIndices.size(); ++idx){  
    for(unsigned jdx=0; jdx<maxPtIndices.size(); ++jdx){  
      if( jdx==idx || maxPtIndices[idx]>maxPtIndices[jdx] ) continue;
        reco::Particle::LorentzVector sum = jets[maxPtIndices[idx]].p4()+jets[maxPtIndices[jdx]].p4();
        if( wDist<0. || wDist>fabs(sum.mass()-wMass_) ){
          wDist=fabs(sum.mass()-wMass_);
          wMassIndices.clear();
          wMassIndices.push_back(maxPtIndices[idx]);
          wMassIndices.push_back(maxPtIndices[jdx]);
        }
    }
  }
  massWBoson_= (jets[wMassIndices[0]].p4()+
                jets[wMassIndices[1]].p4()).mass();
}

Member Data Documentation

bool Calculate::failed_ [private]

indicate failed associations

Definition at line 84 of file TopDQMHelpers.h.

Referenced by massTopQuark(), massWBoson(), and operator()().

double Calculate::massTopQuark_ [private]

cache of top quark mass estimate

Definition at line 92 of file TopDQMHelpers.h.

Referenced by massTopQuark(), and operator()().

double Calculate::massWBoson_ [private]

cache of w boson mass estimate

Definition at line 90 of file TopDQMHelpers.h.

Referenced by massWBoson(), and operator()().

int Calculate::maxNJets_ [private]

max. number of jets to be considered

Definition at line 86 of file TopDQMHelpers.h.

Referenced by operator()().

double Calculate::wMass_ [private]

paramater of the w boson mass

Definition at line 88 of file TopDQMHelpers.h.

Referenced by operator()().