CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC4_patch1/src/RecoTauTag/HLTProducers/src/L2TauJetsMerger.cc

Go to the documentation of this file.
00001 #include "RecoTauTag/HLTProducers/interface/L2TauJetsMerger.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 L2TauJetsMerger::L2TauJetsMerger(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 L2TauJetsMerger::~L2TauJetsMerger(){ }
00027 
00028 void L2TauJetsMerger::produce(edm::Event& iEvent, const edm::EventSetup& iES)
00029 {
00030 
00031  using namespace edm;
00032  using namespace std;
00033  using namespace reco;
00034 
00035  //Getting all the L1Seeds
00036 
00037  
00038  //Getting the Collections of L2ReconstructedJets from L1Seeds
00039  //and removing the collinear jets
00040  myL2L1JetsMap.clear();
00041  CaloJetCollection myTmpJets;
00042  myTmpJets.clear();
00043 
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    for(CaloJetCollection::const_iterator iTau = tauJets->begin();iTau !=tauJets->end();iTau++)
00049      { 
00050      //Create a Map to associate to every Jet its L1SeedId, i.e. 0,1,2 or 3
00051        if(iTau->et() > mEt_Min) {
00052 
00053          //Add the Pdg Id here 
00054          CaloJet myJet = *iTau;
00055          myJet.setPdgId(15);
00056          myTmpJets.push_back(myJet);
00057        }
00058      }
00059    iL1Jet++;
00060  }
00061 
00062  std::auto_ptr<CaloJetCollection> tauL2jets(new CaloJetCollection); 
00063 
00064  //Removing the collinear jets correctly!
00065 
00066  //First sort the jets you have merged
00067  SorterByPt sorter;
00068  std::sort(myTmpJets.begin(),myTmpJets.end(),sorter);
00069  
00070 //Remove Collinear Jets by prefering the highest ones!
00071 
00072    while(myTmpJets.size()>0) {
00073      tauL2jets->push_back(myTmpJets.at(0));
00074      CaloJetCollection tmp;
00075      for(unsigned int i=1 ;i<myTmpJets.size();++i) {
00076        double DR = ROOT::Math::VectorUtil::DeltaR(myTmpJets.at(0).p4(),myTmpJets.at(i).p4());
00077        if(DR>0.1) 
00078          tmp.push_back(myTmpJets.at(i));
00079      }
00080      myTmpJets.swap(tmp);
00081      tmp.clear();
00082    }
00083 
00084 
00085   iEvent.put(tauL2jets);
00086 
00087 }