CMS 3D CMS Logo

RecoTauConstructor.h
Go to the documentation of this file.
1 #ifndef RecoTauTag_RecoTau_RecoTauConstructor_h
2 #define RecoTauTag_RecoTau_RecoTauConstructor_h
3 
4 /*
5  * RecoTauConstructor
6  *
7  * A generalized builder of reco::PFTau objects. Takes a variety of
8  * different collections and converts them to the proper Ref format
9  * needed for PFTau storage. Automatically sets the p4, charge, and
10  * other properties correctly. Optionally, it can determine the
11  * lead track information, and copy the gamma candidates owned by the
12  * reconstructed pi zeros into the appropriate PiZero collection.
13  *
14  * If the gammas are copied from the PiZeros, the four vector will be
15  * built using the PiZeros + Charged Hadrons. This can be different than
16  * the Gammas + Charged Hadrons, as the PiZero may have a mass hypothesis set by
17  * the RecoTauPiZeroProducer
18  *
19  * Note that the p4 of the tau is *always* set as the sum of objects in
20  * signal cone.
21  *
22  * Author: Evan K. Friis, UC Davis
23  *
24  */
25 
34 
35 #include "boost/shared_ptr.hpp"
36 #include <vector>
37 
38 namespace reco { namespace tau {
39 
41  public:
42  enum Region {
45  };
46 
47  enum ParticleType {
52  };
53 
55  RecoTauConstructor(const PFJetRef& jetRef,
57  bool copyGammasFromPiZeros = false,
61 
62  /*
63  * Code to set leading candidates. These are just wrappers about
64  * the embedded taus methods, but with the ability to convert Ptrs
65  * to Refs.
66  */
67 
69  template<typename T> void setleadPFChargedHadrCand(const T& cand) {
70  tau_->setleadPFChargedHadrCand(convertToPtr(cand));
71  }
72 
74  template<typename T> void setleadPFNeutralCand(const T& cand) {
75  tau_->setleadPFNeutralCand(convertToPtr(cand));
76  }
77 
79  template<typename T> void setleadPFCand(const T& cand) {
80  tau_->setleadPFCand(convertToPtr(cand));
81  }
82 
84  void addPFCand(Region region, ParticleType type, const PFCandidateRef& ref, bool skipAddToP4 = false);
85  void addPFCand(Region region, ParticleType type, const PFCandidatePtr& ptr, bool skipAddToP4 = false);
86 
88  void reserve(Region region, ParticleType type, size_t size);
89 
90  // Add a collection of objects to a given collection
91  template<typename InputIterator>
92  void addPFCands(Region region, ParticleType type, const InputIterator& begin, const InputIterator& end)
93  {
94  for(InputIterator iter = begin; iter != end; ++iter) {
95  addPFCand(region, type, convertToPtr(*iter));
96  }
97  }
98 
100  void reserveTauChargedHadron(Region region, size_t size);
101 
104 
106  template<typename InputIterator> void addTauChargedHadrons(Region region, const InputIterator& begin, const InputIterator& end)
107  {
108  for ( InputIterator iter = begin; iter != end; ++iter ) {
109  addTauChargedHadron(region, *iter);
110  }
111  }
112 
114  void reservePiZero(Region region, size_t size);
115 
117  void addPiZero(Region region, const RecoTauPiZero& piZero);
118 
120  template<typename InputIterator> void addPiZeros(Region region, const InputIterator& begin, const InputIterator& end)
121  {
122  for ( InputIterator iter = begin; iter != end; ++iter ) {
123  addPiZero(region, *iter);
124  }
125  }
126 
127  // Build and return the associated tau
128  std::auto_ptr<reco::PFTau> get(bool setupLeadingCandidates=true);
129 
130  // Get the four vector of the signal objects added so far
131  const reco::Candidate::LorentzVector& p4() const { return p4_; }
132 
133  private:
134  typedef std::pair<Region, ParticleType> CollectionKey;
135  typedef std::map<CollectionKey, std::vector<PFCandidatePtr>*> CollectionMap;
136  typedef boost::shared_ptr<std::vector<PFCandidatePtr> > SortedListPtr;
137  typedef std::map<CollectionKey, SortedListPtr> SortedCollectionMap;
138 
140 
146 
147  // Retrieve collection associated to signal/iso and type
148  std::vector<PFCandidatePtr>* getCollection(Region region, ParticleType type);
149  SortedListPtr getSortedCollection(Region region, ParticleType type);
150 
151  // Sort all our collections by PT and copy them into the tau
152  void sortAndCopyIntoTau();
153 
154  // Helper functions for dealing with refs
155  PFCandidatePtr convertToPtr(const PFCandidatePtr& pfPtr) const;
156  PFCandidatePtr convertToPtr(const CandidatePtr& candPtr) const;
157  PFCandidatePtr convertToPtr(const PFCandidateRef& pfRef) const;
158 
160  std::auto_ptr<reco::PFTau> tau_;
161  CollectionMap collections_;
162 
163  // Keep sorted (by descending pt) collections
164  SortedCollectionMap sortedCollections_;
165 
166  // Keep track of the signal cone four vector in case we want it
168 };
169 } } // end reco::tau namespace
170 #endif
size
Write out results.
type
Definition: HCALResponse.h:21
void addPFCand(Region region, ParticleType type, const PFCandidateRef &ref, bool skipAddToP4=false)
Append a PFCandidateRef/Ptr to a given collection.
void addTauChargedHadron(Region region, const PFRecoTauChargedHadron &chargedHadron)
Add a ChargedHadron to the given collection.
reco::Candidate::LorentzVector p4_
void reserve(Region region, ParticleType type, size_t size)
Reserve a set amount of space for a given RefVector.
void addTauChargedHadrons(Region region, const InputIterator &begin, const InputIterator &end)
Add a list of charged hadrons to the input collection.
boost::shared_ptr< std::vector< PFCandidatePtr > > SortedListPtr
std::vector< PFCandidatePtr > * getCollection(Region region, ParticleType type)
std::map< CollectionKey, SortedListPtr > SortedCollectionMap
void addPiZeros(Region region, const InputIterator &begin, const InputIterator &end)
Add a list of pizeros to the input collection.
SortedCollectionMap sortedCollections_
RecoTauConstructor(const PFJetRef &jetRef, const edm::Handle< PFCandidateCollection > &pfCands, bool copyGammasFromPiZeros=false, const StringObjectFunction< reco::PFTau > *signalConeSize=0, double minAbsPhotonSumPt_insideSignalCone=2.5, double minRelPhotonSumPt_insideSignalCone=0., double minAbsPhotonSumPt_outsideSignalCone=1.e+9, double minRelPhotonSumPt_outsideSignalCone=1.e+9)
Constructor with PFCandidate Handle.
const edm::Handle< PFCandidateCollection > & pfCands_
SortedListPtr getSortedCollection(Region region, ParticleType type)
const reco::Candidate::LorentzVector & p4() const
#define end
Definition: vmac.h:39
std::pair< Region, ParticleType > CollectionKey
void setleadPFCand(const T &cand)
Set leading PF candidate.
void addPFCands(Region region, ParticleType type, const InputIterator &begin, const InputIterator &end)
PFCandidatePtr convertToPtr(const PFCandidatePtr &pfPtr) const
std::auto_ptr< reco::PFTau > tau_
void reserveTauChargedHadron(Region region, size_t size)
Reserve a set amount of space for the ChargedHadrons.
void reservePiZero(Region region, size_t size)
Reserve a set amount of space for the PiZeros.
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
void setleadPFChargedHadrCand(const T &cand)
Set leading PFChargedHadron candidate.
fixed size matrix
#define begin
Definition: vmac.h:32
std::map< CollectionKey, std::vector< PFCandidatePtr > * > CollectionMap
void setleadPFNeutralCand(const T &cand)
Set leading PFGamma candidate.
void addPiZero(Region region, const RecoTauPiZero &piZero)
Add a PiZero to the given collection.
const StringObjectFunction< reco::PFTau > * signalConeSize_
long double T