File indexing completed on 2024-04-06 12:18:34
0001
0002
0003
0004
0005
0006
0007
0008
0009 #include "HLTrigger/JetMET/interface/HLTMhtProducer.h"
0010
0011 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0012 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0013 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0014 #include "DataFormats/Common/interface/Handle.h"
0015
0016 #include "DataFormats/METReco/interface/MET.h"
0017 #include "DataFormats/METReco/interface/METFwd.h"
0018
0019
0020 HLTMhtProducer::HLTMhtProducer(const edm::ParameterSet& iConfig)
0021 : usePt_(iConfig.getParameter<bool>("usePt")),
0022 excludePFMuons_(iConfig.getParameter<bool>("excludePFMuons")),
0023 minNJet_(iConfig.getParameter<int>("minNJet")),
0024 minPtJet_(iConfig.getParameter<double>("minPtJet")),
0025 maxEtaJet_(iConfig.getParameter<double>("maxEtaJet")),
0026 jetsLabel_(iConfig.getParameter<edm::InputTag>("jetsLabel")),
0027 pfCandidatesLabel_(iConfig.getParameter<edm::InputTag>("pfCandidatesLabel")) {
0028 m_theJetToken = consumes<reco::CandidateView>(jetsLabel_);
0029 if (pfCandidatesLabel_.label().empty())
0030 excludePFMuons_ = false;
0031 if (excludePFMuons_)
0032 m_thePFCandidateToken = consumes<reco::PFCandidateCollection>(pfCandidatesLabel_);
0033
0034
0035 produces<reco::METCollection>();
0036 }
0037
0038
0039 HLTMhtProducer::~HLTMhtProducer() = default;
0040
0041
0042 void HLTMhtProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0043
0044 edm::ParameterSetDescription desc;
0045 desc.add<bool>("usePt", true);
0046 desc.add<bool>("excludePFMuons", false);
0047 desc.add<int>("minNJet", 0);
0048 desc.add<double>("minPtJet", 0.);
0049 desc.add<double>("maxEtaJet", 999.);
0050 desc.add<edm::InputTag>("jetsLabel", edm::InputTag("hltAntiKT4PFJets"));
0051 desc.add<edm::InputTag>("pfCandidatesLabel", edm::InputTag("hltParticleFlow"));
0052 descriptions.add("hltMhtProducer", desc);
0053 }
0054
0055
0056 void HLTMhtProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0057
0058 std::unique_ptr<reco::METCollection> result(new reco::METCollection());
0059
0060 edm::Handle<reco::CandidateView> jets;
0061 iEvent.getByToken(m_theJetToken, jets);
0062
0063 edm::Handle<reco::PFCandidateCollection> pfCandidates;
0064 if (excludePFMuons_)
0065 iEvent.getByToken(m_thePFCandidateToken, pfCandidates);
0066
0067 int nj = 0;
0068 double sumet = 0., mhx = 0., mhy = 0.;
0069
0070 for (auto const& aJet : *jets) {
0071 double const pt = usePt_ ? aJet.pt() : aJet.et();
0072 double const eta = aJet.eta();
0073 double const phi = aJet.phi();
0074 double const px = usePt_ ? aJet.px() : aJet.et() * cos(phi);
0075 double const py = usePt_ ? aJet.py() : aJet.et() * sin(phi);
0076
0077 if (pt > minPtJet_ && std::abs(eta) < maxEtaJet_) {
0078 mhx -= px;
0079 mhy -= py;
0080 sumet += pt;
0081 ++nj;
0082 }
0083 }
0084
0085 if (excludePFMuons_) {
0086 for (auto const& aCand : *pfCandidates) {
0087 if (std::abs(aCand.pdgId()) == 13) {
0088 mhx += aCand.px();
0089 mhy += aCand.py();
0090 }
0091 }
0092 }
0093
0094 if (nj < minNJet_) {
0095 sumet = 0;
0096 mhx = 0;
0097 mhy = 0;
0098 }
0099
0100 reco::MET::LorentzVector p4(mhx, mhy, 0, sqrt(mhx * mhx + mhy * mhy));
0101 reco::MET::Point vtx(0, 0, 0);
0102 reco::MET mht(sumet, p4, vtx);
0103 result->push_back(mht);
0104
0105
0106 iEvent.put(std::move(result));
0107 }