CMS 3D CMS Logo

Public Member Functions | Protected Attributes

CATopJetHelper Class Reference

#include <CATopJetHelper.h>

Inheritance diagram for CATopJetHelper:
unary_function

List of all members.

Public Member Functions

 CATopJetHelper (double TopMass, double WMass)
reco::CATopJetProperties operator() (reco::Jet const &ihardJet) const

Protected Attributes

double TopMass_
double WMass_

Detailed Description

Definition at line 17 of file CATopJetHelper.h.


Constructor & Destructor Documentation

CATopJetHelper::CATopJetHelper ( double  TopMass,
double  WMass 
) [inline]

Definition at line 20 of file CATopJetHelper.h.

                                               :
  TopMass_(TopMass), WMass_(WMass)
  {}

Member Function Documentation

reco::CATopJetProperties CATopJetHelper::operator() ( reco::Jet const &  ihardJet) const

Definition at line 11 of file CATopJetHelper.cc.

References reco::Jet::getJetConstituents(), reco::LeafCandidate::mass(), reco::CATopJetProperties::minMass, reco::CATopJetProperties::nSubJets, python::multivaluedict::sort(), reco::CATopJetProperties::topMass, reco::CATopJetProperties::wMass, and WMass_.

                                                                                  {
  reco::CATopJetProperties properties;
  // Get subjets
  reco::Jet::Constituents subjets = ihardJet.getJetConstituents();
  properties.nSubJets = subjets.size();  // number of subjets
  properties.topMass = ihardJet.mass();      // jet mass
  properties.wMass = 99999.;                  // best W mass
  properties.minMass = 999999.;            // minimum mass pairing

  // Require at least three subjets in all cases, if not, untagged
  if ( properties.nSubJets >= 3 ) {

    // Take the highest 3 pt subjets for cuts
    sort ( subjets.begin(), subjets.end(), GreaterByPtCandPtr() );
       
    // Now look at the subjets that were formed
    for ( int isub = 0; isub < 2; ++isub ) {

      // Get this subjet
      reco::Jet::Constituent icandJet = subjets[isub];

      // Now look at the "other" subjets than this one, form the minimum invariant mass
      // pairing, as well as the "closest" combination to the W mass
      for ( int jsub = isub + 1; jsub < 3; ++jsub ) {

        // Get the second subjet
        reco::Jet::Constituent jcandJet = subjets[jsub];

        reco::Candidate::LorentzVector wCand = icandJet->p4() + jcandJet->p4();

        // Get the candidate mass
        double imw = wCand.mass();

        // Find the combination closest to the W mass
        if ( fabs( imw - WMass_ ) < fabs(properties.wMass - WMass_) ) {
          properties.wMass = imw;
        }
        // Find the minimum mass pairing. 
        if ( fabs( imw ) < properties.minMass ) {
          properties.minMass = imw;
        }  
      }// end second loop over subjets
    }// end first loop over subjets
  }// endif 3 subjets
 
  if (properties.minMass == 999999){properties.minMass=-1;}

  return properties;
}

Member Data Documentation

double CATopJetHelper::TopMass_ [protected]

Definition at line 27 of file CATopJetHelper.h.

double CATopJetHelper::WMass_ [protected]

Definition at line 28 of file CATopJetHelper.h.

Referenced by operator()().