CMS 3D CMS Logo

HLTJetL1MatchProducer.cc
Go to the documentation of this file.
1 #include <cmath>
2 #include <memory>
3 
12 
13 template <typename T>
15  jetsInput_ = iConfig.template getParameter<edm::InputTag>("jetsInput");
16  L1TauJets_ = iConfig.template getParameter<edm::InputTag>("L1TauJets");
17  L1CenJets_ = iConfig.template getParameter<edm::InputTag>("L1CenJets");
18  L1ForJets_ = iConfig.template getParameter<edm::InputTag>("L1ForJets");
19 
20  // minimum delta-R^2 threshold with sign
21  auto const DeltaR = iConfig.template getParameter<double>("DeltaR");
22  DeltaR2_ = DeltaR * std::abs(DeltaR);
23 
24  typedef std::vector<T> TCollection;
25  m_theJetToken = consumes<TCollection>(jetsInput_);
26  m_theL1TauJetToken = consumes<l1extra::L1JetParticleCollection>(L1TauJets_);
27  m_theL1CenJetToken = consumes<l1extra::L1JetParticleCollection>(L1CenJets_);
28  m_theL1ForJetToken = consumes<l1extra::L1JetParticleCollection>(L1ForJets_);
29  produces<TCollection>();
30 }
31 
32 template <typename T>
34 
35 template <typename T>
37 
38 template <typename T>
41  desc.add<edm::InputTag>("jetsInput", edm::InputTag("hltAntiKT5PFJets"));
42  desc.add<edm::InputTag>("L1TauJets", edm::InputTag("hltL1extraParticles", "Tau"));
43  desc.add<edm::InputTag>("L1CenJets", edm::InputTag("hltL1extraParticles", "Central"));
44  desc.add<edm::InputTag>("L1ForJets", edm::InputTag("hltL1extraParticles", "Forward"));
45  desc.add<double>("DeltaR", 0.5);
47 }
48 
49 template <typename T>
51  auto const& jets = iEvent.get(m_theJetToken);
52 
53  auto result = std::make_unique<std::vector<T>>();
54  result->reserve(jets.size());
55 
56  auto const& l1TauJets = iEvent.get(m_theL1TauJetToken);
57  auto const& l1CenJets = iEvent.get(m_theL1CenJetToken);
58  auto const& l1ForJets = iEvent.get(m_theL1ForJetToken);
59 
60  for (auto const& jet : jets) {
61  bool isMatched = false;
62 
63  for (auto const& l1t_obj : l1TauJets) {
64  if (reco::deltaR2(jet.eta(), jet.phi(), l1t_obj.eta(), l1t_obj.phi()) < DeltaR2_) {
65  isMatched = true;
66  break;
67  }
68  }
69 
70  if (isMatched) {
71  result->emplace_back(jet);
72  continue;
73  }
74 
75  for (auto const& l1t_obj : l1CenJets) {
76  if (reco::deltaR2(jet.eta(), jet.phi(), l1t_obj.eta(), l1t_obj.phi()) < DeltaR2_) {
77  isMatched = true;
78  break;
79  }
80  }
81 
82  if (isMatched) {
83  result->emplace_back(jet);
84  continue;
85  }
86 
87  for (auto const& l1t_obj : l1ForJets) {
88  if (reco::deltaR2(jet.eta(), jet.phi(), l1t_obj.eta(), l1t_obj.phi()) < DeltaR2_) {
89  isMatched = true;
90  break;
91  }
92  }
93 
94  if (isMatched) {
95  result->emplace_back(jet);
96  continue;
97  }
98  }
99 
100  iEvent.put(std::move(result));
101 }
Definition: DeltaR.py:1
std::string defaultModuleLabel()
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
~HLTJetL1MatchProducer() override
int iEvent
Definition: GenABIO.cc:224
HLTJetL1MatchProducer(const edm::ParameterSet &)
void produce(edm::Event &, const edm::EventSetup &) override
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool isMatched(TrackingRecHit const &hit)
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
void add(std::string const &label, ParameterSetDescription const &psetDescription)
def move(src, dest)
Definition: eostools.py:511