CMS 3D CMS Logo

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

Go to the documentation of this file.
00001 #include "RecoTauTag/HLTProducers/interface/HLTTauProducer.h"
00002 #include "DataFormats/TrackReco/interface/Track.h"
00003 //
00004 // class decleration
00005 //
00006 
00007 
00008 HLTTauProducer::HLTTauProducer(const edm::ParameterSet& iConfig)
00009 {
00010   emIsolatedJetsL2_ = iConfig.getParameter<edm::InputTag>("L2EcalIsoJets");
00011   trackIsolatedJetsL25_ = iConfig.getParameter<edm::InputTag>("L25TrackIsoJets");
00012   trackIsolatedJetsL3_ = iConfig.getParameter<edm::InputTag>("L3TrackIsoJets");
00013   matchingCone_ = iConfig.getParameter<double>("MatchingCone");
00014   signalCone_ = iConfig.getParameter<double>("SignalCone");
00015   isolationCone_ = iConfig.getParameter<double>("IsolationCone");
00016   ptMin_ = iConfig.getParameter<double>("MinPtTracks");
00017   produces<reco::HLTTauCollection>();
00018 }
00019 
00020 HLTTauProducer::~HLTTauProducer(){ }
00021 
00022 void HLTTauProducer::produce(edm::Event& iEvent, const edm::EventSetup& iES)
00023 {
00024 
00025   using namespace reco;
00026   using namespace edm;
00027   using namespace std;
00028   
00029 
00030   HLTTauCollection * jetCollection = new HLTTauCollection;
00031  
00032   edm::Handle<L2TauInfoAssociation> tauL2Jets;
00033   iEvent.getByLabel(emIsolatedJetsL2_ , tauL2Jets );
00034 
00035   edm::Handle<IsolatedTauTagInfoCollection> tauL25Jets;
00036   iEvent.getByLabel(trackIsolatedJetsL25_, tauL25Jets );
00037   
00038   edm::Handle<IsolatedTauTagInfoCollection> tauL3Jets;
00039   iEvent.getByLabel(trackIsolatedJetsL3_, tauL3Jets );
00040 
00041   IsolatedTauTagInfoCollection tauL25 = *(tauL25Jets.product());
00042   IsolatedTauTagInfoCollection tauL3 = *(tauL3Jets.product());
00043   
00044   int i=0;
00045   float eta_, phi_, pt_;
00046   int nTracksL25, nTracksL3;
00047     float sumPtTracksL25 = 1000.;
00048     float sumPtTracksL3 = 1000.;
00049     double ptLeadTkL25=0.;
00050     double ptLeadTkL3=0.;
00051  for(L2TauInfoAssociation::const_iterator p = tauL2Jets->begin();p!=tauL2Jets->end();++p)
00052            {
00053              //Retrieve The L2TauIsolationInfo Class from the AssociationMap
00054              const L2TauIsolationInfo l2info = p->val;
00055              //Retrieve the Jet
00056              //      const CaloJet& jet =*(p->key);
00057              
00058 
00059              double emIsol  = l2info.ecalIsolEt();
00060 
00061     JetTracksAssociationRef jetTracks = tauL25[i].jtaRef();
00062     math::XYZVector jetDirL25(jetTracks->first->px(),jetTracks->first->py(),jetTracks->first->pz());   
00063     eta_ = jetDirL25.eta();
00064     phi_ = jetDirL25.phi();
00065     pt_ = jetTracks->first->pt();
00066 
00067     int trackIsolationL25 = (int)tauL25[i].discriminator(jetDirL25,matchingCone_, signalCone_, isolationCone_,1.,ptMin_,0);
00068     const TrackRef leadTkL25 = tauL25[i].leadingSignalTrack(jetDirL25,matchingCone_, 1.);
00069     ptLeadTkL25 = 0.;
00070     nTracksL25 = 1000;
00071     if(!leadTkL25) 
00072       {}else{
00073         ptLeadTkL25 = (*leadTkL25).pt();
00074         const TrackRefVector signalTracks = tauL25[i].tracksInCone((*leadTkL25).momentum(), signalCone_, ptMin_ );
00075         const TrackRefVector isolationTracks = tauL25[i].tracksInCone((*leadTkL25).momentum(), isolationCone_,   ptMin_ );
00076         nTracksL25 = isolationTracks.size() - signalTracks.size();
00077 
00078         for(unsigned int j=0;j<isolationTracks.size();j++)
00079           sumPtTracksL25 = sumPtTracksL25 + isolationTracks[j]->pt();
00080         for(unsigned int j=0;j<signalTracks.size();j++)
00081           sumPtTracksL25 = sumPtTracksL25 - signalTracks[j]->pt();
00082 
00083       }
00084     jetTracks = tauL3[i].jtaRef();
00085     math::XYZVector jetDirL3(jetTracks->first->px(),jetTracks->first->py(),jetTracks->first->pz());     
00086     int  trackIsolationL3 = (int)tauL3[i].discriminator(jetDirL3,matchingCone_, signalCone_, isolationCone_,1.,ptMin_,0);
00087     
00088     const TrackRef leadTkL3 = tauL3[i].leadingSignalTrack(jetDirL3,matchingCone_,1.);
00089     ptLeadTkL3=0.;
00090     nTracksL3 = 1000;
00091     if(!leadTkL3) 
00092       {}else{
00093         ptLeadTkL3 = (*leadTkL3).pt();
00094         const TrackRefVector signalTracks = tauL3[i].tracksInCone((*leadTkL25).momentum(), signalCone_, ptMin_ );
00095         const TrackRefVector isolationTracks = tauL3[i].tracksInCone((*leadTkL25).momentum(), isolationCone_,   ptMin_ );
00096         nTracksL3 = isolationTracks.size() - signalTracks.size();
00097         float sumPtTracksL3 = 0.;
00098         for(unsigned int j=0;j<isolationTracks.size();j++)
00099           sumPtTracksL3 = sumPtTracksL3 + isolationTracks[j]->pt();
00100         for(unsigned int j=0;j<signalTracks.size();j++)
00101           sumPtTracksL3 = sumPtTracksL3 - signalTracks[j]->pt();
00102       }
00103 
00104     HLTTau pippo(eta_,phi_,pt_,emIsol,trackIsolationL25,ptLeadTkL25,trackIsolationL3,ptLeadTkL3);
00105     pippo.setNL25TrackIsolation(nTracksL25);
00106     pippo.setNL3TrackIsolation(nTracksL3);
00107     pippo.setSumPtTracksL25(sumPtTracksL25);
00108     pippo.setSumPtTracksL3(sumPtTracksL3);
00109     pippo.setSeedEcalHitEt(l2info.seedEcalHitEt());
00110     pippo.setEcalClusterShape(l2info.ecalClusterShape());
00111     pippo.setNEcalHits(l2info.nEcalHits());                    
00112     pippo.setHcalIsolEt(l2info.hcalIsolEt());
00113     pippo.setSeedHcalHitEt(l2info.seedHcalHitEt());
00114     pippo.setHcalClusterShape(l2info.hcalClusterShape());
00115     pippo.setNHcalHits(l2info.nHcalHits());
00116     jetCollection->push_back(pippo);
00117       i++;
00118   }
00119   
00120   std::auto_ptr<reco::HLTTauCollection> selectedTaus(jetCollection);
00121   
00122   iEvent.put(selectedTaus);
00123   
00124 
00125 
00126 
00127 }