CMS 3D CMS Logo

HLTMinDPhiMETFilter.cc
Go to the documentation of this file.
1 
10 
14 //#include "DataFormats/HLTReco/interface/TriggerTypeDefs.h"
15 //#include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
17 
18 
19 // Constructor
21  usePt_ (iConfig.getParameter<bool>("usePt")),
22  //excludePFMuons_ (iConfig.getParameter<bool>("excludePFMuons")),
23  triggerType_ (iConfig.getParameter<int>("triggerType")),
24  maxNJets_ (iConfig.getParameter<int>("maxNJets")),
25  minPt_ (iConfig.getParameter<double>("minPt")),
26  maxEta_ (iConfig.getParameter<double>("maxEta")),
27  minDPhi_ (iConfig.getParameter<double>("minDPhi")),
28  metLabel_ (iConfig.getParameter<edm::InputTag>("metLabel")),
29  calometLabel_ (iConfig.getParameter<edm::InputTag>("calometLabel")),
30  jetsLabel_ (iConfig.getParameter<edm::InputTag>("jetsLabel")) {
31  m_theMETToken = consumes<reco::METCollection>(metLabel_);
32  m_theCaloMETToken = consumes<reco::CaloMETCollection>(calometLabel_);
33  m_theJetToken = consumes<reco::JetView>(jetsLabel_);
34 }
35 
36 // Destructor
38 
39 // Fill descriptions
43  desc.add<bool>("usePt", true);
44  //desc.add<bool>("excludePFMuons", false);
45  desc.add<int>("triggerType", trigger::TriggerJet);
46  desc.add<int>("maxNJets", 2);
47  desc.add<double>("minPt", 30.);
48  desc.add<double>("maxEta", 2.5);
49  desc.add<double>("minDPhi", 0.5);
50  desc.add<edm::InputTag>("metLabel", edm::InputTag("hltPFMETProducer"));
51  desc.add<edm::InputTag>("calometLabel", edm::InputTag(""));
52  desc.add<edm::InputTag>("jetsLabel", edm::InputTag("hltAK4PFJetL1FastL2L3Corrected"));
53  descriptions.add("hltMinDPhiMETFilter", desc);
54 }
55 
56 // Make filter decision
58 
59  // The filter object
60  if (saveTags()) filterproduct.addCollectionTag(jetsLabel_);
61 
62  bool usePFMET = (metLabel_.label() != "") || (calometLabel_.label() == "");
63 
66  if (usePFMET) {
67  iEvent.getByToken(m_theMETToken, mets);
68  } else {
69  iEvent.getByToken(m_theMETToken, calomets);
70  }
71 
72  edm::Handle<reco::JetView> jets; // assume to be sorted by pT
73  iEvent.getByToken(m_theJetToken, jets);
74 
75  double minDPhi = 3.141593;
76  int nJets = 0; // nJets counts all jets in the events, not only those that pass pt, eta requirements
77 
78  if (!jets->empty() &&
79  ((usePFMET ? mets->size() : calomets->size()) > 0) ) {
80  double metphi = usePFMET ? mets->front().phi() : calomets->front().phi();
81  for (reco::JetView::const_iterator j = jets->begin(); j != jets->end(); ++j) {
82  if (nJets >= maxNJets_)
83  break;
84 
85  double pt = usePt_ ? j->pt() : j->et();
86  double eta = j->eta();
87  double phi = j->phi();
88  if (pt > minPt_ && std::abs(eta) < maxEta_) {
89  double dPhi = std::abs(reco::deltaPhi(metphi, phi));
90  if (minDPhi > dPhi) {
91  minDPhi = dPhi;
92  }
93 
94  // Not sure what to save, since an event quantity is used
95  //reco::JetBaseRef ref(jets, distance(jets->begin(), j));
96  //filterproduct.addObject(triggerType_, ref);
97  }
98 
99  ++nJets;
100  }
101  }
102 
103  bool accept(minDPhi > minDPhi_);
104 
105  return accept;
106 }
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:22
double minDPhi_
Minium delta phi between a jet and MET.
int maxNJets_
Consider only n leading-pt (or et) jets, n = maxNJets_.
edm::InputTag metLabel_
Input jet, MET collections.
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
bool hltFilter(edm::Event &iEvent, const edm::EventSetup &iSetup, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:30
int iEvent
Definition: GenABIO.cc:224
const_iterator begin() const
~HLTMinDPhiMETFilter() override
vector< PseudoJet > jets
bool empty() const
edm::InputTag calometLabel_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:29
edm::EDGetTokenT< reco::JetView > m_theJetToken
double minPt_
Minimum pt requirement for jets.
edm::EDGetTokenT< reco::CaloMETCollection > m_theCaloMETToken
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
void add(std::string const &label, ParameterSetDescription const &psetDescription)
std::string const & label() const
Definition: InputTag.h:36
bool usePt_
Use pt; otherwise, use et.
bool saveTags() const
Definition: HLTFilter.h:45
HLT enums.
double maxEta_
Maximum (abs) eta requirement for jets.
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
edm::EDGetTokenT< reco::METCollection > m_theMETToken
HLTMinDPhiMETFilter(const edm::ParameterSet &iConfig)
const_iterator end() const