CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
LeptonJetIsolationAngle.cc
Go to the documentation of this file.
1 //
2 // $Id: LeptonJetIsolationAngle.cc,v 1.5 2009/05/26 08:54:22 fabiocos Exp $
3 //
4 
6 
10 
11 #include <vector>
12 
13 
14 using namespace pat;
15 
16 
17 // constructor
19 }
20 
21 
22 // destructor
24 }
25 
26 
27 // calculate the JetIsoA for the lepton object
28 float LeptonJetIsolationAngle::calculate(const Electron & theElectron, const edm::Handle<edm::View<reco::Track> > & trackHandle, const edm::Event & iEvent) {
29  CLHEP::HepLorentzVector theElectronHLV(theElectron.px(), theElectron.py(), theElectron.pz(), theElectron.energy());
30  return this->calculate(theElectronHLV, trackHandle, iEvent);
31 }
32 float LeptonJetIsolationAngle::calculate(const Muon & theMuon, const edm::Handle<edm::View<reco::Track> > & trackHandle, const edm::Event & iEvent) {
33  CLHEP::HepLorentzVector theMuonHLV(theMuon.px(), theMuon.py(), theMuon.pz(), theMuon.energy());
34  return this->calculate(theMuonHLV, trackHandle, iEvent);
35 }
36 
37 
38 // calculate the JetIsoA for the lepton's HLV
39 float LeptonJetIsolationAngle::calculate(const CLHEP::HepLorentzVector & aLepton, const edm::Handle<edm::View<reco::Track> > & trackHandle, const edm::Event & iEvent) {
40  // FIXME: this is an ugly temporary workaround, JetMET+egamma should come up with a better tool
41  // retrieve the jets
43  iEvent.getByLabel("iterativeCone5CaloJets", jetHandle);
44  reco::CaloJetCollection jetColl = *(jetHandle.product());
45  // retrieve the electrons which might be in the jet list
47  iEvent.getByLabel("pixelMatchGsfElectrons", electronsHandle);
48  std::vector<reco::GsfElectron> electrons = *electronsHandle;
49  // determine the set of isolated electrons
50  std::vector<Electron> isoElectrons;
51  for (size_t ie=0; ie<electrons.size(); ie++) {
52  Electron anElectron(electrons[ie]);
53  if (anElectron.pt() > 10 &&
54  trkIsolator_.calculate(anElectron, *trackHandle) < 3.0) {
55  isoElectrons.push_back(electrons[ie]);
56  }
57  }
58  // determine the collections of jets, cleaned from electrons
59  std::vector<reco::CaloJet> theJets;
60  for (reco::CaloJetCollection::const_iterator itJet = jetColl.begin(); itJet != jetColl.end(); itJet++) {
61  float mindr2 = 9999.;
62  for (size_t ie = 0; ie < isoElectrons.size(); ie++) {
63  float dr2 = ::deltaR2(*itJet, isoElectrons[ie]);
64  if (dr2 < mindr2) mindr2 = dr2;
65  }
66  float mindr = std::sqrt(mindr2);
67  // yes, all cuts hardcoded buts, but it's a second-order effect
68  if (itJet->et() > 15 && mindr > 0.3) theJets.push_back(reco::CaloJet(*itJet));
69  }
70  // calculate finally the isolation angle
71  float isoAngle = 1000; // default to some craze impossible number to inhibit compiler warnings
72  for (std::vector<reco::CaloJet>::const_iterator itJet = theJets.begin(); itJet != theJets.end(); itJet++) {
73  float curDR = this->spaceAngle(aLepton, *itJet);
74  if (curDR < isoAngle) isoAngle = curDR;
75  }
76  return isoAngle;
77 }
78 
79 
80 // calculate the angle between two vectors in 3d eucledian space
81 float LeptonJetIsolationAngle::spaceAngle(const CLHEP::HepLorentzVector & aLepton, const reco::CaloJet & aJet) {
82  return acos(sin(aJet.theta()) * cos(aJet.phi()) * sin(aLepton.theta()) * cos(aLepton.phi())
83  + sin(aJet.theta()) * sin(aJet.phi()) * sin(aLepton.theta()) * sin(aLepton.phi())
84  + cos(aJet.theta()) * cos(aLepton.theta()));
85 }
Jets made from CaloTowers.
Definition: CaloJet.h:30
float spaceAngle(const CLHEP::HepLorentzVector &aLepton, const reco::CaloJet &aJet)
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
virtual double energy() const
energy
int iEvent
Definition: GenABIO.cc:243
T sqrt(T t)
Definition: SSEVec.h:46
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
float calculate(const Electron &anElectron, const edm::Handle< edm::View< reco::Track > > &trackHandle, const edm::Event &iEvent)
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
double deltaR2(const Vector1 &v1, const Vector2 &v2)
Definition: VectorUtil.h:78
virtual double theta() const
momentum polar angle
virtual double px() const
x coordinate of momentum vector
virtual double pt() const
transverse momentum
Analysis-level electron class.
Definition: Electron.h:52
T const * product() const
Definition: Handle.h:74
virtual double pz() const
z coordinate of momentum vector
float calculate(const Electron &theElectron, const edm::View< reco::Track > &theTracks, float isoConeElectron=0.3) const
calculate the TrackIsoPt for the lepton object
virtual double phi() const
momentum azimuthal angle
Analysis-level muon class.
Definition: Muon.h:51
virtual double py() const
y coordinate of momentum vector
std::vector< CaloJet > CaloJetCollection
collection of CaloJet objects