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 13 of file MCTruthHelper.h.

Member Function Documentation

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

Definition at line 575 of file MCTruthHelper.h.

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

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

Definition at line 581 of file MCTruthHelper.h.

References funct::abs().

581  {
582  return std::abs(p.pdg_id());
583 }
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 623 of file MCTruthHelper.h.

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

623  {
624  return static_cast<const reco::GenParticle*>(p.daughter(idau));
625 }
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 629 of file MCTruthHelper.h.

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

Definition at line 635 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().

635  {
636  statusFlags.setIsPrompt(isPrompt(p));
643  statusFlags.setIsHardProcess(isHardProcess(p));
644  statusFlags.setFromHardProcess(fromHardProcess(p));
648  statusFlags.setIsFirstCopy(isFirstCopy(p));
649  statusFlags.setIsLastCopy(isLastCopy(p));
651 }
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 537 of file MCTruthHelper.h.

References P.

537  {
538  const P *mo = mother(p);
539  dupCheck_.clear();
540  while (mo && !isDecayedLeptonHadron(*mo)) {
541  dupCheck_.insert(mo);
542  mo = mother(*mo);
543  if (dupCheck_.count(mo)) return 0;
544  }
545  return mo;
546 }
#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 550 of file MCTruthHelper.h.

References P.

550  {
551  const P *mo = mother(p);
552  dupCheck_.clear();
553  while (mo && (absPdgId(*mo)!=abspdgid || !isDecayedLeptonHadron(*mo)) ) {
554  dupCheck_.insert(mo);
555  mo = mother(*mo);
556  if (dupCheck_.count(mo)) return 0;
557  }
558  return mo;
559 }
#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 396 of file MCTruthHelper.h.

References P, and AlCaHLTBitMon_ParallelJobs::p.

396  {
397  const P *pcopy = &p;
398  dupCheck_.clear();
399  while (previousCopy(*pcopy)) {
400  dupCheck_.insert(pcopy);
401  pcopy = previousCopy(*pcopy);
402  if (dupCheck_.count(pcopy)) return 0;
403  }
404  return pcopy;
405 }
#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 305 of file MCTruthHelper.h.

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

Definition at line 337 of file MCTruthHelper.h.

References P, and AlCaHLTBitMon_ParallelJobs::p.

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

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

Definition at line 311 of file MCTruthHelper.h.

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

Definition at line 489 of file MCTruthHelper.h.

References P, and AlCaHLTBitMon_ParallelJobs::p.

489  {
490  //is particle itself is hard process particle
491  if (isHardProcess(p)) return &p;
492 
493  //check if any other copies are hard process particles
494  const P *pcopy = &p;
495  dupCheck_.clear();
496  while (previousCopy(*pcopy)) {
497  dupCheck_.insert(pcopy);
498  pcopy = previousCopy(*pcopy);
499  if (isHardProcess(*pcopy)) return pcopy;
500  if (dupCheck_.count(pcopy)) break;
501  }
502  return 0;
503 }
#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 207 of file MCTruthHelper.h.

207  {
208  return p.status()==2 && (isHadron(p) || absPdgId(p)==13 || absPdgId(p)==15) && isLastCopy(p);
209 }
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 260 of file MCTruthHelper.h.

References P.

260  {
261  const P *um = uniqueMother(p);
262  return um && isHadron(*um) && isDecayedLeptonHadron(*um);
263 }
#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 330 of file MCTruthHelper.h.

References P, and metsig::tau.

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

Definition at line 239 of file MCTruthHelper.h.

References P, and metsig::tau.

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

Definition at line 231 of file MCTruthHelper.h.

References P, and metsig::tau.

231  {
232  const P *tau = findDecayedMother(p,15);
233  const P *dm = findDecayedMother(p);
234  return tau && tau==dm;
235 }
#define P
const P * findDecayedMother(const P &p)
template<typename P>
bool MCTruthHelper< P >::isFirstCopy ( const P p)

Definition at line 365 of file MCTruthHelper.h.

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

Definition at line 267 of file MCTruthHelper.h.

References RecoTau_DiTaus_pt_20-420_cfg::ParticleID, and benchmark_cfg::pdgId.

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

Definition at line 274 of file MCTruthHelper.h.

References P.

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

References P, and metsig::tau.

323  {
324  const P *tau = findDecayedMother(p,15);
325  return tau && fromHardProcessDecayed(*tau);
326 }
#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 371 of file MCTruthHelper.h.

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

Definition at line 377 of file MCTruthHelper.h.

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

Definition at line 247 of file MCTruthHelper.h.

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

Definition at line 193 of file MCTruthHelper.h.

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

Definition at line 212 of file MCTruthHelper.h.

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

Definition at line 202 of file MCTruthHelper.h.

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

Definition at line 253 of file MCTruthHelper.h.

References RPCpg::mu, and P.

253  {
254  const P *mu = findDecayedMother(p,13);
255  return mu && isPrompt(*mu);
256 }
#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 224 of file MCTruthHelper.h.

References P, and metsig::tau.

224  {
225  const P *tau = findDecayedMother(p,15);
226  return tau && isPrompt(*tau);
227 }
#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 218 of file MCTruthHelper.h.

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

Definition at line 409 of file MCTruthHelper.h.

References P, and AlCaHLTBitMon_ParallelJobs::p.

409  {
410  const P *pcopy = &p;
411  dupCheck_.clear();
412  while (nextCopy(*pcopy)) {
413  dupCheck_.insert(pcopy);
414  pcopy = nextCopy(*pcopy);
415  if (dupCheck_.count(pcopy)) return 0;
416  }
417  return pcopy;
418 }
#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 422 of file MCTruthHelper.h.

References P, and benchmark_cfg::pdgId.

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

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

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

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

599  {
600  return static_cast<const reco::GenParticle*>(p.mother(imoth));
601 }
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 605 of file MCTruthHelper.h.

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

Definition at line 522 of file MCTruthHelper.h.

References P, and benchmark_cfg::pdgId.

522  {
523 
524  const unsigned int ndau = numberOfDaughters(p);
525  for (unsigned int idau = 0; idau<ndau; ++idau) {
526  const P *dau = daughter(p,idau);
527  if (pdgId(*dau)==pdgId(p)) {
528  return dau;
529  }
530  }
531 
532  return 0;
533 }
#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 611 of file MCTruthHelper.h.

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

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

Definition at line 617 of file MCTruthHelper.h.

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

Definition at line 587 of file MCTruthHelper.h.

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

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

Definition at line 593 of file MCTruthHelper.h.

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

Definition at line 563 of file MCTruthHelper.h.

References reco::LeafCandidate::pdgId().

Referenced by Particle.Particle::__str__().

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

Definition at line 569 of file MCTruthHelper.h.

Referenced by Particle.Particle::__str__().

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

Definition at line 507 of file MCTruthHelper.h.

References P, and benchmark_cfg::pdgId.

507  {
508 
509  const unsigned int nmoth = numberOfMothers(p);
510  for (unsigned int imoth = 0; imoth<nmoth; ++imoth) {
511  const P *moth = mother(p,imoth);
512  if (pdgId(*moth)==pdgId(p)) {
513  return moth;
514  }
515  }
516 
517  return 0;
518 }
#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 383 of file MCTruthHelper.h.

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

383  {
384  const P *mo = &p;
385  dupCheck_.clear();
386  while (mo && pdgId(*mo)==pdgId(p)) {
387  dupCheck_.insert(mo);
388  mo = mother(*mo);
389  if (dupCheck_.count(mo)) return 0;
390  }
391  return mo;
392 }
#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 182 of file MCTruthHelper.h.