Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:24:07

0001 //
0002 //
0003 
0004 #include "PhysicsTools/PatUtils/interface/TrackerIsolationPt.h"
0005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0006 #include "FWCore/Utilities/interface/Exception.h"
0007 #include "FWCore/Framework/interface/ESHandle.h"
0008 #include "DataFormats/Common/interface/Handle.h"
0009 #include "CLHEP/Vector/LorentzVector.h"
0010 #include "DataFormats/GeometryVector/interface/GlobalVector.h"
0011 #include "FWCore/Utilities/interface/InputTag.h"
0012 #include "DataFormats/Common/interface/View.h"
0013 #include "DataFormats/PatCandidates/interface/Electron.h"
0014 #include "DataFormats/PatCandidates/interface/Muon.h"
0015 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
0016 #include <vector>
0017 
0018 using namespace pat;
0019 
0020 /// constructor
0021 TrackerIsolationPt::TrackerIsolationPt() {}
0022 
0023 /// destructor
0024 TrackerIsolationPt::~TrackerIsolationPt() {}
0025 
0026 /// calculate the TrackIsoPt for the lepton object
0027 float TrackerIsolationPt::calculate(const Electron& theElectron,
0028                                     const edm::View<reco::Track>& theTracks,
0029                                     float isoConeElectron) const {
0030   return this->calculate(*theElectron.gsfTrack(), theTracks, isoConeElectron);
0031 }
0032 
0033 float TrackerIsolationPt::calculate(const Muon& theMuon,
0034                                     const edm::View<reco::Track>& theTracks,
0035                                     float isoConeMuon) const {
0036   return this->calculate(*theMuon.track(), theTracks, isoConeMuon);
0037 }
0038 
0039 /// calculate the TrackIsoPt for the lepton's track
0040 float TrackerIsolationPt::calculate(const reco::Track& theTrack,
0041                                     const edm::View<reco::Track>& theTracks,
0042                                     float isoCone) const {
0043   // initialize some variables
0044   float isoPtLepton = 0;
0045   const reco::Track *closestTrackDRPt = nullptr, *closestTrackDR = nullptr;
0046   float closestDRPt = 10000, closestDR = 10000;
0047   // use all these pointless vector conversions because the momenta from tracks
0048   // are completely unusable; bah, these math-vectors are worthless!
0049   CLHEP::HepLorentzVector lepton(theTrack.px(), theTrack.py(), theTrack.pz(), theTrack.p());
0050   for (edm::View<reco::Track>::const_iterator itTrack = theTracks.begin(); itTrack != theTracks.end(); itTrack++) {
0051     CLHEP::HepLorentzVector track(itTrack->px(), itTrack->py(), itTrack->pz(), itTrack->p());
0052     float dR = lepton.deltaR(track);
0053     if (dR < isoCone) {
0054       isoPtLepton += track.perp();
0055       // find the closest matching track
0056       // FIXME: we could association by hits or chi2 to match
0057       float pRatio = track.perp() / lepton.perp();
0058       if (dR < closestDRPt && pRatio > 0.5 && pRatio < 1.5) {
0059         closestDRPt = dR;
0060         closestTrackDRPt = &*itTrack;
0061       }
0062       if (dR < closestDR) {
0063         closestDR = dR;
0064         closestTrackDR = &*itTrack;
0065       }
0066     }
0067   }
0068   if (closestTrackDRPt) {
0069     GlobalVector closestTrackVector(closestTrackDRPt->px(), closestTrackDRPt->py(), closestTrackDRPt->pz());
0070     isoPtLepton -= closestTrackVector.perp();
0071   } else if (closestTrackDR) {
0072     GlobalVector closestTrackVector(closestTrackDR->px(), closestTrackDR->py(), closestTrackDR->pz());
0073     isoPtLepton -= closestTrackVector.perp();
0074   }
0075   // back to normal sum - S.L. 30/10/2007
0076   if (isoPtLepton < 0)
0077     isoPtLepton = 0;
0078   //  isoPtLepton <= 0.01 ? isoPtLepton = -1 : isoPtLepton = log(isoPtLepton);
0079   return isoPtLepton;
0080 }