File indexing completed on 2024-04-06 12:18:34
0001
0002
0003
0004
0005
0006
0007
0008
0009 #include "HLTrigger/JetMET/interface/HLTMETCleanerUsingJetID.h"
0010 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0011
0012
0013 HLTMETCleanerUsingJetID::HLTMETCleanerUsingJetID(const edm::ParameterSet& iConfig)
0014 : minPt_(iConfig.getParameter<double>("minPt")),
0015 maxEta_(iConfig.getParameter<double>("maxEta")),
0016 metLabel_(iConfig.getParameter<edm::InputTag>("metLabel")),
0017 jetsLabel_(iConfig.getParameter<edm::InputTag>("jetsLabel")),
0018 goodJetsLabel_(iConfig.getParameter<edm::InputTag>("goodJetsLabel")) {
0019 m_theMETToken = consumes<reco::CaloMETCollection>(metLabel_);
0020 m_theJetToken = consumes<reco::CaloJetCollection>(jetsLabel_);
0021 m_theGoodJetToken = consumes<reco::CaloJetCollection>(goodJetsLabel_);
0022
0023
0024 produces<reco::CaloMETCollection>();
0025 }
0026
0027
0028 HLTMETCleanerUsingJetID::~HLTMETCleanerUsingJetID() = default;
0029
0030
0031 void HLTMETCleanerUsingJetID::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0032 edm::ParameterSetDescription desc;
0033 desc.add<double>("minPt", 20.);
0034 desc.add<double>("maxEta", 5.);
0035 desc.add<edm::InputTag>("metLabel", edm::InputTag("hltMet"));
0036 desc.add<edm::InputTag>("jetsLabel", edm::InputTag("hltAntiKT4CaloJets"));
0037 desc.add<edm::InputTag>("goodJetsLabel", edm::InputTag("hltCaloJetIDPassed"));
0038 descriptions.add("hltMETCleanerUsingJetID", desc);
0039 }
0040
0041
0042 void HLTMETCleanerUsingJetID::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0043
0044 std::unique_ptr<reco::CaloMETCollection> result(new reco::CaloMETCollection);
0045
0046 edm::Handle<reco::CaloMETCollection> met;
0047 edm::Handle<reco::CaloJetCollection> jets;
0048 edm::Handle<reco::CaloJetCollection> goodJets;
0049
0050 iEvent.getByToken(m_theMETToken, met);
0051 iEvent.getByToken(m_theJetToken, jets);
0052 iEvent.getByToken(m_theGoodJetToken, goodJets);
0053
0054 double mex_jets = 0.;
0055 double mey_jets = 0.;
0056
0057 if (!jets->empty()) {
0058 for (auto const& j : *jets) {
0059 double pt = j.pt();
0060 double eta = j.eta();
0061 double px = j.px();
0062 double py = j.py();
0063
0064 if (pt > minPt_ && std::abs(eta) < maxEta_) {
0065 mex_jets -= px;
0066 mey_jets -= py;
0067
0068 }
0069 }
0070 }
0071
0072 double mex_goodJets = 0.;
0073 double mey_goodJets = 0.;
0074
0075 if (!goodJets->empty()) {
0076 for (auto const& j : *goodJets) {
0077 double pt = j.pt();
0078 double eta = j.eta();
0079 double px = j.px();
0080 double py = j.py();
0081
0082 if (pt > minPt_ && std::abs(eta) < maxEta_) {
0083 mex_goodJets -= px;
0084 mey_goodJets -= py;
0085
0086 }
0087 }
0088 }
0089
0090 if (!met->empty()) {
0091 double mex_diff = mex_goodJets - mex_jets;
0092 double mey_diff = mey_goodJets - mey_jets;
0093
0094 reco::Candidate::LorentzVector p4_clean(met->front().px() + mex_diff,
0095 mey_diff + met->front().py(),
0096 0,
0097 sqrt((met->front().px() + mex_diff) * (met->front().px() + mex_diff) +
0098 (met->front().py() + mey_diff) * (met->front().py() + mey_diff)));
0099
0100 reco::CaloMET cleanmet = met->front();
0101 cleanmet.setP4(p4_clean);
0102 result->push_back(cleanmet);
0103 }
0104
0105 iEvent.put(std::move(result));
0106 }