CMS 3D CMS Logo

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