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);
24 
25  //is particle prompt and final state
26  bool isPromptFinalState(const P &p);
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);
33 
34  //is particle prompt and decayed
35  bool isPromptDecayed(const P &p);
36 
37  //this particle is a direct or indirect tau decay product
38  bool isTauDecayProduct(const P &p);
39 
40  //this particle is a direct or indirect decay product of a prompt tau
41  bool isPromptTauDecayProduct(const P &p);
42 
43  //this particle is a direct tau decay product
44  bool isDirectTauDecayProduct(const P &p);
45 
46  //this particle is a direct decay product from a prompt tau
47  bool isDirectPromptTauDecayProduct(const P &p);
48 
49  //this particle is a direct or indirect muon decay product
50  bool isMuonDecayProduct(const P &p);
51 
52  //this particle is a direct or indirect decay product of a prompt muon
53  bool isPromptMuonDecayProduct(const P &p);
54 
55  //this particle is a direct decay product from a hadron
56  bool isDirectHadronDecayProduct(const P &p);
57 
58  //is particle a hadron
59  bool isHadron(const P &p);
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);
68 
69  //this particle is the direct descendant of a hard process particle of the same pdg id
70  bool fromHardProcess(const P &p);
71 
72  //this particle is the final state direct descendant of a hard process particle
73  bool fromHardProcessFinalState(const P &p);
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);
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);
82 
83  //this particle is a direct decay product of a tau
84  //from the hard process
85  bool isDirectHardProcessTauDecayProduct(const P &p);
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);
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);
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);
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);
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);
110 
111  //return first copy of particle in chain (may be the particle itself)
112  const P *firstCopy(const P &p);
113 
114  //return last copy of particle in chain (may be the particle itself)
115  const P *lastCopy(const P &p);
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);
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);
122 
123  //return mother copy which is a hard process particle
124  const P *hardProcessMotherCopy(const P &p);
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);
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);
131 
132  //return decayed mother (walk up the chain until found)
133  const P *findDecayedMother(const P &p);
134 
135  //return decayed mother matching requested abs(pdgid) (walk up the chain until found)
136  const P *findDecayedMother(const P &p, int abspdgid);
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);
144 
145  //pdgid
146  int pdgId(const HepMC::GenParticle &p);
147 
148  //abs(pdgid)
149  int absPdgId(const reco::GenParticle &p);
150 
151  //abs(pdgid)
152  int absPdgId(const HepMC::GenParticle &p);
153 
154  //number of mothers
155  unsigned int numberOfMothers(const reco::GenParticle &p);
156 
157  //number of mothers
158  unsigned int numberOfMothers(const HepMC::GenParticle &p);
159 
160  //mother
161  const reco::GenParticle *mother(const reco::GenParticle &p, unsigned int imoth=0);
162 
163  //mother
164  const HepMC::GenParticle *mother(const HepMC::GenParticle &p, unsigned int imoth=0);
165 
166  //number of daughters
167  unsigned int numberOfDaughters(const reco::GenParticle &p);
168 
169  //number of daughters
170  unsigned int numberOfDaughters(const HepMC::GenParticle &p);
171 
172  //ith daughter
173  const reco::GenParticle *daughter(const reco::GenParticle &p, unsigned int idau);
174 
175  //ith daughter
176  const HepMC::GenParticle *daughter(const HepMC::GenParticle &p, unsigned int idau);
177 
179  //Helper function to fill status flags
180  void fillGenStatusFlags(const P &p, reco::GenStatusFlags &statusFlags);
181 
182 protected:
183  std::unordered_set<const P*> dupCheck_;
184 
185 };
186 
187 
188 
190 //implementations
191 
193 template<typename P>
195  //particle from hadron/muon/tau decay -> not prompt
196  //checking all the way up the chain treats properly the radiated photon
197  //case as well
198  return findDecayedMother(p) == 0;
199 }
200 
202 template<typename P>
204  return p.status()==1 && isPrompt(p);
205 }
206 
207 template<typename P>
209  return p.status()==2 && (isHadron(p) || absPdgId(p)==13 || absPdgId(p)==15) && isLastCopy(p);
210 }
211 
212 template<typename P>
214  return isDecayedLeptonHadron(p) && isPrompt(p);
215 }
216 
218 template<typename P>
220  return findDecayedMother(p,15) != 0;
221 }
222 
224 template<typename P>
226  const P *tau = findDecayedMother(p,15);
227  return tau && isPrompt(*tau);
228 }
229 
231 template<typename P>
233  const P *tau = findDecayedMother(p,15);
234  const P *dm = findDecayedMother(p);
235  return tau && tau==dm;
236 }
237 
239 template<typename P>
241  const P *tau = findDecayedMother(p,15);
242  const P *dm = findDecayedMother(p);
243  return tau && tau==dm && isPrompt(*tau);
244 }
245 
247 template<typename P>
249  return findDecayedMother(p,13) != 0;
250 }
251 
253 template<typename P>
255  const P *mu = findDecayedMother(p,13);
256  return mu && isPrompt(*mu);
257 }
258 
260 template<typename P>
262  const P *um = uniqueMother(p);
263  return um && isHadron(*um) && isDecayedLeptonHadron(*um);
264 }
265 
267 template<typename P>
269  HepPDT::ParticleID heppdtid(pdgId(p));
270  return heppdtid.isHadron();
271 }
272 
274 template<typename P>
276 
277  //status 3 in pythia6 means hard process;
278  if (p.status()==3) return true;
279 
280  //hard process codes for pythia8 are 21-29 inclusive (currently 21,22,23,24 are used)
281  if (p.status()>20 && p.status()<30) return true;
282 
283  //if this is a final state or decayed particle,
284  //check if direct mother is a resonance decay in pythia8 but exclude FSR branchings
285  //(In pythia8 if a resonance decay product did not undergo any further branchings
286  //it will be directly stored as status 1 or 2 without any status 23 copy)
287  if (p.status()==1 || p.status()==2) {
288  const P *um = mother(p);
289  if (um) {
290  const P *firstcopy = firstCopy(*um);
291  bool fromResonance = firstcopy && firstcopy->status()==22;
292 
293  const P *umNext = nextCopy(*um);
294  bool fsrBranching = umNext && umNext->status()>50 && umNext->status()<60;
295 
296  if (fromResonance && !fsrBranching) return true;
297  }
298  }
299 
300  return false;
301 
302 }
303 
305 template<typename P>
307  return hardProcessMotherCopy(p) != 0;
308 }
309 
311 template<typename P>
313  return p.status()==1 && fromHardProcess(p);
314 }
315 
317 template<typename P>
319  return isDecayedLeptonHadron(p) && fromHardProcess(p);
320 }
321 
323 template<typename P>
325  const P *tau = findDecayedMother(p,15);
326  return tau && fromHardProcessDecayed(*tau);
327 }
328 
330 template<typename P>
332  const P *tau = findDecayedMother(p,15);
333  const P *dm = findDecayedMother(p);
334  return tau && tau==dm && fromHardProcess(*tau);
335 }
336 
337 template<typename P>
339  //pythia 6 documentation line roughly corresponds to this condition
340  if (p.status()==3) return true;
341 
342  //check hard process mother properties
343  const P *hpc = hardProcessMotherCopy(p);
344  if (!hpc) return false;
345 
346  //for incoming partons in pythia8, more useful information is not
347  //easily available, so take only the incoming parton itself
348  if (hpc->status()==21 && (&p)==hpc) return true;
349 
350  //for intermediate particles in pythia 8, just take the last copy
351  if (hpc->status()==22 && isLastCopy(p)) return true;
352 
353  //for outgoing particles in pythia 8, explicitly find the last copy
354  //before FSR starting from the hardProcess particle, and take only
355  //this one
356  if ( (hpc->status()==23 || hpc->status()==1) && (&p)==lastDaughterCopyBeforeFSR(*hpc) ) return true;
357 
358 
359  //didn't satisfy any of the conditions
360  return false;
361 
362 }
363 
365 template<typename P>
367  return &p == firstCopy(p);
368 }
369 
371 template<typename P>
373  return &p == lastCopy(p);
374 }
375 
377 template<typename P>
379  return &p == lastCopyBeforeFSR(p);
380 }
381 
383 template<typename P>
385  const P *mo = &p;
386  dupCheck_.clear();
387  while (mo && pdgId(*mo)==pdgId(p)) {
388  dupCheck_.insert(mo);
389  mo = mother(*mo);
390  if (dupCheck_.count(mo)) return 0;
391  }
392  return mo;
393 }
394 
396 template<typename P>
398  const P *pcopy = &p;
399  dupCheck_.clear();
400  while (previousCopy(*pcopy)) {
401  dupCheck_.insert(pcopy);
402  pcopy = previousCopy(*pcopy);
403  if (dupCheck_.count(pcopy)) return 0;
404  }
405  return pcopy;
406 }
407 
409 template<typename P>
410 const P *MCTruthHelper<P>::lastCopy(const P &p) {
411  const P *pcopy = &p;
412  dupCheck_.clear();
413  while (nextCopy(*pcopy)) {
414  dupCheck_.insert(pcopy);
415  pcopy = nextCopy(*pcopy);
416  if (dupCheck_.count(pcopy)) return 0;
417  }
418  return pcopy;
419 }
420 
422 template<typename P>
424  //start with first copy and then walk down until there is FSR
425  const P *pcopy = firstCopy(p);
426  if (!pcopy) return 0;
427  bool hasDaughterCopy = true;
428  dupCheck_.clear();
429  while (hasDaughterCopy) {
430  dupCheck_.insert(pcopy);
431  hasDaughterCopy = false;
432  const unsigned int ndau = numberOfDaughters(*pcopy);
433  //look for FSR
434  for (unsigned int idau = 0; idau<ndau; ++idau) {
435  const P *dau = daughter(*pcopy,idau);
436  if (pdgId(*dau)==21 || pdgId(*dau)==22) {
437  //has fsr (or else decayed and is the last copy by construction)
438  return pcopy;
439  }
440  }
441  //look for daughter copy
442  for (unsigned int idau = 0; idau<ndau; ++idau) {
443  const P *dau = daughter(*pcopy,idau);
444  if (pdgId(*dau)==pdgId(p)) {
445  pcopy = dau;
446  hasDaughterCopy = true;
447  break;
448  }
449  }
450  if (dupCheck_.count(pcopy)) return 0;
451  }
452  return pcopy;
453 }
454 
456 template<typename P>
458  //start with this particle and then walk down until there is FSR
459  const P *pcopy = &p;
460  bool hasDaughterCopy = true;
461  dupCheck_.clear();
462  while (hasDaughterCopy) {
463  dupCheck_.insert(pcopy);
464  hasDaughterCopy = false;
465  const unsigned int ndau = numberOfDaughters(*pcopy);
466  //look for FSR
467  for (unsigned int idau = 0; idau<ndau; ++idau) {
468  const P *dau = daughter(*pcopy,idau);
469  if (pdgId(*dau)==21 || pdgId(*dau)==22) {
470  //has fsr (or else decayed and is the last copy by construction)
471  return pcopy;
472  }
473  }
474  //look for daughter copy
475  for (unsigned int idau = 0; idau<ndau; ++idau) {
476  const P *dau = daughter(*pcopy,idau);
477  if (pdgId(*dau)==pdgId(p)) {
478  pcopy = dau;
479  hasDaughterCopy = true;
480  break;
481  }
482  }
483  if (dupCheck_.count(pcopy)) return 0;
484  }
485  return pcopy;
486 }
487 
489 template<typename P>
491  //is particle itself is hard process particle
492  if (isHardProcess(p)) return &p;
493 
494  //check if any other copies are hard process particles
495  const P *pcopy = &p;
496  dupCheck_.clear();
497  while (previousCopy(*pcopy)) {
498  dupCheck_.insert(pcopy);
499  pcopy = previousCopy(*pcopy);
500  if (isHardProcess(*pcopy)) return pcopy;
501  if (dupCheck_.count(pcopy)) break;
502  }
503  return 0;
504 }
505 
507 template<typename P>
509 
510  const unsigned int nmoth = numberOfMothers(p);
511  for (unsigned int imoth = 0; imoth<nmoth; ++imoth) {
512  const P *moth = mother(p,imoth);
513  if (pdgId(*moth)==pdgId(p)) {
514  return moth;
515  }
516  }
517 
518  return 0;
519 }
520 
522 template<typename P>
523 const P *MCTruthHelper<P>::nextCopy(const P &p) {
524 
525  const unsigned int ndau = numberOfDaughters(p);
526  for (unsigned int idau = 0; idau<ndau; ++idau) {
527  const P *dau = daughter(p,idau);
528  if (pdgId(*dau)==pdgId(p)) {
529  return dau;
530  }
531  }
532 
533  return 0;
534 }
535 
537 template<typename P>
539  const P *mo = mother(p);
540  dupCheck_.clear();
541  while (mo && !isDecayedLeptonHadron(*mo)) {
542  dupCheck_.insert(mo);
543  mo = mother(*mo);
544  if (dupCheck_.count(mo)) return 0;
545  }
546  return mo;
547 }
548 
550 template<typename P>
551 const P *MCTruthHelper<P>::findDecayedMother(const P &p, int abspdgid) {
552  const P *mo = mother(p);
553  dupCheck_.clear();
554  while (mo && (absPdgId(*mo)!=abspdgid || !isDecayedLeptonHadron(*mo)) ) {
555  dupCheck_.insert(mo);
556  mo = mother(*mo);
557  if (dupCheck_.count(mo)) return 0;
558  }
559  return mo;
560 }
561 
563 template<typename P>
565  return p.pdgId();
566 }
567 
569 template<typename P>
571  return p.pdg_id();
572 }
573 
575 template<typename P>
577  return std::abs(p.pdgId());
578 }
579 
581 template<typename P>
583  return std::abs(p.pdg_id());
584 }
585 
587 template<typename P>
589  return p.numberOfMothers();
590 }
591 
593 template<typename P>
595  return p.production_vertex() ? p.production_vertex()->particles_in_size() : 0;
596 }
597 
599 template<typename P>
600 const reco::GenParticle *MCTruthHelper<P>::mother(const reco::GenParticle &p, unsigned int imoth) {
601  return static_cast<const reco::GenParticle*>(p.mother(imoth));
602 }
603 
605 template<typename P>
607  return p.production_vertex() && p.production_vertex()->particles_in_size() ? *(p.production_vertex()->particles_in_const_begin() + imoth) : 0;
608 }
609 
611 template<typename P>
613  return p.numberOfDaughters();
614 }
615 
617 template<typename P>
619  return p.end_vertex() ? p.end_vertex()->particles_out_size() : 0;
620 }
621 
623 template<typename P>
625  return static_cast<const reco::GenParticle*>(p.daughter(idau));
626 }
627 
629 template<typename P>
631  return *(p.end_vertex()->particles_out_const_begin() + idau);
632 }
633 
635 template<typename P>
637  statusFlags.setIsPrompt(isPrompt(p));
639  statusFlags.setIsTauDecayProduct(isTauDecayProduct(p));
644  statusFlags.setIsHardProcess(isHardProcess(p));
645  statusFlags.setFromHardProcess(fromHardProcess(p));
649  statusFlags.setIsFirstCopy(isFirstCopy(p));
650  statusFlags.setIsLastCopy(isLastCopy(p));
652 }
653 
654 
655 #endif
const P * lastDaughterCopyBeforeFSR(const P &p)
bool isTauDecayProduct(const P &p)
void fillGenStatusFlags(const P &p, reco::GenStatusFlags &statusFlags)
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)
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) ...
virtual int pdgId() const final
PDG identifier.
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_
std::pair< OmniClusterRef, TrackingParticleRef > P
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)