CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Protected Attributes
MCTruthHelper< P > Class Template Reference

#include <MCTruthHelper.h>

Public Member Functions

int absPdgId (const reco::GenParticle &p)
 
int absPdgId (const HepMC::GenParticle &p)
 
const reco::GenParticledaughter (const reco::GenParticle &p, unsigned int idau)
 
const HepMC::GenParticle * daughter (const HepMC::GenParticle &p, unsigned int idau)
 
void fillGenStatusFlags (const P &p, reco::GenStatusFlags &statusFlags)
 
const PfindDecayedMother (const P &p)
 
const PfindDecayedMother (const P &p, int abspdgid)
 
const PfirstCopy (const P &p)
 
bool fromHardProcess (const P &p)
 
bool fromHardProcessBeforeFSR (const P &p)
 
bool fromHardProcessDecayed (const P &p)
 
bool fromHardProcessFinalState (const P &p)
 
const PhardProcessMotherCopy (const P &p)
 
bool isDecayedLeptonHadron (const P &p)
 
bool isDirectHadronDecayProduct (const P &p)
 
bool isDirectHardProcessTauDecayProduct (const P &p)
 
bool isDirectPromptTauDecayProduct (const P &p)
 
bool isDirectTauDecayProduct (const P &p)
 
bool isFirstCopy (const P &p)
 
bool isHadron (const P &p)
 
bool isHardProcess (const P &p)
 
bool isHardProcessTauDecayProduct (const P &p)
 
bool isLastCopy (const P &p)
 
bool isLastCopyBeforeFSR (const P &p)
 
bool isMuonDecayProduct (const P &p)
 
bool isPrompt (const P &p)
 
bool isPromptDecayed (const P &p)
 
bool isPromptFinalState (const P &p)
 
bool isPromptMuonDecayProduct (const P &p)
 
bool isPromptTauDecayProduct (const P &p)
 
bool isTauDecayProduct (const P &p)
 
const PlastCopy (const P &p)
 
const PlastCopyBeforeFSR (const P &p)
 
const PlastDaughterCopyBeforeFSR (const P &p)
 
const reco::GenParticlemother (const reco::GenParticle &p, unsigned int imoth=0)
 
const HepMC::GenParticle * mother (const HepMC::GenParticle &p, unsigned int imoth=0)
 
const PnextCopy (const P &p)
 
unsigned int numberOfDaughters (const reco::GenParticle &p)
 
unsigned int numberOfDaughters (const HepMC::GenParticle &p)
 
unsigned int numberOfMothers (const reco::GenParticle &p)
 
unsigned int numberOfMothers (const HepMC::GenParticle &p)
 
int pdgId (const reco::GenParticle &p)
 
int pdgId (const HepMC::GenParticle &p)
 
const PpreviousCopy (const P &p)
 
const PuniqueMother (const P &p)
 

Protected Attributes

std::unordered_set< const P * > dupCheck_
 

Detailed Description

template<typename P>
class MCTruthHelper< P >

Definition at line 14 of file MCTruthHelper.h.

Member Function Documentation

template<typename P >
int MCTruthHelper< P >::absPdgId ( const reco::GenParticle p)

Definition at line 576 of file MCTruthHelper.h.

References funct::abs(), and reco::LeafCandidate::pdgId().

576  {
577  return std::abs(p.pdgId());
578 }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
virtual int pdgId() const final
PDG identifier.
template<typename P >
int MCTruthHelper< P >::absPdgId ( const HepMC::GenParticle &  p)

Definition at line 582 of file MCTruthHelper.h.

References funct::abs().

582  {
583  return std::abs(p.pdg_id());
584 }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
template<typename P >
const reco::GenParticle * MCTruthHelper< P >::daughter ( const reco::GenParticle p,
unsigned int  idau 
)

Definition at line 624 of file MCTruthHelper.h.

References reco::CompositeRefCandidateT< D >::daughter().

624  {
625  return static_cast<const reco::GenParticle*>(p.daughter(idau));
626 }
virtual const Candidate * daughter(size_type) const
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
template<typename P >
const HepMC::GenParticle * MCTruthHelper< P >::daughter ( const HepMC::GenParticle &  p,
unsigned int  idau 
)

