CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
RecoTauConstructor.cc
Go to the documentation of this file.
2 
5 
7 
9 
10 #include <boost/foreach.hpp>
11 #include <boost/bind.hpp>
12 
13 namespace reco { namespace tau {
14 
16  bool copyGammasFromPiZeros,
17  const StringObjectFunction<reco::PFTau>* signalConeSize,
18  double minAbsPhotonSumPt_insideSignalCone, double minRelPhotonSumPt_insideSignalCone,
19  double minAbsPhotonSumPt_outsideSignalCone, double minRelPhotonSumPt_outsideSignalCone)
20  : signalConeSize_(signalConeSize),
21  minAbsPhotonSumPt_insideSignalCone_(minAbsPhotonSumPt_insideSignalCone),
22  minRelPhotonSumPt_insideSignalCone_(minRelPhotonSumPt_insideSignalCone),
23  minAbsPhotonSumPt_outsideSignalCone_(minAbsPhotonSumPt_outsideSignalCone),
24  minRelPhotonSumPt_outsideSignalCone_(minRelPhotonSumPt_outsideSignalCone),
25  pfCands_(pfCands)
26 {
27  // Initialize tau
28  tau_.reset(new PFTau());
29 
30  copyGammas_ = copyGammasFromPiZeros;
31  // Initialize our Accessors
32  collections_[std::make_pair(kSignal, kChargedHadron)] =
33  &tau_->selectedSignalPFChargedHadrCands_;
34  collections_[std::make_pair(kSignal, kGamma)] =
35  &tau_->selectedSignalPFGammaCands_;
36  collections_[std::make_pair(kSignal, kNeutralHadron)] =
37  &tau_->selectedSignalPFNeutrHadrCands_;
38  collections_[std::make_pair(kSignal, kAll)] =
39  &tau_->selectedSignalPFCands_;
40 
41  collections_[std::make_pair(kIsolation, kChargedHadron)] =
42  &tau_->selectedIsolationPFChargedHadrCands_;
43  collections_[std::make_pair(kIsolation, kGamma)] =
44  &tau_->selectedIsolationPFGammaCands_;
45  collections_[std::make_pair(kIsolation, kNeutralHadron)] =
46  &tau_->selectedIsolationPFNeutrHadrCands_;
47  collections_[std::make_pair(kIsolation, kAll)] =
48  &tau_->selectedIsolationPFCands_;
49 
50  // Build our temporary sorted collections, since you can't use stl sorts on
51  // RefVectors
52  BOOST_FOREACH(const CollectionMap::value_type& colkey, collections_) {
53  // Build an empty list for each collection
54  sortedCollections_[colkey.first] = SortedListPtr(
55  new SortedListPtr::element_type);
56  }
57 
58  tau_->setjetRef(jet);
59 }
60 
62  LogDebug("TauConstructorAddPFCand") << " region = " << region << ", type = " << type << ": Pt = " << ref->pt() << ", eta = " << ref->eta() << ", phi = " << ref->phi();
63  if ( region == kSignal ) {
64  // Keep track of the four vector of the signal vector products added so far.
65  // If a photon add it if we are not using PiZeros to build the gammas
66  if ( ((type != kGamma) || !copyGammas_) && !skipAddToP4 ) {
67  LogDebug("TauConstructorAddPFCand") << "--> adding PFCand to tauP4." ;
68  p4_ += ref->p4();
69  }
70  }
71  getSortedCollection(region, type)->push_back(edm::refToPtr<PFCandidateCollection>(ref));
72  // Add to global collection
73  getSortedCollection(region, kAll)->push_back(edm::refToPtr<PFCandidateCollection>(ref));
74 }
75 
77  LogDebug("TauConstructorAddPFCand") << " region = " << region << ", type = " << type << ": Pt = " << ptr->pt() << ", eta = " << ptr->eta() << ", phi = " << ptr->phi();
78  if ( region == kSignal ) {
79  // Keep track of the four vector of the signal vector products added so far.
80  // If a photon add it if we are not using PiZeros to build the gammas
81  if ( ((type != kGamma) || !copyGammas_) && !skipAddToP4 ) {
82  LogDebug("TauConstructorAddPFCand") << "--> adding PFCand to tauP4." ;
83  p4_ += ptr->p4();
84  }
85  }
86  getSortedCollection(region, type)->push_back(ptr);
87  // Add to global collection
88  getSortedCollection(region, kAll)->push_back(ptr);
89 }
90 
92  getSortedCollection(region, type)->reserve(size);
93  getCollection(region, type)->reserve(size);
94  // Reserve global collection as well
95  getSortedCollection(region, kAll)->reserve(
96  getSortedCollection(region, kAll)->size() + size);
97  getCollection(region, kAll)->reserve(
98  getCollection(region, kAll)->size() + size);
99 }
100 
102 {
103  if ( region == kSignal ) {
104  tau_->signalTauChargedHadronCandidates_.reserve(size);
105  tau_->selectedSignalPFChargedHadrCands_.reserve(size);
106  } else {
107  tau_->isolationTauChargedHadronCandidates_.reserve(size);
108  tau_->selectedIsolationPFChargedHadrCands_.reserve(size);
109  }
110 }
111 
112 namespace
113 {
114  void checkOverlap(const PFCandidatePtr& neutral, const std::vector<PFCandidatePtr>& pfGammas, bool& isUnique)
115  {
116  LogDebug("TauConstructorCheckOverlap") << " pfGammas: #entries = " << pfGammas.size();
117  for ( std::vector<PFCandidatePtr>::const_iterator pfGamma = pfGammas.begin();
118  pfGamma != pfGammas.end(); ++pfGamma ) {
119  LogDebug("TauConstructorCheckOverlap") << "pfGamma = " << pfGamma->id() << ":" << pfGamma->key();
120  if ( (*pfGamma) == neutral ) isUnique = false;
121  }
122  }
123 
124  void checkOverlap(const PFCandidatePtr& neutral, const std::vector<reco::RecoTauPiZero>& piZeros, bool& isUnique)
125  {
126  LogDebug("TauConstructorCheckOverlap") << " piZeros: #entries = " << piZeros.size();
127  for ( std::vector<reco::RecoTauPiZero>::const_iterator piZero = piZeros.begin();
128  piZero != piZeros.end(); ++piZero ) {
129  size_t numPFGammas = piZero->numberOfDaughters();
130  for ( size_t iPFGamma = 0; iPFGamma < numPFGammas; ++iPFGamma ) {
131  reco::CandidatePtr pfGamma = piZero->daughterPtr(iPFGamma);
132  LogDebug("TauConstructorCheckOverlap") << "pfGamma = " << pfGamma.id() << ":" << pfGamma.key();
133  if ( pfGamma.id() == neutral.id() && pfGamma.key() == neutral.key() ) isUnique = false;
134  }
135  }
136  }
137 }
138 
140 {
141  LogDebug("TauConstructorAddChH") << " region = " << region << ": Pt = " << chargedHadron.pt() << ", eta = " << chargedHadron.eta() << ", phi = " << chargedHadron.phi();
142  // CV: need to make sure that PFGammas merged with ChargedHadrons are not part of PiZeros
143  const std::vector<PFCandidatePtr>& neutrals = chargedHadron.getNeutralPFCandidates();
144  std::vector<PFCandidatePtr> neutrals_cleaned;
145  for ( std::vector<PFCandidatePtr>::const_iterator neutral = neutrals.begin();
146  neutral != neutrals.end(); ++neutral ) {
147  LogDebug("TauConstructorAddChH") << "neutral = " << neutral->id() << ":" << neutral->key();
148  bool isUnique = true;
149  if ( copyGammas_ ) checkOverlap(*neutral, *getSortedCollection(kSignal, kGamma), isUnique);
150  else checkOverlap(*neutral, tau_->signalPiZeroCandidates_, isUnique);
151  if ( region == kIsolation ) {
152  if ( copyGammas_ ) checkOverlap(*neutral, *getSortedCollection(kIsolation, kGamma), isUnique);
153  else checkOverlap(*neutral, tau_->isolationPiZeroCandidates_, isUnique);
154  }
155  LogDebug("TauConstructorAddChH") << "--> isUnique = " << isUnique;
156  if ( isUnique ) neutrals_cleaned.push_back(*neutral);
157  }
158  PFRecoTauChargedHadron chargedHadron_cleaned = chargedHadron;
159  if ( neutrals_cleaned.size() != neutrals.size() ) {
160  chargedHadron_cleaned.neutralPFCandidates_ = neutrals_cleaned;
161  setChargedHadronP4(chargedHadron_cleaned);
162  }
163  if ( region == kSignal ) {
164  tau_->signalTauChargedHadronCandidates_.push_back(chargedHadron_cleaned);
165  p4_ += chargedHadron_cleaned.p4();
166  if ( chargedHadron_cleaned.getChargedPFCandidate().isNonnull() ) {
167  addPFCand(kSignal, kChargedHadron, chargedHadron_cleaned.getChargedPFCandidate(), true);
168  }
169  const std::vector<PFCandidatePtr>& neutrals = chargedHadron_cleaned.getNeutralPFCandidates();
170  for ( std::vector<PFCandidatePtr>::const_iterator neutral = neutrals.begin();
171  neutral != neutrals.end(); ++neutral ) {
172  if ( (*neutral)->particleId() == reco::PFCandidate::gamma ) addPFCand(kSignal, kGamma, *neutral, true);
173  else if ( (*neutral)->particleId() == reco::PFCandidate::h0 ) addPFCand(kSignal, kNeutralHadron, *neutral, true);
174  };
175  } else {
176  tau_->isolationTauChargedHadronCandidates_.push_back(chargedHadron_cleaned);
177  if ( chargedHadron_cleaned.getChargedPFCandidate().isNonnull() ) {
178  if ( chargedHadron_cleaned.getChargedPFCandidate()->particleId() == reco::PFCandidate::h ) addPFCand(kIsolation, kChargedHadron, chargedHadron_cleaned.getChargedPFCandidate());
179  else if ( chargedHadron_cleaned.getChargedPFCandidate()->particleId() == reco::PFCandidate::h0 ) addPFCand(kIsolation, kNeutralHadron, chargedHadron_cleaned.getChargedPFCandidate());
180  }
181  const std::vector<PFCandidatePtr>& neutrals = chargedHadron_cleaned.getNeutralPFCandidates();
182  for ( std::vector<PFCandidatePtr>::const_iterator neutral = neutrals.begin();
183  neutral != neutrals.end(); ++neutral ) {
184  if ( (*neutral)->particleId() == reco::PFCandidate::gamma ) addPFCand(kIsolation, kGamma, *neutral);
185  else if ( (*neutral)->particleId() == reco::PFCandidate::h0 ) addPFCand(kIsolation, kNeutralHadron, *neutral);
186  };
187  }
188 }
189 
191 {
192  if ( region == kSignal ) {
193  tau_->signalPiZeroCandidates_.reserve(size);
194  // If we are building the gammas with the pizeros, resize that
195  // vector as well
196  if ( copyGammas_ ) reserve(kSignal, kGamma, 2*size);
197  } else {
198  tau_->isolationPiZeroCandidates_.reserve(size);
199  if ( copyGammas_ ) reserve(kIsolation, kGamma, 2*size);
200  }
201 }
202 
204 {
205  LogDebug("TauConstructorAddPi0") << " region = " << region << ": Pt = " << piZero.pt() << ", eta = " << piZero.eta() << ", phi = " << piZero.phi();
206  if ( region == kSignal ) {
207  tau_->signalPiZeroCandidates_.push_back(piZero);
208  // Copy the daughter gammas into the gamma collection if desired
209  if ( copyGammas_ ) {
210  // If we are using the pizeros to build the gammas, make sure we update
211  // the four vector correctly.
212  p4_ += piZero.p4();
213  addPFCands(kSignal, kGamma, piZero.daughterPtrVector().begin(),
214  piZero.daughterPtrVector().end());
215  }
216  } else {
217  tau_->isolationPiZeroCandidates_.push_back(piZero);
218  if ( copyGammas_ ) {
219  addPFCands(kIsolation, kGamma, piZero.daughterPtrVector().begin(),
220  piZero.daughterPtrVector().end());
221  }
222  }
223 }
224 
225 std::vector<PFCandidatePtr>*
227  return collections_[std::make_pair(region, type)];
228 }
229 
232  return sortedCollections_[std::make_pair(region, type)];
233 }
234 
235 // Trivial converter needed for polymorphism
237  const PFCandidatePtr& pfPtr) const {
238  return pfPtr;
239 }
240 
241 namespace {
242 // Make sure the two products come from the same EDM source
243 template<typename T1, typename T2>
244 void checkMatchedProductIds(const T1& t1, const T2& t2) {
245  if (t1.id() != t2.id()) {
246  throw cms::Exception("MismatchedPFCandSrc") << "Error: the input tag"
247  << " for the PF candidate collection provided to the RecoTauBuilder "
248  << " does not match the one that was used to build the source jets."
249  << " Please update the pfCandSrc paramters for the PFTau builders.";
250  }
251 }
252 }
253 
255  const PFCandidateRef& pfRef) const {
256  if(pfRef.isNonnull()) {
257  checkMatchedProductIds(pfRef, pfCands_);
258  return PFCandidatePtr(pfCands_, pfRef.key());
259  } else return PFCandidatePtr();
260 }
261 
262 // Convert from a CandidateRef to a Ptr
264  const CandidatePtr& candPtr) const {
265  if(candPtr.isNonnull()) {
266  checkMatchedProductIds(candPtr, pfCands_);
267  return PFCandidatePtr(pfCands_, candPtr.key());
268  } else return PFCandidatePtr();
269 }
270 
271 namespace {
272 template<typename T> bool ptDescending(const T& a, const T& b) {
273  return a.pt() > b.pt();
274 }
275 template<typename T> bool ptDescendingPtr(const T& a, const T& b) {
276  return a->pt() > b->pt();
277 }
278 }
279 
281  // The charged hadrons and pizeros are a special case, as we can sort them in situ
282  std::sort(tau_->signalTauChargedHadronCandidates_.begin(),
283  tau_->signalTauChargedHadronCandidates_.end(),
284  ptDescending<PFRecoTauChargedHadron>);
285  std::sort(tau_->isolationTauChargedHadronCandidates_.begin(),
286  tau_->isolationTauChargedHadronCandidates_.end(),
287  ptDescending<PFRecoTauChargedHadron>);
288  std::sort(tau_->signalPiZeroCandidates_.begin(),
289  tau_->signalPiZeroCandidates_.end(),
290  ptDescending<RecoTauPiZero>);
291  std::sort(tau_->isolationPiZeroCandidates_.begin(),
292  tau_->isolationPiZeroCandidates_.end(),
293  ptDescending<RecoTauPiZero>);
294 
295  // Sort each of our sortable collections, and copy them into the final
296  // tau RefVector.
297  BOOST_FOREACH ( const CollectionMap::value_type& colkey, collections_ ) {
298  SortedListPtr sortedCollection = sortedCollections_[colkey.first];
299  std::sort(sortedCollection->begin(),
300  sortedCollection->end(),
301  ptDescendingPtr<PFCandidatePtr>);
302  // Copy into the real tau collection
303  for ( std::vector<PFCandidatePtr>::const_iterator particle = sortedCollection->begin();
304  particle != sortedCollection->end(); ++particle ) {
305  colkey.second->push_back(*particle);
306  }
307  }
308 }
309 
310 namespace
311 {
312  PFTau::hadronicDecayMode calculateDecayMode(const reco::PFTau& tau, double dRsignalCone,
313  double minAbsPhotonSumPt_insideSignalCone, double minRelPhotonSumPt_insideSignalCone,
314  double minAbsPhotonSumPt_outsideSignalCone, double minRelPhotonSumPt_outsideSignalCone)
315  {
316  unsigned int nCharged = tau.signalTauChargedHadronCandidates().size();
317  // If no tracks exist, this is definitely not a tau!
318  if ( !nCharged ) return PFTau::kNull;
319 
320  unsigned int nPiZeros = 0;
321  const std::vector<RecoTauPiZero>& piZeros = tau.signalPiZeroCandidates();
322  for ( std::vector<RecoTauPiZero>::const_iterator piZero = piZeros.begin();
323  piZero != piZeros.end(); ++piZero ) {
324  double photonSumPt_insideSignalCone = 0.;
325  double photonSumPt_outsideSignalCone = 0.;
326  int numPhotons = piZero->numberOfDaughters();
327  for ( int idxPhoton = 0; idxPhoton < numPhotons; ++idxPhoton ) {
328  const reco::Candidate* photon = piZero->daughter(idxPhoton);
329  double dR = deltaR(photon->p4(), tau.p4());
330  if ( dR < dRsignalCone ) {
331  photonSumPt_insideSignalCone += photon->pt();
332  } else {
333  photonSumPt_outsideSignalCone += photon->pt();
334  }
335  }
336  if ( photonSumPt_insideSignalCone > minAbsPhotonSumPt_insideSignalCone || photonSumPt_insideSignalCone > (minRelPhotonSumPt_insideSignalCone*tau.pt()) ||
337  photonSumPt_outsideSignalCone > minAbsPhotonSumPt_outsideSignalCone || photonSumPt_outsideSignalCone > (minRelPhotonSumPt_outsideSignalCone*tau.pt()) ) ++nPiZeros;
338  }
339 
340  // Find the maximum number of PiZeros our parameterization can hold
341  const unsigned int maxPiZeros = PFTau::kOneProngNPiZero;
342 
343  // Determine our track index
344  unsigned int trackIndex = (nCharged - 1)*(maxPiZeros + 1);
345 
346  // Check if we handle the given number of tracks
347  if ( trackIndex >= PFTau::kRareDecayMode ) return PFTau::kRareDecayMode;
348 
349  if ( nPiZeros > maxPiZeros ) nPiZeros = maxPiZeros;
350  return static_cast<PFTau::hadronicDecayMode>(trackIndex + nPiZeros);
351  }
352 }
353 
354 std::auto_ptr<reco::PFTau> RecoTauConstructor::get(bool setupLeadingObjects)
355 {
356  LogDebug("TauConstructorGet") << "Start getting" ;
357 
358  // Copy the sorted collections into the interal tau refvectors
360 
361  // Setup all the important member variables of the tau
362  // Set charge of tau
363 // tau_->setCharge(
364 // sumPFCandCharge(getCollection(kSignal, kChargedHadron)->begin(),
365 // getCollection(kSignal, kChargedHadron)->end()));
366  // CV: take charge of highest pT charged hadron as charge of tau,
367  // in case tau does not have three reconstructed tracks
368  // (either because tau is reconstructed as 2prong or because PFRecoTauChargedHadron is built from a PFNeutralHadron)
369  unsigned int nCharged = 0;
370  int charge = 0;
371  double leadChargedHadronPt = 0.;
372  int leadChargedHadronCharge = 0;
373  for ( std::vector<PFRecoTauChargedHadron>::const_iterator chargedHadron = tau_->signalTauChargedHadronCandidates_.begin();
374  chargedHadron != tau_->signalTauChargedHadronCandidates_.end(); ++chargedHadron ) {
375  if ( chargedHadron->algoIs(PFRecoTauChargedHadron::kChargedPFCandidate) || chargedHadron->algoIs(PFRecoTauChargedHadron::kTrack) ) {
376  ++nCharged;
377  charge += chargedHadron->charge();
378  if ( chargedHadron->pt() > leadChargedHadronPt ) {
379  leadChargedHadronPt = chargedHadron->pt();
380  leadChargedHadronCharge = chargedHadron->charge();
381  }
382  }
383  }
384  if ( nCharged == 3 ) tau_->setCharge(charge);
385  else tau_->setCharge(leadChargedHadronCharge);
386 
387  // Set PDG id
388  tau_->setPdgId(tau_->charge() < 0 ? 15 : -15);
389 
390  // Set P4
391  tau_->setP4(p4_);
392 
393  // Set Decay Mode
394  double dRsignalCone = ( signalConeSize_ ) ? (*signalConeSize_)(*tau_) : 0.5;
395  tau_->setSignalConeSize(dRsignalCone);
396  PFTau::hadronicDecayMode dm = calculateDecayMode(
397  *tau_,
398  dRsignalCone,
400  tau_->setDecayMode(dm);
401 
402  LogDebug("TauConstructorGet") << "Pt = " << tau_->pt() << ", eta = " << tau_->eta() << ", phi = " << tau_->phi() << ", mass = " << tau_->mass() << ", dm = " << tau_->decayMode() ;
403 
404  // Set charged isolation quantities
405  tau_->setisolationPFChargedHadrCandsPtSum(
406  sumPFCandPt(
409 
410  // Set gamma isolation quantities
411  tau_->setisolationPFGammaCandsEtSum(
412  sumPFCandPt(
415 
416  // Set em fraction
417  tau_->setemFraction(sumPFCandPt(
419  getCollection(kSignal, kGamma)->end()) / tau_->pt());
420 
421  if ( setupLeadingObjects ) {
422  typedef std::vector<PFCandidatePtr>::const_iterator Iter;
423  // Find the highest PT object in the signal cone
424  Iter leadingCand = leadPFCand(
427 
428  if ( leadingCand != getCollection(kSignal, kAll)->end() )
429  tau_->setleadPFCand(*leadingCand);
430 
431  // Hardest charged object in signal cone
432  Iter leadingChargedCand = leadPFCand(
435 
436  if ( leadingChargedCand != getCollection(kSignal, kChargedHadron)->end() )
437  tau_->setleadPFChargedHadrCand(*leadingChargedCand);
438 
439  // Hardest gamma object in signal cone
440  Iter leadingGammaCand = leadPFCand(
443 
444  if(leadingGammaCand != getCollection(kSignal, kGamma)->end())
445  tau_->setleadPFNeutralCand(*leadingGammaCand);
446  }
447  return tau_;
448 }
449 }} // end namespace reco::tau
#define LogDebug(id)
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_
virtual const Candidate * daughter(size_type i) const =0
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
std::auto_ptr< reco::PFTau > get(bool setupLeadingCandidates=true)
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:250
void reserve(Region region, ParticleType type, size_t size)
Reserve a set amount of space for a given RefVector.
double sumPFCandPt(InputIterator begin, InputIterator end)
Sum the pT of a collection of PFCandidates.
key_type key() const
Definition: Ptr.h:169
boost::shared_ptr< std::vector< PFCandidatePtr > > SortedListPtr
std::vector< PFCandidatePtr > * getCollection(Region region, ParticleType type)
virtual double pt() const =0
transverse momentum
std::vector< PFCandidatePtr > pfGammas(const PFJet &jet, bool sort=true)
Extract all pfGammas from a PFJet.
double deltaR(const T1 &t1, const T2 &t2)
Definition: deltaR.h:48
key_type key() const
Accessor for product key.
Definition: Ref.h:266
virtual double eta() const
momentum pseudorapidity
virtual double pt() const
transverse momentum
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_
const std::vector< RecoTauPiZero > & signalPiZeroCandidates() const
Retrieve the association of signal region gamma candidates into candidate PiZeros.
Definition: PFTau.cc:97
SortedListPtr getSortedCollection(Region region, ParticleType type)
#define end
Definition: vmac.h:37
Container::value_type value_type
void addPFCands(Region region, ParticleType type, const InputIterator &begin, const InputIterator &end)
InputIterator leadPFCand(InputIterator begin, InputIterator end)
bool isNonnull() const
Checks for non-null.
Definition: Ptr.h:152
PFCandidatePtr convertToPtr(const PFCandidatePtr &pfPtr) const
unsigned int nCharged(const GenJet &jet)
std::auto_ptr< reco::PFTau > tau_
void reserveTauChargedHadron(Region region, size_t size)
Reserve a set amount of space for the ChargedHadrons.
const std::vector< PFCandidatePtr > & getNeutralPFCandidates() const
references to additional neutral PFCandidates
void reservePiZero(Region region, size_t size)
Reserve a set amount of space for the PiZeros.
double b
Definition: hdecay.h:120
std::vector< PFCandidatePtr > neutralPFCandidates_
ProductID id() const
Accessor for product ID.
Definition: Ptr.h:164
hadronicDecayMode
Definition: PFTau.h:35
#define begin
Definition: vmac.h:30
double a
Definition: hdecay.h:121
const daughters & daughterPtrVector() const
references to daughtes
const std::vector< PFRecoTauChargedHadron > & signalTauChargedHadronCandidates() const
Retrieve the association of signal region PF candidates into candidate PFRecoTauChargedHadrons.
Definition: PFTau.cc:144
void addPiZero(Region region, const RecoTauPiZero &piZero)
Add a PiZero to the given collection.
const StringObjectFunction< reco::PFTau > * signalConeSize_
ProductIndex id() const
Definition: ProductID.h:38
void setChargedHadronP4(reco::PFRecoTauChargedHadron &chargedHadron, double scaleFactor_neutralPFCands=1.0)
long double T
edm::Ptr< PFCandidate > PFCandidatePtr
persistent Ptr to a PFCandidate
virtual double phi() const
momentum azimuthal angle
virtual const LorentzVector & p4() const
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
tuple size
Write out results.
const PFCandidatePtr & getChargedPFCandidate() const
reference to &quot;charged&quot; PFCandidate (either charged PFCandidate or PFNeutralHadron) ...
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector