CMS 3D CMS Logo

List of all members | Public Member Functions
MCTruthHelper< P > Class Template Reference

#include <MCTruthHelper.h>

Public Member Functions

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

Detailed Description

template<typename P>
class MCTruthHelper< P >

Definition at line 13 of file MCTruthHelper.h.

Member Function Documentation

◆ absPdgId() [1/2]

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

Definition at line 581 of file MCTruthHelper.h.

581  {
582  return std::abs(p.pdgId());
583 }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22

◆ absPdgId() [2/2]

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

Definition at line 587 of file MCTruthHelper.h.

587  {
588  return std::abs(p.pdg_id());
589 }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22

◆ daughter() [1/2]

template<typename P >
const reco::GenParticle * MCTruthHelper< P >::daughter ( const reco::GenParticle p,
unsigned int  idau 
) const

Definition at line 631 of file MCTruthHelper.h.

631  {
632  return static_cast<const reco::GenParticle *>(p.daughter(idau));
633 }

◆ daughter() [2/2]

template<typename P >
const HepMC::GenParticle * MCTruthHelper< P >::daughter ( const HepMC::GenParticle &  p,
unsigned int  idau 
) const

Definition at line 637 of file MCTruthHelper.h.

637  {
638  return *(p.end_vertex()->particles_out_const_begin() + idau);
639 }

◆ fillGenStatusFlags()

template<typename P >
void MCTruthHelper< P >::fillGenStatusFlags ( const P p,
reco::GenStatusFlags statusFlags 
) const

Definition at line 643 of file MCTruthHelper.h.

Referenced by GenParticleProducer::convertParticle().

643  {
644  statusFlags.setIsPrompt(isPrompt(p));
645  statusFlags.setIsDecayedLeptonHadron(isDecayedLeptonHadron(p));
646  statusFlags.setIsTauDecayProduct(isTauDecayProduct(p));
647  statusFlags.setIsPromptTauDecayProduct(isPromptTauDecayProduct(p));
648  statusFlags.setIsDirectTauDecayProduct(isDirectTauDecayProduct(p));
649  statusFlags.setIsDirectPromptTauDecayProduct(isDirectPromptTauDecayProduct(p));
650  statusFlags.setIsDirectHadronDecayProduct(isDirectHadronDecayProduct(p));
651  statusFlags.setIsHardProcess(isHardProcess(p));
652  statusFlags.setFromHardProcess(fromHardProcess(p));
653  statusFlags.setIsHardProcessTauDecayProduct(isHardProcessTauDecayProduct(p));
654  statusFlags.setIsDirectHardProcessTauDecayProduct(isDirectHardProcessTauDecayProduct(p));
655  statusFlags.setFromHardProcessBeforeFSR(fromHardProcessBeforeFSR(p));
656  statusFlags.setIsFirstCopy(isFirstCopy(p));
657  statusFlags.setIsLastCopy(isLastCopy(p));
658  statusFlags.setIsLastCopyBeforeFSR(isLastCopyBeforeFSR(p));
659 }
bool isDirectHadronDecayProduct(const P &p) const
bool isLastCopy(const P &p) const
bool isPrompt(const P &p) const
bool isHardProcess(const P &p) const
bool isLastCopyBeforeFSR(const P &p) const
bool isHardProcessTauDecayProduct(const P &p) const
bool isFirstCopy(const P &p) const
bool isPromptTauDecayProduct(const P &p) const
bool isDecayedLeptonHadron(const P &p) const
bool fromHardProcessBeforeFSR(const P &p) const
bool isDirectPromptTauDecayProduct(const P &p) const
bool isDirectTauDecayProduct(const P &p) const
bool isDirectHardProcessTauDecayProduct(const P &p) const
bool isTauDecayProduct(const P &p) const
bool fromHardProcess(const P &p) const

◆ findDecayedMother() [1/2]

template<typename P >
const P * MCTruthHelper< P >::findDecayedMother ( const P p) const

Definition at line 541 of file MCTruthHelper.h.

541  {
542  const P *mo = mother(p);
543  std::unordered_set<const P *> dupCheck;
544  while (mo && !isDecayedLeptonHadron(*mo)) {
545  dupCheck.insert(mo);
546  mo = mother(*mo);
547  if (dupCheck.count(mo))
548  return nullptr;
549  }
550  return mo;
551 }
const reco::GenParticle * mother(const reco::GenParticle &p, unsigned int imoth=0) const
bool isDecayedLeptonHadron(const P &p) const
std::pair< OmniClusterRef, TrackingParticleRef > P

◆ findDecayedMother() [2/2]

template<typename P >
const P * MCTruthHelper< P >::findDecayedMother ( const P p,
int  abspdgid 
) const

Definition at line 555 of file MCTruthHelper.h.

555  {
556  const P *mo = mother(p);
557  std::unordered_set<const P *> dupCheck;
558  while (mo && (absPdgId(*mo) != abspdgid || !isDecayedLeptonHadron(*mo))) {
559  dupCheck.insert(mo);
560  mo = mother(*mo);
561  if (dupCheck.count(mo))
562  return nullptr;
563  }
564  return mo;
565 }
const reco::GenParticle * mother(const reco::GenParticle &p, unsigned int imoth=0) const
int absPdgId(const reco::GenParticle &p) const
bool isDecayedLeptonHadron(const P &p) const
std::pair< OmniClusterRef, TrackingParticleRef > P

◆ firstCopy()

template<typename P >
const P * MCTruthHelper< P >::firstCopy ( const P p) const

Definition at line 394 of file MCTruthHelper.h.

394  {
395  const P *pcopy = &p;
396  std::unordered_set<const P *> dupCheck;
397  while (previousCopy(*pcopy)) {
398  dupCheck.insert(pcopy);
399  pcopy = previousCopy(*pcopy);
400  if (dupCheck.count(pcopy))
401  return nullptr;
402  }
403  return pcopy;
404 }
std::pair< OmniClusterRef, TrackingParticleRef > P
const P * previousCopy(const P &p) const

◆ fromHardProcess()

template<typename P >
bool MCTruthHelper< P >::fromHardProcess ( const P p) const

Definition at line 299 of file MCTruthHelper.h.

299  {
300  return hardProcessMotherCopy(p) != nullptr;
301 }
const P * hardProcessMotherCopy(const P &p) const

◆ fromHardProcessBeforeFSR()

template<typename P >
bool MCTruthHelper< P >::fromHardProcessBeforeFSR ( const P p) const

Definition at line 331 of file MCTruthHelper.h.

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

◆ fromHardProcessDecayed()

template<typename P >
bool MCTruthHelper< P >::fromHardProcessDecayed ( const P p) const

Definition at line 311 of file MCTruthHelper.h.

311  {
313 }
bool isDecayedLeptonHadron(const P &p) const
bool fromHardProcess(const P &p) const

◆ fromHardProcessFinalState()

template<typename P >
bool MCTruthHelper< P >::fromHardProcessFinalState ( const P p) const

Definition at line 305 of file MCTruthHelper.h.

305  {
306  return p.status() == 1 && fromHardProcess(p);
307 }
bool fromHardProcess(const P &p) const

◆ hardProcessMotherCopy()

template<typename P >
const P * MCTruthHelper< P >::hardProcessMotherCopy ( const P p) const

Definition at line 492 of file MCTruthHelper.h.

492  {
493  //is particle itself is hard process particle
494  if (isHardProcess(p))
495  return &p;
496 
497  //check if any other copies are hard process particles
498  const P *pcopy = &p;
499  std::unordered_set<const P *> dupCheck;
500  while (previousCopy(*pcopy)) {
501  dupCheck.insert(pcopy);
502  pcopy = previousCopy(*pcopy);
503  if (isHardProcess(*pcopy))
504  return pcopy;
505  if (dupCheck.count(pcopy))
506  break;
507  }
508  return nullptr;
509 }
bool isHardProcess(const P &p) const
std::pair< OmniClusterRef, TrackingParticleRef > P
const P * previousCopy(const P &p) const

◆ isDecayedLeptonHadron()

template<typename P >
bool MCTruthHelper< P >::isDecayedLeptonHadron ( const P p) const

Definition at line 200 of file MCTruthHelper.h.

200  {
201  return p.status() == 2 && (isHadron(p) || absPdgId(p) == 13 || absPdgId(p) == 15) && isLastCopy(p);
202 }
int absPdgId(const reco::GenParticle &p) const
bool isHadron(const P &p) const
bool isLastCopy(const P &p) const

