CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/RecoTauTag/HLTProducers/src/L2TauJetsProvider.cc

Go to the documentation of this file.
00001 #include "RecoTauTag/HLTProducers/interface/L2TauJetsProvider.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 L2TauJetsProvider::L2TauJetsProvider(const edm::ParameterSet& iConfig)
00016 {
00017   jetSrc = iConfig.getParameter<vtag>("JetSrc");
00018   l1ParticlesTau = iConfig.getParameter<InputTag>("L1ParticlesTau");
00019   l1ParticlesJet = iConfig.getParameter<InputTag>("L1ParticlesJet");
00020   tauTrigger = iConfig.getParameter<InputTag>("L1TauTrigger");
00021   mEt_Min = iConfig.getParameter<double>("EtMin");
00022   
00023   produces<CaloJetCollection>();
00024 }
00025 
00026 L2TauJetsProvider::~L2TauJetsProvider(){ }
00027 
00028 void L2TauJetsProvider::produce(edm::Event& iEvent, const edm::EventSetup& iES)
00029 {
00030 
00031  using namespace edm;
00032  using namespace std;
00033  using namespace reco;
00034  using namespace trigger;
00035  using namespace l1extra;
00036 
00037 
00038  //Getting all the L1Seeds
00039 
00040  
00041  //Getting the Collections of L2ReconstructedJets from L1Seeds
00042  //and removing the collinear jets
00043  myL2L1JetsMap.clear();
00044  int iL1Jet = 0;
00045  for( vtag::const_iterator s = jetSrc.begin(); s != jetSrc.end(); ++ s ) {
00046    edm::Handle<CaloJetCollection> tauJets;
00047    iEvent.getByLabel( * s, tauJets );
00048    CaloJetCollection::const_iterator iTau = tauJets->begin();
00049    if(iTau != tauJets->end()){
00050      //Create a Map to associate to every Jet its L1SeedId, i.e. 0,1,2 or 3
00051      myL2L1JetsMap.insert(std::pair<int, const CaloJet>(iL1Jet, *(iTau)));
00052    }
00053    iL1Jet++;
00054  }
00055  std::auto_ptr<CaloJetCollection> tauL2jets(new CaloJetCollection);
00056 
00057  //Loop over the jetSrc to select the proper jets
00058   
00059   double deltaR = 0.1;
00060   double matchingR = 0.01;
00061   //Loop over the Map to find which jets has fired the trigger
00062   //myL1Tau is the Collection of L1TauCandidates (from 0 to max  4 elements)
00063   //get the list of trigger candidates from the HLTL1SeedGT filter
00064   edm::Handle< L1JetParticleCollection > tauColl ; 
00065   edm::Handle< L1JetParticleCollection > jetColl ; 
00066   
00067   iEvent.getByLabel( l1ParticlesTau, tauColl );
00068   iEvent.getByLabel(l1ParticlesJet, jetColl);
00069 
00070   const L1JetParticleCollection  myL1Tau = *(tauColl.product());
00071 
00072   const L1JetParticleCollection  myL1Jet = *(jetColl.product());  
00073 
00074     L1JetParticleCollection myL1Obj;
00075     myL1Obj.reserve(8);
00076     
00077     for(unsigned int i=0;i<myL1Tau.size();i++)
00078       {
00079         myL1Obj.push_back(myL1Tau[i]);
00080       }
00081     for(unsigned int j=0;j<myL1Jet.size();j++)
00082       {
00083         myL1Obj.push_back(myL1Jet[j]);
00084       }
00085 
00086 
00087     edm::Handle<trigger::TriggerFilterObjectWithRefs> l1TriggeredTaus;
00088     if(iEvent.getByLabel(tauTrigger,l1TriggeredTaus)){
00089     
00090       
00091       tauCandRefVec.clear();
00092       jetCandRefVec.clear();
00093       
00094       l1TriggeredTaus->getObjects( trigger::TriggerL1TauJet,tauCandRefVec);
00095       l1TriggeredTaus->getObjects( trigger::TriggerL1CenJet,jetCandRefVec);
00096       
00097       for( unsigned int iL1Tau=0; iL1Tau <tauCandRefVec.size();iL1Tau++)
00098         {  
00099           for(unsigned int iJet=0;iJet<myL1Obj.size();iJet++)
00100             {
00101               //Find the relative L2TauJets, to see if it has been reconstructed
00102             std::map<int, const reco::CaloJet>::const_iterator myL2itr = myL2L1JetsMap.find(iJet);
00103             if(myL2itr!=myL2L1JetsMap.end()){
00104               //Calculate the DeltaR between L1TauCandidate and L1Tau which fired the trigger
00105               if(&tauCandRefVec[iL1Tau]) 
00106                 deltaR = ROOT::Math::VectorUtil::DeltaR(myL1Obj[iJet].p4().Vect(), (tauCandRefVec[iL1Tau]->p4()).Vect());
00107               if(deltaR < matchingR ) {
00108               //              Getting back from the map the L2TauJet
00109                 const CaloJet myL2TauJet = myL2itr->second;
00110                 if(myL2TauJet.pt() > mEt_Min) tauL2jets->push_back(myL2TauJet);
00111                 myL2L1JetsMap.erase(myL2itr->first);
00112                 break;
00113                 
00114               }
00115             }
00116             
00117           }
00118       }
00119       
00120     for(unsigned int iL1Tau=0; iL1Tau <jetCandRefVec.size();iL1Tau++)
00121       {  
00122         for(unsigned int iJet=0;iJet<myL1Obj.size();iJet++)
00123           {
00124             //Find the relative L2TauJets, to see if it has been reconstructed
00125             std::map<int, const reco::CaloJet>::const_iterator myL2itr = myL2L1JetsMap.find(iJet);
00126             if(myL2itr!=myL2L1JetsMap.end()){
00127               //Calculate the DeltaR between L1TauCandidate and L1Tau which fired the trigger
00128               if(&jetCandRefVec[iL1Tau])
00129                 deltaR = ROOT::Math::VectorUtil::DeltaR(myL1Obj[iJet].p4().Vect(), (jetCandRefVec[iL1Tau]->p4()).Vect());
00130               if(deltaR < matchingR) {
00131                 // Getting back from the map the L2TauJet
00132                 const CaloJet myL2TauJet = myL2itr->second;
00133                 
00134                 if(myL2TauJet.pt() > mEt_Min) tauL2jets->push_back(myL2TauJet);
00135                 myL2L1JetsMap.erase(myL2itr->first);
00136                 break;
00137                 
00138               }
00139             }
00140             
00141           }
00142       }
00143     
00144   }
00145   //  std::cout <<"Size of L2 jets "<<tauL2jets->size()<<std::endl;
00146 
00147   iEvent.put(tauL2jets);
00148 
00149 }