CMS 3D CMS Logo

PFCandWithSuperClusterExtractor.cc
Go to the documentation of this file.
2 
9 
10 using namespace edm;
11 using namespace reco;
12 
14  : thePFCandToken(iC.consumes<PFCandidateCollection>(par.getParameter<edm::InputTag>("inputCandView"))),
15  theDepositLabel(par.getUntrackedParameter<std::string>("DepositLabel")),
16  theVetoSuperClusterMatch(par.getParameter<bool>("SCMatch_Veto")),
17  theMissHitVetoSuperClusterMatch(par.getParameter<bool>("MissHitSCMatch_Veto")),
18  theDiff_r(par.getParameter<double>("Diff_r")),
19  theDiff_z(par.getParameter<double>("Diff_z")),
20  theDR_Max(par.getParameter<double>("DR_Max")),
21  theDR_Veto(par.getParameter<double>("DR_Veto")) {
22  // std::cout << " Loading PFCandWithSuperClusterExtractor " << std::endl;
23 }
24 /*
25 reco::IsoDeposit::Vetos PFCandWithSuperClusterExtractor::vetos(const edm::Event & ev,
26  const edm::EventSetup & evSetup, const reco::Candidate & cand) const
27 {
28  reco::isodeposit::Direction dir(cand.eta(),cand.phi());
29  return reco::IsoDeposit::Vetos(1,veto(dir));
30 }
31 */
32 
35  result.vetoDir = dir;
36  result.dR = theDR_Veto;
37  return result;
38 }
39 
41  const EventSetup& eventSetup,
42  const Photon& cand) const {
43  reco::isodeposit::Direction candDir(cand.eta(), cand.phi());
44  IsoDeposit deposit(candDir);
45  deposit.setVeto(veto(candDir));
46  deposit.addCandEnergy(cand.pt());
47 
49  event.getByToken(thePFCandToken, PFCandH);
50 
51  double eta = cand.eta(), phi = cand.phi();
52  const reco::Particle::Point& vtx = cand.vertex();
53  for (PFCandidateCollection::const_iterator it = PFCandH->begin(), ed = PFCandH->end(); it != ed; ++it) {
54  double dR = deltaR(it->eta(), it->phi(), eta, phi);
55  // veto SC
56  if (theVetoSuperClusterMatch && cand.superCluster().isNonnull() && it->superClusterRef().isNonnull() &&
57  cand.superCluster() == it->superClusterRef())
58  continue;
59  if ((dR < theDR_Max) && (dR > theDR_Veto) && (std::abs(it->vz() - cand.vz()) < theDiff_z) &&
60  ((it->vertex() - vtx).Rho() < theDiff_r)) {
61  // ok
62  reco::isodeposit::Direction dirTrk(it->eta(), it->phi());
63  deposit.addDeposit(dirTrk, it->pt());
64  }
65  }
66 
67  return deposit;
68 }
69 
71  const EventSetup& eventSetup,
72  const GsfElectron& cand) const {
73  reco::isodeposit::Direction candDir(cand.eta(), cand.phi());
74  IsoDeposit deposit(candDir);
75  deposit.setVeto(veto(candDir));
76  deposit.addCandEnergy(cand.pt());
77 
79  event.getByToken(thePFCandToken, PFCandH);
80 
81  double eta = cand.eta(), phi = cand.phi();
82  const reco::Particle::Point& vtx = cand.vertex();
83  for (PFCandidateCollection::const_iterator it = PFCandH->begin(), ed = PFCandH->end(); it != ed; ++it) {
84  double dR = deltaR(it->eta(), it->phi(), eta, phi);
85  // If MissHits>0 (possibly reconstructed as a photon in the PF in this case, kill the the photon if sharing the same SC)
86  if (cand.gsfTrack()->hitPattern().numberOfLostHits(reco::HitPattern::MISSING_INNER_HITS) > 0 &&
87  theMissHitVetoSuperClusterMatch && it->mva_nothing_gamma() > 0.99 && cand.superCluster().isNonnull() &&
88  it->superClusterRef().isNonnull() && cand.superCluster() == it->superClusterRef()) {
89  continue;
90  }
91  if ((dR < theDR_Max) && (dR > theDR_Veto) && (std::abs(it->vz() - cand.vz()) < theDiff_z) &&
92  ((it->vertex() - vtx).Rho() < theDiff_r)) {
93  // ok
94  reco::isodeposit::Direction dirTrk(it->eta(), it->phi());
95  deposit.addDeposit(dirTrk, it->pt());
96  }
97  }
98 
99  return deposit;
100 }
101 
103  const EventSetup& eventSetup,
104  const Track& cand) const {
105  reco::isodeposit::Direction candDir(cand.eta(), cand.phi());
106  IsoDeposit deposit(candDir);
107  deposit.setVeto(veto(candDir));
108  deposit.addCandEnergy(cand.pt());
110  event.getByToken(thePFCandToken, PFCandH);
111 
112  double eta = cand.eta(), phi = cand.phi();
113  const reco::Particle::Point& vtx = cand.vertex();
114  for (PFCandidateCollection::const_iterator it = PFCandH->begin(), ed = PFCandH->end(); it != ed; ++it) {
115  double dR = deltaR(it->eta(), it->phi(), eta, phi);
116 
117  if ((dR < theDR_Max) && (dR > theDR_Veto) && (std::abs(it->vz() - cand.vz()) < theDiff_z) &&
118  ((it->vertex() - vtx).Rho() < theDiff_r)) {
119  // ok
120  reco::isodeposit::Direction dirTrk(it->eta(), it->phi());
121  deposit.addDeposit(dirTrk, it->pt());
122  }
123  }
124 
125  return deposit;
126 }
127 
129  const EventSetup& eventSetup,
130  const PFCandidate& cand) const {
131  reco::isodeposit::Direction candDir(cand.eta(), cand.phi());
132  IsoDeposit deposit(candDir);
133  deposit.setVeto(veto(candDir));
134  deposit.addCandEnergy(cand.pt());
136  event.getByToken(thePFCandToken, PFCandH);
137 
138  double eta = cand.eta(), phi = cand.phi();
139  const reco::Particle::Point& vtx = cand.vertex();
140  for (PFCandidateCollection::const_iterator it = PFCandH->begin(), ed = PFCandH->end(); it != ed; ++it) {
141  // veto SC
142  if (theVetoSuperClusterMatch && cand.superClusterRef().isNonnull() && it->superClusterRef().isNonnull() &&
143  cand.superClusterRef() == it->superClusterRef())
144  continue;
145  double dR = deltaR(it->eta(), it->phi(), eta, phi);
146 
147  if ((dR < theDR_Max) && (dR > theDR_Veto) && (std::abs(it->vz() - cand.vz()) < theDiff_z) &&
148  ((it->vertex() - vtx).Rho() < theDiff_r)) {
149  // ok
150  reco::isodeposit::Direction dirTrk(it->eta(), it->phi());
151  deposit.addDeposit(dirTrk, it->pt());
152  }
153  }
154 
155  return deposit;
156 }
157 
160 
GsfTrackRef gsfTrack() const override
reference to a GsfTrack
Definition: GsfElectron.h:156
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:238
Definition: Photon.py:1
double eta() const final
momentum pseudorapidity
edm::EDGetTokenT< reco::PFCandidateCollection > thePFCandToken
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:614
double pt() const final
transverse momentum
reco::SuperClusterRef superCluster() const override
Ref to SuperCluster.
const Point & vertex() const override
vertex position (overwritten by PF...)
Definition: PFCandidate.cc:602
const Point & vertex() const
reference point on the track. This method is DEPRECATED, please use referencePoint() instead ...
Definition: TrackBase.h:641
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:617
double pt() const
track transverse momentum
Definition: TrackBase.h:602
const Point & vertex() const override
vertex position (overwritten by PF...)
math::XYZPoint Point
point in the space
Definition: Particle.h:25
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
reco::IsoDeposit depositFromObject(const edm::Event &ev, const edm::EventSetup &evSetup, const reco::Photon &cand) const
reco::IsoDeposit::Veto veto(const reco::IsoDeposit::Direction &dir) const
double vz() const override
z coordinate of vertex position
double vz() const
z coordinate of the reference point on track
Definition: TrackBase.h:626
std::vector< reco::PFCandidate > PFCandidateCollection
collection of PFCandidates
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:40
fixed size matrix
HLT enums.
SuperClusterRef superCluster() const override
reference to a SuperCluster
Definition: GsfElectron.h:155
#define DEFINE_EDM_PLUGIN(factory, type, name)
reco::IsoDeposit deposit(const edm::Event &ev, const edm::EventSetup &evSetup, const reco::Track &muon) const override
double phi() const final
momentum azimuthal angle
Definition: event.py:1
double vz() const override
z coordinate of vertex position
Definition: PFCandidate.h:411
reco::SuperClusterRef superClusterRef() const
return a reference to the corresponding SuperCluster if any
Definition: PFCandidate.cc:558