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

Member Function Documentation

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

Definition at line 573 of file MCTruthHelper.h.

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

Referenced by MCTruthHelper< P >::findDecayedMother(), and MCTruthHelper< P >::isDecayedLeptonHadron().

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

Definition at line 579 of file MCTruthHelper.h.

References funct::abs().

579  {
580  return std::abs(p.pdg_id());
581 }
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 
) const

Definition at line 621 of file MCTruthHelper.h.

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

Referenced by MCTruthHelper< P >::lastCopyBeforeFSR(), MCTruthHelper< P >::lastDaughterCopyBeforeFSR(), and MCTruthHelper< P >::nextCopy().

621  {
622  return static_cast<const reco::GenParticle*>(p.daughter(idau));
623 }
const Candidate * daughter(size_type) const override
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 
) const

Definition at line 627 of file MCTruthHelper.h.

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

Definition at line 633 of file MCTruthHelper.h.

References MCTruthHelper< P >::fromHardProcess(), MCTruthHelper< P >::fromHardProcessBeforeFSR(), MCTruthHelper< P >::isDecayedLeptonHadron(), MCTruthHelper< P >::isDirectHadronDecayProduct(), MCTruthHelper< P >::isDirectHardProcessTauDecayProduct(), MCTruthHelper< P >::isDirectPromptTauDecayProduct(), MCTruthHelper< P >::isDirectTauDecayProduct(), MCTruthHelper< P >::isFirstCopy(), MCTruthHelper< P >::isHardProcess(), MCTruthHelper< P >::isHardProcessTauDecayProduct(), MCTruthHelper< P >::isLastCopy(), MCTruthHelper< P >::isLastCopyBeforeFSR(), MCTruthHelper< P >::isPrompt(), MCTruthHelper< P >::isPromptTauDecayProduct(), MCTruthHelper< P >::isTauDecayProduct(), 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().

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

Definition at line 535 of file MCTruthHelper.h.

References MCTruthHelper< P >::isDecayedLeptonHadron(), and MCTruthHelper< P >::mother().

Referenced by MCTruthHelper< P >::isDirectHardProcessTauDecayProduct(), MCTruthHelper< P >::isDirectPromptTauDecayProduct(), MCTruthHelper< P >::isDirectTauDecayProduct(), MCTruthHelper< P >::isHardProcessTauDecayProduct(), MCTruthHelper< P >::isMuonDecayProduct(), MCTruthHelper< P >::isPrompt(), MCTruthHelper< P >::isPromptMuonDecayProduct(), MCTruthHelper< P >::isPromptTauDecayProduct(), and MCTruthHelper< P >::isTauDecayProduct().

535  {
536  const P *mo = mother(p);
537  std::unordered_set<const P*> dupCheck;
538  while (mo && !isDecayedLeptonHadron(*mo)) {
539  dupCheck.insert(mo);
540  mo = mother(*mo);
541  if (dupCheck.count(mo)) return nullptr;
542  }
543  return mo;
544 }
bool isDecayedLeptonHadron(const P &p) const
const reco::GenParticle * mother(const reco::GenParticle &p, unsigned int imoth=0) const
std::pair< OmniClusterRef, TrackingParticleRef > P
template<typename P>
const P * MCTruthHelper< P >::findDecayedMother ( const P p,
int  abspdgid 
) const

Definition at line 548 of file MCTruthHelper.h.

References MCTruthHelper< P >::absPdgId(), MCTruthHelper< P >::isDecayedLeptonHadron(), and MCTruthHelper< P >::mother().

548  {
549  const P *mo = mother(p);
550  std::unordered_set<const P*> dupCheck;
551  while (mo && (absPdgId(*mo)!=abspdgid || !isDecayedLeptonHadron(*mo)) ) {
552  dupCheck.insert(mo);
553  mo = mother(*mo);
554  if (dupCheck.count(mo)) return nullptr;
555  }
556  return mo;
557 }
int absPdgId(const reco::GenParticle &p) const
bool isDecayedLeptonHadron(const P &p) const
const reco::GenParticle * mother(const reco::GenParticle &p, unsigned int imoth=0) const
std::pair< OmniClusterRef, TrackingParticleRef > P
template<typename P>
const P * MCTruthHelper< P >::firstCopy ( const P p) const

Definition at line 394 of file MCTruthHelper.h.

References AlCaHLTBitMon_ParallelJobs::p, and MCTruthHelper< P >::previousCopy().

