CMS 3D CMS Logo

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

#include <CATopJetHelper.h>

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 16 of file CATopJetHelper.h.

Constructor & Destructor Documentation

◆ CATopJetHelper()

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

Member Function Documentation

◆ operator()()

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

Definition at line 9 of file CATopJetHelper.cc.

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

9  {
10  reco::CATopJetProperties properties;
11  // Get subjets
12  reco::Jet::Constituents subjets = ihardJet.getJetConstituents();
13  properties.nSubJets = subjets.size(); // number of subjets
14  properties.topMass = ihardJet.mass(); // jet mass
15  properties.wMass = 99999.; // best W mass
16  properties.minMass = 999999.; // minimum mass pairing
17 
18  // Require at least three subjets in all cases, if not, untagged
19  if (properties.nSubJets >= 3) {
20  // Take the highest 3 pt subjets for cuts
21  sort(subjets.begin(), subjets.end(), GreaterByPtCandPtr());
22 
23  // Now look at the subjets that were formed
24  for (int isub = 0; isub < 2; ++isub) {
25  // Get this subjet
26  reco::Jet::Constituent icandJet = subjets[isub];
27 
28  // Now look at the "other" subjets than this one, form the minimum invariant mass
29  // pairing, as well as the "closest" combination to the W mass
30  for (int jsub = isub + 1; jsub < 3; ++jsub) {
31  // Get the second subjet
32  reco::Jet::Constituent jcandJet = subjets[jsub];
33 
34  reco::Candidate::LorentzVector wCand = icandJet->p4() + jcandJet->p4();
35 
36  // Get the candidate mass
37  double imw = wCand.mass();
38 
39  // Find the combination closest to the W mass
40  if (fabs(imw - WMass_) < fabs(properties.wMass - WMass_)) {
41  properties.wMass = imw;
42  }
43  // Find the minimum mass pairing.
44  if (fabs(imw) < properties.minMass) {
45  properties.minMass = imw;
46  }
47  } // end second loop over subjets
48  } // end first loop over subjets
49  } // endif 3 subjets
50 
51  if (properties.minMass == 999999) {
52  properties.minMass = -1;
53  }
54 
55  return properties;
56 }
std::vector< Constituent > Constituents
Definition: Jet.h:23
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:36

Member Data Documentation

◆ TopMass_

double CATopJetHelper::TopMass_
protected

Definition at line 23 of file CATopJetHelper.h.

◆ WMass_

double CATopJetHelper::WMass_
protected

Definition at line 24 of file CATopJetHelper.h.

Referenced by operator()().