CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
RecoTauBuilderPlugins.h
Go to the documentation of this file.
1 #ifndef RecoTauTag_RecoTau_RecoTauBuilderPlugin_h
2 #define RecoTauTag_RecoTau_RecoTauBuilderPlugin_h
3 
4 /*
5  * RecoTauBuilderPlugins
6  *
7  * Author: Evan K. Friis (UC Davis)
8  *
9  * Classes for building new and modifying existing PFTaus from PFJets and
10  * reconstructed PiZeros.
11  *
12  * RecoTauBuilderPlugin is the base class for any algorithm that constructs
13  * taus. Algorithms should override the abstract function
14  *
15  * std::vector<PFTau> operator()(const PFJet&, const
16  * std::vector<RecoTauPiZero>&) const;
17  *
18  * implementing it such that a list of taus a produced for a given jet and its
19  * associated collection of PiZeros.
20  *
21  * RecoTauModifierPlugin takes an input tau and modifies it.
22  *
23  * Both plugins inherit from RecoTauEventHolderPlugin, which provides the
24  * methods
25  *
26  * const edm::Event* evt() const; const edm::EventSetup* evtSetup()
27  *
28  * to retrieve the current event if necessary.
29  *
30  * $Id $
31  *
32  */
33 
34 #include <boost/ptr_container/ptr_vector.hpp>
35 
37 
43 
46 
47 #include <vector>
48 
49 namespace reco { namespace tau {
50 
51 /* Class that constructs PFTau(s) from a PFJet and its associated PiZeros */
53  public:
54  typedef boost::ptr_vector<reco::PFTau> output_type;
55  typedef std::auto_ptr<output_type> return_type;
56 
59  pfCandSrc_ = pset.getParameter<edm::InputTag>("pfCandSrc");
60  pvSrc_ = pset.getParameter<edm::InputTag>("primaryVertexSrc");
61  };
62 
63  virtual ~RecoTauBuilderPlugin() {}
64 
69  const std::vector<reco::RecoTauPiZero>& piZeros,
70  const std::vector<PFCandidatePtr>& regionalExtras) const = 0;
71 
74  return pfCands_; };
75 
77  const reco::VertexRef& primaryVertex() const { return pv_; }
78 
79  // Hook called by base class at the beginning of each event. Used to update
80  // handle to PFCandidates
81  virtual void beginEvent();
82 
83  private:
86 
88  // Handle to PFCandidates needed to build Refs
90 
91 };
92 
93 /* Class that updates a PFTau's members (i.e. electron variables) */
95  public:
99  // Modify an existing PFTau (i.e. add electron rejection, etc)
100  virtual void operator()(PFTau&) const = 0;
101  virtual void beginEvent() {}
102 };
103 
104 /* Class that returns a double value indicating the quality of a given tau */
106  public:
108  RecoTauEventHolderPlugin(pset){};
110  // Modify an existing PFTau (i.e. add electron rejection, etc)
111  virtual double operator()(const PFTauRef&) const = 0;
112  virtual void beginEvent() {}
113 };
114 } } // end namespace reco::tau
115 
123 
124 #endif
T getParameter(std::string const &) const
const edm::Handle< PFCandidateCollection > & getPFCands() const
Hack to be able to convert Ptrs to Refs.
virtual double operator()(const PFTauRef &) const =0
edm::Handle< PFCandidateCollection > pfCands_
virtual void operator()(PFTau &) const =0
RecoTauCleanerPlugin(const edm::ParameterSet &pset)
boost::ptr_vector< reco::PFTau > output_type
std::auto_ptr< output_type > return_type
RecoTauModifierPlugin(const edm::ParameterSet &pset)
virtual return_type operator()(const reco::PFJetRef &jet, const std::vector< reco::RecoTauPiZero > &piZeros, const std::vector< PFCandidatePtr > &regionalExtras) const =0
RecoTauBuilderPlugin(const edm::ParameterSet &pset)
const reco::VertexRef & primaryVertex() const
Get primary vertex associated to this event.