◆ isDirectHadronDecayProduct()

template<typename P >
bool MCTruthHelper< P >::isDirectHadronDecayProduct ( const P p) const

Definition at line 253 of file MCTruthHelper.h.

253  {
254  const P *um = uniqueMother(p);
255  return um && isHadron(*um) && isDecayedLeptonHadron(*um);
256 }
const P * uniqueMother(const P &p) const
bool isHadron(const P &p) const
bool isDecayedLeptonHadron(const P &p) const
std::pair< OmniClusterRef, TrackingParticleRef > P

◆ isDirectHardProcessTauDecayProduct()

template<typename P >
bool MCTruthHelper< P >::isDirectHardProcessTauDecayProduct ( const P p) const

Definition at line 324 of file MCTruthHelper.h.

324  {
325  const P *tau = findDecayedMother(p, 15);
326  const P *dm = findDecayedMother(p);
327  return tau && tau == dm && fromHardProcess(*tau);
328 }
const P * findDecayedMother(const P &p) const
std::pair< OmniClusterRef, TrackingParticleRef > P
bool fromHardProcess(const P &p) const

◆ isDirectPromptTauDecayProduct()

template<typename P >
bool MCTruthHelper< P >::isDirectPromptTauDecayProduct ( const P p) const

Definition at line 232 of file MCTruthHelper.h.

232  {
233  const P *tau = findDecayedMother(p, 15);
234  const P *dm = findDecayedMother(p);
235  return tau && tau == dm && isPrompt(*tau);
236 }
const P * findDecayedMother(const P &p) const
bool isPrompt(const P &p) const
std::pair< OmniClusterRef, TrackingParticleRef > P

◆ isDirectTauDecayProduct()

template<typename P >
bool MCTruthHelper< P >::isDirectTauDecayProduct ( const P p) const

Definition at line 224 of file MCTruthHelper.h.

224  {
225  const P *tau = findDecayedMother(p, 15);
226  const P *dm = findDecayedMother(p);
227  return tau && tau == dm;
228 }
const P * findDecayedMother(const P &p) const
std::pair< OmniClusterRef, TrackingParticleRef > P

◆ isFirstCopy()

template<typename P >
bool MCTruthHelper< P >::isFirstCopy ( const P p) const

Definition at line 362 of file MCTruthHelper.h.

362  {
363  return &p == firstCopy(p);
364 }
const P * firstCopy(const P &p) const

◆ isHadron()

template<typename P >
bool MCTruthHelper< P >::isHadron ( const P p) const

Definition at line 260 of file MCTruthHelper.h.

260  {
261  HepPDT::ParticleID heppdtid(pdgId(p));
262  return heppdtid.isHadron();
263 }
int pdgId(const reco::GenParticle &p) const

◆ isHardProcess()

template<typename P >
bool MCTruthHelper< P >::isHardProcess ( const P p) const

Definition at line 267 of file MCTruthHelper.h.

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

◆ isHardProcessTauDecayProduct()

template<typename P >
bool MCTruthHelper< P >::isHardProcessTauDecayProduct ( const P p) const

Definition at line 317 of file MCTruthHelper.h.

317  {
318  const P *tau = findDecayedMother(p, 15);
319  return tau && fromHardProcessDecayed(*tau);
320 }
const P * findDecayedMother(const P &p) const
bool fromHardProcessDecayed(const P &p) const
std::pair< OmniClusterRef, TrackingParticleRef > P

◆ isLastCopy()

template<typename P >
bool MCTruthHelper< P >::isLastCopy ( const P p) const

Definition at line 368 of file MCTruthHelper.h.

368  {
369  return &p == lastCopy(p);
370 }
const P * lastCopy(const P &p) const

◆ isLastCopyBeforeFSR()

template<typename P >
bool MCTruthHelper< P >::isLastCopyBeforeFSR ( const P p) const

Definition at line 374 of file MCTruthHelper.h.

374  {
375  return &p == lastCopyBeforeFSR(p);
376 }
const P * lastCopyBeforeFSR(const P &p) const

◆ isMuonDecayProduct()

template<typename P >
bool MCTruthHelper< P >::isMuonDecayProduct ( const P p) const

Definition at line 240 of file MCTruthHelper.h.

240  {
241  return findDecayedMother(p, 13) != 0;
242 }
const P * findDecayedMother(const P &p) const