Referenced by MCTruthHelper< P >::isFirstCopy(), MCTruthHelper< P >::isHardProcess(), and MCTruthHelper< P >::lastCopyBeforeFSR().

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)) return nullptr;
401  }
402  return pcopy;
403 }
const P * previousCopy(const P &p) const
std::pair< OmniClusterRef, TrackingParticleRef > P
template<typename P>
bool MCTruthHelper< P >::fromHardProcess ( const P p) const
template<typename P>
bool MCTruthHelper< P >::fromHardProcessBeforeFSR ( const P p) const

Definition at line 335 of file MCTruthHelper.h.

References MCTruthHelper< P >::hardProcessMotherCopy(), MCTruthHelper< P >::isLastCopy(), MCTruthHelper< P >::lastDaughterCopyBeforeFSR(), and AlCaHLTBitMon_ParallelJobs::p.

Referenced by MCTruthHelper< P >::fillGenStatusFlags().

335  {
336  //pythia 6 documentation line roughly corresponds to this condition
337  if (p.status()==3) return true;
338 
339  //check hard process mother properties
340  const P *hpc = hardProcessMotherCopy(p);
341  if (!hpc) return false;
342 
343  //for incoming partons in pythia8, more useful information is not
344  //easily available, so take only the incoming parton itself
345  if (hpc->status()==21 && (&p)==hpc) return true;
346 
347  //for intermediate particles in pythia 8, just take the last copy
348  if (hpc->status()==22 && isLastCopy(p)) 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) ) return true;
354 
355 
356  //didn't satisfy any of the conditions
357  return false;
358 
359 }
bool isLastCopy(const P &p) const
const P * hardProcessMotherCopy(const P &p) const
const P * lastDaughterCopyBeforeFSR(const P &p) const
std::pair< OmniClusterRef, TrackingParticleRef > P
template<typename P>
bool MCTruthHelper< P >::fromHardProcessDecayed ( const P p) const
template<typename P>
bool MCTruthHelper< P >::fromHardProcessFinalState ( const P p) const

Definition at line 309 of file MCTruthHelper.h.

References MCTruthHelper< P >::fromHardProcess().

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

Definition at line 487 of file MCTruthHelper.h.

References MCTruthHelper< P >::isHardProcess(), AlCaHLTBitMon_ParallelJobs::p, and MCTruthHelper< P >::previousCopy().

Referenced by MCTruthHelper< P >::fromHardProcess(), and MCTruthHelper< P >::fromHardProcessBeforeFSR().

487  {
488  //is particle itself is hard process particle
489  if (isHardProcess(p)) return &p;
490 
491  //check if any other copies are hard process particles
492  const P *pcopy = &p;
493  std::unordered_set<const P*> dupCheck;
494  while (previousCopy(*pcopy)) {
495  dupCheck.insert(pcopy);
496  pcopy = previousCopy(*pcopy);
497  if (isHardProcess(*pcopy)) return pcopy;
498  if (dupCheck.count(pcopy)) break;
499  }
500  return nullptr;
501 }
const P * previousCopy(const P &p) const
bool isHardProcess(const P &p) const
std::pair< OmniClusterRef, TrackingParticleRef > P
template<typename P>
bool MCTruthHelper< P >::isDecayedLeptonHadron ( const P p) const
template<typename P>
bool MCTruthHelper< P >::isDirectHadronDecayProduct ( const P p) const

Definition at line 258 of file MCTruthHelper.h.

References MCTruthHelper< P >::isDecayedLeptonHadron(), MCTruthHelper< P >::isHadron(), and MCTruthHelper< P >::uniqueMother().

Referenced by MCTruthHelper< P >::fillGenStatusFlags().

258  {
259  const P *um = uniqueMother(p);
260  return um && isHadron(*um) && isDecayedLeptonHadron(*um);
261 }
const P * uniqueMother(const P &p) const
bool isDecayedLeptonHadron(const P &p) const
bool isHadron(const P &p) const
std::pair< OmniClusterRef, TrackingParticleRef > P
template<typename P>
bool MCTruthHelper< P >::isDirectHardProcessTauDecayProduct ( const P p) const

Definition at line 328 of file MCTruthHelper.h.

References symbols::dm, MCTruthHelper< P >::findDecayedMother(), MCTruthHelper< P >::fromHardProcess(), and metsig::tau.

Referenced by MCTruthHelper< P >::fillGenStatusFlags().

328  {
329  const P *tau = findDecayedMother(p,15);
330  const P *dm = findDecayedMother(p);
331  return tau && tau==dm && fromHardProcess(*tau);
332 }
bool fromHardProcess(const P &p) const
std::pair< OmniClusterRef, TrackingParticleRef > P
const P * findDecayedMother(const P &p) const
template<typename P>
bool MCTruthHelper< P >::isDirectPromptTauDecayProduct ( const P p) const

Definition at line 237 of file MCTruthHelper.h.

