CMS 3D CMS Logo

HLTJetL1TMatchProducer.cc
Go to the documentation of this file.
1 #include <string>
2 
7 
12 
13 template <typename T>
15  jetsInput_ = iConfig.template getParameter<edm::InputTag>("jetsInput");
16  L1Jets_ = iConfig.template getParameter<edm::InputTag>("L1Jets");
17  DeltaR_ = iConfig.template getParameter<double>("DeltaR");
18 
19  typedef std::vector<T> TCollection;
20  m_theJetToken = consumes<TCollection>(jetsInput_);
21  m_theL1JetToken = consumes<l1t::JetBxCollection>(L1Jets_);
22  produces<TCollection>();
23 }
24 
25 template <typename T>
27 
28 template <typename T>
30 
31 template <typename T>
34  desc.add<edm::InputTag>("jetsInput", edm::InputTag("hltAntiKT5PFJets"));
35  desc.add<edm::InputTag>("L1Jets", edm::InputTag("hltCaloStage2Digis"));
36  desc.add<double>("DeltaR", 0.5);
38 }
39 
40 template <typename T>
42  typedef std::vector<T> TCollection;
43 
45  iEvent.getByToken(m_theJetToken, jets);
46 
47  std::unique_ptr<TCollection> result(new TCollection);
48 
50  iEvent.getByToken(m_theL1JetToken, l1Jets);
51  bool trigger_bx_only = true; // selection of BX not implemented
52 
53  if (l1Jets.isValid()) {
54  typename TCollection::const_iterator jet_iter;
55  for (jet_iter = jets->begin(); jet_iter != jets->end(); ++jet_iter) {
56  bool isMatched = false;
57  for (int ibx = l1Jets->getFirstBX(); ibx <= l1Jets->getLastBX(); ++ibx) {
58  if (trigger_bx_only && (ibx != 0))
59  continue;
60  for (auto it = l1Jets->begin(ibx); it != l1Jets->end(ibx); it++) {
61  if (it->et() == 0)
62  continue; // if you don't care about L1T candidates with zero ET.
63  const double deltaeta = jet_iter->eta() - it->eta();
64  const double deltaphi = deltaPhi(jet_iter->phi(), it->phi());
65  if (sqrt(deltaeta * deltaeta + deltaphi * deltaphi) < DeltaR_)
66  isMatched = true;
67  //cout << "bx: " << ibx << " et: " << it->et() << " eta: " << it->eta() << " phi: " << it->phi() << "\n";
68  }
69  }
70  if (isMatched == true)
71  result->push_back(*jet_iter);
72  } // jet_iter
73  } else {
74  edm::LogWarning("MissingProduct") << "L1Upgrade l1Jets bx collection not found." << std::endl;
75  }
76 
77  iEvent.put(std::move(result));
78 }
int getLastBX() const
int getFirstBX() const
std::string defaultModuleLabel()
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
const_iterator begin(int bx) const
int iEvent
Definition: GenABIO.cc:224
T sqrt(T t)
Definition: SSEVec.h:19
bool isMatched(TrackingRecHit const &hit)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
bool isValid() const
Definition: HandleBase.h:70
const_iterator end(int bx) const
void produce(edm::Event &, const edm::EventSetup &) override
~HLTJetL1TMatchProducer() override
HLTJetL1TMatchProducer(const edm::ParameterSet &)
Log< level::Warning, false > LogWarning
def move(src, dest)
Definition: eostools.py:511