CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 
12 template<typename P>
14 
15 public:
17  //these are robust, generator-independent functions for categorizing
18  //mainly final state particles, but also intermediate hadrons
19  //or radiating leptons
20 
21  //is particle prompt (not from hadron, muon, or tau decay)
22  bool isPrompt(const P &p);
23 
24  //is particle prompt and final state
25  bool isPromptFinalState(const P &p);
26 
27  //is particle a decayed hadron, muon, or tau (does not include resonance decays like W,Z,Higgs,top,etc)
28  //This flag is equivalent to status 2 in the current HepMC standard
29  //but older generators (pythia6, herwig6) predate this and use status 2 also for other intermediate
30  //particles/states
31  bool isDecayedLeptonHadron(const P &p);
32 
33  //is particle prompt and decayed
34  bool isPromptDecayed(const P &p);
35 
36  //this particle is a direct or indirect tau decay product
37  bool isTauDecayProduct(const P &p);
38 
39  //this particle is a direct or indirect decay product of a prompt tau
40  bool isPromptTauDecayProduct(const P &p);
41 
42  //this particle is a direct tau decay product
43  bool isDirectTauDecayProduct(const P &p);
44 
45  //this particle is a direct decay product from a prompt tau
46  bool isDirectPromptTauDecayProduct(const P &p);
47 
48  //this particle is a direct or indirect muon decay product
49  bool isMuonDecayProduct(const P &p);
50 
51  //this particle is a direct or indirect decay product of a prompt muon
52  bool isPromptMuonDecayProduct(const P &p);
53 
54  //this particle is a direct decay product from a hadron
55  bool isDirectHadronDecayProduct(const P &p);
56 
57  //is particle a hadron
58  bool isHadron(const P &p);
59 
61  //these are generator history-dependent functions for tagging particles
62  //associated with the hard process
63  //Currently implemented for Pythia 6 and Pythia 8 status codes and history
64 
65  //this particle is part of the hard process
66  bool isHardProcess(const P &p);
67 
68  //this particle is the direct descendant of a hard process particle of the same pdg id
69  bool fromHardProcess(const P &p);
70 
71  //this particle is the final state direct descendant of a hard process particle
72  bool fromHardProcessFinalState(const P &p);
73 
74  //this particle is the decayed direct descendant of a hard process particle
75  //such as a tau from the hard process
76  bool fromHardProcessDecayed(const P &p);
77 
78  //this particle is a direct or indirect decay product of a tau
79  //from the hard process
80  bool isHardProcessTauDecayProduct(const P &p);
81 
82  //this particle is a direct decay product of a tau
83  //from the hard process
84  bool isDirectHardProcessTauDecayProduct(const P &p);
85 
86  //this particle is the direct descendant of a hard process particle of the same pdg id
87  //For outgoing particles the kinematics are those before QCD or QED FSR
88  //This corresponds roughly to status code 3 in pythia 6
89  bool fromHardProcessBeforeFSR(const P &p);
90 
91  //this particle is the first copy of the particle in the chain with the same pdg id
92  bool isFirstCopy(const P &p);
93 
94  //this particle is the last copy of the particle in the chain with the same pdg id
95  //(and therefore is more likely, but not guaranteed, to carry the final physical momentum)
96  bool isLastCopy(const P &p);
97 
98  //this particle is the last copy of the particle in the chain with the same pdg id
99  //before QED or QCD FSR
100  //(and therefore is more likely, but not guaranteed, to carry the momentum after ISR)
101  //This flag only really makes sense for outgoing particles
102  bool isLastCopyBeforeFSR(const P &p);
103 
105  //These are utility functions used by the above
106 
107  //first mother in chain with a different pdg than the particle
108  const P *uniqueMother(const P &p);
109 
110  //return first copy of particle in chain (may be the particle itself)
111  const P *firstCopy(const P &p);
112 
113  //return last copy of particle in chain (may be the particle itself)
114  const P *lastCopy(const P &p);
115 
116  //return last copy of particle in chain before QED or QCD FSR (may be the particle itself)
117  const P *lastCopyBeforeFSR(const P &p);
118 
119  //return last copy of particle in chain before QED or QCD FSR, starting from itself (may be the particle itself)
120  const P *lastDaughterCopyBeforeFSR(const P &p);
121 
122  //return mother copy which is a hard process particle
123  const P *hardProcessMotherCopy(const P &p);
124 
125  //return previous copy of particle in chain (0 in case this is already the first copy)
126  const P *previousCopy(const P &p);
127 
128  //return next copy of particle in chain (0 in case this is already the last copy)
129  const P *nextCopy(const P &p);
130 
131  //return decayed mother (walk up the chain until found)
132  const P *findDecayedMother(const P &p);
133 
134  //return decayed mother matching requested abs(pdgid) (walk up the chain until found)
135  const P *findDecayedMother(const P &p, int abspdgid);
136 
138  //These are very basic utility functions to implement a common interface for reco::GenParticle
139  //and HepMC::GenParticle
140 
141  //pdgid
142  int pdgId(const reco::GenParticle &p);
143 
144  //pdgid
145  int pdgId(const HepMC::GenParticle &p);
146 
147  //abs(pdgid)
148  int absPdgId(const reco::GenParticle &p);
149 
150  //abs(pdgid)
151  int absPdgId(const HepMC::GenParticle &p);
152 
153  //number of mothers
154  unsigned int numberOfMothers(const reco::GenParticle &p);
155 
156  //number of mothers
157  unsigned int numberOfMothers(const HepMC::GenParticle &p);
158 
159  //mother
160  const reco::GenParticle *mother(const reco::GenParticle &p, unsigned int imoth=0);
161 
162  //mother
163  const HepMC::GenParticle *mother(const HepMC::GenParticle &p, unsigned int imoth=0);
164 
165  //number of daughters
166  unsigned int numberOfDaughters(const reco::GenParticle &p);
167 
168  //number of daughters
169  unsigned int numberOfDaughters(const HepMC::GenParticle &p);
170 
171  //ith daughter
172  const reco::GenParticle *daughter(const reco::GenParticle &p, unsigned int idau);
173 
174  //ith daughter
175  const HepMC::GenParticle *daughter(const HepMC::GenParticle &p, unsigned int idau);
176 
178  //Helper function to fill status flags
179  void fillGenStatusFlags(const P &p, reco::GenStatusFlags &statusFlags);
180 
181 protected:
182  std::unordered_set<const P*> dupCheck_;
183 
184 };
185 
186 
187 
189 //implementations
190 
192 template<typename P>
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 }
199 
201 template<typename P>
203  return p.status()==1 && isPrompt(p);
204 }
205 
206 template<typename P>
208  return p.status()==2 && (isHadron(p) || absPdgId(p)==13 || absPdgId(p)==15) && isLastCopy(p);
209 }
210 
211 template<typename P>
213  return isDecayedLeptonHadron(p) && isPrompt(p);
214 }
215 
217 template<typename P>
219  return findDecayedMother(p,15) != 0;
220 }
221 
223 template<typename P>
225  const P *tau = findDecayedMother(p,15);
226  return tau && isPrompt(*tau);
227 }
228 
230 template<typename P>
232  const P *tau = findDecayedMother(p,15);
233  const P *dm = findDecayedMother(p);
234  return tau && tau==dm;
235 }
236 
238 template<typename P>
240  const P *tau = findDecayedMother(p,15);
241  const P *dm = findDecayedMother(p);
242  return tau && tau==dm && isPrompt(*tau);
243 }
244 
246 template<typename P>
248  return findDecayedMother(p,13) != 0;
249 }
250 
252 template<typename P>
254  const P *mu = findDecayedMother(p,13);
255  return mu && isPrompt(*mu);
256 }
257 
259 template<typename P>
261  const P *um = uniqueMother(p);
262  return um && isHadron(*um) && isDecayedLeptonHadron(*um);
263 }
264 
266 template<typename P>
268  HepPDT::ParticleID heppdtid(pdgId(p));
269  return heppdtid.isHadron();
270 }
271 
273 template<typename P>
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 }
302 
304 template<typename P>
306  return hardProcessMotherCopy(p) != 0;
307 }
308 
310 template<typename P>
312  return p.status()==1 && fromHardProcess(p);
313 }
314 
316 template<typename P>
318  return isDecayedLeptonHadron(p) && fromHardProcess(p);
319 }
320 
322 template<typename P>
324  const P *tau = findDecayedMother(p,15);
325  return tau && fromHardProcessDecayed(*tau);
326 }
327 
329 template<typename P>
331  const P *tau = findDecayedMother(p,15);
332  const P *dm = findDecayedMother(p);
333  return tau && tau==dm && fromHardProcess(*tau);
334 }
335 
336 template<typename P>
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 }
362 
364 template<typename P>
366  return &p == firstCopy(p);
367 }
368 
370 template<typename P>
372  return &p == lastCopy(p);
373 }
374 
376 template<typename P>
378  return &p == lastCopyBeforeFSR(p);
379 }
380 
382 template<typename P>
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 }
393 
395 template<typename P>
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 }
406 
408 template<typename P>
409 const P *MCTruthHelper<P>::lastCopy(const P &p) {
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 }
419 
421 template<typename P>
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 }
453 
455 template<typename P>
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 }
486 
488 template<typename P>
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 }
504 
506 template<typename P>
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 }
519 
521 template<typename P>
522 const P *MCTruthHelper<P>::nextCopy(const P &p) {
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 }
534 
536 template<typename P>
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 }
547 
549 template<typename P>
550 const P *MCTruthHelper<P>::findDecayedMother(const P &p, int abspdgid) {
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 }
560 
562 template<typename P>
564  return p.pdgId();
565 }
566 
568 template<typename P>
570  return p.pdg_id();
571 }
572 
574 template<typename P>
576  return std::abs(p.pdgId());
577 }
578 
580 template<typename P>
582  return std::abs(p.pdg_id());
583 }
584 
586 template<typename P>
588  return p.numberOfMothers();
589 }
590 
592 template<typename P>
594  return p.production_vertex() ? p.production_vertex()->particles_in_size() : 0;
595 }
596 
598 template<typename P>
599 const reco::GenParticle *MCTruthHelper<P>::mother(const reco::GenParticle &p, unsigned int imoth) {
600  return static_cast<const reco::GenParticle*>(p.mother(imoth));
601 }
602 
604 template<typename P>
606  return p.production_vertex() && p.production_vertex()->particles_in_size() ? *(p.production_vertex()->particles_in_const_begin() + imoth) : 0;
607 }
608 
610 template<typename P>
612  return p.numberOfDaughters();
613 }
614 
616 template<typename P>
618  return p.end_vertex() ? p.end_vertex()->particles_out_size() : 0;
619 }
620 
622 template<typename P>
624  return static_cast<const reco::GenParticle*>(p.daughter(idau));
625 }
626 
628 template<typename P>
630  return *(p.end_vertex()->particles_out_const_begin() + idau);
631 }
632 
634 template<typename P>
636  statusFlags.setIsPrompt(isPrompt(p));
637  statusFlags.setIsDecayedLeptonHadron(isDecayedLeptonHadron(p));
638  statusFlags.setIsTauDecayProduct(isTauDecayProduct(p));
639  statusFlags.setIsPromptTauDecayProduct(isPromptTauDecayProduct(p));
640  statusFlags.setIsDirectTauDecayProduct(isDirectTauDecayProduct(p));
641  statusFlags.setIsDirectPromptTauDecayProduct(isDirectPromptTauDecayProduct(p));
642  statusFlags.setIsDirectHadronDecayProduct(isDirectHadronDecayProduct(p));
643  statusFlags.setIsHardProcess(isHardProcess(p));
644  statusFlags.setFromHardProcess(fromHardProcess(p));
645  statusFlags.setIsHardProcessTauDecayProduct(isHardProcessTauDecayProduct(p));
646  statusFlags.setIsDirectHardProcessTauDecayProduct(isDirectHardProcessTauDecayProduct(p));
647  statusFlags.setFromHardProcessBeforeFSR(fromHardProcessBeforeFSR(p));
648  statusFlags.setIsFirstCopy(isFirstCopy(p));
649  statusFlags.setIsLastCopy(isLastCopy(p));
650  statusFlags.setIsLastCopyBeforeFSR(isLastCopyBeforeFSR(p));
651 }
652 
653 
654 #endif
const P * lastDaughterCopyBeforeFSR(const P &p)
bool isTauDecayProduct(const P &p)
void fillGenStatusFlags(const P &p, reco::GenStatusFlags &statusFlags)
virtual int pdgId() const
PDG identifier.
bool fromHardProcess(const P &p)
bool isDirectHadronDecayProduct(const P &p)
const P * hardProcessMotherCopy(const P &p)
bool isPromptMuonDecayProduct(const P &p)
bool fromHardProcessFinalState(const P &p)
#define P
unsigned int numberOfMothers(const reco::GenParticle &p)
bool isDirectPromptTauDecayProduct(const P &p)
const reco::GenParticle * daughter(const reco::GenParticle &p, unsigned int idau)
void setIsLastCopyBeforeFSR(bool b)
void setIsPrompt(bool b)
void setIsTauDecayProduct(bool b)
void setIsDirectPromptTauDecayProduct(bool b)
bool isPromptFinalState(const P &p)
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)
unsigned int numberOfDaughters(const reco::GenParticle &p)
bool isPromptDecayed(const P &p)
bool isDirectHardProcessTauDecayProduct(const P &p)
const P * uniqueMother(const P &p)
const P * firstCopy(const P &p)
const P * findDecayedMother(const P &p)
const P * previousCopy(const P &p)
virtual size_t numberOfMothers() const
number of mothers
bool isFirstCopy(const P &p)
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) ...
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void setIsHardProcessTauDecayProduct(bool b)
const int mu
Definition: Constants.h:22
bool isHardProcessTauDecayProduct(const P &p)
void setFromHardProcess(bool b)
bool isDirectTauDecayProduct(const P &p)
const P * nextCopy(const P &p)
bool isHadron(const P &p)
bool fromHardProcessDecayed(const P &p)
void setIsPromptTauDecayProduct(bool b)
int absPdgId(const reco::GenParticle &p)
bool isPromptTauDecayProduct(const P &p)
void setIsDecayedLeptonHadron(bool b)
void setIsFirstCopy(bool b)
std::unordered_set< const P * > dupCheck_
bool isDecayedLeptonHadron(const P &p)
bool fromHardProcessBeforeFSR(const P &p)
bool isMuonDecayProduct(const P &p)
void setIsDirectTauDecayProduct(bool b)
void setFromHardProcessBeforeFSR(bool b)
bool isLastCopy(const P &p)
const P * lastCopy(const P &p)
void setIsLastCopy(bool b)
const P * lastCopyBeforeFSR(const P &p)
int pdgId(const reco::GenParticle &p)
virtual const Candidate * mother(size_type=0) const
return mother at a given position, i = 0, ... numberOfMothers() - 1 (read only mode) ...
const reco::GenParticle * mother(const reco::GenParticle &p, unsigned int imoth=0)