CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 
33 
34 #include "boost/shared_ptr.hpp"
35 #include <vector>
36 
37 namespace reco { namespace tau {
38 
40  public:
41  enum Region {
44  };
45 
46  enum ParticleType {
51  };
52 
54  RecoTauConstructor(const PFJetRef& jetRef,
56  bool copyGammasFromPiZeros=false);
57 
58  /*
59  * Code to set leading candidates. These are just wrappers about
60  * the embedded taus methods, but with the ability to convert Ptrs
61  * to Refs.
62  */
63 
65  template<typename T> void setleadPFChargedHadrCand(const T& cand) {
66  tau_->setleadPFChargedHadrCand(convertToPtr(cand));
67  }
68 
70  template<typename T> void setleadPFNeutralCand(const T& cand) {
71  tau_->setleadPFNeutralCand(convertToPtr(cand));
72  }
73 
75  template<typename T> void setleadPFCand(const T& cand) {
76  tau_->setleadPFCand(convertToPtr(cand));
77  }
78 
80  void addPFCand(Region region, ParticleType type, const PFCandidateRef& ref, bool skipAddToP4 = false);
81  void addPFCand(Region region, ParticleType type, const PFCandidatePtr& ptr, bool skipAddToP4 = false);
82 
84  void reserve(Region region, ParticleType type, size_t size);
85 
86  // Add a collection of objects to a given collection
87  template<typename InputIterator>
88  void addPFCands(Region region, ParticleType type, const InputIterator& begin, const InputIterator& end)
89  {
90  for(InputIterator iter = begin; iter != end; ++iter) {
91  addPFCand(region, type, convertToPtr(*iter));
92  }
93  }
94 
96  void reserveTauChargedHadron(Region region, size_t size);
97 
99  void addTauChargedHadron(Region region, const PFRecoTauChargedHadron& chargedHadron);
100 
102  template<typename InputIterator> void addTauChargedHadrons(Region region, const InputIterator& begin, const InputIterator& end)
103  {
104  for ( InputIterator iter = begin; iter != end; ++iter ) {
105  addTauChargedHadron(region, *iter);
106  }
107  }
108 
110  void reservePiZero(Region region, size_t size);
111 
113  void addPiZero(Region region, const RecoTauPiZero& piZero);
114 
116  template<typename InputIterator> void addPiZeros(Region region, const InputIterator& begin, const InputIterator& end)
117  {
118  for ( InputIterator iter = begin; iter != end; ++iter ) {
119  addPiZero(region, *iter);
120  }
121  }
122 
123  // Build and return the associated tau
124  std::auto_ptr<reco::PFTau> get(bool setupLeadingCandidates=true);
125 
126  // Get the four vector of the signal objects added so far
127  const reco::Candidate::LorentzVector& p4() const { return p4_; }
128 
129  private:
130  typedef std::pair<Region, ParticleType> CollectionKey;
131  typedef std::map<CollectionKey, std::vector<PFCandidatePtr>*> CollectionMap;
132  typedef boost::shared_ptr<std::vector<PFCandidatePtr> > SortedListPtr;
133  typedef std::map<CollectionKey, SortedListPtr> SortedCollectionMap;
134 
136  // Retrieve collection associated to signal/iso and type
137  std::vector<PFCandidatePtr>* getCollection(Region region, ParticleType type);
139 
140  // Sort all our collections by PT and copy them into the tau
141  void sortAndCopyIntoTau();
142 
143  // Helper functions for dealing with refs
144  PFCandidatePtr convertToPtr(const PFCandidatePtr& pfPtr) const;
145  PFCandidatePtr convertToPtr(const CandidatePtr& candPtr) const;
146  PFCandidatePtr convertToPtr(const PFCandidateRef& pfRef) const;
147 
149  std::auto_ptr<reco::PFTau> tau_;
151 
152  // Keep sorted (by descending pt) collections
154 
155  // Keep track of the signal cone four vector in case we want it
157 };
158 } } // end reco::tau namespace
159 #endif
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_
const edm::Handle< PFCandidateCollection > & pfCands_
SortedListPtr getSortedCollection(Region region, ParticleType type)
ParticleType
const reco::Candidate::LorentzVector & p4() const
#define end
Definition: vmac.h:37
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.
#define begin
Definition: vmac.h:30
RecoTauConstructor(const PFJetRef &jetRef, const edm::Handle< PFCandidateCollection > &pfCands, bool copyGammasFromPiZeros=false)
Constructor with PFCandidate Handle.
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.
long double T
tuple size
Write out results.