CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_1/src/RecoTauTag/RecoTau/interface/RecoTauQualityCuts.h

Go to the documentation of this file.
00001 #ifndef RecoTauTag_RecoTau_RecoTauQualityCuts_h
00002 #define RecoTauTag_RecoTau_RecoTauQualityCuts_h
00003 
00004 /*
00005  * RecoTauQualityCuts
00006  *
00007  * Author: Evan K. Friis
00008  *
00009  * Constructs a number of independent requirements on PFCandidates by building
00010  * binary predicate functions.  These are held in a number of lists of
00011  * functions.  Each of these lists is mapped to a PFCandidate particle type
00012  * (like hadron, gamma, etc).  When a PFCandidate is passed to filter(),
00013  * the correct list is looked up, and the result is the AND of all the predicate
00014  * functions.  See the .cc files for the QCut functions.
00015  *
00016  * Note that for some QCuts, the primary vertex must be updated every event.
00017  * Others require the lead track be defined for each tau before filter(..)
00018  * is called.
00019  *
00020  */
00021 
00022 #include <boost/function.hpp>
00023 #include <boost/foreach.hpp>
00024 
00025 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00026 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
00027 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidateFwd.h"
00028 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00029 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00030 
00031 namespace reco { namespace tau {
00032 
00033 class RecoTauQualityCuts {
00034   public:
00035     // Quality cut types
00036     typedef boost::function<bool (const PFCandidate&)> QCutFunc;
00037     typedef std::vector<QCutFunc> QCutFuncCollection;
00038     typedef std::map<PFCandidate::ParticleType, QCutFuncCollection> QCutFuncMap;
00039 
00040     explicit RecoTauQualityCuts(const edm::ParameterSet &qcuts);
00041 
00043     void setPV(const reco::VertexRef& vtx) const { pv_ = vtx; }
00044 
00046     void setLeadTrack(const reco::PFCandidate& leadCand) const;
00047 
00050     void setLeadTrack(const reco::PFCandidateRef& leadCand) const;
00051 
00053     const QCutFunc& predicate() const { return predicate_; }
00054 
00056     bool filter(const reco::PFCandidate& cand) const {
00057       return predicate_(cand);
00058     }
00059 
00061     template<typename PFCandRefType>
00062     bool filterRef(const PFCandRefType& cand) const { return filter(*cand); }
00063 
00065     template<typename Coll> Coll filterRefs(
00066         const Coll& refcoll, bool invert=false) const {
00067       Coll output;
00068       BOOST_FOREACH(const typename Coll::value_type cand, refcoll) {
00069         if (filterRef(cand)^invert)
00070           output.push_back(cand);
00071       }
00072       return output;
00073     }
00074 
00075   private:
00076     // The current primary vertex
00077     mutable reco::VertexRef pv_;
00078     // The current lead track references
00079     mutable reco::TrackBaseRef leadTrack_;
00080     // A mapping from particle type to a set of QCuts
00081     QCutFuncMap qcuts_;
00082     // Our entire predicate function
00083     QCutFunc predicate_;
00084 };
00085 
00086 // Split an input set of quality cuts into those that need to be inverted
00087 // to select PU (the first member) and those that are general quality cuts.
00088 std::pair<edm::ParameterSet, edm::ParameterSet> factorizePUQCuts(
00089     const edm::ParameterSet& inputSet);
00090 
00091 }}  // end reco::tau:: namespace
00092 #endif