CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/PhysicsTools/PatUtils/src/LeptonJetIsolationAngle.cc

Go to the documentation of this file.
00001 //
00002 // $Id: LeptonJetIsolationAngle.cc,v 1.5 2009/05/26 08:54:22 fabiocos Exp $
00003 //
00004 
00005 #include "PhysicsTools/PatUtils/interface/LeptonJetIsolationAngle.h"
00006 
00007 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00008 #include "FWCore/Utilities/interface/Exception.h"
00009 #include "DataFormats/Math/interface/deltaR.h"
00010 
00011 #include <vector>
00012 
00013 
00014 using namespace pat;
00015 
00016 
00017 // constructor
00018 LeptonJetIsolationAngle::LeptonJetIsolationAngle() {
00019 }
00020 
00021 
00022 // destructor
00023 LeptonJetIsolationAngle::~LeptonJetIsolationAngle() {
00024 }
00025 
00026 
00027 // calculate the JetIsoA for the lepton object
00028 float LeptonJetIsolationAngle::calculate(const Electron & theElectron, const edm::Handle<edm::View<reco::Track> > & trackHandle, const edm::Event & iEvent) {
00029   CLHEP::HepLorentzVector theElectronHLV(theElectron.px(), theElectron.py(), theElectron.pz(), theElectron.energy());
00030   return this->calculate(theElectronHLV, trackHandle, iEvent);
00031 }
00032 float LeptonJetIsolationAngle::calculate(const Muon & theMuon, const edm::Handle<edm::View<reco::Track> > & trackHandle, const edm::Event & iEvent) {
00033   CLHEP::HepLorentzVector theMuonHLV(theMuon.px(), theMuon.py(), theMuon.pz(), theMuon.energy());
00034   return this->calculate(theMuonHLV, trackHandle, iEvent);
00035 }
00036 
00037 
00038 // calculate the JetIsoA for the lepton's HLV
00039 float LeptonJetIsolationAngle::calculate(const CLHEP::HepLorentzVector & aLepton, const edm::Handle<edm::View<reco::Track> > & trackHandle, const edm::Event & iEvent) {
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<reco::GsfElectron> > electronsHandle;
00047   iEvent.getByLabel("pixelMatchGsfElectrons", electronsHandle);
00048   std::vector<reco::GsfElectron> 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 }
00078 
00079 
00080 // calculate the angle between two vectors in 3d eucledian space
00081 float LeptonJetIsolationAngle::spaceAngle(const CLHEP::HepLorentzVector & aLepton, const reco::CaloJet & aJet) {
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 }