File indexing completed on 2023-03-17 11:09:35
0001 #include <string>
0002
0003 #include "HLTrigger/JetMET/interface/HLTJetL1MatchProducer.h"
0004 #include "DataFormats/Common/interface/Handle.h"
0005 #include "FWCore/Framework/interface/ESHandle.h"
0006 #include "DataFormats/Math/interface/deltaPhi.h"
0007
0008 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0009 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0010 #include "HLTrigger/HLTcore/interface/defaultModuleLabel.h"
0011
0012 template <typename T>
0013 HLTJetL1MatchProducer<T>::HLTJetL1MatchProducer(const edm::ParameterSet& iConfig) {
0014 jetsInput_ = iConfig.template getParameter<edm::InputTag>("jetsInput");
0015 L1TauJets_ = iConfig.template getParameter<edm::InputTag>("L1TauJets");
0016 L1CenJets_ = iConfig.template getParameter<edm::InputTag>("L1CenJets");
0017 L1ForJets_ = iConfig.template getParameter<edm::InputTag>("L1ForJets");
0018 DeltaR_ = iConfig.template getParameter<double>("DeltaR");
0019
0020 typedef std::vector<T> TCollection;
0021 m_theJetToken = consumes<TCollection>(jetsInput_);
0022 m_theL1TauJetToken = consumes<l1extra::L1JetParticleCollection>(L1TauJets_);
0023 m_theL1CenJetToken = consumes<l1extra::L1JetParticleCollection>(L1CenJets_);
0024 m_theL1ForJetToken = consumes<l1extra::L1JetParticleCollection>(L1ForJets_);
0025 produces<TCollection>();
0026 }
0027
0028 template <typename T>
0029 void HLTJetL1MatchProducer<T>::beginJob() {}
0030
0031 template <typename T>
0032 HLTJetL1MatchProducer<T>::~HLTJetL1MatchProducer() = default;
0033
0034 template <typename T>
0035 void HLTJetL1MatchProducer<T>::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0036 edm::ParameterSetDescription desc;
0037 desc.add<edm::InputTag>("jetsInput", edm::InputTag("hltAntiKT5PFJets"));
0038 desc.add<edm::InputTag>("L1TauJets", edm::InputTag("hltL1extraParticles", "Tau"));
0039 desc.add<edm::InputTag>("L1CenJets", edm::InputTag("hltL1extraParticles", "Central"));
0040 desc.add<edm::InputTag>("L1ForJets", edm::InputTag("hltL1extraParticles", "Forward"));
0041 desc.add<double>("DeltaR", 0.5);
0042 descriptions.add(defaultModuleLabel<HLTJetL1MatchProducer<T>>(), desc);
0043 }
0044
0045 template <typename T>
0046 void HLTJetL1MatchProducer<T>::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0047 typedef std::vector<T> TCollection;
0048
0049 edm::Handle<TCollection> jets;
0050 iEvent.getByToken(m_theJetToken, jets);
0051
0052 std::unique_ptr<TCollection> result(new TCollection);
0053
0054 edm::Handle<l1extra::L1JetParticleCollection> l1TauJets;
0055 iEvent.getByToken(m_theL1TauJetToken, l1TauJets);
0056
0057 edm::Handle<l1extra::L1JetParticleCollection> l1CenJets;
0058 iEvent.getByToken(m_theL1CenJetToken, l1CenJets);
0059
0060 edm::Handle<l1extra::L1JetParticleCollection> l1ForJets;
0061 iEvent.getByToken(m_theL1ForJetToken, l1ForJets);
0062
0063 typename TCollection::const_iterator jet_iter;
0064 for (jet_iter = jets->begin(); jet_iter != jets->end(); ++jet_iter) {
0065 bool isMatched = false;
0066
0067
0068 for (unsigned int jetc = 0; jetc < l1TauJets->size(); ++jetc) {
0069 const double deltaeta = jet_iter->eta() - (*l1TauJets)[jetc].eta();
0070 const double deltaphi = deltaPhi(jet_iter->phi(), (*l1TauJets)[jetc].phi());
0071
0072 if (sqrt(deltaeta * deltaeta + deltaphi * deltaphi) < DeltaR_)
0073 isMatched = true;
0074 }
0075
0076 for (unsigned int jetc = 0; jetc < l1CenJets->size(); ++jetc) {
0077 const double deltaeta = jet_iter->eta() - (*l1CenJets)[jetc].eta();
0078 const double deltaphi = deltaPhi(jet_iter->phi(), (*l1CenJets)[jetc].phi());
0079 if (sqrt(deltaeta * deltaeta + deltaphi * deltaphi) < DeltaR_)
0080 isMatched = true;
0081 }
0082
0083 for (unsigned int jetc = 0; jetc < l1ForJets->size(); ++jetc) {
0084 const double deltaeta = jet_iter->eta() - (*l1ForJets)[jetc].eta();
0085 const double deltaphi = deltaPhi(jet_iter->phi(), (*l1ForJets)[jetc].phi());
0086 if (sqrt(deltaeta * deltaeta + deltaphi * deltaphi) < DeltaR_)
0087 isMatched = true;
0088 }
0089
0090 if (isMatched == true)
0091 result->push_back(*jet_iter);
0092
0093 }
0094
0095 iEvent.put(std::move(result));
0096 }