CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTJetL1MatchProducer.cc
Go to the documentation of this file.
8 
10 {
11  jetsInput_ = iConfig.getParameter<edm::InputTag>("jetsInput");
12  L1TauJets_ = iConfig.getParameter<edm::InputTag>("L1TauJets");
13  L1CenJets_ = iConfig.getParameter<edm::InputTag>("L1CenJets");
14  L1ForJets_ = iConfig.getParameter<edm::InputTag>("L1ForJets");
15  DeltaR_ = iConfig.getParameter<double>("DeltaR");
16 
17  produces< reco::CaloJetCollection > ();
18 
19 }
20 
22 {
23 
24 }
25 
27 {
28 
29 }
30 
32 {
33 
35  iEvent.getByLabel(jetsInput_, calojets);
36 
37  std::auto_ptr<reco::CaloJetCollection> result (new reco::CaloJetCollection);
38 
39 
41  iEvent.getByLabel(L1TauJets_,l1TauJets);
42 
44  iEvent.getByLabel(L1CenJets_,l1CenJets);
45 
47  iEvent.getByLabel(L1ForJets_,l1ForJets);
48 
49 
50  for (reco::CaloJetCollection::const_iterator calojetc = calojets->begin();
51  calojetc != calojets->end(); ++calojetc) {
52 
53  bool isMatched=false;
54 
55  //std::cout << "FL: l1TauJets.size = " << l1TauJets->size() << std::endl;
56  for (unsigned int jetc=0;jetc<l1TauJets->size();++jetc)
57  {
58  const double deltaeta=calojetc->eta()-(*l1TauJets)[jetc].eta();
59  const double deltaphi=deltaPhi(calojetc->phi(),(*l1TauJets)[jetc].phi());
60  //std::cout << "FL: sqrt(2) = " << sqrt(2) << std::endl;
61  if (sqrt(deltaeta*deltaeta+deltaphi*deltaphi) < DeltaR_) isMatched=true;
62  }
63 
64  for (unsigned int jetc=0;jetc<l1CenJets->size();++jetc)
65  {
66  const double deltaeta=calojetc->eta()-(*l1CenJets)[jetc].eta();
67  const double deltaphi=deltaPhi(calojetc->phi(),(*l1CenJets)[jetc].phi());
68  if (sqrt(deltaeta*deltaeta+deltaphi*deltaphi) < DeltaR_) isMatched=true;
69  }
70 
71  for (unsigned int jetc=0;jetc<l1ForJets->size();++jetc)
72  {
73  const double deltaeta=calojetc->eta()-(*l1ForJets)[jetc].eta();
74  const double deltaphi=deltaPhi(calojetc->phi(),(*l1ForJets)[jetc].phi());
75  if (sqrt(deltaeta*deltaeta+deltaphi*deltaphi) < DeltaR_) isMatched=true;
76  }
77 
78 
79  if (isMatched==true) result->push_back(*calojetc);
80 
81  } // calojetc
82 
83  iEvent.put( result);
84 
85 }
86 
87 
T getParameter(std::string const &) const
double deltaPhi(float phi1, float phi2)
Definition: VectorUtil.h:30
virtual void produce(edm::Event &, const edm::EventSetup &)
int iEvent
Definition: GenABIO.cc:243
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:84
T sqrt(T t)
Definition: SSEVec.h:28
tuple result
Definition: query.py:137
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:359
HLTJetL1MatchProducer(const edm::ParameterSet &)
std::vector< CaloJet > CaloJetCollection
collection of CaloJet objects