CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch13/src/RecoTauTag/HLTProducers/src/L1HLTJetsMatching.cc

Go to the documentation of this file.
00001 #include "RecoTauTag/HLTProducers/interface/L1HLTJetsMatching.h"
00002 #include "Math/GenVector/VectorUtil.h"
00003 #include "DataFormats/L1Trigger/interface/L1JetParticle.h"
00004 #include "DataFormats/L1Trigger/interface/L1JetParticleFwd.h"
00005 #include "DataFormats/HLTReco/interface/TriggerTypeDefs.h"
00006 #include "FWCore/Utilities/interface/EDMException.h"
00007 //
00008 // class decleration
00009 //
00010 using namespace reco;
00011 using namespace std;
00012 using namespace edm;
00013 using namespace l1extra;
00014 
00015 L1HLTJetsMatching::L1HLTJetsMatching(const edm::ParameterSet& iConfig)
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 }
00023 
00024 L1HLTJetsMatching::~L1HLTJetsMatching(){ }
00025 
00026 void L1HLTJetsMatching::produce(edm::Event& iEvent, const edm::EventSetup& iES)
00027 {
00028         
00029         using namespace edm;
00030         using namespace std;
00031         using namespace reco;
00032         using namespace trigger;
00033         using namespace l1extra;
00034         
00035         typedef std::vector<LeafCandidate> LeafCandidateCollection;
00036         
00037         auto_ptr<CaloJetCollection> tauL2jets(new CaloJetCollection);
00038         
00039         double deltaR = 1.0;
00040         double matchingR = 0.5;
00041         //Getting HLT jets to be matched
00042         edm::Handle<edm::View<Candidate> > tauJets;
00043         iEvent.getByLabel( jetSrc, tauJets );
00044 
00045 //              std::cout <<"Size of input jet collection "<<tauJets->size()<<std::endl;
00046                 
00047         edm::Handle<trigger::TriggerFilterObjectWithRefs> l1TriggeredTaus;
00048     iEvent.getByLabel(tauTrigger,l1TriggeredTaus);
00049                 
00050                 
00051                 tauCandRefVec.clear();
00052                 jetCandRefVec.clear();
00053                 
00054                 l1TriggeredTaus->getObjects( trigger::TriggerL1TauJet,tauCandRefVec);
00055                 l1TriggeredTaus->getObjects( trigger::TriggerL1CenJet,jetCandRefVec);
00056                 math::XYZPoint a(0.,0.,0.);
00057                 CaloJet::Specific f;
00058                 
00059                 for( unsigned int iL1Tau=0; iL1Tau <tauCandRefVec.size();iL1Tau++)
00060                 {  
00061                         for(unsigned int iJet=0;iJet<tauJets->size();iJet++)
00062                         {
00063                                 //Find the relative L2TauJets, to see if it has been reconstructed
00064                                 const Candidate &  myJet = (*tauJets)[iJet];
00065                                 deltaR = ROOT::Math::VectorUtil::DeltaR(myJet.p4().Vect(), (tauCandRefVec[iL1Tau]->p4()).Vect());
00066                                 if(deltaR < matchingR ) {
00067                                         //               LeafCandidate myLC(myJet);
00068                                         CaloJet myCaloJet(myJet.p4(),a,f);
00069                                         if(myJet.pt() > mEt_Min) {
00070                                                 //                tauL2LC->push_back(myLC);
00071                                                 tauL2jets->push_back(myCaloJet);
00072                                         }
00073                                         break;
00074                                 }
00075                         }
00076                 }  
00077                 
00078                 for(unsigned int iL1Tau=0; iL1Tau <jetCandRefVec.size();iL1Tau++)
00079                 {  
00080                         for(unsigned int iJet=0;iJet<tauJets->size();iJet++)
00081                         {
00082                                 const Candidate & myJet = (*tauJets)[iJet];
00083                                 //Find the relative L2TauJets, to see if it has been reconstructed
00084                                 deltaR = ROOT::Math::VectorUtil::DeltaR(myJet.p4().Vect(), (jetCandRefVec[iL1Tau]->p4()).Vect());
00085                                 if(deltaR < matchingR ) {
00086                                         //               LeafCandidate myLC(myJet);
00087                                         CaloJet myCaloJet(myJet.p4(),a,f);
00088                                         if(myJet.pt() > mEt_Min) {
00089                                                 //tauL2LC->push_back(myLC);
00090                                                 tauL2jets->push_back(myCaloJet);
00091                                         }
00092                                         break;
00093                                 }
00094                         }
00095                 }
00096         
00097 
00098 //std::cout <<"Size of L1HLT matched jets "<<tauL2jets->size()<<std::endl;
00099 
00100 iEvent.put(tauL2jets);
00101 // iEvent.put(tauL2LC);
00102 }