◆ isPrompt()

template<typename P >
bool MCTruthHelper< P >::isPrompt ( const P p) const

Definition at line 186 of file MCTruthHelper.h.

186  {
187  //particle from hadron/muon/tau decay -> not prompt
188  //checking all the way up the chain treats properly the radiated photon
189  //case as well
190  return findDecayedMother(p) == nullptr;
191 }
const P * findDecayedMother(const P &p) const

◆ isPromptDecayed()

template<typename P >
bool MCTruthHelper< P >::isPromptDecayed ( const P p) const

Definition at line 205 of file MCTruthHelper.h.

205  {
206  return isDecayedLeptonHadron(p) && isPrompt(p);
207 }
bool isPrompt(const P &p) const
bool isDecayedLeptonHadron(const P &p) const

◆ isPromptFinalState()

template<typename P >
bool MCTruthHelper< P >::isPromptFinalState ( const P p) const

Definition at line 195 of file MCTruthHelper.h.

195  {
196  return p.status() == 1 && isPrompt(p);
197 }
bool isPrompt(const P &p) const

◆ isPromptMuonDecayProduct()

template<typename P >
bool MCTruthHelper< P >::isPromptMuonDecayProduct ( const P p) const

Definition at line 246 of file MCTruthHelper.h.

246  {
247  const P *mu = findDecayedMother(p, 13);
248  return mu && isPrompt(*mu);
249 }
const P * findDecayedMother(const P &p) const
bool isPrompt(const P &p) const
std::pair< OmniClusterRef, TrackingParticleRef > P

◆ isPromptTauDecayProduct()

template<typename P >
bool MCTruthHelper< P >::isPromptTauDecayProduct ( const P p) const

Definition at line 217 of file MCTruthHelper.h.

217  {
218  const P *tau = findDecayedMother(p, 15);
219  return tau && isPrompt(*tau);
220 }
const P * findDecayedMother(const P &p) const
bool isPrompt(const P &p) const
std::pair< OmniClusterRef, TrackingParticleRef > P

◆ isTauDecayProduct()

template<typename P >
bool MCTruthHelper< P >::isTauDecayProduct ( const P p) const

Definition at line 211 of file MCTruthHelper.h.

211  {
212  return findDecayedMother(p, 15) != nullptr;
213 }
const P * findDecayedMother(const P &p) const

◆ lastCopy()

template<typename P >
const P * MCTruthHelper< P >::lastCopy ( const P p) const

Definition at line 408 of file MCTruthHelper.h.

408  {
409  const P *pcopy = &p;
410  std::unordered_set<const P *> dupCheck;
411  while (nextCopy(*pcopy)) {
412  dupCheck.insert(pcopy);
413  pcopy = nextCopy(*pcopy);
414  if (dupCheck.count(pcopy))
415  return nullptr;
416  }
417  return pcopy;
418 }
const P * nextCopy(const P &p) const
std::pair< OmniClusterRef, TrackingParticleRef > P

◆ lastCopyBeforeFSR()

template<typename P >
const P * MCTruthHelper< P >::lastCopyBeforeFSR ( const P p) const

Definition at line 422 of file MCTruthHelper.h.

422  {
423  //start with first copy and then walk down until there is FSR
424  const P *pcopy = firstCopy(p);
425  if (!pcopy)
426  return nullptr;
427  bool hasDaughterCopy = true;
428  std::unordered_set<const P *> dupCheck;
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))
451  return nullptr;
452  }
453  return pcopy;
454 }
const reco::GenParticle * daughter(const reco::GenParticle &p, unsigned int idau) const
int pdgId(const reco::GenParticle &p) const
std::pair< OmniClusterRef, TrackingParticleRef > P
unsigned int numberOfDaughters(const reco::GenParticle &p) const
const P * firstCopy(const P &p) const

◆ lastDaughterCopyBeforeFSR()

template<typename P >
const P * MCTruthHelper< P >::lastDaughterCopyBeforeFSR ( const P p) const

Definition at line 458 of file MCTruthHelper.h.

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

◆ mother() [1/2]

template<typename P >
const reco::GenParticle * MCTruthHelper< P >::mother ( const reco::GenParticle p,
unsigned int  imoth = 0 
) const

Definition at line 605 of file MCTruthHelper.h.