Definition at line 630 of file MCTruthHelper.h.

630  {
631  return *(p.end_vertex()->particles_out_const_begin() + idau);
632 }
template<typename P>
void MCTruthHelper< P >::fillGenStatusFlags ( const P p,
reco::GenStatusFlags statusFlags 
)

Definition at line 636 of file MCTruthHelper.h.

References reco::GenStatusFlags::setFromHardProcess(), reco::GenStatusFlags::setFromHardProcessBeforeFSR(), reco::GenStatusFlags::setIsDecayedLeptonHadron(), reco::GenStatusFlags::setIsDirectHadronDecayProduct(), reco::GenStatusFlags::setIsDirectHardProcessTauDecayProduct(), reco::GenStatusFlags::setIsDirectPromptTauDecayProduct(), reco::GenStatusFlags::setIsDirectTauDecayProduct(), reco::GenStatusFlags::setIsFirstCopy(), reco::GenStatusFlags::setIsHardProcess(), reco::GenStatusFlags::setIsHardProcessTauDecayProduct(), reco::GenStatusFlags::setIsLastCopy(), reco::GenStatusFlags::setIsLastCopyBeforeFSR(), reco::GenStatusFlags::setIsPrompt(), reco::GenStatusFlags::setIsPromptTauDecayProduct(), and reco::GenStatusFlags::setIsTauDecayProduct().

Referenced by GenParticleProducer::convertParticle().

636  {
637  statusFlags.setIsPrompt(isPrompt(p));
644  statusFlags.setIsHardProcess(isHardProcess(p));
645  statusFlags.setFromHardProcess(fromHardProcess(p));
649  statusFlags.setIsFirstCopy(isFirstCopy(p));
650  statusFlags.setIsLastCopy(isLastCopy(p));
652 }
bool isTauDecayProduct(const P &p)
bool fromHardProcess(const P &p)
bool isDirectHadronDecayProduct(const P &p)
bool isDirectPromptTauDecayProduct(const P &p)
void setIsLastCopyBeforeFSR(bool b)
void setIsPrompt(bool b)
void setIsTauDecayProduct(bool b)
void setIsDirectPromptTauDecayProduct(bool b)
bool isHardProcess(const P &p)
void setIsDirectHadronDecayProduct(bool b)
bool isPrompt(const P &p)
void setIsDirectHardProcessTauDecayProduct(bool b)
void setIsHardProcess(bool b)
bool isLastCopyBeforeFSR(const P &p)
bool isDirectHardProcessTauDecayProduct(const P &p)
bool isFirstCopy(const P &p)
void setIsHardProcessTauDecayProduct(bool b)
bool isHardProcessTauDecayProduct(const P &p)
void setFromHardProcess(bool b)
bool isDirectTauDecayProduct(const P &p)
void setIsPromptTauDecayProduct(bool b)
bool isPromptTauDecayProduct(const P &p)
void setIsDecayedLeptonHadron(bool b)
void setIsFirstCopy(bool b)
bool isDecayedLeptonHadron(const P &p)
bool fromHardProcessBeforeFSR(const P &p)
void setIsDirectTauDecayProduct(bool b)
void setFromHardProcessBeforeFSR(bool b)
bool isLastCopy(const P &p)
void setIsLastCopy(bool b)
template<typename P>
const P * MCTruthHelper< P >::findDecayedMother ( const P p)

Definition at line 538 of file MCTruthHelper.h.

References P.

538  {
539  const P *mo = mother(p);
540  dupCheck_.clear();
541  while (mo && !isDecayedLeptonHadron(*mo)) {
542  dupCheck_.insert(mo);
543  mo = mother(*mo);
544  if (dupCheck_.count(mo)) return 0;
545  }
546  return mo;
547 }
#define P
std::unordered_set< const P * > dupCheck_
bool isDecayedLeptonHadron(const P &p)
const reco::GenParticle * mother(const reco::GenParticle &p, unsigned int imoth=0)
template<typename P>
const P * MCTruthHelper< P >::findDecayedMother ( const P p,
int  abspdgid 
)

Definition at line 551 of file MCTruthHelper.h.

References P.

