CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
GsfEleDzCut.cc
Go to the documentation of this file.
6 
8 public:
10 
11  result_type operator()(const reco::GsfElectronPtr&) const override final;
12 
13  void setConsumes(edm::ConsumesCollector&) override final;
14  void getEventContent(const edm::EventBase&) override final;
15 
16  double value(const reco::CandidatePtr& cand) const override final;
17 
18  CandidateType candidateType() const override final {
19  return ELECTRON;
20  }
21 
22 private:
25 };
26 
29  "GsfEleDzCut");
30 
33  _dzCutValueEB(c.getParameter<double>("dzCutValueEB")),
34  _dzCutValueEE(c.getParameter<double>("dzCutValueEE")),
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 
56 CutApplicatorBase::result_type
58 operator()(const reco::GsfElectronPtr& cand) const{
59  const float dzCutValue =
60  ( std::abs(cand->superCluster()->position().eta()) < _barrelCutOff ?
62 
63  const reco::VertexCollection& vtxs = *_vtxs;
64  const double dz = ( vtxs.size() ?
65  cand->gsfTrack()->dz(vtxs[0].position()) :
66  cand->gsfTrack()->dz() );
67  return std::abs(dz) < dzCutValue;
68 }
69 
70 double GsfEleDzCut::value(const reco::CandidatePtr& cand) const {
71  reco::GsfElectronPtr ele(cand);
72  const reco::VertexCollection& vtxs = *_vtxs;
73  const double dz = ( vtxs.size() ?
74  ele->gsfTrack()->dz(vtxs[0].position()) :
75  ele->gsfTrack()->dz() );
76  return std::abs(dz);
77 }
T getParameter(std::string const &) const
CandidateType candidateType() const overridefinal
Definition: GsfEleDzCut.cc:18
double value(const reco::CandidatePtr &cand) const overridefinal
Definition: GsfEleDzCut.cc:70
EDGetTokenT< ProductType > mayConsume(edm::InputTag const &tag)
std::unordered_map< std::string, edm::EDGetToken > contentTokens_
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
bool ev
std::unordered_map< std::string, edm::InputTag > contentTags_
const double _barrelCutOff
Definition: GsfEleDzCut.cc:23
result_type operator()(const reco::GsfElectronPtr &) const overridefinal
Definition: GsfEleDzCut.cc:58
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void setConsumes(edm::ConsumesCollector &) overridefinal
Definition: GsfEleDzCut.cc:42
bool isValid() const
Definition: HandleBase.h:75
void getEventContent(const edm::EventBase &) overridefinal
Definition: GsfEleDzCut.cc:49
GsfEleDzCut(const edm::ParameterSet &c)
Definition: GsfEleDzCut.cc:31
const double _dzCutValueEB
Definition: GsfEleDzCut.cc:23
string const
Definition: compareJSON.py:14
const double _dzCutValueEE
Definition: GsfEleDzCut.cc:23
bool getByLabel(InputTag const &, Handle< T > &) const
Definition: EventBase.h:90
static int position[264][3]
Definition: ReadPGInfo.cc:509
edm::Handle< reco::VertexCollection > _vtxs
Definition: GsfEleDzCut.cc:24
#define DEFINE_EDM_PLUGIN(factory, type, name)