CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch1/src/RecoTauTag/HLTProducers/src/L2TauPixelTrackMatch.cc

Go to the documentation of this file.
00001 // $Id: L2TauPixelTrackMatch.cc,v 1.1 2012/02/02 23:15:20 khotilov Exp $
00002 
00003 #include "RecoTauTag/HLTProducers/interface/L2TauPixelTrackMatch.h"
00004 #include "FWCore/Utilities/interface/Exception.h"
00005 #include "DataFormats/Common/interface/Handle.h"
00006 #include "DataFormats/Common/interface/RefToBase.h"
00007 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00008 #include "DataFormats/TrackReco/interface/Track.h"
00009 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00010 #include "DataFormats/JetReco/interface/CaloJet.h"
00011 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00012 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
00013 #include "DataFormats/HLTReco/interface/TriggerTypeDefs.h"
00014 #include "DataFormats/Math/interface/deltaPhi.h"
00015 
00016 
00017 // all these debug printouts will need to be removed at some point
00018 //#define DBG_PRINT(arg) (arg)
00019 #define DBG_PRINT(arg) 
00020 
00021 
00022 L2TauPixelTrackMatch::L2TauPixelTrackMatch(const edm::ParameterSet& conf)
00023 {
00024   m_jetSrc           = conf.getParameter<edm::InputTag>("JetSrc");
00025   m_jetMinPt         = conf.getParameter<double>("JetMinPt");
00026   m_jetMaxEta        = conf.getParameter<double>("JetMaxEta");
00027   //m_jetMinN        = conf.getParameter<int>("JetMinN");
00028   m_trackSrc         = conf.getParameter<edm::InputTag>("TrackSrc");
00029   m_trackMinPt       = conf.getParameter<double>("TrackMinPt");
00030   m_deltaR           = conf.getParameter<double>("deltaR");
00031   m_beamSpotTag      = conf.getParameter<edm::InputTag>("BeamSpotSrc");
00032 
00033   produces<reco::CaloJetCollection>();
00034 }
00035 
00036 
00037 L2TauPixelTrackMatch::~L2TauPixelTrackMatch(){}
00038 
00039 
00040 void L2TauPixelTrackMatch::produce(edm::Event& ev, const edm::EventSetup& es)
00041 {
00042   using namespace std;
00043   using namespace reco;
00044 
00045   // *** Pick up beam spot ***
00046   
00047   // use beam spot for vertex x,y
00048   edm::Handle<BeamSpot> bsHandle;
00049   ev.getByLabel( m_beamSpotTag, bsHandle);
00050   const reco::BeamSpot & bs = *bsHandle;
00051   math::XYZPoint beam_spot(bs.x0(), bs.y0(), bs.z0());
00052   DBG_PRINT(cout<<endl<<"beamspot "<<beam_spot<<endl);
00053 
00054   // *** Pick up pixel tracks ***
00055 
00056   edm::Handle<TrackCollection> tracksHandle;
00057   ev.getByLabel(m_trackSrc,tracksHandle);
00058 
00059   // *** Pick up L2 tau jets that were previously selected by some other filter ***
00060   
00061   // first, get L2 object refs by the label
00062   edm::Handle<trigger::TriggerFilterObjectWithRefs> jetsHandle;
00063   ev.getByLabel(m_jetSrc, jetsHandle);
00064 
00065   // now we can get pre-selected L2 tau jets
00066   std::vector<CaloJetRef> tau_jets;
00067   jetsHandle->getObjects(trigger::TriggerTau, tau_jets);
00068   const size_t n_jets = tau_jets.size();
00069 
00070   // *** Selects interesting tracks ***
00071 
00072   vector<TinyTrack> good_tracks;
00073   for(TrackCollection::const_iterator itrk = tracksHandle->begin(); itrk != tracksHandle->end(); ++itrk)
00074   {
00075     if (itrk->pt() < m_trackMinPt) continue;
00076     if ( std::abs(itrk->eta()) > m_jetMaxEta + m_deltaR ) continue;
00077        
00078     TinyTrack trk;
00079     trk.pt = itrk->pt();
00080     trk.phi = itrk->phi();
00081     trk.eta = itrk->eta();
00082     double dz = itrk->dz(beam_spot);
00083     trk.vtx = math::XYZPoint(bs.x(dz), bs.y(dz), dz);
00084     good_tracks.push_back(trk);
00085   }
00086   DBG_PRINT(cout<<"got "<<good_tracks.size()<<" good tracks;   "<<n_jets<<" tau jets"<<endl);
00087 
00088 
00089   // *** Match tau jets to intertesting tracks  and assign them new vertices ***
00090 
00091   // the new product
00092   std::auto_ptr<CaloJetCollection> new_tau_jets(new CaloJetCollection);
00093   int n_uniq = 0;
00094   if (good_tracks.size()) for (size_t i=0; i < n_jets; ++i)
00095   {
00096     reco::CaloJetRef jet = tau_jets[i];
00097     if ( jet->pt() < m_jetMinPt || std::abs(jet->eta()) > m_jetMaxEta ) continue;
00098 
00099     DBG_PRINT(cout<<i<<" :"<<endl);
00100     
00101     size_t n0 = new_tau_jets->size();
00102     
00103     for(vector<TinyTrack>::const_iterator itrk = good_tracks.begin(); itrk != good_tracks.end(); ++itrk)
00104     {
00105       DBG_PRINT(cout<<"  trk pt,eta,phi,z: "<<itrk->pt<<" "<<itrk->eta<<" "<<itrk->phi<<" "<<itrk->vtx.z()<<" \t\t ");
00106 
00107       math::XYZTLorentzVector new_jet_dir = Jet::physicsP4(itrk->vtx, *jet, itrk->vtx);
00108       float dphi = reco::deltaPhi(new_jet_dir.phi(), itrk->phi);
00109       float deta = new_jet_dir.eta() - itrk->eta;
00110       
00111       DBG_PRINT(cout<<" jet pt,deta,dphi,dr: "<<jet->pt()<<" "<<deta<<" "<<dphi<<" "<<sqrt(dphi*dphi + deta*deta)<<endl);
00112       
00113       if ( dphi*dphi + deta*deta > m_deltaR*m_deltaR ) continue;
00114       
00115       DBG_PRINT(cout<<"  jet-trk match!"<<endl);
00116       
00117       // create a jet copy and assign a new vertex to it
00118       CaloJet new_jet = *jet;
00119       new_jet.setVertex(itrk->vtx);
00120       
00121       new_tau_jets->push_back(new_jet);
00122     }
00123     DBG_PRINT(cout<<"  nmatchedjets "<<new_tau_jets->size() - n0<<endl);
00124     if (new_tau_jets->size() - n0 > 0 ) n_uniq++;
00125     
00127   }
00128   DBG_PRINT(cout<<"n_uniq_matched_jets "<<n_uniq<<endl<<"storing njets "<<new_tau_jets->size()<<endl);
00129   
00130   // store the result
00131   ev.put(new_tau_jets);
00132 }