551  {
552  const P *mo = mother(p);
553  dupCheck_.clear();
554  while (mo && (absPdgId(*mo)!=abspdgid || !isDecayedLeptonHadron(*mo)) ) {
555  dupCheck_.insert(mo);
556  mo = mother(*mo);
557  if (dupCheck_.count(mo)) return 0;
558  }
559  return mo;
560 }
#define P
int absPdgId(const reco::GenParticle &p)
std::unordered_set< const P * > dupCheck_
bool isDecayedLeptonHadron(const P &p)
const reco::GenParticle * mother(const reco::GenParticle &p, unsigned int imoth=0)
template<typename P>
const P * MCTruthHelper< P >::firstCopy ( const P p)

Definition at line 397 of file MCTruthHelper.h.

References P, and AlCaHLTBitMon_ParallelJobs::p.

397  {
398  const P *pcopy = &p;
399  dupCheck_.clear();
400  while (previousCopy(*pcopy)) {
401  dupCheck_.insert(pcopy);
402  pcopy = previousCopy(*pcopy);
403  if (dupCheck_.count(pcopy)) return 0;
404  }
405  return pcopy;
406 }
#define P
const P * previousCopy(const P &p)
std::unordered_set< const P * > dupCheck_
template<typename P>
bool MCTruthHelper< P >::fromHardProcess ( const P p)

Definition at line 306 of file MCTruthHelper.h.

306  {
307  return hardProcessMotherCopy(p) != 0;
308 }
const P * hardProcessMotherCopy(const P &p)
template<typename P>
bool MCTruthHelper< P >::fromHardProcessBeforeFSR ( const P p)

Definition at line 338 of file MCTruthHelper.h.

References P, and AlCaHLTBitMon_ParallelJobs::p.

338  {
339  //pythia 6 documentation line roughly corresponds to this condition
340  if (p.status()==3) return true;
341 
342  //check hard process mother properties
343  const P *hpc = hardProcessMotherCopy(p);
344  if (!hpc) return false;
345 
346  //for incoming partons in pythia8, more useful information is not
347  //easily available, so take only the incoming parton itself
348  if (hpc->status()==21 && (&p)==hpc) return true;
349 
350  //for intermediate particles in pythia 8, just take the last copy
351  if (hpc->status()==22 && isLastCopy(p)) return true;
352 
353  //for outgoing particles in pythia 8, explicitly find the last copy
354  //before FSR starting from the hardProcess particle, and take only
355  //this one
356  if ( (hpc->status()==23 || hpc->status()==1) && (&p)==lastDaughterCopyBeforeFSR(*hpc) ) return true;
357 
358 
359  //didn't satisfy any of the conditions
360  return false;
361 
362 }
const P * lastDaughterCopyBeforeFSR(const P &p)
const P * hardProcessMotherCopy(const P &p)
#define P
bool isLastCopy(const P &p)
template<typename P>
bool MCTruthHelper< P >::fromHardProcessDecayed ( const P p)

Definition at line 318 of file MCTruthHelper.h.

318  {
320 }
bool fromHardProcess(const P &p)
bool isDecayedLeptonHadron(const P &p)
template<typename P>
bool MCTruthHelper< P >::fromHardProcessFinalState ( const P p)

Definition at line 312 of file MCTruthHelper.h.

312  {
313  return p.status()==1 && fromHardProcess(p);
314 }
bool fromHardProcess(const P &p)
template<typename P>
const P * MCTruthHelper< P >::hardProcessMotherCopy ( const P p)

Definition at line 490 of file MCTruthHelper.h.

References P, and AlCaHLTBitMon_ParallelJobs::p.

490  {
491  //is particle itself is hard process particle
492  if (isHardProcess(p)) return &p;
493 
494  //check if any other copies are hard process particles
495  const P *pcopy = &p;
496  dupCheck_.clear();
497  while (previousCopy(*pcopy)) {
498  dupCheck_.insert(pcopy);
499  pcopy = previousCopy(*pcopy);
500  if (isHardProcess(*pcopy)) return pcopy;
501  if (dupCheck_.count(pcopy)) break;
502  }
503  return 0;
504 }
#define P
bool isHardProcess(const P &p)
const P * previousCopy(const P &p)
std::unordered_set< const P * > dupCheck_
template<typename P>
bool MCTruthHelper< P >::isDecayedLeptonHadron ( const P p)

