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 
31 
32 namespace reco { namespace tau {
33 
35 {
36  public:
37  // Quality cut types
38  typedef boost::function<bool (const TrackBaseRef&)> TrackQCutFunc;
39  typedef std::vector<TrackQCutFunc> TrackQCutFuncCollection;
40  typedef boost::function<bool (const PFCandidate&)> CandQCutFunc;
41  typedef std::vector<CandQCutFunc> CandQCutFuncCollection;
42  typedef std::map<PFCandidate::ParticleType, CandQCutFuncCollection> CandQCutFuncMap;
43 
44  explicit RecoTauQualityCuts(const edm::ParameterSet& qcuts);
45 
47  void setPV(const reco::VertexRef& vtx) const { pv_ = vtx; }
48 
50  void setLeadTrack(const reco::TrackRef& leadTrack) const;
51  void setLeadTrack(const reco::PFCandidate& leadCand) const;
52 
55  void setLeadTrack(const reco::PFCandidateRef& leadCand) const;
56 
58  bool filterTrack(const reco::TrackBaseRef& track) const;
59  bool filterTrack(const reco::TrackRef& track) const;
60 
62  template<typename Coll>
63  Coll filterTracks(const Coll& coll, bool invert = false) const
64  {
65  Coll output;
66  BOOST_FOREACH( const typename Coll::value_type track, coll ) {
67  if ( filterTrack(track)^invert ) output.push_back(track);
68  }
69  return output;
70  }
71 
73  bool filterCand(const reco::PFCandidate& cand) const;
74 
76  template<typename PFCandRefType>
77  bool filterCandRef(const PFCandRefType& cand) const { return filterCand(*cand); }
78 
80  template<typename Coll>
81  Coll filterCandRefs(const Coll& refcoll, bool invert = false) const
82  {
83  Coll output;
84  BOOST_FOREACH( const typename Coll::value_type cand, refcoll ) {
85  if ( filterCandRef(cand)^invert ) output.push_back(cand);
86  }
87  return output;
88  }
89 
90  private:
91  template <typename T> bool filterTrack_(const T& trackRef) const;
92  bool filterGammaCand(const reco::PFCandidate& cand) const;
93  bool filterNeutralHadronCand(const reco::PFCandidate& cand) const;
94  bool filterCandByType(const reco::PFCandidate& cand) const;
95 
96  // The current primary vertex
98  // The current lead track references
100 
101  double minTrackPt_;
106  double maxDeltaZ_;
109  double minGammaEt_;
112  bool checkPV_;
113 };
114 
115 // Split an input set of quality cuts into those that need to be inverted
116 // to select PU (the first member) and those that are general quality cuts.
117 std::pair<edm::ParameterSet, edm::ParameterSet> factorizePUQCuts(const edm::ParameterSet& inputSet);
118 
119 }} // end reco::tau:: namespace
120 
121 #endif
bool filterCandByType(const reco::PFCandidate &cand) const
RecoTauQualityCuts(const edm::ParameterSet &qcuts)
Coll filterCandRefs(const Coll &refcoll, bool invert=false) const
Filter a ref vector of PFCandidates.
bool filterTrack_(const T &trackRef) const
bool filterNeutralHadronCand(const reco::PFCandidate &cand) const
std::vector< CandQCutFunc > CandQCutFuncCollection
bool filterGammaCand(const reco::PFCandidate &cand) const
void setPV(const reco::VertexRef &vtx) const
Update the primary vertex.
std::map< PFCandidate::ParticleType, CandQCutFuncCollection > CandQCutFuncMap
std::pair< edm::ParameterSet, edm::ParameterSet > factorizePUQCuts(const edm::ParameterSet &inputSet)
std::vector< TrackQCutFunc > TrackQCutFuncCollection
Container::value_type value_type
boost::function< bool(const TrackBaseRef &)> TrackQCutFunc
boost::function< bool(const PFCandidate &)> CandQCutFunc
bool filterCandRef(const PFCandRefType &cand) const
Filter a PFCandidate held by a smart pointer or Ref.
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:38
void setLeadTrack(const reco::TrackRef &leadTrack) const
Update the leading track.
bool filterCand(const reco::PFCandidate &cand) const
Filter a single PFCandidate.
Coll filterTracks(const Coll &coll, bool invert=false) const
Filter a collection of Tracks.
long double T
bool filterTrack(const reco::TrackBaseRef &track) const
Filter a single Track.