CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/HLTrigger/JetMET/src/HLTJetL1MatchProducer.cc

Go to the documentation of this file.
00001 #include "HLTrigger/JetMET/interface/HLTJetL1MatchProducer.h"
00002 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00003 #include "DataFormats/Common/interface/Handle.h"
00004 #include "FWCore/Framework/interface/ESHandle.h"
00005 #include "DataFormats/L1Trigger/interface/L1JetParticleFwd.h"
00006 #include "DataFormats/L1Trigger/interface/L1JetParticle.h"
00007 #include "DataFormats/Math/interface/deltaPhi.h"
00008 
00009 HLTJetL1MatchProducer::HLTJetL1MatchProducer(const edm::ParameterSet& iConfig)
00010 {
00011   jetsInput_ = iConfig.getParameter<edm::InputTag>("jetsInput");
00012   L1TauJets_ = iConfig.getParameter<edm::InputTag>("L1TauJets");
00013   L1CenJets_ = iConfig.getParameter<edm::InputTag>("L1CenJets");
00014   L1ForJets_ = iConfig.getParameter<edm::InputTag>("L1ForJets");
00015   DeltaR_ = iConfig.getParameter<double>("DeltaR");
00016 
00017   produces< reco::CaloJetCollection > ();
00018 
00019 }
00020 
00021 void HLTJetL1MatchProducer::beginJob()
00022 {
00023 
00024 }
00025 
00026 HLTJetL1MatchProducer::~HLTJetL1MatchProducer()
00027 {
00028 
00029 }
00030 
00031 void HLTJetL1MatchProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00032 {
00033 
00034   edm::Handle<reco::CaloJetCollection> calojets;
00035   iEvent.getByLabel(jetsInput_, calojets);
00036 
00037   std::auto_ptr<reco::CaloJetCollection> result (new reco::CaloJetCollection);
00038 
00039 
00040   edm::Handle<l1extra::L1JetParticleCollection> l1TauJets;
00041   iEvent.getByLabel(L1TauJets_,l1TauJets);
00042 
00043   edm::Handle<l1extra::L1JetParticleCollection> l1CenJets;
00044   iEvent.getByLabel(L1CenJets_,l1CenJets);
00045 
00046   edm::Handle<l1extra::L1JetParticleCollection> l1ForJets;
00047   iEvent.getByLabel(L1ForJets_,l1ForJets);
00048 
00049 
00050   for (reco::CaloJetCollection::const_iterator calojetc = calojets->begin(); 
00051        calojetc != calojets->end(); ++calojetc) {
00052 
00053     bool isMatched=false;
00054 
00055     //std::cout << "FL: l1TauJets.size  = " << l1TauJets->size() << std::endl;
00056     for (unsigned int jetc=0;jetc<l1TauJets->size();++jetc)
00057     {
00058       const double deltaeta=calojetc->eta()-(*l1TauJets)[jetc].eta();
00059       const double deltaphi=deltaPhi(calojetc->phi(),(*l1TauJets)[jetc].phi());
00060       //std::cout << "FL: sqrt(2) = " << sqrt(2) << std::endl;
00061       if (sqrt(deltaeta*deltaeta+deltaphi*deltaphi) < DeltaR_) isMatched=true;
00062     }
00063 
00064     for (unsigned int jetc=0;jetc<l1CenJets->size();++jetc)
00065     {
00066       const double deltaeta=calojetc->eta()-(*l1CenJets)[jetc].eta();
00067       const double deltaphi=deltaPhi(calojetc->phi(),(*l1CenJets)[jetc].phi());
00068       if (sqrt(deltaeta*deltaeta+deltaphi*deltaphi) < DeltaR_) isMatched=true;
00069     }
00070 
00071     for (unsigned int jetc=0;jetc<l1ForJets->size();++jetc)
00072     {
00073       const double deltaeta=calojetc->eta()-(*l1ForJets)[jetc].eta();
00074       const double deltaphi=deltaPhi(calojetc->phi(),(*l1ForJets)[jetc].phi());
00075       if (sqrt(deltaeta*deltaeta+deltaphi*deltaphi) < DeltaR_) isMatched=true;
00076     }
00077 
00078 
00079     if (isMatched==true) result->push_back(*calojetc);
00080 
00081   } // calojetc
00082 
00083   iEvent.put( result);
00084 
00085 }
00086 
00087