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 
8 
15 
18 
21 
22 #include <algorithm>
23 
25  public:
28  bool filter(edm::Event& evt, const edm::EventSetup& es) override;
29  static void fillDescriptions(edm::ConfigurationDescriptions & descriptions);
30  private:
32  double minPt_;
33  bool filter_;
35 };
36 
39  pset.getParameter<double>("minTrackSumPt")) {
40  src_ = pset.getParameter<edm::InputTag>("src");
41  token = consumes<reco::VertexCollection>(src_);
42  filter_ = pset.getParameter<bool>("filter");
43  produces<reco::VertexCollection>();
44 }
45 
46 
48  edm::Event& evt, const edm::EventSetup& es) {
50  evt.getByToken(token, vertices_);
51  auto output = std::make_unique<reco::VertexCollection>();
52  // If there is only one vertex, there are no PU vertices!
53  if (vertices_->size() > 1) {
54  // Copy over all the vertices that have associatd tracks with pt greater
55  // than the threshold. The predicate function is the VertexTrackPtSumFilter
56  // better name: copy_if_not
57  std::remove_copy_if(vertices_->begin()+1, vertices_->end(),
58  std::back_inserter(*output), [this](auto const& vtx){
59  double trackPtSum = 0.;
60  for ( reco::Vertex::trackRef_iterator track = vtx.tracks_begin();
61  track != vtx.tracks_end(); ++track ) {
62  trackPtSum += (*track)->pt();
63  }
64  return trackPtSum > this->minPt_;
65  });
66  }
67  size_t nPUVtx = output->size();
68  evt.put(std::move(output));
69  // If 'filter' is enabled, return whether true if there are PU vertices
70  if (!filter_)
71  return true;
72  else
73  return nPUVtx;
74 }
75 
76 void
78  // recoTauPileUpVertexSelector
80 
81  desc.add<double>("minTrackSumPt");
82  desc.add<edm::InputTag>("src");
83  desc.add<bool>("filter");
84 
85  descriptions.add("recoTauPileUpVertexSelector", desc);
86 }
87 
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
#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)
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector<TrackRef>
Definition: Vertex.h:37
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
def move(src, dest)
Definition: eostools.py:511