File indexing completed on 2024-04-06 12:33:17
0001
0002
0003
0004 #include "FWCore/Framework/interface/Event.h"
0005 #include "FWCore/Framework/interface/MakerMacros.h"
0006 #include "FWCore/Framework/interface/global/EDProducer.h"
0007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0008 #include "FWCore/Utilities/interface/InputTag.h"
0009
0010 #include "CommonTools/Utils/interface/StringCutObjectSelector.h"
0011
0012 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
0013 #include "DataFormats/EgammaCandidates/interface/GsfElectronFwd.h"
0014 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
0015 #include "DataFormats/GsfTrackReco/interface/GsfTrackExtra.h"
0016 #include "DataFormats/GsfTrackReco/interface/GsfTrackExtraFwd.h"
0017 #include "DataFormats/GsfTrackReco/interface/GsfTrackFwd.h"
0018 #include "DataFormats/VertexReco/interface/Vertex.h"
0019 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0020
0021 #include <algorithm>
0022 #include <cmath>
0023 #include <memory>
0024 #include <vector>
0025
0026
0027
0028
0029 class GsfElectronFromPVSelector : public edm::global::EDProducer<> {
0030 public:
0031 explicit GsfElectronFromPVSelector(edm::ParameterSet const&);
0032 void produce(edm::StreamID, edm::Event&, edm::EventSetup const&) const override;
0033
0034 private:
0035 double max_dxy_;
0036 double max_dz_;
0037 edm::EDGetTokenT<std::vector<reco::Vertex>> v_recoVertexToken_;
0038 edm::EDGetTokenT<std::vector<reco::GsfElectron>> v_recoGsfElectronToken_;
0039 };
0040
0041
0042
0043
0044
0045 GsfElectronFromPVSelector::GsfElectronFromPVSelector(edm::ParameterSet const& iConfig)
0046 : max_dxy_{iConfig.getParameter<double>("max_dxy")},
0047 max_dz_{iConfig.getParameter<double>("max_dz")},
0048 v_recoVertexToken_{consumes<std::vector<reco::Vertex>>(iConfig.getParameter<edm::InputTag>("srcVertex"))},
0049 v_recoGsfElectronToken_{
0050 consumes<std::vector<reco::GsfElectron>>(iConfig.getParameter<edm::InputTag>("srcElectron"))} {
0051 produces<std::vector<reco::GsfElectron>>();
0052 }
0053
0054
0055
0056
0057
0058 void GsfElectronFromPVSelector::produce(edm::StreamID, edm::Event& iEvent, edm::EventSetup const&) const {
0059 edm::Handle<std::vector<reco::Vertex>> vertices;
0060 iEvent.getByToken(v_recoVertexToken_, vertices);
0061
0062 edm::Handle<std::vector<reco::GsfElectron>> gsfElectrons;
0063 iEvent.getByToken(v_recoGsfElectronToken_, gsfElectrons);
0064
0065 auto goodGsfElectrons = std::make_unique<std::vector<reco::GsfElectron>>();
0066
0067 if (!vertices->empty() && !gsfElectrons->empty()) {
0068 auto const& pv = vertices->front();
0069 std::copy_if(std::cbegin(*gsfElectrons),
0070 std::cend(*gsfElectrons),
0071 std::back_inserter(*goodGsfElectrons),
0072 [this, &pv](auto const& GsfElectron) {
0073 return std::abs(GsfElectron.gsfTrack()->dxy(pv.position())) < max_dxy_ &&
0074 std::abs(GsfElectron.gsfTrack()->dz(pv.position())) < max_dz_;
0075 });
0076 }
0077
0078 iEvent.put(std::move(goodGsfElectrons));
0079 }
0080
0081 DEFINE_FWK_MODULE(GsfElectronFromPVSelector);