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