CMS 3D CMS Logo

MCTruthHelper.h
Go to the documentation of this file.
1 #ifndef PhysicsTools_HepMCCandAlgos_GenParticlesHelper_h
2 #define PhysicsTools_HepMCCandAlgos_GenParticlesHelper_h
3 
5 #include "HepPDT/ParticleID.hh"
8 
9 
10 #include <iostream>
11 #include <unordered_set>
12 
13 template<typename P>
15 
16 public:
18  //these are robust, generator-independent functions for categorizing
19  //mainly final state particles, but also intermediate hadrons
20  //or radiating leptons
21 
22  //is particle prompt (not from hadron, muon, or tau decay)
23  bool isPrompt(const P &p) const;
24 
25  //is particle prompt and final state
26  bool isPromptFinalState(const P &p) const;
27 
28  //is particle a decayed hadron, muon, or tau (does not include resonance decays like W,Z,Higgs,top,etc)
29  //This flag is equivalent to status 2 in the current HepMC standard
30  //but older generators (pythia6, herwig6) predate this and use status 2 also for other intermediate
31  //particles/states
32  bool isDecayedLeptonHadron(const P &p) const;
33 
34  //is particle prompt and decayed
35  bool isPromptDecayed(const P &p) const;
36 
37  //this particle is a direct or indirect tau decay product
38  bool isTauDecayProduct(const P &p) const;
39 
40  //this particle is a direct or indirect decay product of a prompt tau
41  bool isPromptTauDecayProduct(const P &p) const;
42 
43  //this particle is a direct tau decay product
44  bool isDirectTauDecayProduct(const P &p) const;
45 
46  //this particle is a direct decay product from a prompt tau
47  bool isDirectPromptTauDecayProduct(const P &p) const;
48 
49  //this particle is a direct or indirect muon decay product
50  bool isMuonDecayProduct(const P &p) const;
51 
52  //this particle is a direct or indirect decay product of a prompt muon
53  bool isPromptMuonDecayProduct(const P &p) const;
54 
55  //this particle is a direct decay product from a hadron
56  bool isDirectHadronDecayProduct(const P &p) const;
57 
58  //is particle a hadron
59  bool isHadron(const P &p) const;
60 
62  //these are generator history-dependent functions for tagging particles
63  //associated with the hard process
64  //Currently implemented for Pythia 6 and Pythia 8 status codes and history
65 
66  //this particle is part of the hard process
67  bool isHardProcess(const P &p) const;
68 
69  //this particle is the direct descendant of a hard process particle of the same pdg id
70  bool fromHardProcess(const P &p) const;
71 
72  //this particle is the final state direct descendant of a hard process particle
73  bool fromHardProcessFinalState(const P &p) const;
74 
75  //this particle is the decayed direct descendant of a hard process particle
76  //such as a tau from the hard process
77  bool fromHardProcessDecayed(const P &p) const;
78 
79  //this particle is a direct or indirect decay product of a tau
80  //from the hard process
81  bool isHardProcessTauDecayProduct(const P &p) const;
82 
83  //this particle is a direct decay product of a tau
84  //from the hard process
85  bool isDirectHardProcessTauDecayProduct(const P &p) const;
86 
87  //this particle is the direct descendant of a hard process particle of the same pdg id
88  //For outgoing particles the kinematics are those before QCD or QED FSR
89  //This corresponds roughly to status code 3 in pythia 6
90  bool fromHardProcessBeforeFSR(const P &p) const;
91 
92  //this particle is the first copy of the particle in the chain with the same pdg id
93  bool isFirstCopy(const P &p) const;
94 
95  //this particle is the last copy of the particle in the chain with the same pdg id
96  //(and therefore is more likely, but not guaranteed, to carry the final physical momentum)
97  bool isLastCopy(const P &p) const;
98 
99  //this particle is the last copy of the particle in the chain with the same pdg id
100  //before QED or QCD FSR
101  //(and therefore is more likely, but not guaranteed, to carry the momentum after ISR)
102  //This flag only really makes sense for outgoing particles
103  bool isLastCopyBeforeFSR(const P &p) const;
104 
106  //These are utility functions used by the above
107 
108  //first mother in chain with a different pdg than the particle
109  const P *uniqueMother(const P &p) const;
110 
111  //return first copy of particle in chain (may be the particle itself)
112  const P *firstCopy(const P &p) const;
113 
114  //return last copy of particle in chain (may be the particle itself)
115  const P *lastCopy(const P &p) const;
116 
117  //return last copy of particle in chain before QED or QCD FSR (may be the particle itself)
118  const P *lastCopyBeforeFSR(const P &p) const;
119 
120  //return last copy of particle in chain before QED or QCD FSR, starting from itself (may be the particle itself)
121  const P *lastDaughterCopyBeforeFSR(const P &p) const;
122 
123  //return mother copy which is a hard process particle
124  const P *hardProcessMotherCopy(const P &p) const;
125 
126  //return previous copy of particle in chain (0 in case this is already the first copy)
127  const P *previousCopy(const P &p) const;
128 
129  //return next copy of particle in chain (0 in case this is already the last copy)
130  const P *nextCopy(const P &p) const;
131 
132  //return decayed mother (walk up the chain until found)
133  const P *findDecayedMother(const P &p) const;
134 
135  //return decayed mother matching requested abs(pdgid) (walk up the chain until found)
136  const P *findDecayedMother(const P &p, int abspdgid) const;
137 
139  //These are very basic utility functions to implement a common interface for reco::GenParticle
140  //and HepMC::GenParticle
141 
142  //pdgid
143  int pdgId(const reco::GenParticle &p) const;
144 
145  //pdgid
146  int pdgId(const HepMC::GenParticle &p) const;
147 
148  //abs(pdgid)
149  int absPdgId(const reco::GenParticle &p) const;
150 
151  //abs(pdgid)
152  int absPdgId(const HepMC::GenParticle &p) const;
153 
154  //number of mothers
155  unsigned int numberOfMothers(const reco::GenParticle &p) const;
156 
157  //number of mothers
158  unsigned int numberOfMothers(const HepMC::GenParticle &p) const;
159 
160  //mother
161  const reco::GenParticle *mother(const reco::GenParticle &p, unsigned int imoth=0) const;
162 
163  //mother
164  const HepMC::GenParticle *mother(const HepMC::GenParticle &p, unsigned int imoth=0) const;
165 
166  //number of daughters
167  unsigned int numberOfDaughters(const reco::GenParticle &p) const;
168 
169  //number of daughters
170  unsigned int numberOfDaughters(const HepMC::GenParticle &p) const;
171 
172  //ith daughter
173  const reco::GenParticle *daughter(const reco::GenParticle &p, unsigned int idau) const;
174 
175  //ith daughter
176  const HepMC::GenParticle *daughter(const HepMC::GenParticle &p, unsigned int idau) const;
177 
179  //Helper function to fill status flags
180  void fillGenStatusFlags(const P &p, reco::GenStatusFlags &statusFlags) const;
181 
182 };
183 
184 
185 
187 //implementations
188 
190 template<typename P>
191 bool MCTruthHelper<P>::isPrompt(const P &p) const {
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) == 0;
196 }
197 
199 template<typename P>
201  return p.status()==1 && isPrompt(p);
202 }
203 
204 template<typename P>
206  return p.status()==2 && (isHadron(p) || absPdgId(p)==13 || absPdgId(p)==15) && isLastCopy(p);
207 }
208 
209 template<typename P>
211  return isDecayedLeptonHadron(p) && isPrompt(p);
212 }
213 
215 template<typename P>
217  return findDecayedMother(p,15) != 0;
218 }
219 
221 template<typename P>
223  const P *tau = findDecayedMother(p,15);
224  return tau && isPrompt(*tau);
225 }
226 
228 template<typename P>
230  const P *tau = findDecayedMother(p,15);
231  const P *dm = findDecayedMother(p);
232  return tau && tau==dm;
233 }
234 
236 template<typename P>
238  const P *tau = findDecayedMother(p,15);
239  const P *dm = findDecayedMother(p);
240  return tau && tau==dm && isPrompt(*tau);
241 }
242 
244 template<typename P>
246  return findDecayedMother(p,13) != 0;
247 }
248 
250 template<typename P>
252  const P *mu = findDecayedMother(p,13);
253  return mu && isPrompt(*mu);
254 }
255 
257 template<typename P>
259  const P *um = uniqueMother(p);
260  return um && isHadron(*um) && isDecayedLeptonHadron(*um);
261 }
262 
264 template<typename P>
265 bool MCTruthHelper<P>::isHadron(const P &p) const {
266  HepPDT::ParticleID heppdtid(pdgId(p));
267  return heppdtid.isHadron();
268 }
269 
271 template<typename P>
272 bool MCTruthHelper<P>::isHardProcess(const P &p) const {
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 }
300 
302 template<typename P>
304  return hardProcessMotherCopy(p) != 0;
305 }
306 
308 template<typename P>
310  return p.status()==1 && fromHardProcess(p);
311 }
312 
314 template<typename P>
316  return isDecayedLeptonHadron(p) && fromHardProcess(p);
317 }
318 
320 template<typename P>
322  const P *tau = findDecayedMother(p,15);
323  return tau && fromHardProcessDecayed(*tau);
324 }
325 
327 template<typename P>
329  const P *tau = findDecayedMother(p,15);
330  const P *dm = findDecayedMother(p);
331  return tau && tau==dm && fromHardProcess(*tau);
332 }
333 
334 template<typename P>
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 }
360 
362 template<typename P>
363 bool MCTruthHelper<P>::isFirstCopy(const P &p) const {
364  return &p == firstCopy(p);
365 }
366 
368 template<typename P>
369 bool MCTruthHelper<P>::isLastCopy(const P &p) const {
370  return &p == lastCopy(p);
371 }
372 
374 template<typename P>
376  return &p == lastCopyBeforeFSR(p);
377 }
378 
380 template<typename P>
381 const P *MCTruthHelper<P>::uniqueMother(const P &p) const {
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 0;
388  }
389  return mo;
390 }
391 
393 template<typename P>
394 const P *MCTruthHelper<P>::firstCopy(const P &p) const {
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 0;
401  }
402  return pcopy;
403 }
404 
406 template<typename P>
407 const P *MCTruthHelper<P>::lastCopy(const P &p) const {
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 0;
414  }
415  return pcopy;
416 }
417 
419 template<typename P>
420 const P *MCTruthHelper<P>::lastCopyBeforeFSR(const P &p) const {
421  //start with first copy and then walk down until there is FSR
422  const P *pcopy = firstCopy(p);
423  if (!pcopy) return 0;
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 0;
448  }
449  return pcopy;
450 }
451 
453 template<typename P>
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 0;
481  }
482  return pcopy;
483 }
484 
486 template<typename P>
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 0;
501 }
502 
504 template<typename P>
505 const P *MCTruthHelper<P>::previousCopy(const P &p) const {
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 0;
516 }
517 
519 template<typename P>
520 const P *MCTruthHelper<P>::nextCopy(const P &p) const {
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 0;
531 }
532 
534 template<typename P>
535 const P *MCTruthHelper<P>::findDecayedMother(const P &p) const {
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 0;
542  }
543  return mo;
544 }
545 
547 template<typename P>
548 const P *MCTruthHelper<P>::findDecayedMother(const P &p, int abspdgid) const {
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 0;
555  }
556  return mo;
557 }
558 
560 template<typename P>
562  return p.pdgId();
563 }
564 
566 template<typename P>
568  return p.pdg_id();
569 }
570 
572 template<typename P>
574  return std::abs(p.pdgId());
575 }
576 
578 template<typename P>
580  return std::abs(p.pdg_id());
581 }
582 
584 template<typename P>
586  return p.numberOfMothers();
587 }
588 
590 template<typename P>
592  return p.production_vertex() ? p.production_vertex()->particles_in_size() : 0;
593 }
594 
596 template<typename P>
597 const reco::GenParticle *MCTruthHelper<P>::mother(const reco::GenParticle &p, unsigned int imoth) const {
598  return static_cast<const reco::GenParticle*>(p.mother(imoth));
599 }
600 
602 template<typename P>
603 const HepMC::GenParticle *MCTruthHelper<P>::mother(const HepMC::GenParticle &p, unsigned int imoth) const {
604  return p.production_vertex() && p.production_vertex()->particles_in_size() ? *(p.production_vertex()->particles_in_const_begin() + imoth) : 0;
605 }
606 
608 template<typename P>
610  return p.numberOfDaughters();
611 }
612 
614 template<typename P>
616  return p.end_vertex() ? p.end_vertex()->particles_out_size() : 0;
617 }
618 
620 template<typename P>
621 const reco::GenParticle *MCTruthHelper<P>::daughter(const reco::GenParticle &p, unsigned int idau) const {
622  return static_cast<const reco::GenParticle*>(p.daughter(idau));
623 }
624 
626 template<typename P>
627 const HepMC::GenParticle *MCTruthHelper<P>::daughter(const HepMC::GenParticle &p, unsigned int idau) const {
628  return *(p.end_vertex()->particles_out_const_begin() + idau);
629 }
630 
632 template<typename P>
634  statusFlags.setIsPrompt(isPrompt(p));
636  statusFlags.setIsTauDecayProduct(isTauDecayProduct(p));
641  statusFlags.setIsHardProcess(isHardProcess(p));
642  statusFlags.setFromHardProcess(fromHardProcess(p));
646  statusFlags.setIsFirstCopy(isFirstCopy(p));
647  statusFlags.setIsLastCopy(isLastCopy(p));
649 }
650 
651 
652 #endif
const P * uniqueMother(const P &p) const
bool isPromptDecayed(const P &p) const
bool isDirectHardProcessTauDecayProduct(const P &p) const
const P * previousCopy(const P &p) const
bool isTauDecayProduct(const P &p) const
bool isLastCopy(const P &p) const
const P * nextCopy(const P &p) const
bool isPromptFinalState(const P &p) const
int absPdgId(const reco::GenParticle &p) const
bool isDecayedLeptonHadron(const P &p) const
bool fromHardProcessDecayed(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)
bool fromHardProcessFinalState(const P &p) const
void setIsTauDecayProduct(bool b)
void setIsDirectPromptTauDecayProduct(bool b)
bool isHardProcess(const P &p) const
void fillGenStatusFlags(const P &p, reco::GenStatusFlags &statusFlags) const
void setIsDirectHadronDecayProduct(bool b)
void setIsDirectHardProcessTauDecayProduct(bool b)
void setIsHardProcess(bool b)
unsigned int numberOfMothers(const reco::GenParticle &p) const
virtual size_t numberOfMothers() const
number of mothers
const P * hardProcessMotherCopy(const P &p) const
virtual size_t numberOfDaughters() const
number of daughters
virtual const Candidate * daughter(size_type) const
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
virtual int pdgId() const final
PDG identifier.
const P * lastDaughterCopyBeforeFSR(const P &p) const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void setIsHardProcessTauDecayProduct(bool b)
const int mu
Definition: Constants.h:22
bool isFirstCopy(const P &p) const
bool isMuonDecayProduct(const P &p) const
bool fromHardProcess(const P &p) const
const reco::GenParticle * mother(const reco::GenParticle &p, unsigned int imoth=0) const
unsigned int numberOfDaughters(const reco::GenParticle &p) const
bool fromHardProcessBeforeFSR(const P &p) const
void setFromHardProcess(bool b)
bool isPromptTauDecayProduct(const P &p) const
bool isPromptMuonDecayProduct(const P &p) const
bool isHardProcessTauDecayProduct(const P &p) const
int pdgId(const reco::GenParticle &p) const
const P * lastCopyBeforeFSR(const P &p) const
void setIsPromptTauDecayProduct(bool b)
void setIsDecayedLeptonHadron(bool b)
bool isHadron(const P &p) const
void setIsFirstCopy(bool b)
std::pair< OmniClusterRef, TrackingParticleRef > P
bool isDirectTauDecayProduct(const P &p) const
const P * findDecayedMother(const P &p) const
void setIsDirectTauDecayProduct(bool b)
void setFromHardProcessBeforeFSR(bool b)
const P * lastCopy(const P &p) const
bool isDirectHadronDecayProduct(const P &p) const
void setIsLastCopy(bool b)
virtual const Candidate * mother(size_type=0) const
return mother at a given position, i = 0, ... numberOfMothers() - 1 (read only mode) ...
const P * firstCopy(const P &p) const
const reco::GenParticle * daughter(const reco::GenParticle &p, unsigned int idau) const