CMS 3D CMS Logo

RecoTauPileUpVertexSelector.cc
Go to the documentation of this file.
1 /*
2  * Select reco:Vertices consistent with pileup.
3  *
4  * Author: Evan K. Friis
5  *
6  */
7 
13 
16 
19 
20 #include <algorithm>
21 
23 public:
26  bool filter(edm::Event& evt, const edm::EventSetup& es) override;
27  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
28 
29 private:
31  double minPt_;
32  bool filter_;
34 };
35 
37  : minPt_(pset.getParameter<double>("minTrackSumPt")) {
38  src_ = pset.getParameter<edm::InputTag>("src");
39  token = consumes<reco::VertexCollection>(src_);
40  filter_ = pset.getParameter<bool>("filter");
41  produces<reco::VertexCollection>();
42 }
43 
46  evt.getByToken(token, vertices_);
47  auto output = std::make_unique<reco::VertexCollection>();
48  // If there is only one vertex, there are no PU vertices!
49  if (vertices_->size() > 1) {
50  // Copy over all the vertices that have associatd tracks with pt greater
51  // than the threshold. The predicate function is the VertexTrackPtSumFilter
52  // better name: copy_if_not
53  std::remove_copy_if(vertices_->begin() + 1, vertices_->end(), std::back_inserter(*output), [this](auto const& vtx) {
54  double trackPtSum = 0.;
55  for (reco::Vertex::trackRef_iterator track = vtx.tracks_begin(); track != vtx.tracks_end(); ++track) {
56  trackPtSum += (*track)->pt();
57  }
58  return trackPtSum > this->minPt_;
59  });
60  }
61  size_t nPUVtx = output->size();
62  evt.put(std::move(output));
63  // If 'filter' is enabled, return whether true if there are PU vertices
64  if (!filter_)
65  return true;
66  else
67  return nPUVtx;
68 }
69 
71  // recoTauPileUpVertexSelector
73 
74  desc.add<double>("minTrackSumPt");
75  desc.add<edm::InputTag>("src");
76  desc.add<bool>("filter");
77 
78  descriptions.add("recoTauPileUpVertexSelector", desc);
79 }
80 
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:536
edm::EDGetTokenT< reco::VertexCollection > token
RecoTauPileUpVertexSelector(const edm::ParameterSet &pset)
bool filter(edm::Event &evt, const edm::EventSetup &es) override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void add(std::string const &label, ParameterSetDescription const &psetDescription)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
def move(src, dest)
Definition: eostools.py:511
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector<TrackRef>
Definition: Vertex.h:38