CMS 3D CMS Logo

All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 
14 
17 
20 
21 #include <algorithm>
22 
24 public:
27  bool filter(edm::Event& evt, const edm::EventSetup& es) override;
28  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
29 
30 private:
32  double minPt_;
33  bool filter_;
35 };
36 
38  : minPt_(pset.getParameter<double>("minTrackSumPt")) {
39  src_ = pset.getParameter<edm::InputTag>("src");
40  token = consumes<reco::VertexCollection>(src_);
41  filter_ = pset.getParameter<bool>("filter");
42  produces<reco::VertexCollection>();
43 }
44 
47  evt.getByToken(token, vertices_);
48  auto output = std::make_unique<reco::VertexCollection>();
49  // If there is only one vertex, there are no PU vertices!
50  if (vertices_->size() > 1) {
51  // Copy over all the vertices that have associatd tracks with pt greater
52  // than the threshold. The predicate function is the VertexTrackPtSumFilter
53  // better name: copy_if_not
54  std::remove_copy_if(vertices_->begin() + 1, vertices_->end(), std::back_inserter(*output), [this](auto const& vtx) {
55  double trackPtSum = 0.;
56  for (reco::Vertex::trackRef_iterator track = vtx.tracks_begin(); track != vtx.tracks_end(); ++track) {
57  trackPtSum += (*track)->pt();
58  }
59  return trackPtSum > this->minPt_;
60  });
61  }
62  size_t nPUVtx = output->size();
63  evt.put(std::move(output));
64  // If 'filter' is enabled, return whether true if there are PU vertices
65  if (!filter_)
66  return true;
67  else
68  return nPUVtx;
69 }
70 
72  // recoTauPileUpVertexSelector
74 
75  desc.add<double>("minTrackSumPt");
76  desc.add<edm::InputTag>("src");
77  desc.add<bool>("filter");
78 
79  descriptions.add("recoTauPileUpVertexSelector", desc);
80 }
81 
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::EDGetTokenT< reco::VertexCollection > token
RecoTauPileUpVertexSelector(const edm::ParameterSet &pset)
bool filter(edm::Event &evt, const edm::EventSetup &es) override
ParameterDescriptionBase * add(U const &iLabel, T const &value)
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:37