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:
20 
21 #ifndef __GCCXML__
23 #endif
24 
26  push_back("PV NDOF", params.getParameter<double>("minNdof"));
27  push_back("PV Z", params.getParameter<double>("maxZ"));
28  push_back("PV RHO", params.getParameter<double>("maxRho"));
29  set("PV NDOF");
30  set("PV Z");
31  set("PV RHO");
32 
33  indexNDOF_ = index_type(&bits_, "PV NDOF");
34  indexZ_ = index_type(&bits_, "PV Z");
35  indexRho_ = index_type(&bits_, "PV RHO");
36 
37  if (params.exists("cutsToIgnore"))
38  setIgnoredCuts(params.getParameter<std::vector<std::string> >("cutsToIgnore"));
39 
41  }
42 
43  bool operator()(reco::Vertex const& pv, pat::strbitset& ret) override {
44  if (pv.isFake())
45  return false;
46 
47  if (pv.ndof() >= cut(indexNDOF_, double()) || ignoreCut(indexNDOF_)) {
48  passCut(ret, indexNDOF_);
49  if (fabs(pv.z()) <= cut(indexZ_, double()) || ignoreCut(indexZ_)) {
50  passCut(ret, indexZ_);
51  if (fabs(pv.position().Rho()) <= cut(indexRho_, double()) || ignoreCut(indexRho_)) {
52  passCut(ret, indexRho_);
53  }
54  }
55  }
56 
57  setIgnored(ret);
58 
59  return (bool)ret;
60  }
61 
63 
64 private:
68 };
69 
70 #endif
T getParameter(std::string const &) const
ret
prodAgent to be discontinued
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: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
const Point & position() const
position
Definition: Vertex.h:113
PVObjectSelector(edm::ParameterSet const &params, edm::ConsumesCollector &&iC)
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
bool ignoreCut(std::string const &s) const
ignore the cut at index "s"
Definition: Selector.h:127
def pv(vc)
Definition: MetAnalyzer.py:7
double z() const
z coordinate
Definition: Vertex.h:119
virtual void push_back(std::string const &s)
This is the registration of an individual cut string.
Definition: Selector.h:42
PVObjectSelector(edm::ParameterSet const &params)
Functor that operates on <T>
Definition: Selector.h:22
double ndof() const
Definition: Vertex.h:109
bool isFake() const
Definition: Vertex.h:75
index_type indexRho_
index_type indexNDOF_
pat::strbitset getBitTemplate() const
Get an empty bitset with the proper names.
Definition: Selector.h:168
void setIgnoredCuts(std::vector< std::string > const &bitsToIgnore)
set the bits to ignore from a vector
Definition: Selector.h:131
int cut(index_type const &i, int val) const
Access the int cut values at index "s".
Definition: Selector.h:158