CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
pat::LeptonJetIsolationAngle Class Reference

Calculates a lepton's jet isolation angle. More...

#include "PhysicsTools/PatUtils/interface/LeptonJetIsolationAngle.h"

Public Member Functions

float calculate (const Electron &anElectron, const edm::Handle< edm::View< reco::Track > > &trackHandle, const edm::Event &iEvent)
 
float calculate (const Muon &aMuon, const edm::Handle< edm::View< reco::Track > > &trackHandle, const edm::Event &iEvent)
 
 LeptonJetIsolationAngle ()
 
 ~LeptonJetIsolationAngle ()
 

Private Member Functions

float calculate (const CLHEP::HepLorentzVector &aLepton, const edm::Handle< edm::View< reco::Track > > &trackHandle, const edm::Event &iEvent)
 
float spaceAngle (const CLHEP::HepLorentzVector &aLepton, const reco::CaloJet &aJet)
 

Private Attributes

TrackerIsolationPt trkIsolator_
 

Detailed Description

Calculates a lepton's jet isolation angle.

LeptonJetIsolationAngle calculates an isolation angle w.r.t. a list of given jets as the minimal angle to a jet in Euclidean space, as defined in CMS Note 2006/024

Author
Steven Lowette
Version
Id:
LeptonJetIsolationAngle.h,v 1.3 2008/03/05 14:51:02 fronga Exp

Definition at line 34 of file LeptonJetIsolationAngle.h.

Constructor & Destructor Documentation

LeptonJetIsolationAngle::LeptonJetIsolationAngle ( )

Definition at line 18 of file LeptonJetIsolationAngle.cc.

18  {
19 }
LeptonJetIsolationAngle::~LeptonJetIsolationAngle ( )

Definition at line 23 of file LeptonJetIsolationAngle.cc.

23  {
24 }

Member Function Documentation

float LeptonJetIsolationAngle::calculate ( const Electron anElectron,
const edm::Handle< edm::View< reco::Track > > &  trackHandle,
const edm::Event iEvent 
)

Definition at line 28 of file LeptonJetIsolationAngle.cc.

References reco::LeafCandidate::energy(), reco::LeafCandidate::px(), reco::LeafCandidate::py(), and reco::LeafCandidate::pz().

Referenced by calculate().

28  {
29  CLHEP::HepLorentzVector theElectronHLV(theElectron.px(), theElectron.py(), theElectron.pz(), theElectron.energy());
30  return this->calculate(theElectronHLV, trackHandle, iEvent);
31 }
float calculate(const Electron &anElectron, const edm::Handle< edm::View< reco::Track > > &trackHandle, const edm::Event &iEvent)
float LeptonJetIsolationAngle::calculate ( const Muon aMuon,
const edm::Handle< edm::View< reco::Track > > &  trackHandle,
const edm::Event iEvent 
)

Definition at line 32 of file LeptonJetIsolationAngle.cc.

References calculate(), reco::LeafCandidate::energy(), reco::LeafCandidate::px(), reco::LeafCandidate::py(), and reco::LeafCandidate::pz().

32  {
33  CLHEP::HepLorentzVector theMuonHLV(theMuon.px(), theMuon.py(), theMuon.pz(), theMuon.energy());
34  return this->calculate(theMuonHLV, trackHandle, iEvent);
35 }
float calculate(const Electron &anElectron, const edm::Handle< edm::View< reco::Track > > &trackHandle, const edm::Event &iEvent)
float LeptonJetIsolationAngle::calculate ( const CLHEP::HepLorentzVector &  aLepton,
const edm::Handle< edm::View< reco::Track > > &  trackHandle,
const edm::Event iEvent 
)
private

Definition at line 39 of file LeptonJetIsolationAngle.cc.

References pat::TrackerIsolationPt::calculate(), Geom::deltaR2(), edm::Event::getByLabel(), edm::Handle< T >::product(), reco::LeafCandidate::pt(), spaceAngle(), mathSSE::sqrt(), and trkIsolator_.

