CMS 3D CMS Logo

HLTJetL1MatchProducer.cc
Go to the documentation of this file.
1 #include <string>
2 
7 
11 
12 template<typename T>
14 {
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  DeltaR_ = iConfig.template getParameter<double>("DeltaR");
20 
21  typedef std::vector<T> TCollection;
22  m_theJetToken = consumes<TCollection>(jetsInput_);
23  m_theL1TauJetToken = consumes<l1extra::L1JetParticleCollection>(L1TauJets_);
24  m_theL1CenJetToken = consumes<l1extra::L1JetParticleCollection>(L1CenJets_);
25  m_theL1ForJetToken = consumes<l1extra::L1JetParticleCollection>(L1ForJets_);
26  produces<TCollection> ();
27 
28 }
29 
30 template<typename T>
32 {
33 
34 }
35 
36 template<typename T>
38 
39 template<typename T>
42  desc.add<edm::InputTag>("jetsInput",edm::InputTag("hltAntiKT5PFJets"));
43  desc.add<edm::InputTag>("L1TauJets",edm::InputTag("hltL1extraParticles","Tau"));
44  desc.add<edm::InputTag>("L1CenJets",edm::InputTag("hltL1extraParticles","Central"));
45  desc.add<edm::InputTag>("L1ForJets",edm::InputTag("hltL1extraParticles","Forward"));
46  desc.add<double>("DeltaR",0.5);
47  descriptions.add(defaultModuleLabel<HLTJetL1MatchProducer<T>>(), desc);
48 }
49 
50 template<typename T>
52 {
53 
54  typedef std::vector<T> TCollection;
55 
57  iEvent.getByToken(m_theJetToken, jets);
58 
59  std::unique_ptr<TCollection> result (new TCollection);
60 
61 
63  iEvent.getByToken(m_theL1TauJetToken,l1TauJets);
64 
66  iEvent.getByToken(m_theL1CenJetToken,l1CenJets);
67 
69  iEvent.getByToken(m_theL1ForJetToken,l1ForJets);
70 
71  typename TCollection::const_iterator jet_iter;
72  for (jet_iter = jets->begin(); jet_iter != jets->end(); ++jet_iter) {
73 
74  bool isMatched=false;
75 
76  //std::cout << "FL: l1TauJets.size = " << l1TauJets->size() << std::endl;
77  for (unsigned int jetc=0;jetc<l1TauJets->size();++jetc)
78  {
79  const double deltaeta=jet_iter->eta()-(*l1TauJets)[jetc].eta();
80  const double deltaphi=deltaPhi(jet_iter->phi(),(*l1TauJets)[jetc].phi());
81  //std::cout << "FL: sqrt(2) = " << sqrt(2) << std::endl;
82  if (sqrt(deltaeta*deltaeta+deltaphi*deltaphi) < DeltaR_) isMatched=true;
83  }
84 
85  for (unsigned int jetc=0;jetc<l1CenJets->size();++jetc)
86  {
87  const double deltaeta=jet_iter->eta()-(*l1CenJets)[jetc].eta();
88  const double deltaphi=deltaPhi(jet_iter->phi(),(*l1CenJets)[jetc].phi());
89  if (sqrt(deltaeta*deltaeta+deltaphi*deltaphi) < DeltaR_) isMatched=true;
90  }
91 
92  for (unsigned int jetc=0;jetc<l1ForJets->size();++jetc)
93  {
94  const double deltaeta=jet_iter->eta()-(*l1ForJets)[jetc].eta();
95  const double deltaphi=deltaPhi(jet_iter->phi(),(*l1ForJets)[jetc].phi());
96  if (sqrt(deltaeta*deltaeta+deltaphi*deltaphi) < DeltaR_) isMatched=true;
97  }
98 
99 
100  if (isMatched==true) result->push_back(*jet_iter);
101 
102  } // jet_iter
103 
104  iEvent.put(std::move(result));
105 
106 }
107 
108 
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)
~HLTJetL1MatchProducer() override
int iEvent
Definition: GenABIO.cc:230
T sqrt(T t)
Definition: SSEVec.h:18
vector< PseudoJet > jets
HLTJetL1MatchProducer(const edm::ParameterSet &)
void produce(edm::Event &, const edm::EventSetup &) override
bool isMatched(TrackingRecHit const &hit)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
def move(src, dest)
Definition: eostools.py:510