Definition at line 208 of file MCTruthHelper.h.

208  {
209  return p.status()==2 && (isHadron(p) || absPdgId(p)==13 || absPdgId(p)==15) && isLastCopy(p);
210 }
bool isHadron(const P &p)
int absPdgId(const reco::GenParticle &p)
bool isLastCopy(const P &p)
template<typename P>
bool MCTruthHelper< P >::isDirectHadronDecayProduct ( const P p)

Definition at line 261 of file MCTruthHelper.h.

References P.

261  {
262  const P *um = uniqueMother(p);
263  return um && isHadron(*um) && isDecayedLeptonHadron(*um);
264 }
#define P
const P * uniqueMother(const P &p)
bool isHadron(const P &p)
bool isDecayedLeptonHadron(const P &p)
template<typename P>
bool MCTruthHelper< P >::isDirectHardProcessTauDecayProduct ( const P p)

Definition at line 331 of file MCTruthHelper.h.

References symbols::dm, P, and metsig::tau.

331  {
332  const P *tau = findDecayedMother(p,15);
333  const P *dm = findDecayedMother(p);
334  return tau && tau==dm && fromHardProcess(*tau);
335 }
bool fromHardProcess(const P &p)
#define P
const P * findDecayedMother(const P &p)
tuple dm
Definition: symbols.py:65
template<typename P>
bool MCTruthHelper< P >::isDirectPromptTauDecayProduct ( const P p)

Definition at line 240 of file MCTruthHelper.h.

References symbols::dm, P, and metsig::tau.

240  {
241  const P *tau = findDecayedMother(p,15);
242  const P *dm = findDecayedMother(p);
243  return tau && tau==dm && isPrompt(*tau);
244 }
#define P
bool isPrompt(const P &p)
const P * findDecayedMother(const P &p)
tuple dm
Definition: symbols.py:65
template<typename P>
bool MCTruthHelper< P >::isDirectTauDecayProduct ( const P p)

Definition at line 232 of file MCTruthHelper.h.

References symbols::dm, P, and metsig::tau.

232  {
233  const P *tau = findDecayedMother(p,15);
234  const P *dm = findDecayedMother(p);
235  return tau && tau==dm;
236 }
#define P
const P * findDecayedMother(const P &p)
tuple dm
Definition: symbols.py:65
template<typename P>
bool MCTruthHelper< P >::isFirstCopy ( const P p)

Definition at line 366 of file MCTruthHelper.h.

366  {
367  return &p == firstCopy(p);
368 }
const P * firstCopy(const P &p)
template<typename P>
bool MCTruthHelper< P >::isHadron ( const P p)

Definition at line 268 of file MCTruthHelper.h.

References benchmark_cfg::pdgId.

268  {
269  HepPDT::ParticleID heppdtid(pdgId(p));
270  return heppdtid.isHadron();
271 }
int pdgId(const reco::GenParticle &p)
template<typename P>
bool MCTruthHelper< P >::isHardProcess ( const P p)

Definition at line 275 of file MCTruthHelper.h.

References P.

275  {
276 
277  //status 3 in pythia6 means hard process;
278  if (p.status()==3) return true;
279 
280  //hard process codes for pythia8 are 21-29 inclusive (currently 21,22,23,24 are used)
281  if (p.status()>20 && p.status()<30) return true;
282 
283  //if this is a final state or decayed particle,
284  //check if direct mother is a resonance decay in pythia8 but exclude FSR branchings
285  //(In pythia8 if a resonance decay product did not undergo any further branchings
286  //it will be directly stored as status 1 or 2 without any status 23 copy)
287  if (p.status()==1 || p.status()==2) {
288  const P *um = mother(p);
289  if (um) {
290  const P *firstcopy = firstCopy(*um);
291  bool fromResonance = firstcopy && firstcopy->status()==22;
292 
293  const P *umNext = nextCopy(*um);
294  bool fsrBranching = umNext && umNext->status()>50 && umNext->status()<60;
295 
296  if (fromResonance && !fsrBranching) return true;
297  }
298  }
299 
300  return false;
301 
302 }
#define P
const P * firstCopy(const P &p)
const P * nextCopy(const P &p)
const reco::GenParticle * mother(const reco::GenParticle &p, unsigned int imoth=0)
template<typename P>
bool MCTruthHelper< P >::isHardProcessTauDecayProduct ( const P p)

