CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
RecoTauQualityCuts.h
Go to the documentation of this file.
1 #ifndef RecoTauTag_RecoTau_RecoTauQualityCuts_h
2 #define RecoTauTag_RecoTau_RecoTauQualityCuts_h
3 
4 /*
5  * RecoTauQualityCuts
6  *
7  * Author: Evan K. Friis
8  *
9  * Constructs a number of independent requirements on PFCandidates by building
10  * binary predicate functions. These are held in a number of lists of
11  * functions. Each of these lists is mapped to a PFCandidate particle type
12  * (like hadron, gamma, etc). When a PFCandidate is passed to filter(),
13  * the correct list is looked up, and the result is the AND of all the predicate
14  * functions. See the .cc files for the QCut functions.
15  *
16  * Note that for some QCuts, the primary vertex must be updated every event.
17  * Others require the lead track be defined for each tau before filter(..)
18  * is called.
19  *
20  */
21 
22 #include <boost/function.hpp>
23 #include <boost/foreach.hpp>
24 
30 
31 namespace reco { namespace tau {
32 
34  public:
35  // Quality cut types
36  typedef boost::function<bool (const PFCandidate&)> QCutFunc;
37  typedef std::vector<QCutFunc> QCutFuncCollection;
38  typedef std::map<PFCandidate::ParticleType, QCutFuncCollection> QCutFuncMap;
39 
40  explicit RecoTauQualityCuts(const edm::ParameterSet &qcuts);
41 
43  void setPV(const reco::VertexRef& vtx) const { pv_ = vtx; }
44 
46  void setLeadTrack(const reco::PFCandidate& leadCand) const;
47 
50  void setLeadTrack(const reco::PFCandidateRef& leadCand) const;
51 
53  const QCutFunc& predicate() const { return predicate_; }
54 
56  bool filter(const reco::PFCandidate& cand) const {
57  return predicate_(cand);
58  }
59 
61  template<typename PFCandRefType>
62  bool filterRef(const PFCandRefType& cand) const { return filter(*cand); }
63 
65  template<typename Coll> Coll filterRefs(
66  const Coll& refcoll, bool invert=false) const {
67  Coll output;
68  BOOST_FOREACH(const typename Coll::value_type cand, refcoll) {
69  if (filterRef(cand)^invert)
70  output.push_back(cand);
71  }
72  return output;
73  }
74 
75  private:
76  // The current primary vertex
78  // The current lead track references
80  // A mapping from particle type to a set of QCuts
82  // Our entire predicate function
84 };
85 
86 // Split an input set of quality cuts into those that need to be inverted
87 // to select PU (the first member) and those that are general quality cuts.
88 std::pair<edm::ParameterSet, edm::ParameterSet> factorizePUQCuts(
89  const edm::ParameterSet& inputSet);
90 
91 }} // end reco::tau:: namespace
92 #endif
Coll filterRefs(const Coll &refcoll, bool invert=false) const
Filter a ref vector of PFCandidates.
RecoTauQualityCuts(const edm::ParameterSet &qcuts)
boost::function< bool(const PFCandidate &)> QCutFunc
bool filterRef(const PFCandRefType &cand) const
Filter a PFCandidate held by a smart pointer or Ref.
bool filter(const reco::PFCandidate &cand) const
Filter a single PFCandidate.
void setLeadTrack(const reco::PFCandidate &leadCand) const
Update the leading track.
void setPV(const reco::VertexRef &vtx) const
Update the primary vertex.
const QCutFunc & predicate() const
Get the predicate used to filter.
std::map< PFCandidate::ParticleType, QCutFuncCollection > QCutFuncMap
std::pair< edm::ParameterSet, edm::ParameterSet > factorizePUQCuts(const edm::ParameterSet &inputSet)
Container::value_type value_type
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:33
std::vector< QCutFunc > QCutFuncCollection