#include <PhysicsTools/PatUtils/interface/LeptonJetIsolationAngle.h>
Public Member Functions | |
float | calculate (const Muon &aMuon, const edm::Handle< edm::View< reco::Track > > &trackHandle, const edm::Event &iEvent) |
float | calculate (const Electron &anElectron, const edm::Handle< edm::View< reco::Track > > &trackHandle, const edm::Event &iEvent) |
LeptonJetIsolationAngle () | |
~LeptonJetIsolationAngle () | |
Private Member Functions | |
float | calculate (const HepLorentzVector &aLepton, const edm::Handle< edm::View< reco::Track > > &trackHandle, const edm::Event &iEvent) |
float | spaceAngle (const HepLorentzVector &aLepton, const reco::CaloJet &aJet) |
Private Attributes | |
TrackerIsolationPt | trkIsolator_ |
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
Definition at line 34 of file LeptonJetIsolationAngle.h.
LeptonJetIsolationAngle::LeptonJetIsolationAngle | ( | ) |
LeptonJetIsolationAngle::~LeptonJetIsolationAngle | ( | ) |
float LeptonJetIsolationAngle::calculate | ( | const 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(), deltaR2(), edm::Event::getByLabel(), jetColl, edm::Handle< T >::product(), spaceAngle(), funct::sqrt(), theJets, and trkIsolator_.
00039 { 00040 // FIXME: this is an ugly temporary workaround, JetMET+egamma should come up with a better tool 00041 // retrieve the jets 00042 edm::Handle<reco::CaloJetCollection> jetHandle; 00043 iEvent.getByLabel("iterativeCone5CaloJets", jetHandle); 00044 reco::CaloJetCollection jetColl = *(jetHandle.product()); 00045 // retrieve the electrons which might be in the jet list 00046 edm::Handle<std::vector<ElectronType> > electronsHandle; 00047 iEvent.getByLabel("pixelMatchGsfElectrons", electronsHandle); 00048 std::vector<ElectronType> electrons = *electronsHandle; 00049 // determine the set of isolated electrons 00050 std::vector<Electron> isoElectrons; 00051 for (size_t ie=0; ie<electrons.size(); ie++) { 00052 Electron anElectron(electrons[ie]); 00053 if (anElectron.pt() > 10 && 00054 trkIsolator_.calculate(anElectron, *trackHandle) < 3.0) { 00055 isoElectrons.push_back(electrons[ie]); 00056 } 00057 } 00058 // determine the collections of jets, cleaned from electrons 00059 std::vector<reco::CaloJet> theJets; 00060 for (reco::CaloJetCollection::const_iterator itJet = jetColl.begin(); itJet != jetColl.end(); itJet++) { 00061 float mindr2 = 9999.; 00062 for (size_t ie = 0; ie < isoElectrons.size(); ie++) { 00063 float dr2 = ::deltaR2(*itJet, isoElectrons[ie]); 00064 if (dr2 < mindr2) mindr2 = dr2; 00065 } 00066 float mindr = std::sqrt(mindr2); 00067 // yes, all cuts hardcoded buts, but it's a second-order effect 00068 if (itJet->et() > 15 && mindr > 0.3) theJets.push_back(reco::CaloJet(*itJet)); 00069 } 00070 // calculate finally the isolation angle 00071 float isoAngle = 1000; // default to some craze impossible number to inhibit compiler warnings 00072 for (std::vector<reco::CaloJet>::const_iterator itJet = theJets.begin(); itJet != theJets.end(); itJet++) { 00073 float curDR = this->spaceAngle(aLepton, *itJet); 00074 if (curDR < isoAngle) isoAngle = curDR; 00075 } 00076 return isoAngle; 00077 }
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().
00032 { 00033 HepLorentzVector theMuonHLV(theMuon.px(), theMuon.py(), theMuon.pz(), theMuon.energy()); 00034 return this->calculate(theMuonHLV, trackHandle, iEvent); 00035 }
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.
Referenced by calculate().
00028 { 00029 HepLorentzVector theElectronHLV(theElectron.px(), theElectron.py(), theElectron.pz(), theElectron.energy()); 00030 return this->calculate(theElectronHLV, trackHandle, iEvent); 00031 }
float LeptonJetIsolationAngle::spaceAngle | ( | const HepLorentzVector & | aLepton, | |
const reco::CaloJet & | aJet | |||
) | [private] |
Definition at line 81 of file LeptonJetIsolationAngle.cc.
References funct::cos(), reco::Particle::phi(), funct::sin(), and reco::Particle::theta().
Referenced by calculate().
00081 { 00082 return acos(sin(aJet.theta()) * cos(aJet.phi()) * sin(aLepton.theta()) * cos(aLepton.phi()) 00083 + sin(aJet.theta()) * sin(aJet.phi()) * sin(aLepton.theta()) * sin(aLepton.phi()) 00084 + cos(aJet.theta()) * cos(aLepton.theta())); 00085 }