CMS 3D CMS Logo

HLTJetL1MatchProducer.cc
Go to the documentation of this file.
1 #include <string>
2 
7 
11 
12 template <typename T>
14  jetsInput_ = iConfig.template getParameter<edm::InputTag>("jetsInput");
15  L1TauJets_ = iConfig.template getParameter<edm::InputTag>("L1TauJets");
16  L1CenJets_ = iConfig.template getParameter<edm::InputTag>("L1CenJets");
17  L1ForJets_ = iConfig.template getParameter<edm::InputTag>("L1ForJets");
18  DeltaR_ = iConfig.template getParameter<double>("DeltaR");
19 
20  typedef std::vector<T> TCollection;
21  m_theJetToken = consumes<TCollection>(jetsInput_);
22  m_theL1TauJetToken = consumes<l1extra::L1JetParticleCollection>(L1TauJets_);
23  m_theL1CenJetToken = consumes<l1extra::L1JetParticleCollection>(L1CenJets_);
24  m_theL1ForJetToken = consumes<l1extra::L1JetParticleCollection>(L1ForJets_);
25  produces<TCollection>();
26 }
27 
28 template <typename T>
30 
31 template <typename T>
33 
34 template <typename T>
37  desc.add<edm::InputTag>("jetsInput", edm::InputTag("hltAntiKT5PFJets"));
38  desc.add<edm::InputTag>("L1TauJets", edm::InputTag("hltL1extraParticles", "Tau"));
39  desc.add<edm::InputTag>("L1CenJets", edm::InputTag("hltL1extraParticles", "Central"));
40  desc.add<edm::InputTag>("L1ForJets", edm::InputTag("hltL1extraParticles", "Forward"));
41  desc.add<double>("DeltaR", 0.5);
43 }
44 
45 template <typename T>
47  typedef std::vector<T> TCollection;
48 
50  iEvent.getByToken(m_theJetToken, jets);
51 
52  std::unique_ptr<TCollection> result(new TCollection);
53 
55  iEvent.getByToken(m_theL1TauJetToken, l1TauJets);
56 
58  iEvent.getByToken(m_theL1CenJetToken, l1CenJets);
59 
61  iEvent.getByToken(m_theL1ForJetToken, l1ForJets);
62 
63  typename TCollection::const_iterator jet_iter;
64  for (jet_iter = jets->begin(); jet_iter != jets->end(); ++jet_iter) {
65  bool isMatched = false;
66 
67  //std::cout << "FL: l1TauJets.size = " << l1TauJets->size() << std::endl;
68  for (unsigned int jetc = 0; jetc < l1TauJets->size(); ++jetc) {
69  const double deltaeta = jet_iter->eta() - (*l1TauJets)[jetc].eta();
70  const double deltaphi = deltaPhi(jet_iter->phi(), (*l1TauJets)[jetc].phi());
71  //std::cout << "FL: sqrt(2) = " << sqrt(2) << std::endl;
72  if (sqrt(deltaeta * deltaeta + deltaphi * deltaphi) < DeltaR_)
73  isMatched = true;
74  }
75 
76  for (unsigned int jetc = 0; jetc < l1CenJets->size(); ++jetc) {
77  const double deltaeta = jet_iter->eta() - (*l1CenJets)[jetc].eta();
78  const double deltaphi = deltaPhi(jet_iter->phi(), (*l1CenJets)[jetc].phi());
79  if (sqrt(deltaeta * deltaeta + deltaphi * deltaphi) < DeltaR_)
80  isMatched = true;
81  }
82 
83  for (unsigned int jetc = 0; jetc < l1ForJets->size(); ++jetc) {
84  const double deltaeta = jet_iter->eta() - (*l1ForJets)[jetc].eta();
85  const double deltaphi = deltaPhi(jet_iter->phi(), (*l1ForJets)[jetc].phi());
86  if (sqrt(deltaeta * deltaeta + deltaphi * deltaphi) < DeltaR_)
87  isMatched = true;
88  }
89 
90  if (isMatched == true)
91  result->push_back(*jet_iter);
92 
93  } // jet_iter
94 
95  iEvent.put(std::move(result));
96 }
std::string defaultModuleLabel()
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
~HLTJetL1MatchProducer() override
int iEvent
Definition: GenABIO.cc:224
T sqrt(T t)
Definition: SSEVec.h:19
HLTJetL1MatchProducer(const edm::ParameterSet &)
void produce(edm::Event &, const edm::EventSetup &) override
bool isMatched(TrackingRecHit const &hit)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
def move(src, dest)
Definition: eostools.py:511