CMS 3D CMS Logo

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::CaloJetCollectionjetToken_
 
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 30 of file LeptonJetIsolationAngle.h.

Constructor & Destructor Documentation

◆ LeptonJetIsolationAngle()

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

Definition at line 12 of file LeptonJetIsolationAngle.cc.

13  : jetToken_(iC.consumes<reco::CaloJetCollection>(edm::InputTag("iterativeCone5CaloJets"))),
14  electronsToken_(iC.consumes<std::vector<reco::GsfElectron> >(edm::InputTag("pixelMatchGsfElectrons"))) {}

◆ ~LeptonJetIsolationAngle()

LeptonJetIsolationAngle::~LeptonJetIsolationAngle ( )

Definition at line 17 of file LeptonJetIsolationAngle.cc.

17 {}

Member Function Documentation

◆ calculate() [1/3]

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

Definition at line 34 of file LeptonJetIsolationAngle.cc.

36  {
37  // FIXME: this is an ugly temporary workaround, JetMET+egamma should come up with a better tool
38  // retrieve the jets
40  iEvent.getByToken(jetToken_, jetHandle);
41  reco::CaloJetCollection jetColl = *(jetHandle.product());
42  // retrieve the electrons which might be in the jet list
44  iEvent.getByToken(electronsToken_, electronsHandle);
45  std::vector<reco::GsfElectron> electrons = *electronsHandle;
46  // determine the set of isolated electrons
47  std::vector<Electron> isoElectrons;
48  for (size_t ie = 0; ie < electrons.size(); ie++) {
49  Electron anElectron(electrons[ie]);
50  if (anElectron.pt() > 10 && trkIsolator_.calculate(anElectron, *trackHandle) < 3.0) {
51  isoElectrons.push_back(electrons[ie]);
52  }
53  }
54  // determine the collections of jets, cleaned from electrons
55  std::vector<reco::CaloJet> theJets;
56  for (reco::CaloJetCollection::const_iterator itJet = jetColl.begin(); itJet != jetColl.end(); itJet++) {
57  float mindr2 = 9999.;
58  for (size_t ie = 0; ie < isoElectrons.size(); ie++) {
59  float dr2 = ::deltaR2(*itJet, isoElectrons[ie]);
60  if (dr2 < mindr2)
61  mindr2 = dr2;
62  }
63  float mindr = std::sqrt(mindr2);
64  // yes, all cuts hardcoded buts, but it's a second-order effect
65  if (itJet->et() > 15 && mindr > 0.3)
66  theJets.push_back(reco::CaloJet(*itJet));
67  }
68  // calculate finally the isolation angle
69  float isoAngle = 1000; // default to some craze impossible number to inhibit compiler warnings
70  for (std::vector<reco::CaloJet>::const_iterator itJet = theJets.begin(); itJet != theJets.end(); itJet++) {
71  float curDR = this->spaceAngle(aLepton, *itJet);
72  if (curDR < isoAngle)
73  isoAngle = curDR;
74  }
75  return isoAngle;
76 }

References pat::TrackerIsolationPt::calculate(), HLTMuonOfflineAnalyzer_cfi::deltaR2, pwdgSkimBPark_cfi::electrons, electronsToken_, iEvent, jetToken_, MTVHistoProducerAlgoForTrackerBlock_cfi::mindr, edm::Handle< T >::product(), reco::LeafCandidate::pt(), spaceAngle(), mathSSE::sqrt(), and trkIsolator_.

◆ calculate() [2/3]

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

Definition at line 20 of file LeptonJetIsolationAngle.cc.

22  {
23  CLHEP::HepLorentzVector theElectronHLV(theElectron.px(), theElectron.py(), theElectron.pz(), theElectron.energy());
24  return this->calculate(theElectronHLV, trackHandle, iEvent);
25 }

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

Referenced by calculate().

◆ calculate() [3/3]

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

Definition at line 26 of file LeptonJetIsolationAngle.cc.

