CMS 3D CMS Logo

TrackerIsolationPt.cc
Go to the documentation of this file.
1 //
2 //
3 
9 #include "CLHEP/Vector/LorentzVector.h"
16 #include <vector>
17 
18 using namespace pat;
19 
22 
25 
27 float TrackerIsolationPt::calculate(const Electron& theElectron,
28  const edm::View<reco::Track>& theTracks,
29  float isoConeElectron) const {
30  return this->calculate(*theElectron.gsfTrack(), theTracks, isoConeElectron);
31 }
32 
33 float TrackerIsolationPt::calculate(const Muon& theMuon,
34  const edm::View<reco::Track>& theTracks,
35  float isoConeMuon) const {
36  return this->calculate(*theMuon.track(), theTracks, isoConeMuon);
37 }
38 
41  const edm::View<reco::Track>& theTracks,
42  float isoCone) const {
43  // initialize some variables
44  float isoPtLepton = 0;
45  const reco::Track *closestTrackDRPt = nullptr, *closestTrackDR = nullptr;
46  float closestDRPt = 10000, closestDR = 10000;
47  // use all these pointless vector conversions because the momenta from tracks
48  // are completely unusable; bah, these math-vectors are worthless!
49  CLHEP::HepLorentzVector lepton(theTrack.px(), theTrack.py(), theTrack.pz(), theTrack.p());
50  for (edm::View<reco::Track>::const_iterator itTrack = theTracks.begin(); itTrack != theTracks.end(); itTrack++) {
51  CLHEP::HepLorentzVector track(itTrack->px(), itTrack->py(), itTrack->pz(), itTrack->p());
52  float dR = lepton.deltaR(track);
53  if (dR < isoCone) {
54  isoPtLepton += track.perp();
55  // find the closest matching track
56  // FIXME: we could association by hits or chi2 to match
57  float pRatio = track.perp() / lepton.perp();
58  if (dR < closestDRPt && pRatio > 0.5 && pRatio < 1.5) {
59  closestDRPt = dR;
60  closestTrackDRPt = &*itTrack;
61  }
62  if (dR < closestDR) {
63  closestDR = dR;
64  closestTrackDR = &*itTrack;
65  }
66  }
67  }
68  if (closestTrackDRPt) {
69  GlobalVector closestTrackVector(closestTrackDRPt->px(), closestTrackDRPt->py(), closestTrackDRPt->pz());
70  isoPtLepton -= closestTrackVector.perp();
71  } else if (closestTrackDR) {
72  GlobalVector closestTrackVector(closestTrackDR->px(), closestTrackDR->py(), closestTrackDR->pz());
73  isoPtLepton -= closestTrackVector.perp();
74  }
75  // back to normal sum - S.L. 30/10/2007
76  if (isoPtLepton < 0)
77  isoPtLepton = 0;
78  // isoPtLepton <= 0.01 ? isoPtLepton = -1 : isoPtLepton = log(isoPtLepton);
79  return isoPtLepton;
80 }
T perp() const
Definition: PV3DBase.h:69
double px() const
x coordinate of momentum vector
Definition: TrackBase.h:640
double p() const
momentum vector magnitude
Definition: TrackBase.h:631
double py() const
y coordinate of momentum vector
Definition: TrackBase.h:643
reco::GsfTrackRef gsfTrack() const override
override the reco::GsfElectron::gsfTrack method, to access the internal storage of the supercluster ...
virtual ~TrackerIsolationPt()
destructor
Definition: HeavyIon.h:7
float calculate(const Electron &theElectron, const edm::View< reco::Track > &theTracks, float isoConeElectron=0.3) const
calculate the TrackIsoPt for the lepton object
Definition: Muon.py:1
reco::TrackRef track() const override
reference to Track reconstructed in the tracker only (reimplemented from reco::Muon) ...
double pz() const
z coordinate of momentum vector
Definition: TrackBase.h:646
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
const_iterator begin() const
const_iterator end() const