CMS 3D CMS Logo

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>

List of all members.

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.4 2009/05/26 08:54:21 fabiocos Exp

Definition at line 34 of file LeptonJetIsolationAngle.h.


Constructor & Destructor Documentation

LeptonJetIsolationAngle::LeptonJetIsolationAngle ( )

Definition at line 18 of file LeptonJetIsolationAngle.cc.

                                                 {
}
LeptonJetIsolationAngle::~LeptonJetIsolationAngle ( )

Definition at line 23 of file LeptonJetIsolationAngle.cc.

                                                  {
}

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().

                                                                                                                                                  {
  CLHEP::HepLorentzVector theElectronHLV(theElectron.px(), theElectron.py(), theElectron.pz(), theElectron.energy());
  return this->calculate(theElectronHLV, trackHandle, 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().

                                                                                                                                          {
  CLHEP::HepLorentzVector theMuonHLV(theMuon.px(), theMuon.py(), theMuon.pz(), theMuon.energy());
  return this->calculate(theMuonHLV, trackHandle, 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(), HI_PhotonSkim_cff::electrons, edm::Event::getByLabel(), edm::Handle< T >::product(), reco::LeafCandidate::pt(), spaceAngle(), mathSSE::sqrt(), and trkIsolator_.

                                                                                                                                                           {
  // FIXME: this is an ugly temporary workaround, JetMET+egamma should come up with a better tool
  // retrieve the jets
  edm::Handle<reco::CaloJetCollection> jetHandle;
  iEvent.getByLabel("iterativeCone5CaloJets", jetHandle);
  reco::CaloJetCollection jetColl = *(jetHandle.product());
  // retrieve the electrons which might be in the jet list
  edm::Handle<std::vector<reco::GsfElectron> > electronsHandle;
  iEvent.getByLabel("pixelMatchGsfElectrons", electronsHandle);
  std::vector<reco::GsfElectron> electrons = *electronsHandle;
  // determine the set of isolated electrons
  std::vector<Electron> isoElectrons;
  for (size_t ie=0; ie<electrons.size(); ie++) {
    Electron anElectron(electrons[ie]);
    if (anElectron.pt() > 10 &&
        trkIsolator_.calculate(anElectron, *trackHandle) < 3.0) {
      isoElectrons.push_back(electrons[ie]);
    }
  }
  // determine the collections of jets, cleaned from electrons
  std::vector<reco::CaloJet> theJets;
  for (reco::CaloJetCollection::const_iterator itJet = jetColl.begin(); itJet != jetColl.end(); itJet++) {
    float mindr2 = 9999.;
    for (size_t ie = 0; ie < isoElectrons.size(); ie++) {
      float dr2 = ::deltaR2(*itJet, isoElectrons[ie]);
      if (dr2 < mindr2) mindr2 = dr2;
    }
    float mindr = std::sqrt(mindr2);
    // yes, all cuts hardcoded buts, but it's a second-order effect
    if (itJet->et() > 15 && mindr > 0.3) theJets.push_back(reco::CaloJet(*itJet));
  }
  // calculate finally the isolation angle
  float isoAngle = 1000; // default to some craze impossible number to inhibit compiler warnings
  for (std::vector<reco::CaloJet>::const_iterator itJet = theJets.begin(); itJet != theJets.end(); itJet++) {
    float curDR = this->spaceAngle(aLepton, *itJet);
    if (curDR < isoAngle) isoAngle = curDR;
  }
  return isoAngle;
}
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().

                                                                                                         {
  return acos(sin(aJet.theta()) * cos(aJet.phi()) * sin(aLepton.theta()) * cos(aLepton.phi())
            + sin(aJet.theta()) * sin(aJet.phi()) * sin(aLepton.theta()) * sin(aLepton.phi())
            + cos(aJet.theta()) * cos(aLepton.theta()));
}

Member Data Documentation

Definition at line 51 of file LeptonJetIsolationAngle.h.

Referenced by calculate().