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 
14 // make a selector for this selection
15 class PVSelector : public Selector<edm::EventBase> {
16 public:
17 
19 
20  PVSelector( edm::ParameterSet const & params ) :
21  pvSrc_ (params.getParameter<edm::InputTag>("pvSrc") ),
22  pvSel_ (params)
23  {
24  push_back("NPV", params.getParameter<int>("NPV") );
25  set("NPV");
27  indexNPV_ = index_type(&bits_, "NPV");
28  }
29 
30 #ifndef __GCCXML__
32  PVSelector(params)
33  {
34  pvSrcToken_ = iC.consumes<std::vector<reco::Vertex> >(pvSrc_);
35  }
36 #endif
37 
38  bool operator() ( edm::EventBase const & event, pat::strbitset & ret ) override {
39  ret.set(false);
40  event.getByLabel(pvSrc_, h_primVtx);
41 
42  // check if there is a good primary vertex
43 
44  if ( h_primVtx->empty() ) return false;
45 
46  // Loop over PV's and count those that pass
47  int npv = 0;
48  int _ntotal = 0;
49  mvSelPvs.clear();
50  for ( std::vector<reco::Vertex>::const_iterator ibegin = h_primVtx->begin(),
51  iend = h_primVtx->end(), i = ibegin; i != iend; ++i ) {
52  reco::Vertex const & pv = *i;
53  bool ipass = pvSel_(pv);
54  if ( ipass ) {
55  ++npv;
56  mvSelPvs.push_back(edm::Ptr<reco::Vertex>(h_primVtx,_ntotal));
57  }
58  ++_ntotal;
59  }
60 
61  // cache npv
62  mNpv = npv;
63 
64  // Set the strbitset
65  if ( npv >= cut(indexNPV_, int() ) || ignoreCut(indexNPV_) ) {
66  passCut(ret, indexNPV_);
67  }
68 
69  // Check if there is anything to ignore
70  setIgnored(ret);
71 
72  // Return status
73  bool pass = (bool)ret;
74  return pass;
75  }
76 
78 
80 
81 
82 
83  // get NPV from the last check
84  int GetNpv(void){return mNpv;}
85 
86 
87 
88  std::vector<edm::Ptr<reco::Vertex> > const &
89  GetSelectedPvs() const { return mvSelPvs; }
90 
91 
92 
93 private:
95 #ifndef __GCCXML__
97 #endif
100  std::vector<edm::Ptr<reco::Vertex> > mvSelPvs; // selected vertices
102  int mNpv; // cache number of PVs
103 };
104 
105 #endif
std::vector< edm::Ptr< reco::Vertex > > mvSelPvs
Definition: PVSelector.h:100
T getParameter(std::string const &) const
index_type indexNPV_
Definition: PVSelector.h:101
edm::Handle< std::vector< reco::Vertex > > const & vertices() const
Definition: PVSelector.h:79
pat::strbitset::index_type index_type
Definition: Selector.h:29
void setIgnored(pat::strbitset &ret)
set ignored bits
Definition: Selector.h:224
pat::strbitset retInternal_
internal ret if users don&#39;t care about return bits
Definition: Selector.h:287
PVSelector(edm::ParameterSet const &params, edm::ConsumesCollector &&iC)
Definition: PVSelector.h:31
pat::strbitset bits_
the bitset indexed by strings
Definition: Selector.h:286
void passCut(pat::strbitset &ret, std::string const &s)
Passing cuts.
Definition: Selector.h:176
bool ignoreCut(std::string const &s) const
ignore the cut at index "s"
Definition: Selector.h:159
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:46
Functor that operates on <T>
Definition: Selector.h:24
edm::InputTag pvSrc_
Definition: PVSelector.h:94
strbitset & set(bool val=true)
set method of all bits
Definition: strbitset.h:144
std::vector< edm::Ptr< reco::Vertex > > const & GetSelectedPvs() const
Definition: PVSelector.h:89
pat::strbitset getBitTemplate() const
Get an empty bitset with the proper names.
Definition: Selector.h:212
PVObjectSelector pvSel_
Definition: PVSelector.h:98
edm::EDGetTokenT< std::vector< reco::Vertex > > pvSrcToken_
Definition: PVSelector.h:96
HLT enums.
int GetNpv(void)
Definition: PVSelector.h:84
bool operator()(edm::EventBase const &event, pat::strbitset &ret) override
Definition: PVSelector.h:38
edm::Handle< std::vector< reco::Vertex > > h_primVtx
Definition: PVSelector.h:99
PVSelector(edm::ParameterSet const &params)
Definition: PVSelector.h:20
Definition: event.py:1
int cut(index_type const &i, int val) const
Access the int cut values at index "s".
Definition: Selector.h:194