1 #ifndef PhysicsTools_HepMCCandAlgos_GenParticlesHelper_h 2 #define PhysicsTools_HepMCCandAlgos_GenParticlesHelper_h 5 #include "HepPDT/ParticleID.hh" 10 #include <unordered_set> 185 template <
typename P>
190 return findDecayedMother(
p) ==
nullptr;
194 template <
typename P>
196 return p.status() == 1 && isPrompt(
p);
199 template <
typename P>
201 return p.status() == 2 && (isHadron(
p) || absPdgId(
p) == 13 || absPdgId(
p) == 15) && isLastCopy(
p);
204 template <
typename P>
206 return isDecayedLeptonHadron(
p) && isPrompt(
p);
210 template <
typename P>
212 return findDecayedMother(
p, 15) !=
nullptr;
216 template <
typename P>
218 const P *
tau = findDecayedMother(
p, 15);
219 return tau && isPrompt(*
tau);
223 template <
typename P>
225 const P *
tau = findDecayedMother(
p, 15);
226 const P *
dm = findDecayedMother(
p);
231 template <
typename P>
233 const P *
tau = findDecayedMother(
p, 15);
234 const P *
dm = findDecayedMother(
p);
239 template <
typename P>
241 return findDecayedMother(
p, 13) != 0;
245 template <
typename P>
247 const P *
mu = findDecayedMother(
p, 13);
248 return mu && isPrompt(*
mu);
252 template <
typename P>
254 const P *um = uniqueMother(
p);
255 return um && isHadron(*um) && isDecayedLeptonHadron(*um);
259 template <
typename P>
262 return heppdtid.isHadron();
266 template <
typename P>
273 if (
p.status() > 20 &&
p.status() < 30)
280 if (
p.status() == 1 ||
p.status() == 2) {
281 const P *um = mother(
p);
283 const P *firstcopy = firstCopy(*um);
284 bool fromResonance = firstcopy && firstcopy->status() == 22;
286 const P *umNext = nextCopy(*um);
287 bool fsrBranching = umNext && umNext->status() > 50 && umNext->status() < 60;
289 if (fromResonance && !fsrBranching)
298 template <
typename P>
300 return hardProcessMotherCopy(
p) !=
nullptr;
304 template <
typename P>
306 return p.status() == 1 && fromHardProcess(
p);
310 template <
typename P>
312 return isDecayedLeptonHadron(
p) && fromHardProcess(
p);
316 template <
typename P>
318 const P *
tau = findDecayedMother(
p, 15);
319 return tau && fromHardProcessDecayed(*
tau);
323 template <
typename P>
325 const P *
tau = findDecayedMother(
p, 15);
326 const P *
dm = findDecayedMother(
p);
330 template <
typename P>
337 const P *hpc = hardProcessMotherCopy(
p);
343 if (hpc->status() == 21 && (&
p) == hpc)
347 if (hpc->status() == 22 && isLastCopy(
p))
353 if ((hpc->status() == 23 || hpc->status() == 1) && (&
p) == lastDaughterCopyBeforeFSR(*hpc))
361 template <
typename P>
363 return &
p == firstCopy(
p);
367 template <
typename P>
369 return &
p == lastCopy(
p);
373 template <
typename P>
375 return &
p == lastCopyBeforeFSR(
p);
379 template <
typename P>
382 std::unordered_set<const P *> dupCheck;
386 if (dupCheck.count(mo))
393 template <
typename 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))
407 template <
typename 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))
421 template <
typename P>
424 const P *pcopy = firstCopy(
p);
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);
434 for (
unsigned int idau = 0; idau < ndau; ++idau) {
435 const P *dau = daughter(*pcopy, idau);
442 for (
unsigned int idau = 0; idau < ndau; ++idau) {
443 const P *dau = daughter(*pcopy, idau);
446 hasDaughterCopy =
true;
450 if (dupCheck.count(pcopy))
457 template <
typename 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);
468 for (
unsigned int idau = 0; idau < ndau; ++idau) {
469 const P *dau = daughter(*pcopy, idau);
476 for (
unsigned int idau = 0; idau < ndau; ++idau) {
477 const P *dau = daughter(*pcopy, idau);
480 hasDaughterCopy =
true;
484 if (dupCheck.count(pcopy))
491 template <
typename P>
494 if (isHardProcess(
p))
499 std::unordered_set<const P *> dupCheck;
500 while (previousCopy(*pcopy)) {
501 dupCheck.insert(pcopy);
502 pcopy = previousCopy(*pcopy);
503 if (isHardProcess(*pcopy))
505 if (dupCheck.count(pcopy))
512 template <
typename P>
514 const unsigned int nmoth = numberOfMothers(
p);
515 for (
unsigned int imoth = 0; imoth < nmoth; ++imoth) {
516 const P *moth = mother(
p, imoth);
526 template <
typename P>
528 const unsigned int ndau = numberOfDaughters(
p);
529 for (
unsigned int idau = 0; idau < ndau; ++idau) {
530 const P *dau = daughter(
p, idau);
540 template <
typename P>
542 const P *mo = mother(
p);
543 std::unordered_set<const P *> dupCheck;
544 while (mo && !isDecayedLeptonHadron(*mo)) {
547 if (dupCheck.count(mo))
554 template <
typename P>
556 const P *mo = mother(
p);
557 std::unordered_set<const P *> dupCheck;
558 while (mo && (absPdgId(*mo) != abspdgid || !isDecayedLeptonHadron(*mo))) {
561 if (dupCheck.count(mo))
568 template <
typename P>
574 template <
typename P>
580 template <
typename P>
586 template <
typename P>
592 template <
typename P>
594 return p.numberOfMothers();
598 template <
typename P>
600 return p.production_vertex() ?
p.production_vertex()->particles_in_size() : 0;
604 template <
typename P>
610 template <
typename P>
612 return p.production_vertex() &&
p.production_vertex()->particles_in_size()
613 ? *(
p.production_vertex()->particles_in_const_begin() + imoth)
618 template <
typename P>
620 return p.numberOfDaughters();
624 template <
typename P>
626 return p.end_vertex() ?
p.end_vertex()->particles_out_size() : 0;
630 template <
typename P>
636 template <
typename P>
638 return *(
p.end_vertex()->particles_out_const_begin() + idau);
642 template <
typename P>
645 statusFlags.setIsDecayedLeptonHadron(isDecayedLeptonHadron(
p));
647 statusFlags.setIsPromptTauDecayProduct(isPromptTauDecayProduct(
p));
648 statusFlags.setIsDirectTauDecayProduct(isDirectTauDecayProduct(
p));
649 statusFlags.setIsDirectPromptTauDecayProduct(isDirectPromptTauDecayProduct(
p));
650 statusFlags.setIsDirectHadronDecayProduct(isDirectHadronDecayProduct(
p));
653 statusFlags.setIsHardProcessTauDecayProduct(isHardProcessTauDecayProduct(
p));
654 statusFlags.setIsDirectHardProcessTauDecayProduct(isDirectHardProcessTauDecayProduct(
p));
655 statusFlags.setFromHardProcessBeforeFSR(fromHardProcessBeforeFSR(
p));
658 statusFlags.setIsLastCopyBeforeFSR(isLastCopyBeforeFSR(
p));
const reco::GenParticle * mother(const reco::GenParticle &p, unsigned int imoth=0) const
const P * uniqueMother(const P &p) const
int absPdgId(const reco::GenParticle &p) const
bool isHadron(const P &p) const
const P * findDecayedMother(const P &p) const
bool isPromptMuonDecayProduct(const P &p) const
const reco::GenParticle * daughter(const reco::GenParticle &p, unsigned int idau) const
const P * lastCopy(const P &p) const
const P * lastCopyBeforeFSR(const P &p) const
bool isDirectHadronDecayProduct(const P &p) const
const P * hardProcessMotherCopy(const P &p) const
bool isLastCopy(const P &p) const
bool isPrompt(const P &p) const
bool isHardProcess(const P &p) const
bool isPromptFinalState(const P &p) const
bool isPromptDecayed(const P &p) const
const P * nextCopy(const P &p) const
bool isLastCopyBeforeFSR(const P &p) const
bool fromHardProcessDecayed(const P &p) const
bool isHardProcessTauDecayProduct(const P &p) const
Abs< T >::type abs(const T &t)
const P * lastDaughterCopyBeforeFSR(const P &p) const
bool isFirstCopy(const P &p) const
int pdgId(const reco::GenParticle &p) const
bool isPromptTauDecayProduct(const P &p) const
bool isDecayedLeptonHadron(const P &p) const
std::pair< OmniClusterRef, TrackingParticleRef > P
const P * previousCopy(const P &p) const
bool fromHardProcessFinalState(const P &p) const
unsigned int numberOfMothers(const reco::GenParticle &p) const
bool isMuonDecayProduct(const P &p) const
bool fromHardProcessBeforeFSR(const P &p) const
void fillGenStatusFlags(const P &p, reco::GenStatusFlags &statusFlags) const
bool isDirectPromptTauDecayProduct(const P &p) const
bool isDirectTauDecayProduct(const P &p) const
unsigned int numberOfDaughters(const reco::GenParticle &p) const
bool isDirectHardProcessTauDecayProduct(const P &p) const
bool isTauDecayProduct(const P &p) const
bool fromHardProcess(const P &p) const
const P * firstCopy(const P &p) const