CMS 3D CMS Logo

PixelVertexCollectionTrimmer.cc
Go to the documentation of this file.
1 // Original Author: Riccardo Manzoni
2 #include <vector>
3 #include <memory>
4 #include <numeric>
5 #include <algorithm>
6 
18 
20 public:
22 
23  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
24 
25 private:
26  void produce(edm::Event&, const edm::EventSetup&) override;
27 
29  uint const maxVtx_;
30  double const fractionSumPt2_;
31  double const minSumPt2_;
32 
33  std::unique_ptr<PVClusterComparer> pvComparer_;
34 };
35 
37  : vtxToken_(consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("src"))),
38  maxVtx_(iConfig.getParameter<uint>("maxVtx")),
39  fractionSumPt2_(iConfig.getParameter<double>("fractionSumPt2")),
40  minSumPt2_(iConfig.getParameter<double>("minSumPt2")) {
41  if (fractionSumPt2_ > 1)
42  throw cms::Exception("PixelVertexConfiguration") << "value of \"fractionSumPt2\" is larger than 1.";
43 
44  auto const& pvComparerPSet = iConfig.getParameterSet("PVcomparer");
45  auto const track_pt_min = pvComparerPSet.getParameter<double>("track_pt_min");
46  auto const track_pt_max = pvComparerPSet.getParameter<double>("track_pt_max");
47  auto const track_chi2_max = pvComparerPSet.getParameter<double>("track_chi2_max");
48  auto const track_prob_min = pvComparerPSet.getParameter<double>("track_prob_min");
49 
51  throw cms::Exception("PixelVertexConfiguration")
52  << "PVcomparer.track_pt_min (" << track_pt_min << ") >= PVcomparer.track_pt_max (" << track_pt_max
53  << ") : PVClusterComparer will use pT=" << track_pt_max << " for all selected tracks.";
54 
55  pvComparer_ = std::make_unique<PVClusterComparer>(track_pt_min, track_pt_max, track_chi2_max, track_prob_min);
56 
57  produces<reco::VertexCollection>();
58 }
59 
61  auto vtxs_trim = std::make_unique<reco::VertexCollection>();
62 
63  auto const& vtxs = iEvent.get(vtxToken_);
64 
65  if (vtxs.empty())
66  edm::LogWarning("PixelVertexInput") << "Input collection of vertices is empty. Output collection will be empty.";
67  else {
68  std::vector<double> foms(vtxs.size());
69  for (size_t idx = 0; idx < vtxs.size(); ++idx)
70  foms[idx] = pvComparer_->pTSquaredSum(vtxs[idx]);
71 
72  std::vector<size_t> sortIdxs(vtxs.size());
73  std::iota(sortIdxs.begin(), sortIdxs.end(), 0);
74  std::sort(sortIdxs.begin(), sortIdxs.end(), [&](size_t const i1, size_t const i2) { return foms[i1] > foms[i2]; });
75 
76  auto const minFOM_fromFrac = foms[sortIdxs.front()] * fractionSumPt2_;
77 
78  vtxs_trim->reserve(std::min((size_t)maxVtx_, vtxs.size()));
79  for (auto const idx : sortIdxs) {
80  if (vtxs_trim->size() >= maxVtx_)
81  break;
82  if (foms[idx] >= minFOM_fromFrac and foms[idx] > minSumPt2_)
83  vtxs_trim->emplace_back(vtxs[idx]);
84  }
85 
86  if (vtxs_trim->empty())
87  edm::LogInfo("PixelVertexOutput") << "Output collection is empty.";
88  }
89 
90  iEvent.put(std::move(vtxs_trim));
91 }
92 
95  desc.add<edm::InputTag>("src", edm::InputTag(""))->setComment("input (pixel) vertex collection");
96  desc.add<uint>("maxVtx", 100)->setComment("max output collection size (number of accepted vertices)");
97  desc.add<double>("fractionSumPt2", 0.3)->setComment("threshold on sumPt2 fraction of the leading vertex");
98  desc.add<double>("minSumPt2", 0.)->setComment("min sumPt2");
99  edm::ParameterSetDescription PVcomparerPSet;
100  PVcomparerPSet.add<double>("track_pt_min", 1.0)->setComment("min track p_T");
101  PVcomparerPSet.add<double>("track_pt_max", 10.0)->setComment("max track p_T");
102  PVcomparerPSet.add<double>("track_chi2_max", 99999.)->setComment("max track chi2");
103  PVcomparerPSet.add<double>("track_prob_min", -1.)->setComment("min track prob");
104  desc.add<edm::ParameterSetDescription>("PVcomparer", PVcomparerPSet)
105  ->setComment("from RecoTracker/PixelVertexFinding/python/PVClusterComparer_cfi.py");
106  descriptions.add("hltPixelVertexCollectionTrimmer", desc);
107 }
108 
void setComment(std::string const &value)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
ParameterSet const & getParameterSet(std::string const &) const
std::vector< Vertex > VertexCollection
Definition: Vertex.h:31
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void produce(edm::Event &, const edm::EventSetup &) override
void add(std::string const &label, ParameterSetDescription const &psetDescription)
PixelVertexCollectionTrimmer(const edm::ParameterSet &)
fixed size matrix
HLT enums.
edm::EDGetTokenT< reco::VertexCollection > const vtxToken_
std::unique_ptr< PVClusterComparer > pvComparer_
def move(src, dest)
Definition: eostools.py:511