CMS 3D CMS Logo

GsfEleDxyCut.cc
Go to the documentation of this file.
6 
8 public:
10 
11  result_type operator()(const reco::GsfElectronPtr&) const final;
12 
14  void getEventContent(const edm::EventBase&) final;
15 
16  double value(const reco::CandidatePtr& cand) const final;
17 
19 
20 private:
23 };
24 
26 
29  _dxyCutValueEB(c.getParameter<double>("dxyCutValueEB")),
30  _dxyCutValueEE(c.getParameter<double>("dxyCutValueEE")),
31  _barrelCutOff(c.getParameter<double>("barrelCutOff")) {
32  edm::InputTag vertextag = c.getParameter<edm::InputTag>("vertexSrc");
33  edm::InputTag vertextagMiniAOD = c.getParameter<edm::InputTag>("vertexSrcMiniAOD");
34  contentTags_.emplace("vertices", vertextag);
35  contentTags_.emplace("verticesMiniAOD", vertextagMiniAOD);
36 }
37 
39  auto vtcs = cc.mayConsume<reco::VertexCollection>(contentTags_["vertices"]);
40  auto vtcsMiniAOD = cc.mayConsume<reco::VertexCollection>(contentTags_["verticesMiniAOD"]);
41  contentTokens_.emplace("vertices", vtcs);
42  contentTokens_.emplace("verticesMiniAOD", vtcsMiniAOD);
43 }
44 
46  // First try AOD, then go to miniAOD. Use the same Handle since collection class is the same.
47  ev.getByLabel(contentTags_["vertices"], _vtxs);
48  if (!_vtxs.isValid())
49  ev.getByLabel(contentTags_["verticesMiniAOD"], _vtxs);
50 }
51 
53  const float dxyCutValue =
54  (std::abs(cand->superCluster()->position().eta()) < _barrelCutOff ? _dxyCutValueEB : _dxyCutValueEE);
55 
56  const reco::VertexCollection& vtxs = *_vtxs;
57  const double dxy = (!vtxs.empty() ? cand->gsfTrack()->dxy(vtxs[0].position()) : cand->gsfTrack()->dxy());
58  return std::abs(dxy) < dxyCutValue;
59 }
60 
63  const reco::VertexCollection& vtxs = *_vtxs;
64  const double dxy = (!vtxs.empty() ? ele->gsfTrack()->dxy(vtxs[0].position()) : ele->gsfTrack()->dxy());
65  return std::abs(dxy);
66 }
std::unordered_map< std::string, edm::InputTag > contentTags_
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
const double _dxyCutValueEB
Definition: GsfEleDxyCut.cc:21
result_type operator()(const reco::GsfElectronPtr &) const final
Definition: GsfEleDxyCut.cc:52
GsfTrackRef gsfTrack() const override
reference to a GsfTrack
Definition: GsfElectron.h:156
double value(const reco::CandidatePtr &cand) const final
Definition: GsfEleDxyCut.cc:61
edm::Handle< reco::VertexCollection > _vtxs
Definition: GsfEleDxyCut.cc:22
const double _dxyCutValueEE
Definition: GsfEleDxyCut.cc:21
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
CandidateType candidateType() const final
Definition: GsfEleDxyCut.cc:18
void setConsumes(edm::ConsumesCollector &) final
Definition: GsfEleDxyCut.cc:38
void getEventContent(const edm::EventBase &) final
Definition: GsfEleDxyCut.cc:45
GsfEleDxyCut(const edm::ParameterSet &c)
Definition: GsfEleDxyCut.cc:27
std::unordered_map< std::string, edm::EDGetToken > contentTokens_
const double _barrelCutOff
Definition: GsfEleDxyCut.cc:21
bool isValid() const
Definition: HandleBase.h:70
static int position[264][3]
Definition: ReadPGInfo.cc:289
#define DEFINE_EDM_PLUGIN(factory, type, name)