Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:18:31

0001 #ifndef HLTrigger_JetMET_HLTJetsMatchedToFilteredJetsProducer_h
0002 #define HLTrigger_JetMET_HLTJetsMatchedToFilteredJetsProducer_h
0003 
0004 #include "FWCore/Framework/interface/Frameworkfwd.h"
0005 #include "FWCore/Framework/interface/global/EDProducer.h"
0006 #include "FWCore/Framework/interface/Event.h"
0007 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0008 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0009 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0010 #include "FWCore/Utilities/interface/Exception.h"
0011 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
0012 #include "DataFormats/HLTReco/interface/TriggerTypeDefs.h"
0013 #include "DataFormats/Math/interface/deltaR.h"
0014 
0015 #include <vector>
0016 #include <memory>
0017 #include <utility>
0018 
0019 template <typename TriggerJetsType, typename TriggerJetsRefType>
0020 class HLTJetsMatchedToFilteredJetsProducer : public edm::global::EDProducer<> {
0021 public:
0022   explicit HLTJetsMatchedToFilteredJetsProducer(edm::ParameterSet const& iConfig);
0023 
0024   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0025 
0026   void produce(edm::StreamID streamID, edm::Event& iEvent, edm::EventSetup const& iSetup) const override;
0027 
0028 private:
0029   using InputJetCollection = edm::View<TriggerJetsType>;
0030   using OutputJetCollection = std::vector<TriggerJetsType>;
0031 
0032   // collection of input jets
0033   edm::EDGetTokenT<InputJetCollection> const jetsToken_;
0034   // collection of TriggerFilterObjectWithRefs containing TriggerObjects holding refs to Jets stored by an HLTFilter
0035   edm::EDGetTokenT<trigger::TriggerFilterObjectWithRefs> const triggerJetsToken_;
0036   // TriggerType of Jets produced by previous HLTFilter
0037   int const triggerJetsType_;
0038   // maximum Delta-R and Delta-R^2 distances between matched jets
0039   double const maxDeltaR_, maxDeltaR2_;
0040 };
0041 
0042 template <typename TriggerJetsType, typename TriggerJetsRefType>
0043 HLTJetsMatchedToFilteredJetsProducer<TriggerJetsType, TriggerJetsRefType>::HLTJetsMatchedToFilteredJetsProducer(
0044     edm::ParameterSet const& iConfig)
0045     : jetsToken_(consumes(iConfig.getParameter<edm::InputTag>("src"))),
0046       triggerJetsToken_(consumes(iConfig.getParameter<edm::InputTag>("triggerJetsFilter"))),
0047       triggerJetsType_(iConfig.getParameter<int>("triggerJetsType")),
0048       maxDeltaR_(iConfig.getParameter<double>("maxDeltaR")),
0049       maxDeltaR2_(maxDeltaR_ * maxDeltaR_) {
0050   if (maxDeltaR_ <= 0.) {
0051     throw cms::Exception("HLTJetsMatchedToFilteredJetsProducerConfiguration")
0052         << "invalid value for parameter \"maxDeltaR\" (must be > 0): " << maxDeltaR_;
0053   }
0054 
0055   produces<OutputJetCollection>();
0056 }
0057 
0058 template <typename TriggerJetsType, typename TriggerJetsRefType>
0059 void HLTJetsMatchedToFilteredJetsProducer<TriggerJetsType, TriggerJetsRefType>::fillDescriptions(
0060     edm::ConfigurationDescriptions& descriptions) {
0061   edm::ParameterSetDescription desc;
0062   desc.add<edm::InputTag>("src", edm::InputTag("hltJets"));
0063   desc.add<edm::InputTag>("triggerJetsFilter", edm::InputTag("hltCaloJetsFiltered"));
0064   desc.add<int>("triggerJetsType", trigger::TriggerJet);
0065   desc.add<double>("maxDeltaR", 0.5);
0066   descriptions.addWithDefaultLabel(desc);
0067 }
0068 
0069 template <typename TriggerJetsType, typename TriggerJetsRefType>
0070 void HLTJetsMatchedToFilteredJetsProducer<TriggerJetsType, TriggerJetsRefType>::produce(edm::StreamID,
0071                                                                                         edm::Event& iEvent,
0072                                                                                         edm::EventSetup const&) const {
0073   auto const& jets = iEvent.get(jetsToken_);
0074 
0075   std::vector<TriggerJetsRefType> triggerJetsRefVec;
0076   iEvent.get(triggerJetsToken_).getObjects(triggerJetsType_, triggerJetsRefVec);
0077 
0078   auto outJets = std::make_unique<OutputJetCollection>();
0079   outJets->reserve(jets.size());
0080 
0081   for (auto const& jet_i : jets) {
0082     for (auto const& jetRef_j : triggerJetsRefVec) {
0083       if (reco::deltaR2(jet_i.p4(), jetRef_j->p4()) < maxDeltaR2_) {
0084         outJets->emplace_back(jet_i);
0085         break;
0086       }
0087     }
0088   }
0089 
0090   iEvent.put(std::move(outJets));
0091 }
0092 
0093 #endif