References symbols::dm, MCTruthHelper< P >::findDecayedMother(), MCTruthHelper< P >::isPrompt(), and metsig::tau.

Referenced by MCTruthHelper< P >::fillGenStatusFlags().

237  {
238  const P *tau = findDecayedMother(p,15);
239  const P *dm = findDecayedMother(p);
240  return tau && tau==dm && isPrompt(*tau);
241 }
bool isPrompt(const P &p) const
std::pair< OmniClusterRef, TrackingParticleRef > P
const P * findDecayedMother(const P &p) const
template<typename P>
bool MCTruthHelper< P >::isDirectTauDecayProduct ( const P p) const

Definition at line 229 of file MCTruthHelper.h.

References symbols::dm, MCTruthHelper< P >::findDecayedMother(), and metsig::tau.

Referenced by MCTruthHelper< P >::fillGenStatusFlags().

229  {
230  const P *tau = findDecayedMother(p,15);
231  const P *dm = findDecayedMother(p);
232  return tau && tau==dm;
233 }
std::pair< OmniClusterRef, TrackingParticleRef > P
const P * findDecayedMother(const P &p) const
template<typename P>
bool MCTruthHelper< P >::isFirstCopy ( const P p) const

Definition at line 363 of file MCTruthHelper.h.

References MCTruthHelper< P >::firstCopy().

Referenced by MCTruthHelper< P >::fillGenStatusFlags().

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

Definition at line 272 of file MCTruthHelper.h.

References MCTruthHelper< P >::firstCopy(), MCTruthHelper< P >::mother(), and MCTruthHelper< P >::nextCopy().

Referenced by MCTruthHelper< P >::fillGenStatusFlags(), and MCTruthHelper< P >::hardProcessMotherCopy().

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

Definition at line 321 of file MCTruthHelper.h.

References MCTruthHelper< P >::findDecayedMother(), MCTruthHelper< P >::fromHardProcessDecayed(), and metsig::tau.

Referenced by MCTruthHelper< P >::fillGenStatusFlags().

321  {
322  const P *tau = findDecayedMother(p,15);
323  return tau && fromHardProcessDecayed(*tau);
324 }
bool fromHardProcessDecayed(const P &p) const
std::pair< OmniClusterRef, TrackingParticleRef > P
const P * findDecayedMother(const P &p) const
template<typename P>
bool MCTruthHelper< P >::isLastCopy ( const P p) const
template<typename P>
bool MCTruthHelper< P >::isLastCopyBeforeFSR ( const P p) const

Definition at line 375 of file MCTruthHelper.h.

References MCTruthHelper< P >::lastCopyBeforeFSR().

Referenced by MCTruthHelper< P >::fillGenStatusFlags().

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

Definition at line 245 of file MCTruthHelper.h.

References MCTruthHelper< P >::findDecayedMother().

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

Definition at line 191 of file MCTruthHelper.h.

References MCTruthHelper< P >::findDecayedMother().

Referenced by MCTruthHelper< P >::fillGenStatusFlags(), MCTruthHelper< P >::isDirectPromptTauDecayProduct(), MCTruthHelper< P >::isPromptDecayed(), MCTruthHelper< P >::isPromptFinalState(), MCTruthHelper< P >::isPromptMuonDecayProduct(), and MCTruthHelper< P >::isPromptTauDecayProduct().

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

Definition at line 210 of file MCTruthHelper.h.

References MCTruthHelper< P >::isDecayedLeptonHadron(), and MCTruthHelper< P >::isPrompt().

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

Definition at line 200 of file MCTruthHelper.h.

References MCTruthHelper< P >::isPrompt().

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

Definition at line 251 of file MCTruthHelper.h.

References MCTruthHelper< P >::findDecayedMother(), MCTruthHelper< P >::isPrompt(), and RPCpg::mu.

251  {
252  const P *mu = findDecayedMother(p,13);
253  return mu && isPrompt(*mu);
254 }
bool isPrompt(const P &p) const
const int mu
Definition: Constants.h:22
std::pair< OmniClusterRef, TrackingParticleRef > P
const P * findDecayedMother(const P &p) const
template<typename P>
bool MCTruthHelper< P >::isPromptTauDecayProduct ( const P p) const

Definition at line 222 of file MCTruthHelper.h.

References MCTruthHelper< P >::findDecayedMother(), MCTruthHelper< P >::isPrompt(), and metsig::tau.

Referenced by MCTruthHelper< P >::fillGenStatusFlags().

222  {
223  const P *tau = findDecayedMother(p,15);
224  return tau && isPrompt(*tau);
225 }
bool isPrompt(const P &p) const
std::pair< OmniClusterRef, TrackingParticleRef > P
const P * findDecayedMother(const P &p) const
template<typename P>
bool MCTruthHelper< P >::isTauDecayProduct ( const P p) const