28  {
29  CLHEP::HepLorentzVector theMuonHLV(theMuon.px(), theMuon.py(), theMuon.pz(), theMuon.energy());
30  return this->calculate(theMuonHLV, trackHandle, iEvent);
31 }

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

◆ spaceAngle()

float LeptonJetIsolationAngle::spaceAngle ( const CLHEP::HepLorentzVector &  aLepton,
const reco::CaloJet aJet 
)
private

Definition at line 79 of file LeptonJetIsolationAngle.cc.

79  {
80  return acos(sin(aJet.theta()) * cos(aJet.phi()) * sin(aLepton.theta()) * cos(aLepton.phi()) +
81  sin(aJet.theta()) * sin(aJet.phi()) * sin(aLepton.theta()) * sin(aLepton.phi()) +
82  cos(aJet.theta()) * cos(aLepton.theta()));
83 }

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

Referenced by calculate().

Member Data Documentation

◆ electronsToken_

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

Definition at line 51 of file LeptonJetIsolationAngle.h.

Referenced by calculate().

◆ jetToken_

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

Definition at line 50 of file LeptonJetIsolationAngle.h.

Referenced by calculate().

◆ trkIsolator_

TrackerIsolationPt pat::LeptonJetIsolationAngle::trkIsolator_
private

Definition at line 49 of file LeptonJetIsolationAngle.h.

Referenced by calculate().

pat::LeptonJetIsolationAngle::trkIsolator_
TrackerIsolationPt trkIsolator_
Definition: LeptonJetIsolationAngle.h:49
reco::CaloJet
Jets made from CaloTowers.
Definition: CaloJet.h:27
edm::Handle::product
T const * product() const
Definition: Handle.h:70
pat::LeptonJetIsolationAngle::electronsToken_
edm::EDGetTokenT< std::vector< reco::GsfElectron > > electronsToken_
Definition: LeptonJetIsolationAngle.h:51
Electron
Definition: Electron.py:1
pat::LeptonJetIsolationAngle::calculate
float calculate(const Electron &anElectron, const edm::Handle< edm::View< reco::Track > > &trackHandle, const edm::Event &iEvent)
Definition: LeptonJetIsolationAngle.cc:20
pat::TrackerIsolationPt::calculate
float calculate(const Electron &theElectron, const edm::View< reco::Track > &theTracks, float isoConeElectron=0.3) const
calculate the TrackIsoPt for the lepton object
Definition: TrackerIsolationPt.cc:27
edm::Handle< reco::CaloJetCollection >
funct::sin
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
pat::LeptonJetIsolationAngle::spaceAngle
float spaceAngle(const CLHEP::HepLorentzVector &aLepton, const reco::CaloJet &aJet)
Definition: LeptonJetIsolationAngle.cc:79
funct::cos
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
pat::LeptonJetIsolationAngle::jetToken_
edm::EDGetTokenT< reco::CaloJetCollection > jetToken_
Definition: LeptonJetIsolationAngle.h:50
reco::LeafCandidate::theta
double theta() const final
momentum polar angle
Definition: LeafCandidate.h:150
edm::ConsumesCollector::consumes
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
Definition: ConsumesCollector.h:55
iEvent
int iEvent
Definition: GenABIO.cc:224
reco::CaloJetCollection
std::vector< CaloJet > CaloJetCollection
collection of CaloJet objects
Definition: CaloJetCollection.h:15
reco::LeafCandidate::phi
double phi() const final
momentum azimuthal angle
Definition: LeafCandidate.h:148
pwdgSkimBPark_cfi.electrons
electrons
Definition: pwdgSkimBPark_cfi.py:6
HLTMuonOfflineAnalyzer_cfi.deltaR2
deltaR2
Definition: HLTMuonOfflineAnalyzer_cfi.py:105
MTVHistoProducerAlgoForTrackerBlock_cfi.mindr
mindr
Definition: MTVHistoProducerAlgoForTrackerBlock_cfi.py:79
edm::InputTag
Definition: InputTag.h:15