File indexing completed on 2024-04-06 12:25:06
0001 #include "PhysicsTools/SelectorUtils/interface/CutApplicatorWithEventContentBase.h"
0002 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
0003 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0004 #include "DataFormats/VertexReco/interface/Vertex.h"
0005 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
0006
0007 class GsfEleDxyCut : public CutApplicatorWithEventContentBase {
0008 public:
0009 GsfEleDxyCut(const edm::ParameterSet& c);
0010
0011 result_type operator()(const reco::GsfElectronPtr&) const final;
0012
0013 void setConsumes(edm::ConsumesCollector&) final;
0014 void getEventContent(const edm::EventBase&) final;
0015
0016 double value(const reco::CandidatePtr& cand) const final;
0017
0018 CandidateType candidateType() const final { return ELECTRON; }
0019
0020 private:
0021 const double _dxyCutValueEB, _dxyCutValueEE, _barrelCutOff;
0022 edm::Handle<reco::VertexCollection> _vtxs;
0023 };
0024
0025 DEFINE_EDM_PLUGIN(CutApplicatorFactory, GsfEleDxyCut, "GsfEleDxyCut");
0026
0027 GsfEleDxyCut::GsfEleDxyCut(const edm::ParameterSet& c)
0028 : CutApplicatorWithEventContentBase(c),
0029 _dxyCutValueEB(c.getParameter<double>("dxyCutValueEB")),
0030 _dxyCutValueEE(c.getParameter<double>("dxyCutValueEE")),
0031 _barrelCutOff(c.getParameter<double>("barrelCutOff")) {
0032 edm::InputTag vertextag = c.getParameter<edm::InputTag>("vertexSrc");
0033 edm::InputTag vertextagMiniAOD = c.getParameter<edm::InputTag>("vertexSrcMiniAOD");
0034 contentTags_.emplace("vertices", vertextag);
0035 contentTags_.emplace("verticesMiniAOD", vertextagMiniAOD);
0036 }
0037
0038 void GsfEleDxyCut::setConsumes(edm::ConsumesCollector& cc) {
0039 auto vtcs = cc.mayConsume<reco::VertexCollection>(contentTags_["vertices"]);
0040 auto vtcsMiniAOD = cc.mayConsume<reco::VertexCollection>(contentTags_["verticesMiniAOD"]);
0041 contentTokens_.emplace("vertices", vtcs);
0042 contentTokens_.emplace("verticesMiniAOD", vtcsMiniAOD);
0043 }
0044
0045 void GsfEleDxyCut::getEventContent(const edm::EventBase& ev) {
0046
0047 ev.getByLabel(contentTags_["vertices"], _vtxs);
0048 if (!_vtxs.isValid())
0049 ev.getByLabel(contentTags_["verticesMiniAOD"], _vtxs);
0050 }
0051
0052 CutApplicatorBase::result_type GsfEleDxyCut::operator()(const reco::GsfElectronPtr& cand) const {
0053 const float dxyCutValue =
0054 (std::abs(cand->superCluster()->position().eta()) < _barrelCutOff ? _dxyCutValueEB : _dxyCutValueEE);
0055
0056 const reco::VertexCollection& vtxs = *_vtxs;
0057 const double dxy = (!vtxs.empty() ? cand->gsfTrack()->dxy(vtxs[0].position()) : cand->gsfTrack()->dxy());
0058 return std::abs(dxy) < dxyCutValue;
0059 }
0060
0061 double GsfEleDxyCut::value(const reco::CandidatePtr& cand) const {
0062 reco::GsfElectronPtr ele(cand);
0063 const reco::VertexCollection& vtxs = *_vtxs;
0064 const double dxy = (!vtxs.empty() ? ele->gsfTrack()->dxy(vtxs[0].position()) : ele->gsfTrack()->dxy());
0065 return std::abs(dxy);
0066 }