File indexing completed on 2024-04-06 12:33:18
0001
0002
0003
0004 #include "DataFormats/VertexReco/interface/Vertex.h"
0005 #include "FWCore/Framework/interface/Event.h"
0006 #include "FWCore/Framework/interface/MakerMacros.h"
0007 #include "FWCore/Framework/interface/global/EDProducer.h"
0008 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0009 #include "FWCore/Utilities/interface/InputTag.h"
0010
0011 #include <memory>
0012 #include <vector>
0013 #include <sstream>
0014
0015
0016
0017
0018 class bestPVselector : public edm::global::EDProducer<> {
0019 public:
0020 explicit bestPVselector(edm::ParameterSet const& iConfig);
0021 void produce(edm::StreamID, edm::Event&, edm::EventSetup const&) const override;
0022
0023 private:
0024 edm::EDGetTokenT<std::vector<reco::Vertex>> src_;
0025 };
0026
0027
0028
0029
0030
0031
0032 bestPVselector::bestPVselector(edm::ParameterSet const& iConfig)
0033 : src_{consumes<std::vector<reco::Vertex>>(iConfig.getParameter<edm::InputTag>("src"))} {
0034 produces<std::vector<reco::Vertex>>();
0035 }
0036
0037
0038
0039
0040
0041
0042 void bestPVselector::produce(edm::StreamID, edm::Event& iEvent, edm::EventSetup const&) const {
0043 edm::Handle<std::vector<reco::Vertex>> vertices;
0044 iEvent.getByToken(src_, vertices);
0045
0046 auto theBestPV = std::make_unique<std::vector<reco::Vertex>>();
0047
0048 if (!vertices->empty()) {
0049 auto sumSquarePt = [](auto const& pv) { return pv.p4().pt() * pv.p4().pt(); };
0050 auto bestPV =
0051 std::max_element(std::cbegin(*vertices), std::cend(*vertices), [sumSquarePt](auto const& v1, auto const& v2) {
0052 return sumSquarePt(v1) < sumSquarePt(v2);
0053 });
0054 theBestPV->push_back(*bestPV);
0055 }
0056 iEvent.put(std::move(theBestPV));
0057 }
0058
0059 using HighestSumP4PrimaryVertexSelector = bestPVselector;
0060 DEFINE_FWK_MODULE(HighestSumP4PrimaryVertexSelector);