File indexing completed on 2024-09-07 04:38:05
0001
0002 #include <vector>
0003 #include <memory>
0004 #include <numeric>
0005 #include <algorithm>
0006
0007 #include "DataFormats/VertexReco/interface/Vertex.h"
0008 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0009 #include "FWCore/Framework/interface/Frameworkfwd.h"
0010 #include "FWCore/Framework/interface/stream/EDProducer.h"
0011 #include "FWCore/Framework/interface/Event.h"
0012 #include "FWCore/Framework/interface/MakerMacros.h"
0013 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0014 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0015 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0016 #include "FWCore/Utilities/interface/Exception.h"
0017 #include "RecoVertex/PixelVertexFinding/interface/PVClusterComparer.h"
0018
0019 class PixelVertexCollectionTrimmer : public edm::stream::EDProducer<> {
0020 public:
0021 explicit PixelVertexCollectionTrimmer(const edm::ParameterSet&);
0022
0023 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0024
0025 private:
0026 void produce(edm::Event&, const edm::EventSetup&) override;
0027
0028 edm::EDGetTokenT<reco::VertexCollection> const vtxToken_;
0029 uint const maxVtx_;
0030 double const fractionSumPt2_;
0031 double const minSumPt2_;
0032
0033 std::unique_ptr<PVClusterComparer> pvComparer_;
0034 };
0035
0036 PixelVertexCollectionTrimmer::PixelVertexCollectionTrimmer(const edm::ParameterSet& iConfig)
0037 : vtxToken_(consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("src"))),
0038 maxVtx_(iConfig.getParameter<uint>("maxVtx")),
0039 fractionSumPt2_(iConfig.getParameter<double>("fractionSumPt2")),
0040 minSumPt2_(iConfig.getParameter<double>("minSumPt2")) {
0041 if (fractionSumPt2_ > 1)
0042 throw cms::Exception("PixelVertexConfiguration") << "value of \"fractionSumPt2\" is larger than 1.";
0043
0044 auto const& pvComparerPSet = iConfig.getParameterSet("PVcomparer");
0045 auto const track_pt_min = pvComparerPSet.getParameter<double>("track_pt_min");
0046 auto const track_pt_max = pvComparerPSet.getParameter<double>("track_pt_max");
0047 auto const track_chi2_max = pvComparerPSet.getParameter<double>("track_chi2_max");
0048 auto const track_prob_min = pvComparerPSet.getParameter<double>("track_prob_min");
0049
0050 if (track_pt_min >= track_pt_max)
0051 throw cms::Exception("PixelVertexConfiguration")
0052 << "PVcomparer.track_pt_min (" << track_pt_min << ") >= PVcomparer.track_pt_max (" << track_pt_max
0053 << ") : PVClusterComparer will use pT=" << track_pt_max << " for all selected tracks.";
0054
0055 pvComparer_ = std::make_unique<PVClusterComparer>(track_pt_min, track_pt_max, track_chi2_max, track_prob_min);
0056
0057 produces<reco::VertexCollection>();
0058 }
0059
0060 void PixelVertexCollectionTrimmer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0061 auto vtxs_trim = std::make_unique<reco::VertexCollection>();
0062
0063 auto const& vtxs = iEvent.get(vtxToken_);
0064
0065 if (vtxs.empty())
0066 edm::LogWarning("PixelVertexInput") << "Input collection of vertices is empty. Output collection will be empty.";
0067 else {
0068 std::vector<double> foms(vtxs.size());
0069 for (size_t idx = 0; idx < vtxs.size(); ++idx)
0070 foms[idx] = pvComparer_->pTSquaredSum(vtxs[idx]);
0071
0072 std::vector<size_t> sortIdxs(vtxs.size());
0073 std::iota(sortIdxs.begin(), sortIdxs.end(), 0);
0074 std::sort(sortIdxs.begin(), sortIdxs.end(), [&](size_t const i1, size_t const i2) { return foms[i1] > foms[i2]; });
0075
0076 auto const minFOM_fromFrac = foms[sortIdxs.front()] * fractionSumPt2_;
0077
0078 vtxs_trim->reserve(std::min((size_t)maxVtx_, vtxs.size()));
0079 for (auto const idx : sortIdxs) {
0080 if (vtxs_trim->size() >= maxVtx_)
0081 break;
0082 if (foms[idx] >= minFOM_fromFrac and foms[idx] > minSumPt2_)
0083 vtxs_trim->emplace_back(vtxs[idx]);
0084 }
0085
0086 if (vtxs_trim->empty())
0087 edm::LogInfo("PixelVertexOutput") << "Output collection is empty.";
0088 }
0089
0090 iEvent.put(std::move(vtxs_trim));
0091 }
0092
0093 void PixelVertexCollectionTrimmer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0094 edm::ParameterSetDescription desc;
0095 desc.add<edm::InputTag>("src", edm::InputTag(""))->setComment("input (pixel) vertex collection");
0096 desc.add<uint>("maxVtx", 100)->setComment("max output collection size (number of accepted vertices)");
0097 desc.add<double>("fractionSumPt2", 0.3)->setComment("threshold on sumPt2 fraction of the leading vertex");
0098 desc.add<double>("minSumPt2", 0.)->setComment("min sumPt2");
0099 edm::ParameterSetDescription PVcomparerPSet;
0100 PVcomparerPSet.add<double>("track_pt_min", 1.0)->setComment("min track p_T");
0101 PVcomparerPSet.add<double>("track_pt_max", 10.0)->setComment("max track p_T");
0102 PVcomparerPSet.add<double>("track_chi2_max", 99999.)->setComment("max track chi2");
0103 PVcomparerPSet.add<double>("track_prob_min", -1.)->setComment("min track prob");
0104 desc.add<edm::ParameterSetDescription>("PVcomparer", PVcomparerPSet)
0105 ->setComment("from RecoVertex/PixelVertexFinding/python/PVClusterComparer_cfi.py");
0106 descriptions.add("hltPixelVertexCollectionTrimmer", desc);
0107 }
0108
0109 DEFINE_FWK_MODULE(PixelVertexCollectionTrimmer);