CMS 3D CMS Logo

PVSelector.h
Go to the documentation of this file.
1 #ifndef Analysis_AnalysisFilters_interface_PVSelector_h
2 #define Analysis_AnalysisFilters_interface_PVSelector_h
3 
5 #ifndef __GCCXML__
7 #endif
9 
12 
13 // make a selector for this selection
14 class PVSelector : public Selector<edm::EventBase> {
15 public:
17 
18  PVSelector(edm::ParameterSet const& params) : pvSrc_(params.getParameter<edm::InputTag>("pvSrc")), pvSel_(params) {
19  push_back("NPV", params.getParameter<int>("NPV"));
20  set("NPV");
22  indexNPV_ = index_type(&bits_, "NPV");
23  }
24 
25 #ifndef __GCCXML__
27  pvSrcToken_ = iC.consumes<std::vector<reco::Vertex> >(pvSrc_);
28  }
29 #endif
30 
31  bool operator()(edm::EventBase const& event, pat::strbitset& ret) override {
32  ret.set(false);
33  event.getByLabel(pvSrc_, h_primVtx);
34 
35  // check if there is a good primary vertex
36 
37  if (h_primVtx->empty())
38  return false;
39 
40  // Loop over PV's and count those that pass
41  int npv = 0;
42  int _ntotal = 0;
43  mvSelPvs.clear();
44  for (std::vector<reco::Vertex>::const_iterator ibegin = h_primVtx->begin(), iend = h_primVtx->end(), i = ibegin;
45  i != iend;
46  ++i) {
47  reco::Vertex const& pv = *i;
48  bool ipass = pvSel_(pv);
49  if (ipass) {
50  ++npv;
51  mvSelPvs.push_back(edm::Ptr<reco::Vertex>(h_primVtx, _ntotal));
52  }
53  ++_ntotal;
54  }
55 
56  // cache npv
57  mNpv = npv;
58 
59  // Set the strbitset
60  if (npv >= cut(indexNPV_, int()) || ignoreCut(indexNPV_)) {
62  }
63 
64  // Check if there is anything to ignore
65  setIgnored(ret);
66 
67  // Return status
68  bool pass = (bool)ret;
69  return pass;
70  }
71 
73 
75 
76  // get NPV from the last check
77  int GetNpv(void) { return mNpv; }
78 
79  std::vector<edm::Ptr<reco::Vertex> > const& GetSelectedPvs() const { return mvSelPvs; }
80 
81 private:
83 #ifndef __GCCXML__
85 #endif
88  std::vector<edm::Ptr<reco::Vertex> > mvSelPvs; // selected vertices
90  int mNpv; // cache number of PVs
91 };
92 
93 #endif
bool ignoreCut(std::string const &s) const
ignore the cut at index "s"
Definition: Selector.h:127
std::vector< edm::Ptr< reco::Vertex > > mvSelPvs
Definition: PVSelector.h:88
index_type indexNPV_
Definition: PVSelector.h:89
ret
prodAgent to be discontinued
pat::strbitset::index_type index_type
Definition: Selector.h:25
void setIgnored(pat::strbitset &ret)
set ignored bits
Definition: Selector.h:181
pat::strbitset retInternal_
internal ret if users don&#39;t care about return bits
Definition: Selector.h:242
PVSelector(edm::ParameterSet const &params, edm::ConsumesCollector &&iC)
Definition: PVSelector.h:26
pat::strbitset bits_
the bitset indexed by strings
Definition: Selector.h:241
void passCut(pat::strbitset &ret, std::string const &s)
Passing cuts.
Definition: Selector.h:142
def pv(vc)
Definition: MetAnalyzer.py:7
virtual void push_back(std::string const &s)
This is the registration of an individual cut string.
Definition: Selector.h:42
Functor that operates on <T>
Definition: Selector.h:22
pat::strbitset getBitTemplate() const
Get an empty bitset with the proper names.
Definition: Selector.h:168
edm::InputTag pvSrc_
Definition: PVSelector.h:82
int cut(index_type const &i, int val) const
Access the int cut values at index "s".
Definition: Selector.h:158
std::vector< edm::Ptr< reco::Vertex > > const & GetSelectedPvs() const
Definition: PVSelector.h:79
edm::Handle< std::vector< reco::Vertex > > const & vertices() const
Definition: PVSelector.h:74
PVObjectSelector pvSel_
Definition: PVSelector.h:86
edm::EDGetTokenT< std::vector< reco::Vertex > > pvSrcToken_
Definition: PVSelector.h:84
HLT enums.
int GetNpv(void)
Definition: PVSelector.h:77
bool operator()(edm::EventBase const &event, pat::strbitset &ret) override
This provides the interface for base classes to select objects.
Definition: PVSelector.h:31
edm::Handle< std::vector< reco::Vertex > > h_primVtx
Definition: PVSelector.h:87
PVSelector(edm::ParameterSet const &params)
Definition: PVSelector.h:18
Definition: event.py:1