CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2/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/CaloJet.h"
00003 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00004 #include "DataFormats/JetReco/interface/PFJet.h"
00005 #include "DataFormats/JetReco/interface/PFJetCollection.h"
00006 #include "DataFormats/JetReco/interface/TrackJetCollection.h"
00007 #include "DataFormats/JetReco/interface/BasicJetCollection.h"
00008 #include "DataFormats/Common/interface/Handle.h"
00009 #include "FWCore/Framework/interface/ESHandle.h"
00010 #include "DataFormats/L1Trigger/interface/L1JetParticleFwd.h"
00011 #include "DataFormats/L1Trigger/interface/L1JetParticle.h"
00012 #include "DataFormats/Math/interface/deltaPhi.h"
00013 
00014 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
00015 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
00016 
00017 #include<typeinfo>
00018 #include<string>
00019 
00020 template<typename T>
00021 HLTJetL1MatchProducer<T>::HLTJetL1MatchProducer(const edm::ParameterSet& iConfig)
00022 {
00023   jetsInput_ = iConfig.template getParameter<edm::InputTag>("jetsInput");
00024   L1TauJets_ = iConfig.template getParameter<edm::InputTag>("L1TauJets");
00025   L1CenJets_ = iConfig.template getParameter<edm::InputTag>("L1CenJets");
00026   L1ForJets_ = iConfig.template getParameter<edm::InputTag>("L1ForJets");
00027   DeltaR_ = iConfig.template getParameter<double>("DeltaR");
00028 
00029   typedef std::vector<T> TCollection;
00030   produces<TCollection> ();
00031 
00032 }
00033 
00034 template<typename T>
00035 void HLTJetL1MatchProducer<T>::beginJob()
00036 {
00037 
00038 }
00039 
00040 template<typename T>
00041 HLTJetL1MatchProducer<T>::~HLTJetL1MatchProducer()
00042 {
00043 
00044 }
00045 
00046 template<typename T>
00047 void HLTJetL1MatchProducer<T>::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
00048   edm::ParameterSetDescription desc;
00049   desc.add<edm::InputTag>("jetsInput",edm::InputTag("hltAntiKT5PFJets"));
00050   desc.add<edm::InputTag>("L1TauJets",edm::InputTag("hltL1extraParticles","Tau"));
00051   desc.add<edm::InputTag>("L1CenJets",edm::InputTag("hltL1extraParticles","Central"));
00052   desc.add<edm::InputTag>("L1ForJets",edm::InputTag("hltL1extraParticles","Forward"));
00053   desc.add<double>("DeltaR",0.5);
00054   descriptions.add(std::string("hlt")+std::string(typeid(HLTJetL1MatchProducer<T>).name()),desc);
00055 }
00056 
00057 template<typename T>
00058 void HLTJetL1MatchProducer<T>::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00059 {
00060 
00061   typedef std::vector<T> TCollection;
00062 
00063   edm::Handle<TCollection> jets;
00064   iEvent.getByLabel(jetsInput_, jets);
00065 
00066   std::auto_ptr<TCollection> result (new TCollection);
00067 
00068 
00069   edm::Handle<l1extra::L1JetParticleCollection> l1TauJets;
00070   iEvent.getByLabel(L1TauJets_,l1TauJets);
00071 
00072   edm::Handle<l1extra::L1JetParticleCollection> l1CenJets;
00073   iEvent.getByLabel(L1CenJets_,l1CenJets);
00074 
00075   edm::Handle<l1extra::L1JetParticleCollection> l1ForJets;
00076   iEvent.getByLabel(L1ForJets_,l1ForJets);
00077 
00078   typename TCollection::const_iterator jet_iter;
00079   for (jet_iter = jets->begin(); jet_iter != jets->end(); ++jet_iter) {
00080 
00081     bool isMatched=false;
00082 
00083     //std::cout << "FL: l1TauJets.size  = " << l1TauJets->size() << std::endl;
00084     for (unsigned int jetc=0;jetc<l1TauJets->size();++jetc)
00085     {
00086       const double deltaeta=jet_iter->eta()-(*l1TauJets)[jetc].eta();
00087       const double deltaphi=deltaPhi(jet_iter->phi(),(*l1TauJets)[jetc].phi());
00088       //std::cout << "FL: sqrt(2) = " << sqrt(2) << std::endl;
00089       if (sqrt(deltaeta*deltaeta+deltaphi*deltaphi) < DeltaR_) isMatched=true;
00090     }
00091 
00092     for (unsigned int jetc=0;jetc<l1CenJets->size();++jetc)
00093     {
00094       const double deltaeta=jet_iter->eta()-(*l1CenJets)[jetc].eta();
00095       const double deltaphi=deltaPhi(jet_iter->phi(),(*l1CenJets)[jetc].phi());
00096       if (sqrt(deltaeta*deltaeta+deltaphi*deltaphi) < DeltaR_) isMatched=true;
00097     }
00098 
00099     for (unsigned int jetc=0;jetc<l1ForJets->size();++jetc)
00100     {
00101       const double deltaeta=jet_iter->eta()-(*l1ForJets)[jetc].eta();
00102       const double deltaphi=deltaPhi(jet_iter->phi(),(*l1ForJets)[jetc].phi());
00103       if (sqrt(deltaeta*deltaeta+deltaphi*deltaphi) < DeltaR_) isMatched=true;
00104     }
00105 
00106 
00107     if (isMatched==true) result->push_back(*jet_iter);
00108 
00109   } // jet_iter
00110 
00111   iEvent.put( result);
00112 
00113 }
00114 
00115