605  {
606  return static_cast<const reco::GenParticle *>(p.mother(imoth));
607 }

◆ mother() [2/2]

template<typename P >
const HepMC::GenParticle * MCTruthHelper< P >::mother ( const HepMC::GenParticle &  p,
unsigned int  imoth = 0 
) const

Definition at line 611 of file MCTruthHelper.h.

611  {
612  return p.production_vertex() && p.production_vertex()->particles_in_size()
613  ? *(p.production_vertex()->particles_in_const_begin() + imoth)
614  : nullptr;
615 }

◆ nextCopy()

template<typename P >
const P * MCTruthHelper< P >::nextCopy ( const P p) const

Definition at line 527 of file MCTruthHelper.h.

527  {
528  const unsigned int ndau = numberOfDaughters(p);
529  for (unsigned int idau = 0; idau < ndau; ++idau) {
530  const P *dau = daughter(p, idau);
531  if (pdgId(*dau) == pdgId(p)) {
532  return dau;
533  }
534  }
535 
536  return nullptr;
537 }
const reco::GenParticle * daughter(const reco::GenParticle &p, unsigned int idau) const
int pdgId(const reco::GenParticle &p) const
std::pair< OmniClusterRef, TrackingParticleRef > P
unsigned int numberOfDaughters(const reco::GenParticle &p) const

◆ numberOfDaughters() [1/2]

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

Definition at line 619 of file MCTruthHelper.h.

619  {
620  return p.numberOfDaughters();
621 }

◆ numberOfDaughters() [2/2]

template<typename P >
unsigned int MCTruthHelper< P >::numberOfDaughters ( const HepMC::GenParticle &  p) const

Definition at line 625 of file MCTruthHelper.h.

625  {
626  return p.end_vertex() ? p.end_vertex()->particles_out_size() : 0;
627 }

◆ numberOfMothers() [1/2]

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

Definition at line 593 of file MCTruthHelper.h.

593  {
594  return p.numberOfMothers();
595 }

◆ numberOfMothers() [2/2]

template<typename P >
unsigned int MCTruthHelper< P >::numberOfMothers ( const HepMC::GenParticle &  p) const

Definition at line 599 of file MCTruthHelper.h.

599  {
600  return p.production_vertex() ? p.production_vertex()->particles_in_size() : 0;
601 }

◆ pdgId() [1/2]

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

Definition at line 569 of file MCTruthHelper.h.

Referenced by Particle.Particle::__str__().

569  {
570  return p.pdgId();
571 }

◆ pdgId() [2/2]

template<typename P >
int MCTruthHelper< P >::pdgId ( const HepMC::GenParticle &  p) const

Definition at line 575 of file MCTruthHelper.h.

Referenced by Particle.Particle::__str__().

575  {
576  return p.pdg_id();
577 }

◆ previousCopy()

template<typename P >
const P * MCTruthHelper< P >::previousCopy ( const P p) const

Definition at line 513 of file MCTruthHelper.h.

513  {
514  const unsigned int nmoth = numberOfMothers(p);
515  for (unsigned int imoth = 0; imoth < nmoth; ++imoth) {
516  const P *moth = mother(p, imoth);
517  if (pdgId(*moth) == pdgId(p)) {
518  return moth;
519  }
520  }
521 
522  return nullptr;
523 }
const reco::GenParticle * mother(const reco::GenParticle &p, unsigned int imoth=0) const
int pdgId(const reco::GenParticle &p) const
std::pair< OmniClusterRef, TrackingParticleRef > P
unsigned int numberOfMothers(const reco::GenParticle &p) const

◆ uniqueMother()

template<typename P >
const P * MCTruthHelper< P >::uniqueMother ( const P p) const

Definition at line 380 of file MCTruthHelper.h.

380  {
381  const P *mo = &p;
382  std::unordered_set<const P *> dupCheck;
383  while (mo && pdgId(*mo) == pdgId(p)) {
384  dupCheck.insert(mo);
385  mo = mother(*mo);
386  if (dupCheck.count(mo))
387  return nullptr;
388  }
389  return mo;
390 }
const reco::GenParticle * mother(const reco::GenParticle &p, unsigned int imoth=0) const
int pdgId(const reco::GenParticle &p) const
std::pair< OmniClusterRef, TrackingParticleRef > P