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 (edm::ConsumesCollector &&iC)
 
 ~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

edm::EDGetTokenT< std::vector
< reco::GsfElectron > > 
electronsToken_
 
edm::EDGetTokenT
< reco::CaloJetCollection
jetToken_
 
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

Definition at line 33 of file LeptonJetIsolationAngle.h.

Constructor & Destructor Documentation

LeptonJetIsolationAngle::LeptonJetIsolationAngle ( edm::ConsumesCollector &&  iC)

Definition at line 14 of file LeptonJetIsolationAngle.cc.

15 : jetToken_( iC.consumes< reco::CaloJetCollection >( edm::InputTag( "iterativeCone5CaloJets" ) ) )
16 , electronsToken_( iC.consumes< std::vector<reco::GsfElectron > >( edm::InputTag( "pixelMatchGsfElectrons" ) ) )
17 {
18 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
edm::EDGetTokenT< reco::CaloJetCollection > jetToken_
edm::EDGetTokenT< std::vector< reco::GsfElectron > > electronsToken_
std::vector< CaloJet > CaloJetCollection
collection of CaloJet objects
LeptonJetIsolationAngle::~LeptonJetIsolationAngle ( )

Definition at line 22 of file LeptonJetIsolationAngle.cc.

22  {
23 }

Member Function Documentation

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

Definition at line 27 of file LeptonJetIsolationAngle.cc.

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

Referenced by calculate().

27  {
28  CLHEP::HepLorentzVector theElectronHLV(theElectron.px(), theElectron.py(), theElectron.pz(), theElectron.energy());
29  return this->calculate(theElectronHLV, trackHandle, iEvent);
30 }
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 31 of file LeptonJetIsolationAngle.cc.

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

31  {
32  CLHEP::HepLorentzVector theMuonHLV(theMuon.px(), theMuon.py(), theMuon.pz(), theMuon.energy());
33  return this->calculate(theMuonHLV, trackHandle, iEvent);
34 }
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 38 of file LeptonJetIsolationAngle.cc.

References pat::TrackerIsolationPt::calculate(), reco::deltaR2(), HI_PhotonSkim_cff::electrons, electronsToken_, edm::Event::getByToken(), jetToken_, edm::Handle< T >::product(), reco::LeafCandidate::pt(), spaceAngle(), mathSSE::sqrt(), and trkIsolator_.

38  {
39  // FIXME: this is an ugly temporary workaround, JetMET+egamma should come up with a better tool
40  // retrieve the jets
42  iEvent.getByToken(jetToken_, jetHandle);
43  reco::CaloJetCollection jetColl = *(jetHandle.product());
44  // retrieve the electrons which might be in the jet list
46  iEvent.getByToken(electronsToken_, electronsHandle);
47  std::vector<reco::GsfElectron> electrons = *electronsHandle;
48  // determine the set of isolated electrons
49  std::vector<Electron> isoElectrons;
50  for (size_t ie=0; ie<electrons.size(); ie++) {
51  Electron anElectron(electrons[ie]);
52  if (anElectron.pt() > 10 &&
53  trkIsolator_.calculate(anElectron, *trackHandle) < 3.0) {
54  isoElectrons.push_back(electrons[ie]);
55  }
56  }
57  // determine the collections of jets, cleaned from electrons
58  std::vector<reco::CaloJet> theJets;
59  for (reco::CaloJetCollection::const_iterator itJet = jetColl.begin(); itJet != jetColl.end(); itJet++) {
60  float mindr2 = 9999.;
61  for (size_t ie = 0; ie < isoElectrons.size(); ie++) {
62  float dr2 = ::deltaR2(*itJet, isoElectrons[ie]);
63  if (dr2 < mindr2) mindr2 = dr2;
64  }
65  float mindr = std::sqrt(mindr2);
66  // yes, all cuts hardcoded buts, but it's a second-order effect
67  if (itJet->et() > 15 && mindr > 0.3) theJets.push_back(reco::CaloJet(*itJet));
68  }
69  // calculate finally the isolation angle
70  float isoAngle = 1000; // default to some craze impossible number to inhibit compiler warnings
71  for (std::vector<reco::CaloJet>::const_iterator itJet = theJets.begin(); itJet != theJets.end(); itJet++) {
72  float curDR = this->spaceAngle(aLepton, *itJet);
73  if (curDR < isoAngle) isoAngle = curDR;
74  }
75  return isoAngle;
76 }
Jets made from CaloTowers.
Definition: CaloJet.h:29
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:446
float spaceAngle(const CLHEP::HepLorentzVector &aLepton, const reco::CaloJet &aJet)
T sqrt(T t)
Definition: SSEVec.h:48
edm::EDGetTokenT< reco::CaloJetCollection > jetToken_
double deltaR2(const T1 &t1, const T2 &t2)
Definition: deltaR.h:36
edm::EDGetTokenT< std::vector< reco::GsfElectron > > electronsToken_
T const * product() const
Definition: Handle.h:81
Analysis-level electron class.
Definition: Electron.h:52
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 80 of file LeptonJetIsolationAngle.cc.

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

Referenced by calculate().

80  {
81  return acos(sin(aJet.theta()) * cos(aJet.phi()) * sin(aLepton.theta()) * cos(aLepton.phi())
82  + sin(aJet.theta()) * sin(aJet.phi()) * sin(aLepton.theta()) * sin(aLepton.phi())
83  + cos(aJet.theta()) * cos(aLepton.theta()));
84 }
virtual float phi() const
momentum azimuthal angle
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

Member Data Documentation

edm::EDGetTokenT<std::vector<reco::GsfElectron> > pat::LeptonJetIsolationAngle::electronsToken_
private

Definition at line 52 of file LeptonJetIsolationAngle.h.

Referenced by calculate().

edm::EDGetTokenT<reco::CaloJetCollection> pat::LeptonJetIsolationAngle::jetToken_
private

Definition at line 51 of file LeptonJetIsolationAngle.h.

Referenced by calculate().

TrackerIsolationPt pat::LeptonJetIsolationAngle::trkIsolator_
private

Definition at line 50 of file LeptonJetIsolationAngle.h.

Referenced by calculate().