CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Protected Attributes
CATopJetHelper Class Reference

#include <CATopJetHelper.h>

Inheritance diagram for CATopJetHelper:

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.

20  :
21  TopMass_(TopMass), WMass_(WMass)
22  {}

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_.

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

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()().