39  {
40  // FIXME: this is an ugly temporary workaround, JetMET+egamma should come up with a better tool
41  // retrieve the jets
43  iEvent.getByLabel("iterativeCone5CaloJets", jetHandle);
44  reco::CaloJetCollection jetColl = *(jetHandle.product());
45  // retrieve the electrons which might be in the jet list
47  iEvent.getByLabel("pixelMatchGsfElectrons", electronsHandle);
48  std::vector<reco::GsfElectron> electrons = *electronsHandle;
49  // determine the set of isolated electrons
50  std::vector<Electron> isoElectrons;
51  for (size_t ie=0; ie<electrons.size(); ie++) {
52  Electron anElectron(electrons[ie]);
53  if (anElectron.pt() > 10 &&
54  trkIsolator_.calculate(anElectron, *trackHandle) < 3.0) {
55  isoElectrons.push_back(electrons[ie]);
56  }
57  }
58  // determine the collections of jets, cleaned from electrons
59  std::vector<reco::CaloJet> theJets;
60  for (reco::CaloJetCollection::const_iterator itJet = jetColl.begin(); itJet != jetColl.end(); itJet++) {
61  float mindr2 = 9999.;
62  for (size_t ie = 0; ie < isoElectrons.size(); ie++) {
63  float dr2 = ::deltaR2(*itJet, isoElectrons[ie]);
64  if (dr2 < mindr2) mindr2 = dr2;
65  }
66  float mindr = std::sqrt(mindr2);
67  // yes, all cuts hardcoded buts, but it's a second-order effect
68  if (itJet->et() > 15 && mindr > 0.3) theJets.push_back(reco::CaloJet(*itJet));
69  }
70  // calculate finally the isolation angle
71  float isoAngle = 1000; // default to some craze impossible number to inhibit compiler warnings
72  for (std::vector<reco::CaloJet>::const_iterator itJet = theJets.begin(); itJet != theJets.end(); itJet++) {
73  float curDR = this->spaceAngle(aLepton, *itJet);
74  if (curDR < isoAngle) isoAngle = curDR;
75  }
76  return isoAngle;
77 }
Jets made from CaloTowers.
Definition: CaloJet.h:30
float spaceAngle(const CLHEP::HepLorentzVector &aLepton, const reco::CaloJet &aJet)
T sqrt(T t)
Definition: SSEVec.h:28
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:359
double deltaR2(const Vector1 &v1, const Vector2 &v2)
Definition: VectorUtil.h:78
Analysis-level electron class.
Definition: Electron.h:46
T const * product() const
Definition: Handle.h:74
float calculate(const Electron &theElectron, const edm::View< reco::Track > &theTracks, float isoConeElectron=0.3) const
calculate the TrackIsoPt for the lepton object
std::vector< CaloJet > CaloJetCollection
collection of CaloJet objects
float LeptonJetIsolationAngle::spaceAngle ( const CLHEP::HepLorentzVector &  aLepton,
const reco::CaloJet aJet 
)
private

Definition at line 81 of file LeptonJetIsolationAngle.cc.

References funct::cos(), reco::LeafCandidate::phi(), funct::sin(), and reco::LeafCandidate::theta().

Referenced by calculate().

81  {
82  return acos(sin(aJet.theta()) * cos(aJet.phi()) * sin(aLepton.theta()) * cos(aLepton.phi())
83  + sin(aJet.theta()) * sin(aJet.phi()) * sin(aLepton.theta()) * sin(aLepton.phi())
84  + cos(aJet.theta()) * cos(aLepton.theta()));
85 }
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
virtual double theta() const
momentum polar angle
virtual double phi() const
momentum azimuthal angle

Member Data Documentation

TrackerIsolationPt pat::LeptonJetIsolationAngle::trkIsolator_
private

Definition at line 51 of file LeptonJetIsolationAngle.h.

Referenced by calculate().