File indexing completed on 2024-04-06 12:18:31
0001 #ifndef HLTrigger_JetMET_plugins_HLTJetHFCleaner_h
0002 #define HLTrigger_JetMET_plugins_HLTJetHFCleaner_h
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033
0034
0035
0036 #include <memory>
0037 #include <vector>
0038
0039 #include "FWCore/Framework/interface/Frameworkfwd.h"
0040 #include "FWCore/Framework/interface/global/EDProducer.h"
0041 #include "FWCore/Framework/interface/Event.h"
0042 #include "FWCore/Framework/interface/MakerMacros.h"
0043 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0044 #include "FWCore/Utilities/interface/StreamID.h"
0045 #include "FWCore/Utilities/interface/Exception.h"
0046 #include "DataFormats/Common/interface/Ref.h"
0047 #include "DataFormats/JetReco/interface/Jet.h"
0048 #include "DataFormats/Math/interface/deltaPhi.h"
0049 #include "DataFormats/METReco/interface/MET.h"
0050
0051 template <typename JetType>
0052 class HLTJetHFCleaner : public edm::global::EDProducer<> {
0053 public:
0054 typedef std::vector<JetType> JetCollection;
0055 typedef edm::Ref<JetCollection> JetRef;
0056 explicit HLTJetHFCleaner(const edm::ParameterSet&);
0057 ~HLTJetHFCleaner() override = default;
0058
0059 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0060
0061 private:
0062 void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
0063
0064 const edm::EDGetTokenT<std::vector<JetType>> jetToken_;
0065 const edm::EDGetTokenT<edm::View<reco::MET>> metToken_;
0066
0067 const edm::EDGetTokenT<edm::ValueMap<float>> sigmaEtaEtaToken_;
0068 const edm::EDGetTokenT<edm::ValueMap<float>> sigmaPhiPhiToken_;
0069 const edm::EDGetTokenT<edm::ValueMap<int>> centralEtaStripSizeToken_;
0070
0071
0072 const float jetPtMin_;
0073 const float dphiJetMetMin_;
0074
0075
0076 const float jetEtaMin_;
0077 const float jetEtaMax_;
0078
0079 const float sigmaEtaPhiDiffMax_;
0080 const float centralEtaStripSizeMax_;
0081 const float cornerCutSigmaEtaEta_;
0082 const float cornerCutSigmaPhiPhi_;
0083 const bool applySigmaEtaPhiCornerCut_;
0084 const bool applySigmaEtaPhiCut_;
0085 const bool applyStripSizeCut_;
0086 };
0087
0088 template <typename JetType>
0089 HLTJetHFCleaner<JetType>::HLTJetHFCleaner(const edm::ParameterSet& iConfig)
0090 : jetToken_(consumes(iConfig.getParameter<edm::InputTag>("jets"))),
0091 metToken_(consumes(iConfig.getParameter<edm::InputTag>("mets"))),
0092 sigmaEtaEtaToken_(consumes(iConfig.getParameter<edm::InputTag>("sigmaEtaEta"))),
0093 sigmaPhiPhiToken_(consumes(iConfig.getParameter<edm::InputTag>("sigmaPhiPhi"))),
0094 centralEtaStripSizeToken_(consumes(iConfig.getParameter<edm::InputTag>("centralEtaStripSize"))),
0095 jetPtMin_(iConfig.getParameter<double>("jetPtMin")),
0096 dphiJetMetMin_(iConfig.getParameter<double>("dphiJetMetMin")),
0097 jetEtaMin_(iConfig.getParameter<double>("jetEtaMin")),
0098 jetEtaMax_(iConfig.getParameter<double>("jetEtaMax")),
0099 sigmaEtaPhiDiffMax_(iConfig.getParameter<double>("sigmaEtaPhiDiffMax")),
0100 centralEtaStripSizeMax_(iConfig.getParameter<int>("centralEtaStripSizeMax")),
0101 cornerCutSigmaEtaEta_(iConfig.getParameter<double>("cornerCutSigmaEtaEta")),
0102 cornerCutSigmaPhiPhi_(iConfig.getParameter<double>("cornerCutSigmaPhiPhi")),
0103 applySigmaEtaPhiCornerCut_(iConfig.getParameter<bool>("applySigmaEtaPhiCornerCut")),
0104 applySigmaEtaPhiCut_(iConfig.getParameter<bool>("applySigmaEtaPhiCut")),
0105 applyStripSizeCut_(iConfig.getParameter<bool>("applyStripSizeCut")) {
0106 produces<JetCollection>();
0107 }
0108
0109 template <typename JetType>
0110 void HLTJetHFCleaner<JetType>::produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup&) const {
0111 auto const jets = iEvent.getHandle(jetToken_);
0112 auto const& mets = iEvent.get(metToken_);
0113
0114
0115 auto cleaned_jets = std::make_unique<JetCollection>();
0116 cleaned_jets->reserve(jets->size());
0117
0118
0119 if (mets.size() != 1) {
0120 throw cms::Exception("ComparisonFailure") << "Size of reco::MET collection is different from 1";
0121 }
0122 auto const met_phi = mets.front().phi();
0123
0124
0125
0126
0127 auto const& sigmaEtaEtas = iEvent.get(sigmaEtaEtaToken_);
0128 auto const& sigmaPhiPhis = iEvent.get(sigmaPhiPhiToken_);
0129 auto const& centralEtaStripSizes = iEvent.get(centralEtaStripSizeToken_);
0130
0131
0132 for (uint ijet = 0; ijet < jets->size(); ijet++) {
0133 auto const jet = (*jets)[ijet];
0134 auto const withinEtaRange = (std::abs(jet.eta()) < jetEtaMax_) and (std::abs(jet.eta()) > jetEtaMin_);
0135
0136 if (std::abs(reco::deltaPhi(jet.phi(), met_phi)) <= dphiJetMetMin_) {
0137 cleaned_jets->emplace_back(jet);
0138 continue;
0139 }
0140
0141
0142 if ((jet.pt() < jetPtMin_) or (!withinEtaRange)) {
0143 cleaned_jets->emplace_back(jet);
0144 continue;
0145 }
0146
0147 JetRef const jetRef(jets, ijet);
0148
0149
0150
0151
0152 if (applySigmaEtaPhiCut_) {
0153
0154 auto const sigmaEtaEta = sigmaEtaEtas[jetRef];
0155 auto const sigmaPhiPhi = sigmaPhiPhis[jetRef];
0156
0157
0158
0159 if ((sigmaEtaEta < 0) and (sigmaPhiPhi < 0)) {
0160 cleaned_jets->emplace_back(jet);
0161 continue;
0162 }
0163
0164 auto passSigmaEtaPhiCut = (sigmaEtaEta - sigmaPhiPhi) <= sigmaEtaPhiDiffMax_;
0165
0166 if (applySigmaEtaPhiCornerCut_) {
0167 auto const inCorner = (sigmaEtaEta < cornerCutSigmaEtaEta_) and (sigmaPhiPhi < cornerCutSigmaPhiPhi_);
0168 passSigmaEtaPhiCut &= !inCorner;
0169 }
0170 if (!passSigmaEtaPhiCut) {
0171 continue;
0172 }
0173 }
0174
0175
0176 if (applyStripSizeCut_ and (centralEtaStripSizes[jetRef] > centralEtaStripSizeMax_)) {
0177 continue;
0178 }
0179
0180 cleaned_jets->emplace_back(jet);
0181 }
0182
0183 iEvent.put(std::move(cleaned_jets));
0184 }
0185
0186
0187 template <typename JetType>
0188 void HLTJetHFCleaner<JetType>::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0189 edm::ParameterSetDescription desc;
0190 desc.add<edm::InputTag>("jets", edm::InputTag("hltAK4PFJetsTightIDCorrected"))->setComment("Input jet collection.");
0191 desc.add<edm::InputTag>("mets", edm::InputTag("hltMet"))->setComment("Input MET collection.");
0192 desc.add<edm::InputTag>("sigmaEtaEta", edm::InputTag("hltHFJetShowerShape", "sigmaEtaEta"))
0193 ->setComment("Input collection which stores the sigmaEtaEta values per jet.");
0194 desc.add<edm::InputTag>("sigmaPhiPhi", edm::InputTag("hltHFJetShowerShape", "sigmaPhiPhi"))
0195 ->setComment("Input collection which stores the sigmaPhiPhis values per jet.");
0196 desc.add<edm::InputTag>("centralEtaStripSize", edm::InputTag("hltHFJetShowerShape", "centralEtaStripSize"))
0197 ->setComment("Input collection which stores the central strip size values per jet.");
0198 desc.add<double>("jetPtMin", 100)->setComment("The minimum pt value for a jet such that the cuts will be applied.");
0199 desc.add<double>("dphiJetMetMin", 2.5)
0200 ->setComment(
0201 "The minimum value for deltaPhi between jet and MET, such that the cuts will be applied to a given jet.");
0202 desc.add<double>("jetEtaMin", 2.9)->setComment("Minimum value of jet |eta| for which the cuts will be applied.");
0203 desc.add<double>("jetEtaMax", 5.0)->setComment("Maximum value of jet |eta| for which the cuts will be applied.");
0204 desc.add<double>("sigmaEtaPhiDiffMax", 0.05)
0205 ->setComment("Determines the threshold in the following cut: sigmaEtaEta-sigmaPhiPhi <= X");
0206 desc.add<double>("cornerCutSigmaEtaEta", 0.02)
0207 ->setComment(
0208 "Corner cut value for sigmaEtaEta. Jets will be cut if both sigmaEtaEta and sigmaPhiPhi are lower than the "
0209 "corner value.");
0210 desc.add<double>("cornerCutSigmaPhiPhi", 0.02)
0211 ->setComment(
0212 "Corner cut value for sigmaPhiPhi. Jets will be cut if both sigmaEtaEta and sigmaPhiPhi are lower than the "
0213 "corner value.");
0214 desc.add<int>("centralEtaStripSizeMax", 2)
0215 ->setComment("Determines the threshold in the following cut: centralEtaStripSize <= X");
0216 desc.add<bool>("applySigmaEtaPhiCornerCut", true)->setComment("Boolean specifying whether to apply the corner cut.");
0217 desc.add<bool>("applySigmaEtaPhiCut", true)
0218 ->setComment("Boolean specifying whether to apply the sigmaEtaEta-sigmaPhiPhi cut.");
0219 desc.add<bool>("applyStripSizeCut", true)
0220 ->setComment("Boolean specifying whether to apply the centralEtaStripSize cut.");
0221 descriptions.addWithDefaultLabel(desc);
0222 }
0223
0224 #endif