Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:27:46

0001 #ifndef RecoTauTag_RecoTau_RecoTauQualityCuts_h
0002 #define RecoTauTag_RecoTau_RecoTauQualityCuts_h
0003 
0004 /*
0005  * RecoTauQualityCuts
0006  *
0007  * Author: Evan K. Friis
0008  *
0009  * Constructs a number of independent requirements on Candidates by building
0010  * binary predicate functions.  These are held in a number of lists of
0011  * functions.  Each of these lists is mapped to a Candidate particle type
0012  * (like hadron, gamma, etc).  When a Candidate is passed to filter(),
0013  * the correct list is looked up, and the result is the AND of all the predicate
0014  * functions.  See the .cc files for the QCut functions.
0015  *
0016  * Note that for some QCuts, the primary vertex must be updated every event.
0017  * Others require the lead track be defined for each tau before filter(..)
0018  * is called.
0019  *
0020  */
0021 
0022 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0023 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0024 #include "DataFormats/Candidate/interface/Candidate.h"
0025 #include "DataFormats/Candidate/interface/CandidateFwd.h"
0026 #include "DataFormats/TrackReco/interface/Track.h"
0027 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0028 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0029 
0030 #include <functional>
0031 
0032 namespace reco {
0033   namespace tau {
0034 
0035     class RecoTauQualityCuts {
0036     public:
0037       // Quality cut types
0038       typedef std::function<bool(const TrackBaseRef&)> TrackQCutFunc;
0039       typedef std::vector<TrackQCutFunc> TrackQCutFuncCollection;
0040       typedef std::function<bool(const Candidate&)> CandQCutFunc;
0041       typedef std::vector<CandQCutFunc> CandQCutFuncCollection;
0042       typedef std::map<int, CandQCutFuncCollection> CandQCutFuncMap;
0043 
0044       explicit RecoTauQualityCuts(const edm::ParameterSet& qcuts);
0045 
0046       /// Update the primary vertex
0047       void setPV(const reco::VertexRef& vtx) { pv_ = vtx; }
0048 
0049       /// Update the leading track
0050       void setLeadTrack(const reco::Track& leadTrack);
0051       void setLeadTrack(const reco::Candidate& leadCand);
0052 
0053       /// Update the leading track (using reference)
0054       /// If null, this will set the lead track ref null.
0055       void setLeadTrack(const reco::CandidateRef& leadCand);
0056 
0057       /// Filter a single Track
0058       bool filterTrack(const reco::TrackBaseRef& track) const;
0059       bool filterTrack(const reco::TrackRef& track) const;
0060       bool filterTrack(const reco::Track& track) const;
0061       /// or a single charged candidate
0062       bool filterChargedCand(const reco::Candidate& cand) const;
0063 
0064       /// Filter a collection of Tracks
0065       template <typename Coll>
0066       Coll filterTracks(const Coll& coll, bool invert = false) const {
0067         Coll output;
0068         for (auto const& track : coll) {
0069           if (filterTrack(track) ^ invert)
0070             output.push_back(track);
0071         }
0072         return output;
0073       }
0074 
0075       /// Filter a single Candidate
0076       bool filterCand(const reco::Candidate& cand) const;
0077 
0078       /// Filter a Candidate held by a smart pointer or Ref
0079       template <typename CandRefType>
0080       bool filterCandRef(const CandRefType& cand) const {
0081         return filterCand(*cand);
0082       }
0083 
0084       /// Filter a ref vector of Candidates
0085       template <typename Coll>
0086       Coll filterCandRefs(const Coll& refcoll, bool invert = false) const {
0087         Coll output;
0088         for (auto const& cand : refcoll) {
0089           if (filterCandRef(cand) ^ invert)
0090             output.push_back(cand);
0091         }
0092         return output;
0093       }
0094 
0095       /// Declare all parameters read from python config file
0096       static void fillDescriptions(edm::ParameterSetDescription& descriptions);
0097 
0098     private:
0099       bool filterTrack_(const reco::Track* track) const;
0100       bool filterGammaCand(const reco::Candidate& cand) const;
0101       bool filterNeutralHadronCand(const reco::Candidate& cand) const;
0102       bool filterCandByType(const reco::Candidate& cand) const;
0103 
0104       // The current primary vertex
0105       reco::VertexRef pv_;
0106       // The current lead track references
0107       const reco::Track* leadTrack_;
0108 
0109       double minTrackPt_;
0110       double maxTrackChi2_;
0111       int minTrackPixelHits_;
0112       int minTrackHits_;
0113       double maxTransverseImpactParameter_;
0114       double maxDeltaZ_;
0115       double maxDeltaZToLeadTrack_;
0116       double minTrackVertexWeight_;
0117       double minGammaEt_;
0118       double minNeutralHadronEt_;
0119       bool checkHitPattern_;
0120       bool checkPV_;
0121     };
0122 
0123     // Split an input set of quality cuts into those that need to be inverted
0124     // to select PU (the first member) and those that are general quality cuts.
0125     std::pair<edm::ParameterSet, edm::ParameterSet> factorizePUQCuts(const edm::ParameterSet& inputSet);
0126 
0127   }  // namespace tau
0128 }  // namespace reco
0129 
0130 #endif