CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
RecoTauBuilderCombinatoricPlugin.cc
Go to the documentation of this file.
1 #include <vector>
2 
5 
9 
16 
17 namespace reco { namespace tau {
18 
19 typedef std::vector<reco::PFRecoTauChargedHadron> ChargedHadronList;
21 typedef std::vector<RecoTauPiZero> PiZeroList;
23 
25 {
26  public:
29 
31  const reco::PFJetRef&,
32  const std::vector<reco::PFRecoTauChargedHadron>&,
33  const std::vector<RecoTauPiZero>&,
34  const std::vector<PFCandidatePtr>&) const override;
35 
36  private:
38 
40 
41  struct decayModeInfo
42  {
43  uint32_t maxPiZeros_;
44  uint32_t maxPFCHs_;
45  uint32_t nCharged_;
46  uint32_t nPiZeros_;
47  };
48  std::vector<decayModeInfo> decayModesToBuild_;
49 
51 };
52 
54  : RecoTauBuilderPlugin(pset, std::move(iC)),
55  qcuts_(pset.getParameterSet("qualityCuts").getParameterSet("signalQualityCuts")),
56  isolationConeSize_(pset.getParameter<double>("isolationConeSize"))
57 {
58  typedef std::vector<edm::ParameterSet> VPSet;
59  const VPSet& decayModes = pset.getParameter<VPSet>("decayModes");
60  for ( VPSet::const_iterator decayMode = decayModes.begin();
61  decayMode != decayModes.end(); ++decayMode ) {
63  info.nCharged_ = decayMode->getParameter<uint32_t>("nCharged");
64  info.nPiZeros_ = decayMode->getParameter<uint32_t>("nPiZeros");
65  info.maxPFCHs_ = decayMode->getParameter<uint32_t>("maxTracks");
66  info.maxPiZeros_ = decayMode->getParameter<uint32_t>("maxPiZeros");
67  decayModesToBuild_.push_back(info);
68  }
69 
70  verbosity_ = ( pset.exists("verbosity") ) ?
71  pset.getParameter<int>("verbosity") : 0;
72 }
73 
74 // define template specialization for cross-cleaning
75 namespace xclean
76 {
77  template<>
79  {
80  // Get the list of objects we need to clean
81  for ( ChargedHadronCombo::combo_iterator chargedHadron = chargedHadronsBegin; chargedHadron != chargedHadronsEnd; ++chargedHadron ) {
82  // CV: Remove PFGammas that are merged into TauChargedHadrons from isolation PiZeros, but not from signal PiZeros.
83  // The overlap between PFGammas contained in signal PiZeros and merged into TauChargedHadrons
84  // is resolved by RecoTauConstructor::addTauChargedHadron,
85  // which gives preference to PFGammas that are within PiZeros and removes those PFGammas from TauChargedHadrons.
86  if ( mode_ == kRemoveChargedDaughterOverlaps ) {
87  if ( chargedHadron->getChargedPFCandidate().isNonnull() ) toRemove_.insert(reco::CandidatePtr(chargedHadron->getChargedPFCandidate()));
88  } else if ( mode_ == kRemoveChargedAndNeutralDaughterOverlaps ) {
89  const reco::CompositePtrCandidate::daughters& daughters = chargedHadron->daughterPtrVector();
90  for ( reco::CompositePtrCandidate::daughters::const_iterator daughter = daughters.begin();
91  daughter != daughters.end(); ++daughter ) {
92  toRemove_.insert(reco::CandidatePtr(*daughter));
93  }
94  } else assert(0);
95  }
96  }
97 
98  template<>
99  inline void CrossCleanPtrs<PiZeroList::const_iterator>::initialize(const PiZeroList::const_iterator& piZerosBegin, const PiZeroList::const_iterator& piZerosEnd)
100  {
101  BOOST_FOREACH( const PFCandidatePtr &ptr, flattenPiZeros(piZerosBegin, piZerosEnd) ) {
102  toRemove_.insert(CandidatePtr(ptr));
103  }
104  }
105 
106  template<>
108  {
109  //std::cout << "<CrossCleanPtrs<ChargedHadronCombo>::initialize>:" << std::endl;
110  for ( ChargedHadronCombo::combo_iterator chargedHadron = chargedHadronsBegin; chargedHadron != chargedHadronsEnd; ++chargedHadron ) {
111  const reco::CompositePtrCandidate::daughters& daughters = chargedHadron->daughterPtrVector();
112  for ( reco::CompositePtrCandidate::daughters::const_iterator daughter = daughters.begin();
113  daughter != daughters.end(); ++daughter ) {
114  //std::cout << " adding PFCandidate = " << daughter->id() << ":" << daughter->key() << std::endl;
115  toRemove_.insert(reco::CandidatePtr(*daughter));
116  }
117  }
118  }
119 
120  template<>
121  inline void CrossCleanPtrs<ChargedHadronList::const_iterator>::initialize(const ChargedHadronList::const_iterator& chargedHadronsBegin, const ChargedHadronList::const_iterator& chargedHadronsEnd)
122  {
123  //std::cout << "<CrossCleanPtrs<ChargedHadronList>::initialize>:" << std::endl;
124  for ( ChargedHadronList::const_iterator chargedHadron = chargedHadronsBegin; chargedHadron != chargedHadronsEnd; ++chargedHadron ) {
125  const reco::CompositePtrCandidate::daughters& daughters = chargedHadron->daughterPtrVector();
126  for ( reco::CompositePtrCandidate::daughters::const_iterator daughter = daughters.begin();
127  daughter != daughters.end(); ++daughter ) {
128  //std::cout << " adding PFCandidate = " << daughter->id() << ":" << daughter->key() << std::endl;
129  toRemove_.insert(reco::CandidatePtr(*daughter));
130  }
131  }
132  }
133 }
134 
135 namespace
136 {
137  // auxiliary class for sorting pizeros by descending transverse momentum
138  class SortPi0sDescendingPt
139  {
140  public:
141  bool operator()(const RecoTauPiZero& a, const RecoTauPiZero& b) const
142  {
143  return a.pt() > b.pt();
144  }
145  };
146 
147  std::string getPFCandidateType(reco::PFCandidate::ParticleType pfCandidateType)
148  {
149  if ( pfCandidateType == reco::PFCandidate::X ) return "undefined";
150  else if ( pfCandidateType == reco::PFCandidate::h ) return "PFChargedHadron";
151  else if ( pfCandidateType == reco::PFCandidate::e ) return "PFElectron";
152  else if ( pfCandidateType == reco::PFCandidate::mu ) return "PFMuon";
153  else if ( pfCandidateType == reco::PFCandidate::gamma ) return "PFGamma";
154  else if ( pfCandidateType == reco::PFCandidate::h0 ) return "PFNeutralHadron";
155  else if ( pfCandidateType == reco::PFCandidate::h_HF ) return "HF_had";
156  else if ( pfCandidateType == reco::PFCandidate::egamma_HF ) return "HF_em";
157  else assert(0);
158  }
159 }
160 
163  const reco::PFJetRef& jet,
164  const std::vector<reco::PFRecoTauChargedHadron>& chargedHadrons,
165  const std::vector<RecoTauPiZero>& piZeros,
166  const std::vector<PFCandidatePtr>& regionalExtras) const
167 {
168  if ( verbosity_ ) {
169  std::cout << "<RecoTauBuilderCombinatoricPlugin::operator()>:" << std::endl;
170  std::cout << " processing jet: Pt = " << jet->pt() << ", eta = " << jet->eta() << ", phi = " << jet->eta() << ","
171  << " mass = " << jet->mass() << ", area = " << jet->jetArea() << std::endl;
172  }
173 
174  // Define output.
176 
177  reco::VertexRef primaryVertexRef = primaryVertex(jet);
178 
179  // Update the primary vertex used by the quality cuts. The PV is supplied by
180  // the base class.
181  qcuts_.setPV(primaryVertexRef);
182 
183  typedef std::vector<PFCandidatePtr> PFCandPtrs;
184 
185  if ( verbosity_ ) {
186  std::cout << "#chargedHadrons = " << chargedHadrons.size() << std::endl;
187  int idx = 0;
188  for ( ChargedHadronList::const_iterator chargedHadron = chargedHadrons.begin();
189  chargedHadron != chargedHadrons.end(); ++chargedHadron ) {
190  std::cout << "chargedHadron #" << idx << ":" << std::endl;
191  chargedHadron->print(std::cout);
192  ++idx;
193  }
194  std::cout << "#piZeros = " << piZeros.size() << std::endl;
195  idx = 0;
196  for ( PiZeroList::const_iterator piZero = piZeros.begin();
197  piZero != piZeros.end(); ++piZero ) {
198  std::cout << "piZero #" << idx << ": Pt = " << piZero->pt() << ", eta = " << piZero->eta() << ", phi = " << piZero->phi() << std::endl;
199  size_t numDaughters = piZero->numberOfDaughters();
200  for ( size_t iDaughter = 0; iDaughter < numDaughters; ++iDaughter ) {
201  const reco::PFCandidate* daughter = dynamic_cast<const reco::PFCandidate*>(piZero->daughterPtr(iDaughter).get());
202  std::cout << " daughter #" << iDaughter << " (" << getPFCandidateType(daughter->particleId()) << "):"
203  << " Pt = " << daughter->pt() << ", eta = " << daughter->eta() << ", phi = " << daughter->phi() << std::endl;
204  }
205  ++idx;
206  }
207  }
208 
209  PFCandPtrs pfchs = qcuts_.filterCandRefs(pfChargedCands(*jet));
210  PFCandPtrs pfnhs = qcuts_.filterCandRefs(pfCandidates(*jet, reco::PFCandidate::h0));
211 
214  PFCandPtrs regionalJunk = qcuts_.filterCandRefs(regionalExtras);
215 
216  // Loop over the decay modes we want to build
217  for ( std::vector<decayModeInfo>::const_iterator decayMode = decayModesToBuild_.begin();
218  decayMode != decayModesToBuild_.end(); ++decayMode ) {
219  // Find how many piZeros are in this decay mode
220  size_t piZerosToBuild = decayMode->nPiZeros_;
221  // Find how many tracks are in this decay mode
222  size_t tracksToBuild = decayMode->nCharged_;
223  if ( verbosity_ ) {
224  std::cout << "piZerosToBuild = " << piZerosToBuild << std::endl;
225  std::cout << "tracksToBuild = " << tracksToBuild << std::endl;
226  std::cout << "#chargedHadrons = " << chargedHadrons.size() << std::endl;
227  }
228 
229  // Skip decay mode if jet doesn't have the multiplicity to support it
230  if ( chargedHadrons.size() < tracksToBuild ) continue;
231 
232  // Find the start and end of potential signal tracks
233  ChargedHadronList::const_iterator chargedHadron_begin = chargedHadrons.begin();
234  ChargedHadronList::const_iterator chargedHadron_end = chargedHadrons.end();
235  chargedHadron_end = takeNElements(chargedHadron_begin, chargedHadron_end, decayMode->maxPFCHs_);
236 
237  // Build our track combo generator
238  ChargedHadronCombo trackCombos(chargedHadron_begin, chargedHadron_end, tracksToBuild);
239 
240  PFCandPtrs::iterator pfch_end = pfchs.end();
241  pfch_end = takeNElements(pfchs.begin(), pfch_end, decayMode->maxPFCHs_);
242 
243  //-------------------------------------------------------
244  // Begin combinatoric loop for this decay mode
245  //-------------------------------------------------------
246 
247  // Loop over the different combinations of tracks
248  for ( ChargedHadronCombo::iterator trackCombo = trackCombos.begin();
249  trackCombo != trackCombos.end(); ++trackCombo ) {
251  trackCombo->combo_begin(), trackCombo->combo_end(),
253 
254  PiZeroList cleanSignalPiZeros = signalPiZeroXCleaner(piZeros);
255 
256  // CV: sort collection of cross-cleaned pi0s by descending Pt
257  std::sort(cleanSignalPiZeros.begin(), cleanSignalPiZeros.end(), SortPi0sDescendingPt());
258 
259  // Skip decay mode if we don't have enough remaining clean pizeros to
260  // build it.
261  if ( cleanSignalPiZeros.size() < piZerosToBuild ) continue;
262 
263  // Find the start and end of potential signal tracks
264  PiZeroList::iterator signalPiZero_begin = cleanSignalPiZeros.begin();
265  PiZeroList::iterator signalPiZero_end = cleanSignalPiZeros.end();
266  signalPiZero_end = takeNElements(signalPiZero_begin, signalPiZero_end, decayMode->maxPiZeros_);
267 
268  // Build our piZero combo generator
269  PiZeroCombo piZeroCombos(signalPiZero_begin, signalPiZero_end, piZerosToBuild);
270  // Loop over the different combinations of PiZeros
271  for ( PiZeroCombo::iterator piZeroCombo = piZeroCombos.begin();
272  piZeroCombo != piZeroCombos.end(); ++piZeroCombo ) {
273  // Output tau
274  RecoTauConstructor tau(jet, getPFCands(), true);
275  // Reserve space in our collections
276  tau.reserve(
278  RecoTauConstructor::kChargedHadron, tracksToBuild);
279  tau.reserve(
281  RecoTauConstructor::kGamma, 2*piZerosToBuild); // k-factor = 2
282  tau.reservePiZero(RecoTauConstructor::kSignal, piZerosToBuild);
283 
285  trackCombo->combo_begin(), trackCombo->combo_end(),
287 
288  PiZeroList precleanedIsolationPiZeros = isolationPiZeroXCleaner(piZeros);
289  std::set<reco::CandidatePtr> toRemove;
290  for ( PiZeroCombo::combo_iterator signalPiZero = piZeroCombo->combo_begin();
291  signalPiZero != piZeroCombo->combo_end(); ++signalPiZero ) {
292  toRemove.insert(signalPiZero->daughterPtrVector().begin(), signalPiZero->daughterPtrVector().end());
293  }
294  PiZeroList cleanIsolationPiZeros;
295  BOOST_FOREACH( const RecoTauPiZero& precleanedPiZero, precleanedIsolationPiZeros ) {
296  std::set<reco::CandidatePtr> toCheck(precleanedPiZero.daughterPtrVector().begin(), precleanedPiZero.daughterPtrVector().end());
297  std::vector<reco::CandidatePtr> cleanDaughters;
298  std::set_difference(toCheck.begin(), toCheck.end(), toRemove.begin(), toRemove.end(), std::back_inserter(cleanDaughters));
299  // CV: piZero is signal piZero if at least one daughter overlaps
300  if ( cleanDaughters.size() == precleanedPiZero.daughterPtrVector().size() ) {
301  cleanIsolationPiZeros.push_back(precleanedPiZero);
302  }
303  }
304 
305  // FIXME - are all these reserves okay? will they get propagated to the
306  // dataformat size if they are wrong?
307  tau.reserve(
309  RecoTauConstructor::kChargedHadron, chargedHadrons.size() - tracksToBuild);
310  tau.reserve(
313  (piZeros.size() - piZerosToBuild)*2);
314  tau.reservePiZero(
316  (piZeros.size() - piZerosToBuild));
317 
318  // Get signal PiZero constituents and add them to the tau.
319  // The sub-gammas are automatically added.
320  tau.addPiZeros(
322  piZeroCombo->combo_begin(), piZeroCombo->combo_end());
323 
324  // Set signal and isolation components for charged hadrons, after
325  // converting them to a PFCandidateRefVector
326  //
327  // NOTE: signal ChargedHadrons need to be added **after** signal PiZeros
328  // to avoid double-counting PFGammas as part of PiZero and merged with ChargedHadron
329  //
332  trackCombo->combo_begin(), trackCombo->combo_end());
333 
334  // Now build isolation collections
335  // Load our isolation tools
336  using namespace reco::tau::cone;
337  PFCandPtrDRFilter isolationConeFilter(tau.p4(), -0.1, isolationConeSize_);
338 
339  // Cross cleaning predicate: Remove any PFCandidatePtrs that are contained within existing ChargedHadrons or PiZeros.
340  // The predicate will return false for any object that overlaps with chargedHadrons or cleanPiZeros.
341  // 1.) to select charged PFCandidates within jet that are not signalPFChargedHadrons
342  typedef xclean::CrossCleanPtrs<ChargedHadronCombo::combo_iterator> pfChargedHadronXCleanerType;
343  pfChargedHadronXCleanerType pfChargedHadronXCleaner_comboChargedHadrons(trackCombo->combo_begin(), trackCombo->combo_end());
344  // And this cleaning filter predicate with our Iso cone filter
345  xclean::PredicateAND<PFCandPtrDRFilter, pfChargedHadronXCleanerType> pfCandFilter_comboChargedHadrons(isolationConeFilter, pfChargedHadronXCleaner_comboChargedHadrons);
346  // 2.) to select neutral PFCandidates within jet
347  xclean::CrossCleanPtrs<ChargedHadronList::const_iterator> pfChargedHadronXCleaner_allChargedHadrons(chargedHadrons.begin(), chargedHadrons.end());
348  xclean::CrossCleanPtrs<PiZeroList::const_iterator> piZeroXCleaner(cleanIsolationPiZeros.begin(), cleanIsolationPiZeros.end());
350  pfCandXCleanerType pfCandXCleaner_allChargedHadrons(pfChargedHadronXCleaner_allChargedHadrons, piZeroXCleaner);
351  // And this cleaning filter predicate with our Iso cone filter
352  xclean::PredicateAND<PFCandPtrDRFilter, pfCandXCleanerType> pfCandFilter_allChargedHadrons(isolationConeFilter, pfCandXCleaner_allChargedHadrons);
353 
354  ChargedHadronDRFilter isolationConeFilterChargedHadron(tau.p4(), -0.1, isolationConeSize_);
355  PiZeroDRFilter isolationConeFilterPiZero(tau.p4(), -0.1, isolationConeSize_);
356 
357  // Additionally make predicates to select the different PF object types
358  // of the regional junk objects to add
360  PFCandPtrDRFilter> RegionalJunkConeAndIdFilter;
361 
363  pfchCandSelector(reco::PFCandidate::h);
365  pfgammaCandSelector(reco::PFCandidate::gamma);
367  pfnhCandSelector(reco::PFCandidate::h0);
368 
369  RegionalJunkConeAndIdFilter pfChargedJunk(
370  pfchCandSelector, // select charged stuff from junk
371  isolationConeFilter); // only take those in iso cone
372 
373  RegionalJunkConeAndIdFilter pfGammaJunk(
374  pfgammaCandSelector, // select gammas from junk
375  isolationConeFilter); // only take those in iso cone
376 
377  RegionalJunkConeAndIdFilter pfNeutralJunk(
378  pfnhCandSelector, // select neutral stuff from junk
379  isolationConeFilter); // select stuff in iso cone
380 
381  tau.addPiZeros(
383  boost::make_filter_iterator(
384  isolationConeFilterPiZero,
385  cleanIsolationPiZeros.begin(), cleanIsolationPiZeros.end()),
386  boost::make_filter_iterator(
387  isolationConeFilterPiZero,
388  cleanIsolationPiZeros.end(), cleanIsolationPiZeros.end()));
389 
390  // Filter the isolation candidates in a DR cone
391  //
392  // NOTE: isolation ChargedHadrons need to be added **after** signal and isolation PiZeros
393  // to avoid double-counting PFGammas as part of PiZero and merged with ChargedHadron
394  //
395  if ( verbosity_ >= 2 ) {
396  std::cout << "adding isolation PFChargedHadrons from trackCombo:" << std::endl;
397  }
400  boost::make_filter_iterator(
401  isolationConeFilterChargedHadron,
402  trackCombo->remainder_begin(), trackCombo->remainder_end()),
403  boost::make_filter_iterator(
404  isolationConeFilterChargedHadron,
405  trackCombo->remainder_end(), trackCombo->remainder_end()));
406 
407  // Add all the candidates that weren't included in the combinatoric
408  // generation
409  if ( verbosity_ >= 2 ) {
410  std::cout << "adding isolation PFChargedHadrons not considered in trackCombo:" << std::endl;
411  }
412  tau.addPFCands(
414  boost::make_filter_iterator(
415  pfCandFilter_comboChargedHadrons,
416  pfch_end, pfchs.end()),
417  boost::make_filter_iterator(
418  pfCandFilter_comboChargedHadrons,
419  pfchs.end(), pfchs.end()));
420  // Add all charged candidates that are in the iso cone but weren't in the
421  // original PFJet
422  if ( verbosity_ >= 2 ) {
423  std::cout << "adding isolation PFChargedHadrons from 'regional junk':" << std::endl;
424  }
425  tau.addPFCands(
427  boost::make_filter_iterator(
428  pfChargedJunk, regionalJunk.begin(), regionalJunk.end()),
429  boost::make_filter_iterator(
430  pfChargedJunk, regionalJunk.end(), regionalJunk.end()));
431 
432  // Add all gammas that are in the iso cone but weren't in the
433  // orginal PFJet
434  tau.addPFCands(
436  boost::make_filter_iterator(
437  pfGammaJunk, regionalJunk.begin(), regionalJunk.end()),
438  boost::make_filter_iterator(
439  pfGammaJunk, regionalJunk.end(), regionalJunk.end()));
440 
441  // Add all the neutral hadron candidates to the isolation collection
442  tau.addPFCands(
444  boost::make_filter_iterator(
445  pfCandFilter_allChargedHadrons,
446  pfnhs.begin(), pfnhs.end()),
447  boost::make_filter_iterator(
448  pfCandFilter_allChargedHadrons,
449  pfnhs.end(), pfnhs.end()));
450  // Add all the neutral hadrons from the region collection that are in
451  // the iso cone to the tau
452  tau.addPFCands(
454  boost::make_filter_iterator(
455  pfNeutralJunk, regionalJunk.begin(), regionalJunk.end()),
456  boost::make_filter_iterator(
457  pfNeutralJunk, regionalJunk.end(), regionalJunk.end()));
458 
459  std::auto_ptr<reco::PFTau> tauPtr = tau.get(true);
460 
461  if ( primaryVertexRef.isNonnull() ) tauPtr->setVertex(primaryVertexRef->position());
462 
463  output.push_back(tauPtr);
464  }
465  }
466  }
467 
468  return output.release();
469 }
470 
471 }} // end namespace reco::tau
472 
476  "RecoTauBuilderCombinatoricPlugin");
T getParameter(std::string const &) const
const edm::Handle< PFCandidateCollection > & getPFCands() const
Hack to be able to convert Ptrs to Refs.
InputIterator takeNElements(const InputIterator &begin, const InputIterator &end, size_t N)
static const TGPicture * info(bool iBackgroundIsBlack)
ParticleType
particle types
Definition: PFCandidate.h:44
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.
void addTauChargedHadrons(Region region, const InputIterator &begin, const InputIterator &end)
Add a list of charged hadrons to the input collection.
virtual float pt() const
transverse momentum
CombinatoricIterator< T > end()
Coll filterCandRefs(const Coll &refcoll, bool invert=false) const
Filter a ref vector of PFCandidates.
std::vector< reco::PFRecoTauChargedHadron > ChargedHadronList
virtual float phi() const
momentum azimuthal angle
ParameterSet const & getParameterSet(ParameterSetID const &id)
std::vector< reco::PFCandidatePtr > PFCandPtrs
bool exists(std::string const &parameterName) const
checks if a parameter exists
boost::ptr_vector< reco::PFTau > output_type
std::auto_ptr< output_type > return_type
return_type operator()(const reco::PFJetRef &, const std::vector< reco::PFRecoTauChargedHadron > &, const std::vector< RecoTauPiZero > &, const std::vector< PFCandidatePtr > &) const override
RecoTauBuilderCombinatoricPlugin(const edm::ParameterSet &pset, edm::ConsumesCollector &&iC)
tau::CombinatoricGenerator< PiZeroList > PiZeroCombo
std::vector< PFCandidatePtr > pfCandidates(const PFJet &jet, int particleId, bool sort=true)
void addPiZeros(Region region, const InputIterator &begin, const InputIterator &end)
Add a list of pizeros to the input collection.
reco::VertexRef primaryVertex(const reco::PFJetRef &jet) const
Get primary vertex associated to this jet.
void setPV(const reco::VertexRef &vtx) const
Update the primary vertex.
virtual float eta() const
momentum pseudorapidity
const reco::Candidate::LorentzVector & p4() const
void addPFCands(Region region, ParticleType type, const InputIterator &begin, const InputIterator &end)
void initialize(const PtrIter &chargedHadronsBegin, const PtrIter &chargedHadronsEnd)
tau::CombinatoricGenerator< ChargedHadronList > ChargedHadronCombo
std::vector< PFCandidatePtr > flattenPiZeros(const std::vector< RecoTauPiZero >::const_iterator &, const std::vector< RecoTauPiZero >::const_iterator &)
Flatten a list of pi zeros into a list of there constituent PFCandidates.
std::vector< CandidatePtr > daughters
collection of references to daughters
edm::Ptr< Candidate > CandidatePtr
persistent reference to an object in a collection of Candidate objects
Definition: CandidateFwd.h:25
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
void reservePiZero(Region region, size_t size)
Reserve a set amount of space for the PiZeros.
double b
Definition: hdecay.h:120
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:39
std::vector< PFCandidatePtr > pfChargedCands(const PFJet &jet, bool sort=true)
Extract all non-neutral candidates from a PFJet.
double a
Definition: hdecay.h:121
const daughters & daughterPtrVector() const
references to daughtes
tuple cout
Definition: gather_cfg.py:121
#define DEFINE_EDM_PLUGIN(factory, type, name)
virtual ParticleType particleId() const
Definition: PFCandidate.h:369
void initialize(const PtrIter &particlesBegin, const PtrIter &particlesEnd)
std::vector< RecoTauPiZero > PiZeroList
iterator::value_type::ValueIter combo_iterator
Transform a pizero to remove given candidates.
CombinatoricIterator< T > begin()