CMS 3D CMS Logo

All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
GsfEleConversionVetoCut.cc
Go to the documentation of this file.
6 
7 
9 public:
11 
12  result_type operator()(const reco::GsfElectronPtr&) const override final;
13 
14  void setConsumes(edm::ConsumesCollector&) override final;
15  void getEventContent(const edm::EventBase&) override final;
16 
17  CandidateType candidateType() const override final {
18  return ELECTRON;
19  }
20 
21 private:
24 };
25 
28  "GsfEleConversionVetoCut");
29 
32 
33  edm::InputTag conversiontag = c.getParameter<edm::InputTag>("conversionSrc");
34  contentTags_.emplace("conversions",conversiontag);
35 
36  edm::InputTag conversiontagMiniAOD = c.getParameter<edm::InputTag>("conversionSrcMiniAOD");
37  contentTags_.emplace("conversionsMiniAOD",conversiontagMiniAOD);
38 
39  edm::InputTag beamspottag = c.getParameter<edm::InputTag>("beamspotSrc");
40  contentTags_.emplace("beamspot",beamspottag);
41 }
42 
44  auto convs = cc.mayConsume<reco::ConversionCollection>(contentTags_["conversions"]);
45  auto convsMiniAOD = cc.mayConsume<reco::ConversionCollection>(contentTags_["conversionsMiniAOD"]);
46  auto thebs = cc.consumes<reco::BeamSpot>(contentTags_["beamspot"]);
47  contentTokens_.emplace("conversions",convs);
48  contentTokens_.emplace("conversionsMiniAOD",convsMiniAOD);
49  contentTokens_.emplace("beamspot",thebs);
50 }
51 
53 
54  // First try AOD, then go to miniAOD. Use the same Handle since collection class is the same.
55  ev.getByLabel(contentTags_["conversions"],_convs);
56  if (!_convs.isValid())
57  ev.getByLabel(contentTags_["conversionsMiniAOD"],_convs);
58 
59  ev.getByLabel(contentTags_["beamspot"],_thebs);
60 }
61 
62 CutApplicatorBase::result_type
64 operator()(const reco::GsfElectronPtr& cand) const{
65  if( _thebs.isValid() && _convs.isValid() ) {
67  _thebs->position());
68  } else {
69  edm::LogWarning("GsfEleConversionVetoCut")
70  << "Couldn't find a necessary collection, returning true!";
71  }
72  return true;
73 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
EDGetTokenT< ProductType > mayConsume(edm::InputTag const &tag)
edm::Handle< reco::BeamSpot > _thebs
CandidateType candidateType() const overridefinal
void setConsumes(edm::ConsumesCollector &) overridefinal
std::unordered_map< std::string, edm::EDGetToken > contentTokens_
bool ev
std::unordered_map< std::string, edm::InputTag > contentTags_
std::vector< Conversion > ConversionCollection
collectin of Conversion objects
Definition: ConversionFwd.h:9
edm::Handle< reco::ConversionCollection > _convs
GsfEleConversionVetoCut(const edm::ParameterSet &c)
bool isValid() const
Definition: HandleBase.h:75
static bool hasMatchedConversion(const reco::GsfElectron &ele, const edm::Handle< reco::ConversionCollection > &convCol, const math::XYZPoint &beamspot, bool allowCkfMatch=true, float lxyMin=2.0, float probMin=1e-6, unsigned int nHitsBeforeVtxMax=0)
void getEventContent(const edm::EventBase &) overridefinal
result_type operator()(const reco::GsfElectronPtr &) const overridefinal
string const
Definition: compareJSON.py:14
bool getByLabel(InputTag const &, Handle< T > &) const
Definition: EventBase.h:92
#define DEFINE_EDM_PLUGIN(factory, type, name)