File indexing completed on 2024-04-06 12:18:34
0001
0002
0003
0004
0005
0006
0007
0008
0009 #include "HLTrigger/JetMET/interface/HLTPFJetIDProducer.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
0017 HLTPFJetIDProducer::HLTPFJetIDProducer(const edm::ParameterSet& iConfig)
0018 : minPt_(iConfig.getParameter<double>("minPt")),
0019 maxEta_(iConfig.getParameter<double>("maxEta")),
0020 CHF_(iConfig.getParameter<double>("CHF")),
0021 NHF_(iConfig.getParameter<double>("NHF")),
0022 CEF_(iConfig.getParameter<double>("CEF")),
0023 NEF_(iConfig.getParameter<double>("NEF")),
0024 maxCF_(iConfig.getParameter<double>("maxCF")),
0025 NCH_(iConfig.getParameter<int>("NCH")),
0026 NTOT_(iConfig.getParameter<int>("NTOT")),
0027 inputTag_(iConfig.getParameter<edm::InputTag>("jetsInput")) {
0028 m_thePFJetToken = consumes<reco::PFJetCollection>(inputTag_);
0029
0030
0031 produces<reco::PFJetCollection>();
0032 }
0033
0034
0035 HLTPFJetIDProducer::~HLTPFJetIDProducer() = default;
0036
0037
0038 void HLTPFJetIDProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0039 edm::ParameterSetDescription desc;
0040 desc.add<double>("minPt", 20.);
0041 desc.add<double>("maxEta", 1e99);
0042 desc.add<double>("CHF", -99.);
0043 desc.add<double>("NHF", 99.);
0044 desc.add<double>("CEF", 99.);
0045 desc.add<double>("NEF", 99.);
0046 desc.add<double>("maxCF", 99.);
0047 desc.add<int>("NCH", 0);
0048 desc.add<int>("NTOT", 0);
0049 desc.add<edm::InputTag>("jetsInput", edm::InputTag("hltAntiKT4PFJets"));
0050 descriptions.add("hltPFJetIDProducer", desc);
0051 }
0052
0053
0054 void HLTPFJetIDProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0055
0056 std::unique_ptr<reco::PFJetCollection> result(new reco::PFJetCollection());
0057
0058 edm::Handle<reco::PFJetCollection> pfjets;
0059 iEvent.getByToken(m_thePFJetToken, pfjets);
0060
0061 for (auto const& j : *pfjets) {
0062 bool pass = false;
0063 double pt = j.pt();
0064 double eta = j.eta();
0065 double abseta = std::abs(eta);
0066
0067 if (!(pt > 0.))
0068 continue;
0069
0070 if (pt < minPt_) {
0071 pass = true;
0072
0073 } else if (abseta >= maxEta_) {
0074 pass = true;
0075
0076 } else {
0077 double chf = j.chargedHadronEnergyFraction();
0078
0079 double nhf = j.neutralHadronEnergyFraction();
0080 double cef = j.chargedEmEnergyFraction();
0081 double nef = j.neutralEmEnergyFraction();
0082 double cftot = chf + cef + j.chargedMuEnergyFraction();
0083 int nch = j.chargedMultiplicity();
0084 int ntot = j.numberOfDaughters();
0085
0086 pass = true;
0087 pass = pass && (ntot > NTOT_);
0088 pass = pass && (nef < NEF_);
0089 pass = pass && (nhf < NHF_ || abseta >= 2.4);
0090 pass = pass && (cef < CEF_ || abseta >= 2.4);
0091 pass = pass && (chf > CHF_ || abseta >= 2.4);
0092 pass = pass && (nch > NCH_ || abseta >= 2.4);
0093 pass = pass && (cftot < maxCF_ || abseta >= 2.4);
0094 }
0095
0096 if (pass)
0097 result->push_back(j);
0098 }
0099
0100
0101 iEvent.put(std::move(result));
0102 }