Definition at line 324 of file MCTruthHelper.h.

References P, and metsig::tau.

324  {
325  const P *tau = findDecayedMother(p,15);
326  return tau && fromHardProcessDecayed(*tau);
327 }
#define P
const P * findDecayedMother(const P &p)
bool fromHardProcessDecayed(const P &p)
template<typename P>
bool MCTruthHelper< P >::isLastCopy ( const P p)

Definition at line 372 of file MCTruthHelper.h.

372  {
373  return &p == lastCopy(p);
374 }
const P * lastCopy(const P &p)
template<typename P>
bool MCTruthHelper< P >::isLastCopyBeforeFSR ( const P p)

Definition at line 378 of file MCTruthHelper.h.

378  {
379  return &p == lastCopyBeforeFSR(p);
380 }
const P * lastCopyBeforeFSR(const P &p)
template<typename P>
bool MCTruthHelper< P >::isMuonDecayProduct ( const P p)

Definition at line 248 of file MCTruthHelper.h.

248  {
249  return findDecayedMother(p,13) != 0;
250 }
const P * findDecayedMother(const P &p)
template<typename P>
bool MCTruthHelper< P >::isPrompt ( const P p)

Definition at line 194 of file MCTruthHelper.h.

194  {
195  //particle from hadron/muon/tau decay -> not prompt
196  //checking all the way up the chain treats properly the radiated photon
197  //case as well
198  return findDecayedMother(p) == 0;
199 }
const P * findDecayedMother(const P &p)
template<typename P>
bool MCTruthHelper< P >::isPromptDecayed ( const P p)

Definition at line 213 of file MCTruthHelper.h.

213  {
214  return isDecayedLeptonHadron(p) && isPrompt(p);
215 }
bool isPrompt(const P &p)
bool isDecayedLeptonHadron(const P &p)
template<typename P>
bool MCTruthHelper< P >::isPromptFinalState ( const P p)

Definition at line 203 of file MCTruthHelper.h.

203  {
204  return p.status()==1 && isPrompt(p);
205 }
bool isPrompt(const P &p)
template<typename P>
bool MCTruthHelper< P >::isPromptMuonDecayProduct ( const P p)

Definition at line 254 of file MCTruthHelper.h.

References RPCpg::mu, and P.

254  {
255  const P *mu = findDecayedMother(p,13);
256  return mu && isPrompt(*mu);
257 }
#define P
bool isPrompt(const P &p)
const P * findDecayedMother(const P &p)
const int mu
Definition: Constants.h:22
template<typename P>
bool MCTruthHelper< P >::isPromptTauDecayProduct ( const P p)

Definition at line 225 of file MCTruthHelper.h.

References P, and metsig::tau.

225  {
226  const P *tau = findDecayedMother(p,15);
227  return tau && isPrompt(*tau);
228 }
#define P
bool isPrompt(const P &p)
const P * findDecayedMother(const P &p)
template<typename P>
bool MCTruthHelper< P >::isTauDecayProduct ( const P p)

Definition at line 219 of file MCTruthHelper.h.

219  {
220  return findDecayedMother(p,15) != 0;
221 }
const P * findDecayedMother(const P &p)
template<typename P>
const P * MCTruthHelper< P >::lastCopy ( const P p)

Definition at line 410 of file MCTruthHelper.h.

References P, and AlCaHLTBitMon_ParallelJobs::p.

410  {
411  const P *pcopy = &p;
412  dupCheck_.clear();
413  while (nextCopy(*pcopy)) {
414  dupCheck_.insert(pcopy);
415  pcopy = nextCopy(*pcopy);
416  if (dupCheck_.count(pcopy)) return 0;
417  }
418  return pcopy;
419 }
#define P
const P * nextCopy(const P &p)
std::unordered_set< const P * > dupCheck_
template<typename P>
const P * MCTruthHelper< P >::lastCopyBeforeFSR ( const P p)

