File indexing completed on 2024-04-06 12:18:34
0001
0002
0003
0004
0005
0006
0007
0008
0009 #include "HLTrigger/JetMET/interface/HLTMinDPhiMETFilter.h"
0010
0011 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0012 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0013 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0014
0015
0016 #include "DataFormats/Math/interface/deltaPhi.h"
0017
0018
0019 HLTMinDPhiMETFilter::HLTMinDPhiMETFilter(const edm::ParameterSet& iConfig)
0020 : HLTFilter(iConfig),
0021 usePt_(iConfig.getParameter<bool>("usePt")),
0022
0023 triggerType_(iConfig.getParameter<int>("triggerType")),
0024 maxNJets_(iConfig.getParameter<int>("maxNJets")),
0025 minPt_(iConfig.getParameter<double>("minPt")),
0026 maxEta_(iConfig.getParameter<double>("maxEta")),
0027 minDPhi_(iConfig.getParameter<double>("minDPhi")),
0028 metLabel_(iConfig.getParameter<edm::InputTag>("metLabel")),
0029 calometLabel_(iConfig.getParameter<edm::InputTag>("calometLabel")),
0030 jetsLabel_(iConfig.getParameter<edm::InputTag>("jetsLabel")) {
0031 m_theMETToken = consumes<reco::METCollection>(metLabel_);
0032 m_theCaloMETToken = consumes<reco::CaloMETCollection>(calometLabel_);
0033 m_theJetToken = consumes<reco::JetView>(jetsLabel_);
0034 }
0035
0036
0037 HLTMinDPhiMETFilter::~HLTMinDPhiMETFilter() = default;
0038
0039
0040 void HLTMinDPhiMETFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0041 edm::ParameterSetDescription desc;
0042 makeHLTFilterDescription(desc);
0043 desc.add<bool>("usePt", true);
0044
0045 desc.add<int>("triggerType", trigger::TriggerJet);
0046 desc.add<int>("maxNJets", 2);
0047 desc.add<double>("minPt", 30.);
0048 desc.add<double>("maxEta", 2.5);
0049 desc.add<double>("minDPhi", 0.5);
0050 desc.add<edm::InputTag>("metLabel", edm::InputTag("hltPFMETProducer"));
0051 desc.add<edm::InputTag>("calometLabel", edm::InputTag(""));
0052 desc.add<edm::InputTag>("jetsLabel", edm::InputTag("hltAK4PFJetL1FastL2L3Corrected"));
0053 descriptions.add("hltMinDPhiMETFilter", desc);
0054 }
0055
0056
0057 bool HLTMinDPhiMETFilter::hltFilter(edm::Event& iEvent,
0058 const edm::EventSetup& iSetup,
0059 trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0060
0061 if (saveTags())
0062 filterproduct.addCollectionTag(jetsLabel_);
0063
0064 bool usePFMET = (!metLabel_.label().empty()) || (calometLabel_.label().empty());
0065
0066 edm::Handle<reco::METCollection> mets;
0067 edm::Handle<reco::CaloMETCollection> calomets;
0068 if (usePFMET) {
0069 iEvent.getByToken(m_theMETToken, mets);
0070 } else {
0071 iEvent.getByToken(m_theMETToken, calomets);
0072 }
0073
0074 edm::Handle<reco::JetView> jets;
0075 iEvent.getByToken(m_theJetToken, jets);
0076
0077 double minDPhi = 3.141593;
0078 int nJets = 0;
0079
0080 if (!jets->empty() && ((usePFMET ? mets->size() : calomets->size()) > 0)) {
0081 double metphi = usePFMET ? mets->front().phi() : calomets->front().phi();
0082 for (reco::JetView::const_iterator j = jets->begin(); j != jets->end(); ++j) {
0083 if (nJets >= maxNJets_)
0084 break;
0085
0086 double pt = usePt_ ? j->pt() : j->et();
0087 double eta = j->eta();
0088 double phi = j->phi();
0089 if (pt > minPt_ && std::abs(eta) < maxEta_) {
0090 double dPhi = std::abs(reco::deltaPhi(metphi, phi));
0091 if (minDPhi > dPhi) {
0092 minDPhi = dPhi;
0093 }
0094
0095
0096
0097
0098 }
0099
0100 ++nJets;
0101 }
0102 }
0103
0104 bool accept(minDPhi > minDPhi_);
0105
0106 return accept;
0107 }