CMS 3D CMS Logo

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 Candidates 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 Candidate particle type
12  * (like hadron, gamma, etc). When a Candidate 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 
29 
30 #include <functional>
31 
32 namespace reco {
33  namespace tau {
34 
36  public:
37  // Quality cut types
38  typedef std::function<bool(const TrackBaseRef&)> TrackQCutFunc;
39  typedef std::vector<TrackQCutFunc> TrackQCutFuncCollection;
40  typedef std::function<bool(const Candidate&)> CandQCutFunc;
41  typedef std::vector<CandQCutFunc> CandQCutFuncCollection;
42  typedef std::map<int, CandQCutFuncCollection> CandQCutFuncMap;
43 
44  explicit RecoTauQualityCuts(const edm::ParameterSet& qcuts);
45 
47  void setPV(const reco::VertexRef& vtx) { pv_ = vtx; }
48 
50  void setLeadTrack(const reco::Track& leadTrack);
52 
56 
58  bool filterTrack(const reco::TrackBaseRef& track) const;
59  bool filterTrack(const reco::TrackRef& track) const;
60  bool filterTrack(const reco::Track& track) const;
62  bool filterChargedCand(const reco::Candidate& cand) const;
63 
65  template <typename Coll>
66  Coll filterTracks(const Coll& coll, bool invert = false) const {
67  Coll output;
68  for (auto const& track : coll) {
69  if (filterTrack(track) ^ invert)
70  output.push_back(track);
71  }
72  return output;
73  }
74 
76  bool filterCand(const reco::Candidate& cand) const;
77 
79  template <typename CandRefType>
80  bool filterCandRef(const CandRefType& cand) const {
81  return filterCand(*cand);
82  }
83 
85  template <typename Coll>
86  Coll filterCandRefs(const Coll& refcoll, bool invert = false) const {
87  Coll output;
88  for (auto const& cand : refcoll) {
89  if (filterCandRef(cand) ^ invert)
90  output.push_back(cand);
91  }
92  return output;
93  }
94 
96  static void fillDescriptions(edm::ParameterSetDescription& descriptions);
97 
98  private:
99  bool filterTrack_(const reco::Track* track) const;
100  bool filterGammaCand(const reco::Candidate& cand) const;
101  bool filterNeutralHadronCand(const reco::Candidate& cand) const;
102  bool filterCandByType(const reco::Candidate& cand) const;
103 
104  // The current primary vertex
106  // The current lead track references
108 
109  double minTrackPt_;
114  double maxDeltaZ_;
117  double minGammaEt_;
120  bool checkPV_;
121  };
122 
123  // Split an input set of quality cuts into those that need to be inverted
124  // to select PU (the first member) and those that are general quality cuts.
125  std::pair<edm::ParameterSet, edm::ParameterSet> factorizePUQCuts(const edm::ParameterSet& inputSet);
126 
127  } // namespace tau
128 } // namespace reco
129 
130 #endif
std::map< int, CandQCutFuncCollection > CandQCutFuncMap
bool filterTrack_(const reco::Track *track) const
RecoTauQualityCuts(const edm::ParameterSet &qcuts)
static void fillDescriptions(edm::ParameterSetDescription &descriptions)
Declare all parameters read from python config file.
bool filterCandByType(const reco::Candidate &cand) const
bool filterGammaCand(const reco::Candidate &cand) const
bool filterCandRef(const CandRefType &cand) const
Filter a Candidate held by a smart pointer or Ref.
bool filterNeutralHadronCand(const reco::Candidate &cand) const
bool filterTrack(const reco::TrackBaseRef &track) const
Filter a single Track.
InputIterator leadCand(InputIterator begin, InputIterator end)
reco::CandidatePtr CandRefType
Coll filterTracks(const Coll &coll, bool invert=false) const
Filter a collection of Tracks.
std::vector< CandQCutFunc > CandQCutFuncCollection
std::pair< edm::ParameterSet, edm::ParameterSet > factorizePUQCuts(const edm::ParameterSet &inputSet)
std::vector< TrackQCutFunc > TrackQCutFuncCollection
bool filterChargedCand(const reco::Candidate &cand) const
or a single charged candidate
Coll filterCandRefs(const Coll &refcoll, bool invert=false) const
Filter a ref vector of Candidates.
fixed size matrix
bool filterCand(const reco::Candidate &cand) const
Filter a single Candidate.
std::function< bool(const TrackBaseRef &)> TrackQCutFunc
void setLeadTrack(const reco::Track &leadTrack)
Update the leading track.
void setPV(const reco::VertexRef &vtx)
Update the primary vertex.
std::function< bool(const Candidate &)> CandQCutFunc