Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 11:09:35

0001 #include <string>
0002 
0003 #include "HLTrigger/JetMET/interface/HLTJetL1TMatchProducer.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/MessageLogger/interface/MessageLogger.h"
0009 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0010 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0011 #include "HLTrigger/HLTcore/interface/defaultModuleLabel.h"
0012 
0013 template <typename T>
0014 HLTJetL1TMatchProducer<T>::HLTJetL1TMatchProducer(const edm::ParameterSet& iConfig) {
0015   jetsInput_ = iConfig.template getParameter<edm::InputTag>("jetsInput");
0016   L1Jets_ = iConfig.template getParameter<edm::InputTag>("L1Jets");
0017   DeltaR_ = iConfig.template getParameter<double>("DeltaR");
0018 
0019   typedef std::vector<T> TCollection;
0020   m_theJetToken = consumes<TCollection>(jetsInput_);
0021   m_theL1JetToken = consumes<l1t::JetBxCollection>(L1Jets_);
0022   produces<TCollection>();
0023 }
0024 
0025 template <typename T>
0026 void HLTJetL1TMatchProducer<T>::beginJob() {}
0027 
0028 template <typename T>
0029 HLTJetL1TMatchProducer<T>::~HLTJetL1TMatchProducer() = default;
0030 
0031 template <typename T>
0032 void HLTJetL1TMatchProducer<T>::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0033   edm::ParameterSetDescription desc;
0034   desc.add<edm::InputTag>("jetsInput", edm::InputTag("hltAntiKT5PFJets"));
0035   desc.add<edm::InputTag>("L1Jets", edm::InputTag("hltCaloStage2Digis"));
0036   desc.add<double>("DeltaR", 0.5);
0037   descriptions.add(defaultModuleLabel<HLTJetL1TMatchProducer<T>>(), desc);
0038 }
0039 
0040 template <typename T>
0041 void HLTJetL1TMatchProducer<T>::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0042   typedef std::vector<T> TCollection;
0043 
0044   edm::Handle<TCollection> jets;
0045   iEvent.getByToken(m_theJetToken, jets);
0046 
0047   std::unique_ptr<TCollection> result(new TCollection);
0048 
0049   edm::Handle<l1t::JetBxCollection> l1Jets;
0050   iEvent.getByToken(m_theL1JetToken, l1Jets);
0051   bool trigger_bx_only = true;  // selection of BX not implemented
0052 
0053   if (l1Jets.isValid()) {
0054     typename TCollection::const_iterator jet_iter;
0055     for (jet_iter = jets->begin(); jet_iter != jets->end(); ++jet_iter) {
0056       bool isMatched = false;
0057       for (int ibx = l1Jets->getFirstBX(); ibx <= l1Jets->getLastBX(); ++ibx) {
0058         if (trigger_bx_only && (ibx != 0))
0059           continue;
0060         for (auto it = l1Jets->begin(ibx); it != l1Jets->end(ibx); it++) {
0061           if (it->et() == 0)
0062             continue;  // if you don't care about L1T candidates with zero ET.
0063           const double deltaeta = jet_iter->eta() - it->eta();
0064           const double deltaphi = deltaPhi(jet_iter->phi(), it->phi());
0065           if (sqrt(deltaeta * deltaeta + deltaphi * deltaphi) < DeltaR_)
0066             isMatched = true;
0067           //cout << "bx:  " << ibx << "  et:  "  << it->et() << "  eta:  "  << it->eta() << "  phi:  "  << it->phi() << "\n";
0068         }
0069       }
0070       if (isMatched == true)
0071         result->push_back(*jet_iter);
0072     }  // jet_iter
0073   } else {
0074     edm::LogWarning("MissingProduct") << "L1Upgrade l1Jets bx collection not found." << std::endl;
0075   }
0076 
0077   iEvent.put(std::move(result));
0078 }