Definition at line 423 of file MCTruthHelper.h.

References P, and benchmark_cfg::pdgId.

423  {
424  //start with first copy and then walk down until there is FSR
425  const P *pcopy = firstCopy(p);
426  if (!pcopy) return 0;
427  bool hasDaughterCopy = true;
428  dupCheck_.clear();
429  while (hasDaughterCopy) {
430  dupCheck_.insert(pcopy);
431  hasDaughterCopy = false;
432  const unsigned int ndau = numberOfDaughters(*pcopy);
433  //look for FSR
434  for (unsigned int idau = 0; idau<ndau; ++idau) {
435  const P *dau = daughter(*pcopy,idau);
436  if (pdgId(*dau)==21 || pdgId(*dau)==22) {
437  //has fsr (or else decayed and is the last copy by construction)
438  return pcopy;
439  }
440  }
441  //look for daughter copy
442  for (unsigned int idau = 0; idau<ndau; ++idau) {
443  const P *dau = daughter(*pcopy,idau);
444  if (pdgId(*dau)==pdgId(p)) {
445  pcopy = dau;
446  hasDaughterCopy = true;
447  break;
448  }
449  }
450  if (dupCheck_.count(pcopy)) return 0;
451  }
452  return pcopy;
453 }
#define P
const reco::GenParticle * daughter(const reco::GenParticle &p, unsigned int idau)
unsigned int numberOfDaughters(const reco::GenParticle &p)
const P * firstCopy(const P &p)
std::unordered_set< const P * > dupCheck_
int pdgId(const reco::GenParticle &p)
template<typename P>
const P * MCTruthHelper< P >::lastDaughterCopyBeforeFSR ( const P p)

Definition at line 457 of file MCTruthHelper.h.

References P, AlCaHLTBitMon_ParallelJobs::p, and benchmark_cfg::pdgId.

457  {
458  //start with this particle and then walk down until there is FSR
459  const P *pcopy = &p;
460  bool hasDaughterCopy = true;
461  dupCheck_.clear();
462  while (hasDaughterCopy) {
463  dupCheck_.insert(pcopy);
464  hasDaughterCopy = false;
465  const unsigned int ndau = numberOfDaughters(*pcopy);
466  //look for FSR
467  for (unsigned int idau = 0; idau<ndau; ++idau) {
468  const P *dau = daughter(*pcopy,idau);
469  if (pdgId(*dau)==21 || pdgId(*dau)==22) {
470  //has fsr (or else decayed and is the last copy by construction)
471  return pcopy;
472  }
473  }
474  //look for daughter copy
475  for (unsigned int idau = 0; idau<ndau; ++idau) {
476  const P *dau = daughter(*pcopy,idau);
477  if (pdgId(*dau)==pdgId(p)) {
478  pcopy = dau;
479  hasDaughterCopy = true;
480  break;
481  }
482  }
483  if (dupCheck_.count(pcopy)) return 0;
484  }
485  return pcopy;
486 }
#define P
const reco::GenParticle * daughter(const reco::GenParticle &p, unsigned int idau)
unsigned int numberOfDaughters(const reco::GenParticle &p)
std::unordered_set< const P * > dupCheck_
int pdgId(const reco::GenParticle &p)
template<typename P >
const reco::GenParticle * MCTruthHelper< P >::mother ( const reco::GenParticle p,
unsigned int  imoth = 0 
)

Definition at line 600 of file MCTruthHelper.h.

References reco::CompositeRefCandidateT< D >::mother().

600  {
601  return static_cast<const reco::GenParticle*>(p.mother(imoth));
602 }
virtual const Candidate * mother(size_type=0) const
return mother at a given position, i = 0, ... numberOfMothers() - 1 (read only mode) ...
template<typename P >
const HepMC::GenParticle * MCTruthHelper< P >::mother ( const HepMC::GenParticle &  p,
unsigned int  imoth = 0 
)

Definition at line 606 of file MCTruthHelper.h.

606  {
607  return p.production_vertex() && p.production_vertex()->particles_in_size() ? *(p.production_vertex()->particles_in_const_begin() + imoth) : 0;
608 }
template<typename P>
const P * MCTruthHelper< P >::nextCopy ( const P p)