Definition at line 216 of file MCTruthHelper.h.

References MCTruthHelper< P >::findDecayedMother().

Referenced by MCTruthHelper< P >::fillGenStatusFlags().

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

Definition at line 407 of file MCTruthHelper.h.

References MCTruthHelper< P >::nextCopy(), and AlCaHLTBitMon_ParallelJobs::p.

Referenced by MCTruthHelper< P >::isLastCopy().

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

Definition at line 420 of file MCTruthHelper.h.

References MCTruthHelper< P >::daughter(), MCTruthHelper< P >::firstCopy(), MCTruthHelper< P >::numberOfDaughters(), and MCTruthHelper< P >::pdgId().

Referenced by MCTruthHelper< P >::isLastCopyBeforeFSR().

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

Definition at line 454 of file MCTruthHelper.h.

References MCTruthHelper< P >::daughter(), MCTruthHelper< P >::numberOfDaughters(), AlCaHLTBitMon_ParallelJobs::p, and MCTruthHelper< P >::pdgId().

Referenced by MCTruthHelper< P >::fromHardProcessBeforeFSR().

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

Definition at line 597 of file MCTruthHelper.h.

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

Referenced by MCTruthHelper< P >::findDecayedMother(), MCTruthHelper< P >::isHardProcess(), MCTruthHelper< P >::previousCopy(), and MCTruthHelper< P >::uniqueMother().

597  {
598  return static_cast<const reco::GenParticle*>(p.mother(imoth));
599 }
const Candidate * mother(size_type=0) const override
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 
) const

Definition at line 603 of file MCTruthHelper.h.

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

Definition at line 520 of file MCTruthHelper.h.

References MCTruthHelper< P >::daughter(), MCTruthHelper< P >::numberOfDaughters(), and MCTruthHelper< P >::pdgId().

Referenced by MCTruthHelper< P >::isHardProcess(), and MCTruthHelper< P >::lastCopy().

520  {
521 
522  const unsigned int ndau = numberOfDaughters(p);
523  for (unsigned int idau = 0; idau<ndau; ++idau) {
524  const P *dau = daughter(p,idau);
525  if (pdgId(*dau)==pdgId(p)) {
526  return dau;
527  }
528  }
529 
530  return nullptr;
531 }
unsigned int numberOfDaughters(const reco::GenParticle &p) const
int pdgId(const reco::GenParticle &p) const
std::pair< OmniClusterRef, TrackingParticleRef > P
const reco::GenParticle * daughter(const reco::GenParticle &p, unsigned int idau) const
template<typename P >
unsigned int MCTruthHelper< P >::numberOfDaughters ( const reco::GenParticle p) const
template<typename P >
unsigned int MCTruthHelper< P >::numberOfDaughters ( const HepMC::GenParticle &  p) const

Definition at line 615 of file MCTruthHelper.h.

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

Definition at line 585 of file MCTruthHelper.h.

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

Referenced by MCTruthHelper< P >::previousCopy().

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

Definition at line 591 of file MCTruthHelper.h.

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

Definition at line 567 of file MCTruthHelper.h.

Referenced by Particle.Particle::__str__().

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

Definition at line 505 of file MCTruthHelper.h.

References MCTruthHelper< P >::mother(), MCTruthHelper< P >::numberOfMothers(), and MCTruthHelper< P >::pdgId().

Referenced by MCTruthHelper< P >::firstCopy(), and MCTruthHelper< P >::hardProcessMotherCopy().

505  {
506 
507  const unsigned int nmoth = numberOfMothers(p);
508  for (unsigned int imoth = 0; imoth<nmoth; ++imoth) {
509  const P *moth = mother(p,imoth);
510  if (pdgId(*moth)==pdgId(p)) {
511  return moth;
512  }
513  }
514 
515  return nullptr;
516 }
unsigned int numberOfMothers(const reco::GenParticle &p) const
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
template<typename P>
const P * MCTruthHelper< P >::uniqueMother ( const P p) const

Definition at line 381 of file MCTruthHelper.h.

References MCTruthHelper< P >::mother(), AlCaHLTBitMon_ParallelJobs::p, and MCTruthHelper< P >::pdgId().

Referenced by MCTruthHelper< P >::isDirectHadronDecayProduct().

381  {
382  const P *mo = &p;
383  std::unordered_set<const P*> dupCheck;
384  while (mo && pdgId(*mo)==pdgId(p)) {
385  dupCheck.insert(mo);
386  mo = mother(*mo);
387  if (dupCheck.count(mo)) 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