CMS 3D CMS Logo

HLTCATopTagFilter.h
Go to the documentation of this file.
1 #ifndef HLTCATopTagFilter_h
2 #define HLTCATopTagFilter_h
3 
4 // system include files
5 #include <memory>
6 #include <vector>
7 #include <sstream>
8 
9 
10 // user include files
26 #include <Math/VectorUtil.h>
27 
29  public:
30 
32  TopMass_(TopMass)
33  {}
34 
35  reco::CATopJetProperties operator()( reco::Jet const & ihardJet ) const;
36 
37  protected:
38  double TopMass_;
39 
40 };
41 
42 
43 
45  bool operator()( const edm::Ptr<reco::Candidate> & t1, const edm::Ptr<reco::Candidate> & t2 ) const {
46  return t1->pt() > t2->pt();
47  }
48 };
49 
50 
51 
52 //
53 // class declaration
54 //
55 
56 class HLTCATopTagFilter : public HLTFilter {
57  public:
58  explicit HLTCATopTagFilter(const edm::ParameterSet&);
59  ~HLTCATopTagFilter() override;
60  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
61  bool hltFilter( edm::Event&, const edm::EventSetup&, trigger::TriggerFilterObjectWithRefs & filterobject) const override;
62 
63  private:
64  // ----------member data ---------------------------
65 
70  double TopMass_;
71  double minTopMass_;
72  double maxTopMass_;
73  double minMinMass_;
74 
75 
76 };
77 
79  reco::CATopJetProperties properties;
80  // Get subjets
81  reco::Jet::Constituents subjets = ihardJet.getJetConstituents();
82  properties.nSubJets = subjets.size(); // number of subjets
83  properties.topMass = ihardJet.mass(); // jet mass
84  properties.minMass = 999999.; // minimum mass pairing
85 
86  // Require at least three subjets in all cases, if not, untagged
87  if ( properties.nSubJets >= 3 ) {
88 
89  // Take the highest 3 pt subjets for cuts
90  sort ( subjets.begin(), subjets.end(), GreaterByPtCandPtrUser() );
91 
92  // Now look at the subjets that were formed
93  for ( int isub = 0; isub < 2; ++isub ) {
94 
95  // Get this subjet
96  reco::Jet::Constituent icandJet = subjets[isub];
97 
98  // Now look at the "other" subjets than this one, form the minimum invariant mass
99  // pairing, as well as the "closest" combination to the W mass
100  for ( int jsub = isub + 1; jsub < 3; ++jsub ) {
101 
102  // Get the second subjet
103  reco::Jet::Constituent jcandJet = subjets[jsub];
104 
105  reco::Candidate::LorentzVector wCand = icandJet->p4() + jcandJet->p4();
106 
107  // Get the candidate mass
108  double imw = wCand.mass();
109 
110  // Find the minimum mass pairing.
111  if ( fabs( imw ) < properties.minMass ) {
112  properties.minMass = imw;
113  }
114  }// end second loop over subjets
115  }// end first loop over subjets
116  }// endif 3 subjets
117 
118  if (properties.minMass == 999999) properties.minMass=-1;
119 
120  return properties;
121 }
122 
123 #endif
Base class for all types of Jets.
Definition: Jet.h:20
reco::CATopJetProperties operator()(reco::Jet const &ihardJet) const
std::vector< Constituent > Constituents
Definition: Jet.h:23
virtual Constituents getJetConstituents() const
list of constituents
edm::InputTag src_
const edm::EDGetTokenT< reco::BasicJetCollection > inputToken_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
virtual double pt() const =0
transverse momentum
CATopJetHelperUser(double TopMass)
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
const edm::EDGetTokenT< reco::PFJetCollection > inputPFToken_
bool operator()(const edm::Ptr< reco::Candidate > &t1, const edm::Ptr< reco::Candidate > &t2) const
edm::InputTag pfsrc_
double mass() const final
mass