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 
28 
29 #include <functional>
30 
31 namespace reco {
32  namespace tau {
33 
35  public:
36  // Quality cut types
37  typedef std::function<bool(const TrackBaseRef&)> TrackQCutFunc;
38  typedef std::vector<TrackQCutFunc> TrackQCutFuncCollection;
39  typedef std::function<bool(const Candidate&)> CandQCutFunc;
40  typedef std::vector<CandQCutFunc> CandQCutFuncCollection;
41  typedef std::map<int, CandQCutFuncCollection> CandQCutFuncMap;
42 
43  explicit RecoTauQualityCuts(const edm::ParameterSet& qcuts);
44 
46  void setPV(const reco::VertexRef& vtx) { pv_ = vtx; }
47 
49  void setLeadTrack(const reco::Track& leadTrack);
51 
54  void setLeadTrack(const reco::CandidateRef& leadCand);
55 
57  bool filterTrack(const reco::TrackBaseRef& track) const;
58  bool filterTrack(const reco::TrackRef& track) const;
59  bool filterTrack(const reco::Track& track) const;
61  bool filterChargedCand(const reco::Candidate& cand) const;
62 
64  template <typename Coll>
65  Coll filterTracks(const Coll& coll, bool invert = false) const {
66  Coll output;
67  for (auto const& track : coll) {
68  if (filterTrack(track) ^ invert)
69  output.push_back(track);
70  }
71  return output;
72  }
73 
75  bool filterCand(const reco::Candidate& cand) const;
76 
78  template <typename CandRefType>
79  bool filterCandRef(const CandRefType& cand) const {
80  return filterCand(*cand);
81  }
82 
84  template <typename Coll>
85  Coll filterCandRefs(const Coll& refcoll, bool invert = false) const {
86  Coll output;
87  for (auto const& cand : refcoll) {
88  if (filterCandRef(cand) ^ invert)
89  output.push_back(cand);
90  }
91  return output;
92  }
93 
94  private:
95  bool filterTrack_(const reco::Track* track) const;
96  bool filterGammaCand(const reco::Candidate& cand) const;
97  bool filterNeutralHadronCand(const reco::Candidate& cand) const;
98  bool filterCandByType(const reco::Candidate& cand) const;
99 
100  // The current primary vertex
102  // The current lead track references
104 
105  double minTrackPt_;
110  double maxDeltaZ_;
113  double minGammaEt_;
116  bool checkPV_;
117  };
118 
119  // Split an input set of quality cuts into those that need to be inverted
120  // to select PU (the first member) and those that are general quality cuts.
121  std::pair<edm::ParameterSet, edm::ParameterSet> factorizePUQCuts(const edm::ParameterSet& inputSet);
122 
123  } // namespace tau
124 } // namespace reco
125 
126 #endif
std::map< int, CandQCutFuncCollection > CandQCutFuncMap
RecoTauQualityCuts(const edm::ParameterSet &qcuts)
Coll filterCandRefs(const Coll &refcoll, bool invert=false) const
Filter a ref vector of Candidates.
bool filterCandRef(const CandRefType &cand) const
Filter a Candidate held by a smart pointer or Ref.
bool filterCand(const reco::Candidate &cand) const
Filter a single Candidate.
InputIterator leadCand(InputIterator begin, InputIterator end)
reco::CandidatePtr CandRefType
std::vector< CandQCutFunc > CandQCutFuncCollection
bool filterCandByType(const reco::Candidate &cand) const
std::pair< edm::ParameterSet, edm::ParameterSet > factorizePUQCuts(const edm::ParameterSet &inputSet)
std::vector< TrackQCutFunc > TrackQCutFuncCollection
bool filterTrack_(const reco::Track *track) const
JetCorrectorParametersCollection coll
Definition: classes.h:10
bool filterChargedCand(const reco::Candidate &cand) const
or a single charged candidate
bool filterNeutralHadronCand(const reco::Candidate &cand) const
fixed size matrix
std::function< bool(const TrackBaseRef &)> TrackQCutFunc
bool filterGammaCand(const reco::Candidate &cand) const
void setLeadTrack(const reco::Track &leadTrack)
Update the leading track.
Coll filterTracks(const Coll &coll, bool invert=false) const
Filter a collection of Tracks.
void setPV(const reco::VertexRef &vtx)
Update the primary vertex.
bool filterTrack(const reco::TrackBaseRef &track) const
Filter a single Track.
std::function< bool(const Candidate &)> CandQCutFunc