CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/PhysicsTools/PatUtils/src/TrackerIsolationPt.cc

Go to the documentation of this file.
00001 //
00002 // $Id: TrackerIsolationPt.cc,v 1.6 2010/02/11 00:13:22 wmtan Exp $
00003 //
00004 
00005 #include "PhysicsTools/PatUtils/interface/TrackerIsolationPt.h"
00006 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00007 #include "FWCore/Utilities/interface/Exception.h"
00008 #include "FWCore/Framework/interface/ESHandle.h"
00009 #include "DataFormats/Common/interface/Handle.h"
00010 #include "CLHEP/Vector/LorentzVector.h"
00011 #include "DataFormats/GeometryVector/interface/GlobalVector.h"
00012 #include "FWCore/Utilities/interface/InputTag.h"
00013 #include "DataFormats/Common/interface/View.h"
00014 #include "DataFormats/PatCandidates/interface/Electron.h"
00015 #include "DataFormats/PatCandidates/interface/Muon.h"
00016 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
00017 #include <vector>
00018 
00019 using namespace pat;
00020 
00022 TrackerIsolationPt::TrackerIsolationPt() {
00023 }
00024 
00026 TrackerIsolationPt::~TrackerIsolationPt() {
00027 }
00028 
00030 float TrackerIsolationPt::calculate(const Electron & theElectron, const edm::View<reco::Track> & theTracks, float isoConeElectron) const {
00031   return this->calculate(*theElectron.gsfTrack(), theTracks, isoConeElectron);
00032 }
00033 
00034 float TrackerIsolationPt::calculate(const Muon & theMuon, const edm::View<reco::Track> & theTracks, float isoConeMuon) const {
00035   return this->calculate(*theMuon.track(), theTracks, isoConeMuon);
00036 }
00037 
00039 float TrackerIsolationPt::calculate(const reco::Track & theTrack, const edm::View<reco::Track> & theTracks, float isoCone) const {
00040   // initialize some variables
00041   float isoPtLepton = 0;
00042   const reco::Track * closestTrackDRPt = 0, * closestTrackDR = 0;
00043   float closestDRPt = 10000, closestDR = 10000;
00044   // use all these pointless vector conversions because the momenta from tracks
00045   // are completely unusable; bah, these math-vectors are worthless!
00046   CLHEP::HepLorentzVector lepton(theTrack.px(), theTrack.py(), theTrack.pz(), theTrack.p());
00047   for (edm::View<reco::Track>::const_iterator itTrack = theTracks.begin(); itTrack != theTracks.end(); itTrack++) {
00048     CLHEP::HepLorentzVector track(itTrack->px(), itTrack->py(), itTrack->pz(), itTrack->p());
00049     float dR = lepton.deltaR(track);
00050     if (dR < isoCone) {
00051       isoPtLepton += track.perp();
00052       // find the closest matching track
00053       // FIXME: we could association by hits or chi2 to match
00054       float pRatio = track.perp()/lepton.perp();
00055       if (dR < closestDRPt && pRatio > 0.5 && pRatio < 1.5) {
00056         closestDRPt = dR;
00057         closestTrackDRPt = &*itTrack;
00058       }
00059       if (dR < closestDR) {
00060         closestDR = dR;
00061         closestTrackDR = &*itTrack;
00062       }
00063     }
00064   }
00065   if (closestTrackDRPt) {
00066     GlobalVector closestTrackVector(closestTrackDRPt->px(), closestTrackDRPt->py(), closestTrackDRPt->pz());
00067     isoPtLepton -= closestTrackVector.perp();
00068   } else if (closestTrackDR) {
00069     GlobalVector closestTrackVector(closestTrackDR->px(), closestTrackDR->py(), closestTrackDR->pz());
00070     isoPtLepton -= closestTrackVector.perp();
00071   }
00072   // back to normal sum - S.L. 30/10/2007
00073   if (isoPtLepton<0) isoPtLepton = 0;
00074   //  isoPtLepton <= 0.01 ? isoPtLepton = -1 : isoPtLepton = log(isoPtLepton);
00075   return isoPtLepton;
00076 }
00077