CMS 3D CMS Logo

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  double value(const reco::CandidatePtr& cand) const override final;
18 
19  CandidateType candidateType() const override final {
20  return ELECTRON;
21  }
22 
23 private:
26 };
27 
30  "GsfEleConversionVetoCut");
31 
34 
35  edm::InputTag conversiontag = c.getParameter<edm::InputTag>("conversionSrc");
36  contentTags_.emplace("conversions",conversiontag);
37 
38  edm::InputTag conversiontagMiniAOD = c.getParameter<edm::InputTag>("conversionSrcMiniAOD");
39  contentTags_.emplace("conversionsMiniAOD",conversiontagMiniAOD);
40 
41  edm::InputTag beamspottag = c.getParameter<edm::InputTag>("beamspotSrc");
42  contentTags_.emplace("beamspot",beamspottag);
43 }
44 
46  auto convs = cc.mayConsume<reco::ConversionCollection>(contentTags_["conversions"]);
47  auto convsMiniAOD = cc.mayConsume<reco::ConversionCollection>(contentTags_["conversionsMiniAOD"]);
48  auto thebs = cc.consumes<reco::BeamSpot>(contentTags_["beamspot"]);
49  contentTokens_.emplace("conversions",convs);
50  contentTokens_.emplace("conversionsMiniAOD",convsMiniAOD);
51  contentTokens_.emplace("beamspot",thebs);
52 }
53 
55 
56  // First try AOD, then go to miniAOD. Use the same Handle since collection class is the same.
57  ev.getByLabel(contentTags_["conversions"],_convs);
58  if (!_convs.isValid())
59  ev.getByLabel(contentTags_["conversionsMiniAOD"],_convs);
60 
61  ev.getByLabel(contentTags_["beamspot"],_thebs);
62 }
63 
64 CutApplicatorBase::result_type
67  if( _thebs.isValid() && _convs.isValid() ) {
69  _thebs->position());
70  } else {
71  edm::LogWarning("GsfEleConversionVetoCut")
72  << "Couldn't find a necessary collection, returning true!";
73  }
74  return true;
75 }
76 
78  reco::GsfElectronPtr ele(cand);
79  if( _thebs.isValid() && _convs.isValid() ) {
81  _thebs->position());
82  } else {
83  edm::LogWarning("GsfEleConversionVetoCut")
84  << "Couldn't find a necessary collection, returning true!";
85  return true;
86  }
87 }
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
std::unordered_map< std::string, edm::EDGetToken > contentTokens_
void getEventContent(const edm::EventBase &) override final
bool ev
std::unordered_map< std::string, edm::InputTag > contentTags_
std::vector< Conversion > ConversionCollection
collectin of Conversion objects
Definition: ConversionFwd.h:9
double value(const reco::CandidatePtr &cand) const override final
result_type operator()(const reco::GsfElectronPtr &) const override final
edm::Handle< reco::ConversionCollection > _convs
GsfEleConversionVetoCut(const edm::ParameterSet &c)
CandidateType candidateType() const override final
bool isValid() const
Definition: HandleBase.h:74
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)
bool getByLabel(InputTag const &, Handle< T > &) const
Definition: EventBase.h:89
void setConsumes(edm::ConsumesCollector &) override final
const Point & position() const
position
Definition: BeamSpot.h:62
#define DEFINE_EDM_PLUGIN(factory, type, name)