CMS 3D CMS Logo

HLTCAWZTagFilter.h
Go to the documentation of this file.
1 #ifndef HLTCAWZTagFilter_h
2 #define HLTCAWZTagFilter_h
3 
4 // system include files
5 #include <memory>
6 #include <vector>
7 #include <sstream>
8 
9 // user include files
23 #include <Math/VectorUtil.h>
24 
26  public:
27 
28  CAWZJetHelperUser(double massdropcut) :
29  massdropcut_(massdropcut)
30  {}
31 
32  reco::CATopJetProperties operator()( reco::Jet const & ihardJet ) const;
33 
34  protected:
35  double massdropcut_;
36 
37 };
38 
39 //
40 // class declaration
41 //
42 
43 class HLTCAWZTagFilter : public HLTFilter {
44  public:
45  explicit HLTCAWZTagFilter(const edm::ParameterSet&);
46  ~HLTCAWZTagFilter() override;
47  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
48  bool hltFilter( edm::Event&, const edm::EventSetup&, trigger::TriggerFilterObjectWithRefs & filterobject) const override;
49 
50 
51  private:
52  // ----------member data ---------------------------
53 
58  double minWMass_;
59  double maxWMass_;
60  double massdropcut_;
61 
62 };
63 
65  reco::CATopJetProperties properties;
66  // Get subjets
67  reco::Jet::Constituents subjets = ihardJet.getJetConstituents();
68  properties.nSubJets = subjets.size(); // number of subjets
69  properties.wMass = 999999.; // best W mass
70  properties.topMass = 999999.;
71  properties.minMass = -1;
72 
73  if (properties.nSubJets == 2) {
74 
75  sort ( subjets.begin(), subjets.end(), [](auto const& t1, auto const& t2){ return t1->pt() > t2->pt(); } );
76 
77  reco::Jet::Constituent icandJet = subjets[0];
78 
79  reco::Candidate::LorentzVector isubJet = icandJet->p4();
80  double imass = isubJet.mass();
81  double imw = ihardJet.mass();
82 
83  if (imass/imw < massdropcut_) {
84  // Get the candidate mass
85  properties.wMass = imw;
86  }
87  }
88 
89  return properties;
90 }
91 
92 #endif
edm::InputTag src_
CAWZJetHelperUser(double massdropcut)
Base class for all types of Jets.
Definition: Jet.h:20
std::vector< Constituent > Constituents
Definition: Jet.h:23
virtual Constituents getJetConstituents() const
list of constituents
reco::CATopJetProperties operator()(reco::Jet const &ihardJet) const
const edm::EDGetTokenT< reco::PFJetCollection > inputPFToken_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
const edm::EDGetTokenT< reco::BasicJetCollection > inputToken_
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
double mass() const final
mass
edm::InputTag pfsrc_