File indexing completed on 2024-04-06 12:18:34
0001 #include <cmath>
0002 #include <memory>
0003
0004 #include "DataFormats/Common/interface/Handle.h"
0005 #include "DataFormats/Math/interface/deltaR.h"
0006 #include "FWCore/Framework/interface/ESHandle.h"
0007 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0008 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0009 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0010 #include "HLTrigger/HLTcore/interface/defaultModuleLabel.h"
0011 #include "HLTrigger/JetMET/interface/HLTJetL1MatchProducer.h"
0012
0013 template <typename T>
0014 HLTJetL1MatchProducer<T>::HLTJetL1MatchProducer(const edm::ParameterSet& iConfig) {
0015 jetsInput_ = iConfig.template getParameter<edm::InputTag>("jetsInput");
0016 L1TauJets_ = iConfig.template getParameter<edm::InputTag>("L1TauJets");
0017 L1CenJets_ = iConfig.template getParameter<edm::InputTag>("L1CenJets");
0018 L1ForJets_ = iConfig.template getParameter<edm::InputTag>("L1ForJets");
0019
0020
0021 auto const DeltaR = iConfig.template getParameter<double>("DeltaR");
0022 DeltaR2_ = DeltaR * std::abs(DeltaR);
0023
0024 typedef std::vector<T> TCollection;
0025 m_theJetToken = consumes<TCollection>(jetsInput_);
0026 m_theL1TauJetToken = consumes<l1extra::L1JetParticleCollection>(L1TauJets_);
0027 m_theL1CenJetToken = consumes<l1extra::L1JetParticleCollection>(L1CenJets_);
0028 m_theL1ForJetToken = consumes<l1extra::L1JetParticleCollection>(L1ForJets_);
0029 produces<TCollection>();
0030 }
0031
0032 template <typename T>
0033 void HLTJetL1MatchProducer<T>::beginJob() {}
0034
0035 template <typename T>
0036 HLTJetL1MatchProducer<T>::~HLTJetL1MatchProducer() = default;
0037
0038 template <typename T>
0039 void HLTJetL1MatchProducer<T>::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0040 edm::ParameterSetDescription desc;
0041 desc.add<edm::InputTag>("jetsInput", edm::InputTag("hltAntiKT5PFJets"));
0042 desc.add<edm::InputTag>("L1TauJets", edm::InputTag("hltL1extraParticles", "Tau"));
0043 desc.add<edm::InputTag>("L1CenJets", edm::InputTag("hltL1extraParticles", "Central"));
0044 desc.add<edm::InputTag>("L1ForJets", edm::InputTag("hltL1extraParticles", "Forward"));
0045 desc.add<double>("DeltaR", 0.5);
0046 descriptions.add(defaultModuleLabel<HLTJetL1MatchProducer<T>>(), desc);
0047 }
0048
0049 template <typename T>
0050 void HLTJetL1MatchProducer<T>::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0051 auto const& jets = iEvent.get(m_theJetToken);
0052
0053 auto result = std::make_unique<std::vector<T>>();
0054 result->reserve(jets.size());
0055
0056 auto const& l1TauJets = iEvent.get(m_theL1TauJetToken);
0057 auto const& l1CenJets = iEvent.get(m_theL1CenJetToken);
0058 auto const& l1ForJets = iEvent.get(m_theL1ForJetToken);
0059
0060 for (auto const& jet : jets) {
0061 bool isMatched = false;
0062
0063 for (auto const& l1t_obj : l1TauJets) {
0064 if (reco::deltaR2(jet.eta(), jet.phi(), l1t_obj.eta(), l1t_obj.phi()) < DeltaR2_) {
0065 isMatched = true;
0066 break;
0067 }
0068 }
0069
0070 if (isMatched) {
0071 result->emplace_back(jet);
0072 continue;
0073 }
0074
0075 for (auto const& l1t_obj : l1CenJets) {
0076 if (reco::deltaR2(jet.eta(), jet.phi(), l1t_obj.eta(), l1t_obj.phi()) < DeltaR2_) {
0077 isMatched = true;
0078 break;
0079 }
0080 }
0081
0082 if (isMatched) {
0083 result->emplace_back(jet);
0084 continue;
0085 }
0086
0087 for (auto const& l1t_obj : l1ForJets) {
0088 if (reco::deltaR2(jet.eta(), jet.phi(), l1t_obj.eta(), l1t_obj.phi()) < DeltaR2_) {
0089 isMatched = true;
0090 break;
0091 }
0092 }
0093
0094 if (isMatched) {
0095 result->emplace_back(jet);
0096 continue;
0097 }
0098 }
0099
0100 iEvent.put(std::move(result));
0101 }