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 
18  CandidateType candidateType() const final {
19  return ELECTRON;
20  }
21 
22 private:
25 };
26 
29  "GsfEleDxyCut");
30 
33  _dxyCutValueEB(c.getParameter<double>("dxyCutValueEB")),
34  _dxyCutValueEE(c.getParameter<double>("dxyCutValueEE")),
35  _barrelCutOff(c.getParameter<double>("barrelCutOff")) {
36  edm::InputTag vertextag = c.getParameter<edm::InputTag>("vertexSrc");
37  edm::InputTag vertextagMiniAOD = c.getParameter<edm::InputTag>("vertexSrcMiniAOD");
38  contentTags_.emplace("vertices",vertextag);
39  contentTags_.emplace("verticesMiniAOD",vertextagMiniAOD);
40 }
41 
43  auto vtcs = cc.mayConsume<reco::VertexCollection>(contentTags_["vertices"]);
44  auto vtcsMiniAOD = cc.mayConsume<reco::VertexCollection>(contentTags_["verticesMiniAOD"]);
45  contentTokens_.emplace("vertices",vtcs);
46  contentTokens_.emplace("verticesMiniAOD",vtcsMiniAOD);
47 }
48 
50  // First try AOD, then go to miniAOD. Use the same Handle since collection class is the same.
51  ev.getByLabel(contentTags_["vertices"],_vtxs);
52  if (!_vtxs.isValid())
53  ev.getByLabel(contentTags_["verticesMiniAOD"],_vtxs);
54 }
55 
59  const float dxyCutValue =
60  ( std::abs(cand->superCluster()->position().eta()) < _barrelCutOff ?
62 
63  const reco::VertexCollection& vtxs = *_vtxs;
64  const double dxy = ( !vtxs.empty() ?
65  cand->gsfTrack()->dxy(vtxs[0].position()) :
66  cand->gsfTrack()->dxy() );
67  return std::abs(dxy) < dxyCutValue;
68 }
69 
71  reco::GsfElectronPtr ele(cand);
72  const reco::VertexCollection& vtxs = *_vtxs;
73  const double dxy = ( !vtxs.empty() ?
74  ele->gsfTrack()->dxy(vtxs[0].position()) :
75  ele->gsfTrack()->dxy() );
76  return std::abs(dxy);
77 }
T getParameter(std::string const &) const
GsfTrackRef gsfTrack() const override
reference to a GsfTrack
Definition: GsfElectron.h:186
EDGetTokenT< ProductType > mayConsume(edm::InputTag const &tag)
std::unordered_map< std::string, edm::EDGetToken > contentTokens_
CandidateType candidateType() const final
Definition: GsfEleDxyCut.cc:18
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
bool ev
const double _dxyCutValueEB
Definition: GsfEleDxyCut.cc:23
std::unordered_map< std::string, edm::InputTag > contentTags_
edm::Handle< reco::VertexCollection > _vtxs
Definition: GsfEleDxyCut.cc:24
const double _dxyCutValueEE
Definition: GsfEleDxyCut.cc:23
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void setConsumes(edm::ConsumesCollector &) final
Definition: GsfEleDxyCut.cc:42
bool isValid() const
Definition: HandleBase.h:74
void getEventContent(const edm::EventBase &) final
Definition: GsfEleDxyCut.cc:49
GsfEleDxyCut(const edm::ParameterSet &c)
Definition: GsfEleDxyCut.cc:31
const double _barrelCutOff
Definition: GsfEleDxyCut.cc:23
bool getByLabel(InputTag const &, Handle< T > &) const
Definition: EventBase.h:92
static int position[264][3]
Definition: ReadPGInfo.cc:509
SuperClusterRef superCluster() const override
reference to a SuperCluster
Definition: GsfElectron.h:185
#define DEFINE_EDM_PLUGIN(factory, type, name)
result_type operator()(const reco::GsfElectronPtr &) const final
Definition: GsfEleDxyCut.cc:58
double value(const reco::CandidatePtr &cand) const final
Definition: GsfEleDxyCut.cc:70