CMS 3D CMS Logo

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