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 <vector>
36 
37 namespace reco {
38  namespace tau {
39 
41  public:
43 
45 
47  RecoTauConstructor(const JetBaseRef& jetRef,
48  const edm::Handle<edm::View<reco::Candidate> >& pfCands,
49  bool copyGammasFromPiZeros = false,
55 
56  /*
57  * Code to set leading candidates. These are just wrappers about
58  * the embedded taus methods, but with the ability to convert Ptrs
59  * to Refs.
60  */
61 
63  template <typename T>
64  void setleadChargedHadrCand(const T& cand) {
65  tau_->setleadChargedHadrCand(convertToPtr(cand));
66  }
67 
69  template <typename T>
70  void setleadNeutralCand(const T& cand) {
71  tau_->setleadNeutralCand(convertToPtr(cand));
72  }
73 
75  template <typename T>
76  void setleadCand(const T& cand) {
77  tau_->setleadCand(convertToPtr(cand));
78  }
79 
81  void addPFCand(Region region, ParticleType type, const CandidatePtr& 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  for (InputIterator iter = begin; iter != end; ++iter) {
91  }
92  }
93 
96 
99 
101  template <typename InputIterator>
102  void addTauChargedHadrons(Region region, const InputIterator& begin, const InputIterator& end) {
103  for (InputIterator iter = begin; iter != end; ++iter) {
104  addTauChargedHadron(region, *iter);
105  }
106  }
107 
109  void reservePiZero(Region region, size_t size);
110 
112  void addPiZero(Region region, const RecoTauPiZero& piZero);
113 
115  template <typename InputIterator>
116  void addPiZeros(Region region, const InputIterator& begin, const InputIterator& end) {
117  for (InputIterator iter = begin; iter != end; ++iter) {
118  addPiZero(region, *iter);
119  }
120  }
121 
122  // Build and return the associated tau
123  std::unique_ptr<reco::PFTau> get(bool setupLeadingCandidates = true);
124 
125  // Get the four vector of the signal objects added so far
126  const reco::Candidate::LorentzVector& p4() const { return p4_; }
127 
128  private:
129  typedef std::pair<Region, ParticleType> CollectionKey;
130  typedef std::map<CollectionKey, std::vector<CandidatePtr>*> CollectionMap;
131  typedef std::shared_ptr<std::vector<CandidatePtr> > SortedListPtr;
132  typedef std::map<CollectionKey, SortedListPtr> SortedCollectionMap;
133 
135 
141 
142  // Retrieve collection associated to signal/iso and type
143  std::vector<CandidatePtr>* getCollection(Region region, ParticleType type);
145 
146  // Sort all our collections by PT and copy them into the tau
147  void sortAndCopyIntoTau();
148 
149  // Helper functions for dealing with refs
150  CandidatePtr convertToPtr(const PFCandidatePtr& pfPtr) const;
151  CandidatePtr convertToPtr(const CandidatePtr& candPtr) const;
152 
154  std::unique_ptr<reco::PFTau> tau_;
156 
157  // Keep sorted (by descending pt) collections
159 
160  // Keep track of the signal cone four vector in case we want it
162  };
163  } // namespace tau
164 } // namespace reco
165 #endif
void addPFCand(Region region, ParticleType type, const CandidatePtr &ptr, bool skipAddToP4=false)
Append a PFCandidateRef/Ptr to a given collection.
size
Write out results.
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.
void setleadChargedHadrCand(const T &cand)
Set leading PFChargedHadron candidate.
std::map< CollectionKey, SortedListPtr > SortedCollectionMap
const reco::Candidate::LorentzVector & p4() const
std::vector< CandidatePtr > * getCollection(Region region, ParticleType type)
void setleadNeutralCand(const T &cand)
Set leading PFGamma candidate.
std::unique_ptr< reco::PFTau > tau_
void addPiZeros(Region region, const InputIterator &begin, const InputIterator &end)
Add a list of pizeros to the input collection.
const edm::Handle< edm::View< reco::Candidate > > & pfCands_
SortedCollectionMap sortedCollections_
SortedListPtr getSortedCollection(Region region, ParticleType type)
std::pair< Region, ParticleType > CollectionKey
void addPFCands(Region region, ParticleType type, const InputIterator &begin, const InputIterator &end)
RecoTauConstructor(const JetBaseRef &jetRef, const edm::Handle< edm::View< reco::Candidate > > &pfCands, bool copyGammasFromPiZeros=false, const StringObjectFunction< reco::PFTau > *signalConeSize=nullptr, 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.
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:36
std::shared_ptr< std::vector< CandidatePtr > > SortedListPtr
CandidatePtr convertToPtr(const PFCandidatePtr &pfPtr) const
fixed size matrix
void setleadCand(const T &cand)
Set leading PF candidate.
void addPiZero(Region region, const RecoTauPiZero &piZero)
Add a PiZero to the given collection.
const StringObjectFunction< reco::PFTau > * signalConeSize_
std::map< CollectionKey, std::vector< CandidatePtr > * > CollectionMap
long double T