![]() |
![]() |
#include <RecoTauTag/HLTProducers/interface/L1HLTJetsMatching.h>
Public Member Functions | |
L1HLTJetsMatching (const edm::ParameterSet &) | |
virtual void | produce (edm::Event &, const edm::EventSetup &) |
~L1HLTJetsMatching () | |
Private Attributes | |
std::vector < l1extra::L1JetParticleRef > | jetCandRefVec |
edm::InputTag | jetSrc |
double | mEt_Min |
std::vector < l1extra::L1JetParticleRef > | objL1CandRefVec |
l1extra::L1JetParticleRef | tauCandRef |
std::vector < l1extra::L1JetParticleRef > | tauCandRefVec |
edm::InputTag | tauTrigger |
Definition at line 21 of file L1HLTJetsMatching.h.
L1HLTJetsMatching::L1HLTJetsMatching | ( | const edm::ParameterSet & | iConfig | ) | [explicit] |
Definition at line 15 of file L1HLTJetsMatching.cc.
References edm::ParameterSet::getParameter(), jetSrc, mEt_Min, and tauTrigger.
00016 { 00017 jetSrc = iConfig.getParameter<InputTag>("JetSrc"); 00018 tauTrigger = iConfig.getParameter<InputTag>("L1TauTrigger"); 00019 mEt_Min = iConfig.getParameter<double>("EtMin"); 00020 00021 produces<CaloJetCollection>(); 00022 }
L1HLTJetsMatching::~L1HLTJetsMatching | ( | ) |
void L1HLTJetsMatching::produce | ( | edm::Event & | iEvent, | |
const edm::EventSetup & | iES | |||
) | [virtual] |
Implements edm::EDProducer.
Definition at line 26 of file L1HLTJetsMatching.cc.
References deltaR(), edm::Event::getByLabel(), jetCandRefVec, jetSrc, mEt_Min, p4, edm::Handle< T >::product(), edm::Event::put(), HcalSimpleRecAlgoImpl::reco(), std, tauCandRefVec, tauTrigger, trigger::TriggerL1CenJet, and trigger::TriggerL1TauJet.
00027 { 00028 00029 using namespace edm; 00030 using namespace std; 00031 using namespace reco; 00032 using namespace trigger; 00033 using namespace l1extra; 00034 00035 00036 auto_ptr<CaloJetCollection> tauL2jets(new CaloJetCollection); 00037 00038 double deltaR = 1.0; 00039 double matchingR = 0.5; 00040 //Getting HLT jets to be matched 00041 edm::Handle<CaloJetCollection> tauJets; 00042 if(iEvent.getByLabel( jetSrc, tauJets )) 00043 { 00044 const CaloJetCollection myJets = *(tauJets.product()); 00045 00046 Handle<trigger::TriggerFilterObjectWithRefs> l1TriggeredTaus; 00047 if(iEvent.getByLabel(tauTrigger,l1TriggeredTaus)){ 00048 00049 00050 tauCandRefVec.clear(); 00051 jetCandRefVec.clear(); 00052 00053 l1TriggeredTaus->getObjects( trigger::TriggerL1TauJet,tauCandRefVec); 00054 l1TriggeredTaus->getObjects( trigger::TriggerL1CenJet,jetCandRefVec); 00055 00056 for( unsigned int iL1Tau=0; iL1Tau <tauCandRefVec.size();iL1Tau++) 00057 { 00058 for(unsigned int iJet=0;iJet<myJets.size();iJet++) 00059 { 00060 //Find the relative L2TauJets, to see if it has been reconstructed 00061 deltaR = ROOT::Math::VectorUtil::DeltaR(myJets[iJet].p4().Vect(), (tauCandRefVec[iL1Tau]->p4()).Vect()); 00062 if(deltaR < matchingR ) { 00063 if(myJets[iJet].pt() > mEt_Min) tauL2jets->push_back(myJets[iJet]); 00064 break; 00065 } 00066 } 00067 } 00068 00069 for(unsigned int iL1Tau=0; iL1Tau <jetCandRefVec.size();iL1Tau++) 00070 { 00071 for(unsigned int iJet=0;iJet<myJets.size();iJet++) 00072 { 00073 //Find the relative L2TauJets, to see if it has been reconstructed 00074 deltaR = ROOT::Math::VectorUtil::DeltaR(myJets[iJet].p4().Vect(), (jetCandRefVec[iL1Tau]->p4()).Vect()); 00075 if(deltaR < matchingR ) { 00076 if(myJets[iJet].pt() > mEt_Min) tauL2jets->push_back(myJets[iJet]); 00077 break; 00078 } 00079 } 00080 } 00081 } 00082 00083 } 00084 00085 // cout <<"Size of L2 jets "<<tauL2jets->size()<<endl; 00086 00087 iEvent.put(tauL2jets); 00088 00089 }
std::vector<l1extra::L1JetParticleRef> L1HLTJetsMatching::jetCandRefVec [private] |
edm::InputTag L1HLTJetsMatching::jetSrc [private] |
Definition at line 33 of file L1HLTJetsMatching.h.
Referenced by L1HLTJetsMatching(), and produce().
double L1HLTJetsMatching::mEt_Min [private] |
Definition at line 35 of file L1HLTJetsMatching.h.
Referenced by L1HLTJetsMatching(), and produce().
std::vector<l1extra::L1JetParticleRef> L1HLTJetsMatching::objL1CandRefVec [private] |
Definition at line 30 of file L1HLTJetsMatching.h.
Definition at line 31 of file L1HLTJetsMatching.h.
std::vector<l1extra::L1JetParticleRef> L1HLTJetsMatching::tauCandRefVec [private] |
edm::InputTag L1HLTJetsMatching::tauTrigger [private] |
Definition at line 34 of file L1HLTJetsMatching.h.
Referenced by L1HLTJetsMatching(), and produce().