Definition at line 523 of file MCTruthHelper.h.

References P, and benchmark_cfg::pdgId.

523  {
524 
525  const unsigned int ndau = numberOfDaughters(p);
526  for (unsigned int idau = 0; idau<ndau; ++idau) {
527  const P *dau = daughter(p,idau);
528  if (pdgId(*dau)==pdgId(p)) {
529  return dau;
530  }
531  }
532 
533  return 0;
534 }
#define P
const reco::GenParticle * daughter(const reco::GenParticle &p, unsigned int idau)
unsigned int numberOfDaughters(const reco::GenParticle &p)
int pdgId(const reco::GenParticle &p)
template<typename P >
unsigned int MCTruthHelper< P >::numberOfDaughters ( const reco::GenParticle p)

Definition at line 612 of file MCTruthHelper.h.

References reco::CompositeRefCandidateT< D >::numberOfDaughters().

612  {
613  return p.numberOfDaughters();
614 }
virtual size_t numberOfDaughters() const
number of daughters
template<typename P >
unsigned int MCTruthHelper< P >::numberOfDaughters ( const HepMC::GenParticle &  p)

Definition at line 618 of file MCTruthHelper.h.

618  {
619  return p.end_vertex() ? p.end_vertex()->particles_out_size() : 0;
620 }
template<typename P >
unsigned int MCTruthHelper< P >::numberOfMothers ( const reco::GenParticle p)

Definition at line 588 of file MCTruthHelper.h.

References reco::CompositeRefCandidateT< D >::numberOfMothers().

588  {
589  return p.numberOfMothers();
590 }
virtual size_t numberOfMothers() const
number of mothers
template<typename P >
unsigned int MCTruthHelper< P >::numberOfMothers ( const HepMC::GenParticle &  p)

Definition at line 594 of file MCTruthHelper.h.

594  {
595  return p.production_vertex() ? p.production_vertex()->particles_in_size() : 0;
596 }
template<typename P >
int MCTruthHelper< P >::pdgId ( const reco::GenParticle p)

Definition at line 564 of file MCTruthHelper.h.

References reco::LeafCandidate::pdgId().

Referenced by Particle.Particle::__str__().

564  {
565  return p.pdgId();
566 }
virtual int pdgId() const final
PDG identifier.
template<typename P >
int MCTruthHelper< P >::pdgId ( const HepMC::GenParticle &  p)

Definition at line 570 of file MCTruthHelper.h.

Referenced by Particle.Particle::__str__().

570  {
571  return p.pdg_id();
572 }
template<typename P>
const P * MCTruthHelper< P >::previousCopy ( const P p)

Definition at line 508 of file MCTruthHelper.h.

References P, and benchmark_cfg::pdgId.

508  {
509 
510  const unsigned int nmoth = numberOfMothers(p);
511  for (unsigned int imoth = 0; imoth<nmoth; ++imoth) {
512  const P *moth = mother(p,imoth);
513  if (pdgId(*moth)==pdgId(p)) {
514  return moth;
515  }
516  }
517 
518  return 0;
519 }
#define P
unsigned int numberOfMothers(const reco::GenParticle &p)
int pdgId(const reco::GenParticle &p)
const reco::GenParticle * mother(const reco::GenParticle &p, unsigned int imoth=0)
template<typename P>
const P * MCTruthHelper< P >::uniqueMother ( const P p)

Definition at line 384 of file MCTruthHelper.h.

References P, AlCaHLTBitMon_ParallelJobs::p, and benchmark_cfg::pdgId.

384  {
385  const P *mo = &p;
386  dupCheck_.clear();
387  while (mo && pdgId(*mo)==pdgId(p)) {
388  dupCheck_.insert(mo);
389  mo = mother(*mo);
390  if (dupCheck_.count(mo)) return 0;
391  }
392  return mo;
393 }
#define P
std::unordered_set< const P * > dupCheck_
int pdgId(const reco::GenParticle &p)
const reco::GenParticle * mother(const reco::GenParticle &p, unsigned int imoth=0)

Member Data Documentation

template<typename P>
std::unordered_set<const P*> MCTruthHelper< P >::dupCheck_
protected

Definition at line 183 of file MCTruthHelper.h.