CMS 3D CMS Logo

PVObjectSelector.h
Go to the documentation of this file.
1 #ifndef Analysis_AnalysisFilters_interface_PVObjectSelector_h
2 #define Analysis_AnalysisFilters_interface_PVObjectSelector_h
3 
4 #ifndef __GCCXML__
6 #endif
9 
12 
13 #include <vector>
14 #include <string>
15 
16 // make a selector for this selection
17 class PVObjectSelector : public Selector<reco::Vertex> {
18 public:
19 
21 
22 #ifndef __GCCXML__
24  PVObjectSelector( params )
25  {}
26 #endif
27 
28  PVObjectSelector( edm::ParameterSet const & params ) {
29  push_back("PV NDOF", params.getParameter<double>("minNdof") );
30  push_back("PV Z", params.getParameter<double>("maxZ") );
31  push_back("PV RHO", params.getParameter<double>("maxRho") );
32  set("PV NDOF");
33  set("PV Z");
34  set("PV RHO");
35 
36  indexNDOF_ = index_type (&bits_, "PV NDOF");
37  indexZ_ = index_type (&bits_, "PV Z");
38  indexRho_ = index_type (&bits_, "PV RHO");
39 
40  if ( params.exists("cutsToIgnore") )
41  setIgnoredCuts( params.getParameter<std::vector<std::string> >("cutsToIgnore") );
42 
44  }
45 
46  bool operator() ( reco::Vertex const & pv, pat::strbitset & ret ) override {
47  if ( pv.isFake() ) return false;
48 
49  if ( pv.ndof() >= cut(indexNDOF_, double() )
50  || ignoreCut(indexNDOF_) ) {
51  passCut(ret, indexNDOF_ );
52  if ( fabs(pv.z()) <= cut(indexZ_, double())
53  || ignoreCut(indexZ_) ) {
54  passCut(ret, indexZ_ );
55  if ( fabs(pv.position().Rho()) <= cut(indexRho_, double() )
56  || ignoreCut(indexRho_) ) {
57  passCut( ret, indexRho_);
58  }
59  }
60  }
61 
62  setIgnored(ret);
63 
64  return (bool)ret;
65  }
66 
68 
69 private:
73 };
74 
75 #endif
T getParameter(std::string const &) const
bool operator()(reco::Vertex const &pv, pat::strbitset &ret) override
bool exists(std::string const &parameterName) const
checks if a parameter exists
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
const Point & position() const
position
Definition: Vertex.h:109
PVObjectSelector(edm::ParameterSet const &params, edm::ConsumesCollector &&iC)
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:6
double z() const
z coordinate
Definition: Vertex.h:115
virtual void push_back(std::string const &s)
This is the registration of an individual cut string.
Definition: Selector.h:46
PVObjectSelector(edm::ParameterSet const &params)
Functor that operates on <T>
Definition: Selector.h:24
double ndof() const
Definition: Vertex.h:105
bool isFake() const
Definition: Vertex.h:72
index_type indexRho_
index_type indexNDOF_
pat::strbitset getBitTemplate() const
Get an empty bitset with the proper names.
Definition: Selector.h:212
void setIgnoredCuts(std::vector< std::string > const &bitsToIgnore)
set the bits to ignore from a vector
Definition: Selector.h:167
int cut(index_type const &i, int val) const
Access the int cut values at index "s".
Definition: Selector.h:194