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