ATLAS Offline Software
AtlasPID.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3 */
4 #ifndef TRUTHUTILS_ATLASPID_H
5 #define TRUTHUTILS_ATLASPID_H
6 #include <vector>
7 #include <cmath>
8 #include <algorithm>
9 #include <array>
10 #include <cstdlib>
16 class DecodedPID: public std::pair<int,std::vector<int>> {
17 public:
18  inline DecodedPID(const int& p){
19  this->first=p;
20  this->second.reserve(10);
21  int ap = std::abs(p);
22  for(; ap; ap/=10) this->second.push_back( ap%10 );
23  std::reverse(this->second.begin(), this->second.end());
24  }
25  inline DecodedPID shift(const size_t n) const { return DecodedPID(this->first%int(std::pow(10,ndigits()-n)));}
26  inline const int& operator()(const size_t n) const { return this->second.at(n);}
27  inline const int& last() const { return this->second.back();}
28  inline const int& pid() const { return this->first;}
29  inline int max_digit(const int m,const int n) const { return *std::max_element(second.rbegin() + m, second.rbegin() + n);}
30  inline int min_digit(const int m,const int n) const { return *std::min_element(second.rbegin() + m, second.rbegin() + n);}
31  inline size_t ndigits() const { return this->second.size();}
32 };
33 
34 static const int TABLESIZE = 100;
35 static const std::array<int,TABLESIZE> triple_charge = {
36  +0, -1, +2, -1, +2, -1, +2, -1, +2, +0,
37  +0, -3, +0, -3, +0, -3, +0, -3, +0, +0,
38  +0, +0, +0, +0, +3, +0, +0, +0, +0, +0,
39  +0, +0, +0, +0, +3, +0, +0, +3, +6, +0,
40  +0, +0, -1, +0, +0, +0, +0, +0, +0, +0,
41  +0, +0, +0, +0, +0, +0, +0, +0, +0, +0,
42  +0, +0, +0, +0, +0, +0, +0, +0, +0, +0,
43  +0, +0, +0, +0, +0, +0, +0, +0, +0, +0,
44  +0, +0, +0, +0, +0, +0, +0, +0, +0, +0,
45  +0, +0, +0, +0, +0, +0, +0, +0, +0, +0
46 };
47 static const std::array<int,TABLESIZE> double_spin = {
48  +0, +1, +1, +1, +1, +1, +1, +1, +1, +0,
49  +0, +1, +1, +1, +1, +1, +1, +1, +1, +0,
50  +2, +2, +2, +2, +2, +0, +0, +0, +0, +0,
51  +0, +0, +2, +2, +2, +0, +0, +0, +0, +4,
52  +0, +0, -1, +0, +0, +0, +0, +0, +0, +0,
53  +0, +0, +1, +2, +0, +2, +0, +0, +0, +0,
54  +0, +0, +0, +0, +0, +0, +0, +0, +0, +0,
55  +0, +0, +0, +0, +0, +0, +0, +0, +0, +0,
56  +0, +0, +0, +0, +0, +0, +0, +0, +0, +0,
57  +0, +0, +0, +0, +0, +0, +0, +0, +0, +0
58 };
59 
60 
61 static const int DQUARK = 1;
62 static const int UQUARK = 2;
63 static const int SQUARK = 3;
64 static const int CQUARK = 4;
65 static const int BQUARK = 5;
66 static const int TQUARK = 6;
67 static const int BPRIME = 7; // 4th Generation quark
68 static const int TPRIME = 8; // 4th Generation quark
69 static const int QUARK_LIMIT = BPRIME; // Quark pdg_ids less than this are considered in (R-)Hadrons and Diquarks
70 
71 static const int ELECTRON = 11;
72 static const int POSITRON = -ELECTRON;
73 static const int NU_E = 12;
74 static const int MUON = 13;
75 static const int NU_MU = 14;
76 static const int TAU = 15;
77 static const int NU_TAU = 16;
78 static const int LPRIME = 17; // 4th Generation charged lepton
79 static const int NUPRIME = 18; // 4th Generation neutrino
80 
81 static const int GLUON = 21;
82 // APID: 9 rather than 21 is used to denote a gluon/gluino in composite states. (From PDG 11g)
83 static const int COMPOSITEGLUON = 9;
84 static const int PHOTON = 22;
85 static const int Z0BOSON = 23;
86 static const int WPLUSBOSON = 24;
87 static const int HIGGSBOSON = 25;
88 static const int ZPRIME = 32; // Z′/Z^0_2
89 static const int ZDBLPRIME = 33; // Z′′/Z^0_3
90 static const int WPLUSPRIME = 34; // W ′/W^+_2
91 static const int HIGGS2 = 35; // H^0/H^0_2 FIXME Any better ideas?
92 static const int HIGGS3 = 36; // A^0/H^0_3 FIXME Any better ideas?
93 static const int HIGGSPLUS = 37; // H^+
94 static const int HIGGSPLUSPLUS = 38; // H^++
95 static const int GRAVITON = 39;
96 static const int HIGGS4 = 40; // a^0/H^0_4 FIXME Any better ideas?
97 static const int LEPTOQUARK = 42;
98 
103 static const int DARKPHOTON = 60000;
104 static const int MAVTOP = 60001;
105 
106 static const int PIPLUS = 211;
107 static const int PIMINUS = -PIPLUS;
108 static const int PI0 = 111;
109 static const int K0L = 130;
110 
111 static const int K0S = 310;
112 static const int K0 = 311;
113 static const int KPLUS = 321;
114 static const int DPLUS = 411;
115 static const int DSTAR = 413;
116 static const int D0 = 421;
117 static const int DSPLUS = 431;
118 static const int JPSI = 443;
119 static const int B0 = 511;
120 static const int BCPLUS = 541;
121 static const int PROTON = 2212;
122 static const int NEUTRON = 2112;
123 static const int LAMBDA0 = 3122;
124 static const int LAMBDACPLUS = 4122;
125 static const int LAMBDAB0 = 5122;
126 static const int PSI2S = 20443;
127 
133 static const int RH_NU_E = 9900012;
134 static const int RH_NU_MU = 9900014;
135 static const int RH_NU_TAU = 9900016;
136 static const int WBOSON_LRSM = 9900024;
137 
138 static const int LEAD = 1000822080;
139 static const int OXYGEN = 1000080160;
140 static const int NEON = 1000100200;
141 static const int HELIUM = 1000020040;
142 
146 static const int POMERON = 990;
147 static const int ODDERON = 9990;
148 static const int REGGEON = 110;
149 
155 static const int GEANTINOPLUS = 998;
156 static const int GEANTINO0 = 999;
157 
158 
164 template<class T> inline bool isQuark(const T& p) {return isQuark(p->pdg_id());}
165 template<> inline bool isQuark(const int& p) { return p != 0 && (std::abs(p) <= TPRIME || std::abs(p) == MAVTOP);}
166 template<> inline bool isQuark(const DecodedPID& p){ return isQuark(p.pid()); }
167 
168 // APID: the fourth generation quarks are not standard model quarks
169 template<class T> inline bool isSMQuark(const T& p) {return isSMQuark(p->pdg_id());}
170 template<> inline bool isSMQuark(const int& p) { return p != 0 && std::abs(p) <= TQUARK;}
171 template<> inline bool isSMQuark(const DecodedPID& p){ return isSMQuark(p.pid()); }
172 
173 template<class T> inline bool isStrange(const T& p) {return isStrange(p->pdg_id());}
174 template<> inline bool isStrange(const int& p){ return std::abs(p) == SQUARK;}
175 
176 template<class T> inline bool isCharm(const T& p){return isCharm(p->pdg_id());}
177 template<> inline bool isCharm(const int& p){ return std::abs(p) == CQUARK;}
178 
179 template<class T> inline bool isBottom(const T& p){return isBottom(p->pdg_id());}
180 template<> inline bool isBottom(const int& p){ return std::abs(p) == BQUARK;}
181 
182 template<class T> inline bool isTop(const T& p){return isTop(p->pdg_id());}
183 template<> inline bool isTop(const int& p){ return std::abs(p) == TQUARK;}
184 
186 template<class T> inline bool isLepton(const T& p){return isLepton(p->pdg_id());}
187 template<> inline bool isLepton(const int& p){ auto sp = std::abs(p); return sp >= ELECTRON && sp <= NUPRIME; }
188 template<> inline bool isLepton(const DecodedPID& p){ return isLepton(p.pid()); }
189 
191 template<class T> inline bool isSMLepton(const T& p){return isSMLepton(p->pdg_id());}
192 template<> inline bool isSMLepton(const int& p){ auto sp = std::abs(p); return sp >= ELECTRON && sp <= NU_TAU; }
193 template<> inline bool isSMLepton(const DecodedPID& p){ return isSMLepton(p.pid()); }
194 
196 template<class T> inline bool isChLepton(const T& p){return isChLepton(p->pdg_id());}
197 template<> inline bool isChLepton(const int& p){ auto sp = std::abs(p); return sp >= ELECTRON && sp <= LPRIME && sp%2 == 1; }
198 
199 template<class T> inline bool isElectron(const T& p){return isElectron(p->pdg_id());}
200 template<> inline bool isElectron(const int& p){ return std::abs(p) == ELECTRON;}
201 
202 template<class T> inline bool isMuon(const T& p){return isMuon(p->pdg_id());}
203 template<> inline bool isMuon(const int& p){ return std::abs(p) == MUON;}
204 
205 template<class T> inline bool isTau(const T& p){return isTau(p->pdg_id());}
206 template<> inline bool isTau(const int& p){ return std::abs(p) == TAU;}
207 
209 template<class T> inline bool isNeutrino(const T& p){return isNeutrino(p->pdg_id());}
210 template<> inline bool isNeutrino(const int& p){ auto sp = std::abs(p); return sp == NU_E || sp == NU_MU || sp == NU_TAU || sp == NUPRIME; }
211 
212 template<class T> inline bool isSMNeutrino(const T& p){return isSMNeutrino(p->pdg_id());}
213 template<> inline bool isSMNeutrino(const int& p){ auto sp = std::abs(p); return sp == NU_E || sp == NU_MU || sp == NU_TAU; }
214 
217 template<class T> inline bool isFourthGeneration(const T& p){return isFourthGeneration(p->pdg_id());}
218 template<> inline bool isFourthGeneration(const int& p) {return std::abs(p) == BPRIME || std::abs(p) == TPRIME || std::abs(p) == LPRIME || std::abs(p) == NUPRIME;}
219 
224 template<class T> inline bool isDiquark(const T& p){return isDiquark(p->pdg_id());}
225 template<> inline bool isDiquark(const DecodedPID& p){
226  if ( p.ndigits() == 4 && p(0) >= p(1) && p(1) !=0 && p(2) == 0 && (p.last() == 1 || p.last() == 3)
227  && p.max_digit(2,4) < QUARK_LIMIT
228  ) return true;
229  return false;
230 }
231 template<> inline bool isDiquark(const int& p){ auto value_digits = DecodedPID(p); return isDiquark(value_digits);}
232 
241 template<class T> inline bool isMeson(const T& p){return isMeson(p->pdg_id());}
242 template<> inline bool isMeson(const DecodedPID& p){
243  if (p.ndigits() < 3 ) return false;
244  if (p.ndigits() == 7 && (p(0) == 1 || p(0) == 2)) return false; // APID don't match SUSY particles
245  if (std::abs(p.pid()) == K0S) return true;
246  if (std::abs(p.pid()) == K0L) return true;
247  if (std::abs(p.pid()) == K0) return true;
248  if (p.last() % 2 != 1 ) return false;
249  if (p.max_digit(1,3) >= QUARK_LIMIT ) return false; // Ignore pdg_ids which would describe states including fourth generation quarks
250  if (p.min_digit(1,3) == 0 ) return false;
251  if (*(p.second.rbegin() + 2) < *(p.second.rbegin() + 1) ) return false; // Quark ordering (nq2 >= nq3)
252  if (*(p.second.rbegin() + 2) == *(p.second.rbegin() + 1) && p.pid() < 0 ) return false; // Illegal antiparticle check (nq2 == nq3)
253  if (p.ndigits() == 3 ) return true;
254  if (*(p.second.rbegin() + 3) != 0 ) return false; // Only two quarks! (nq1 == 0)
255  if (p.ndigits() == 5 && p(0) == 1 ) return true;
256  if (p.ndigits() == 5 && p(0) == 2 && p.last() > 1 ) return true;
257  if (p.ndigits() == 5 && p(0) == 3 && p.last() > 1 ) return true;
258  if (p.ndigits() == 6 && p.last() % 2 == 1 ) return true;
259  if (p.ndigits() == 7 && p(0) == 9 && p(1) == 0 ) return true;
260 
261  return false;
262 }
263 template<> inline bool isMeson(const int& p){ auto value_digits = DecodedPID(p); return isMeson(value_digits);}
264 
270 template<class T> inline bool isQuarkonium(const T& p){return isQuarkonium(p->pdg_id());}
271 template<> inline bool isQuarkonium(const DecodedPID& p) {
272  if (!isMeson(p)) return false; //< all quarkonia are mesons
273  return (*(p.second.rbegin() + 2) > SQUARK && p.last() > 0 && *(p.second.rbegin() + 1) == *(p.second.rbegin() + 2));
274 }
275 template<> inline bool isQuarkonium(const int& p){ auto value_digits = DecodedPID(p); return isQuarkonium(value_digits);}
276 
279 template<class T> inline bool isBaryon(const T& p){return isBaryon(p->pdg_id());}
280 template<> inline bool isBaryon(const DecodedPID& p){
281  if (p.ndigits() < 4 ) return false;
282  if (p.max_digit(1,4) >= QUARK_LIMIT ) return false; // Ignore pdg_ids which would describe states including fourth generation quarks
283  if (p.min_digit(1,4) == 0) return false; // Ignore pdg_ids with zero for nq1, nq2, nq3
284  if (p.ndigits() == 4 && (p.last() == 2 || p.last() == 4|| p.last() == 6|| p.last() == 8) ) return true;
285 
286  if (p.ndigits() == 5 && p(0) == 1 && (p.last() == 2 || p.last() == 4) ) return true;
287  if (p.ndigits() == 5 && p(0) == 3 && (p.last() == 2 || p.last() == 4) ) return true;
288 
289  if (p.ndigits() == 6 ) {
290  if (p(0) == 1 && p(1) == 0 && p.last() == 2 ) return true;
291  if (p(0) == 1 && p(1) == 0 && p.last() == 4 ) return true;
292  if (p(0) == 1 && p(1) == 0 && p.last() == 6 ) return true;
293  if (p(0) == 1 && p(1) == 1 && p.last() == 2 ) return true;
294  if (p(0) == 1 && p(1) == 2 && p.last() == 4 ) return true;
295 
296  if (p(0) == 2 && p(1) == 0 && p.last() == 2 ) return true;
297  if (p(0) == 2 && p(1) == 0 && p.last() == 4 ) return true;
298  if (p(0) == 2 && p(1) == 0 && p.last() == 6 ) return true;
299  if (p(0) == 2 && p(1) == 0 && p.last() == 8 ) return true;
300  if (p(0) == 2 && p(1) == 1 && p.last() == 2 ) return true;
301  }
302 
303  if (p.ndigits() == 5 ) {
304  if (p(0) == 2 && p.last() == 2 ) return true;
305  if (p(0) == 2 && p.last() == 4 ) return true;
306  if (p(0) == 2 && p.last() == 6 ) return true;
307  if (p(0) == 5 && p.last() == 2 ) return true;
308  if (p(0) == 1 && p.last() == 6 ) return true;
309  if (p(0) == 4 && p.last() == 2 ) return true;
310  }
311  return false;
312 }
313 template<> inline bool isBaryon(const int& p){ auto value_digits = DecodedPID(p); return isBaryon(value_digits);}
314 
322 template<class T> inline bool isTetraquark(const T& p){return isTetraquark(p->pdg_id());}
323 template<> inline bool isTetraquark(const DecodedPID& p){
324  return (p.ndigits() == 9 && p(0) == 1 && p(5) == 0 &&
325  p.max_digit(1,3) < QUARK_LIMIT && p.min_digit(1,3) > 0 && // ignore 4th generation quarks for nq3 and nq4
326  p.max_digit(4,6) < QUARK_LIMIT && p.min_digit(4,6) > 0 && // ignore 4th generation quarks for nq1 and nq2
327  ( p(3) >= p(4) && p(6) >= p(7) ) && ( ( p(3) > p(6) ) || ( p(3) == p(6) && (p(4) >= p(7))))
328  );
329 }
330 template<> inline bool isTetraquark(const int& p){ auto value_digits = DecodedPID(p); return isTetraquark(value_digits);}
331 
338 // APID: states with fourth generation quarks are not pentaquarks
339 template<class T> inline bool isPentaquark(const T& p){return isPentaquark(p->pdg_id());}
340 template<> inline bool isPentaquark(const DecodedPID& p){
341  return (p.ndigits() == 9 && p(0) == 1 &&
342  p.max_digit(1,6) < QUARK_LIMIT && p.min_digit(1,6) > 0 && // ignore 4th generation (anti-)quarks
343  ( p(3) >= p(4) && p(4) >= p(5) && p(5) >= p(6)) );
344 }
345 template<> inline bool isPentaquark(const int& p){ auto value_digits = DecodedPID(p); return isPentaquark(value_digits);}
346 
347 // APID Mesons, Baryons, Tetraquarks and Pentaquarks are Hadrons
348 template<class T> inline bool isHadron(const T& p){return isHadron(p->pdg_id());}
349 template<> inline bool isHadron(const DecodedPID& p){ return isMeson(p) || isBaryon(p) || isTetraquark(p) || isPentaquark(p); }
350 template<> inline bool isHadron(const int& p){ auto value_digits = DecodedPID(p); return isHadron(value_digits);}
351 
352 
356 template<class T> inline bool isTrajectory(const T& p){return isTrajectory(p->pdg_id());}
357 template<> inline bool isTrajectory(const int& p){ return std::abs(p) == POMERON || std::abs(p) == ODDERON || std::abs(p) == REGGEON; }
358 
359 
366 template<class T> inline bool isBoson(const T& p){return isBoson(p->pdg_id());}
367 template<> inline bool isBoson(const int& p){ auto sp = std::abs(p); return sp > 20 && sp < 41; }
368 template<> inline bool isBoson(const DecodedPID& p){ return isBoson(p.pid()); }
369 
370 template<class T> inline bool isGluon(const T& p){return isGluon(p->pdg_id());}
371 template<> inline bool isGluon(const int& p){ return p == GLUON; }
372 
373 template<class T> inline bool isPhoton(const T& p){return isPhoton(p->pdg_id());}
374 template<> inline bool isPhoton(const int& p){ return p == PHOTON; }
375 
376 template<class T> inline bool isZ(const T& p){return isZ(p->pdg_id());}
377 template<> inline bool isZ(const int& p){ return p == Z0BOSON; }
378 
379 template<class T> inline bool isW(const T& p){return isW(p->pdg_id());}
380 template<> inline bool isW(const int& p){ return std::abs(p) == WPLUSBOSON; }
381 
383 template<class T> inline bool isHeavyBoson(const T& p){return isHeavyBoson(p->pdg_id());}
384 template<> inline bool isHeavyBoson(const int& p){ return p == ZPRIME || p == ZDBLPRIME || std::abs(p) == WPLUSPRIME; }
385 
387 template<class T> inline bool isHiggs(const T& p){return isHiggs(p->pdg_id());}
388 template<> inline bool isHiggs(const int& p){ return p == HIGGSBOSON; }
389 
391 template<class T> inline bool isMSSMHiggs(const T& p){return isMSSMHiggs(p->pdg_id());}
392 template<> inline bool isMSSMHiggs(const int& p){ return p == HIGGS2 || p == HIGGS3 || std::abs(p) == HIGGSPLUS; }
393 
394 template<class T> inline bool isGraviton(const T& p) {return isGraviton(p->pdg_id());}
395 template<> inline bool isGraviton(const int& p){ return p == GRAVITON; }
396 
397 template<class T> inline bool isResonance(const T& p) { return isZ(p) || isW(p) || isHiggs(p) || isTop(p); } // APID: not including t' (pdg_id=8), Z', Z'' and W'+ or BSM Higgs bosons
398 
405 template<class T> inline bool isLeptoQuark(const T& p){return isLeptoQuark(p->pdg_id());}
406 template<> inline bool isLeptoQuark(const int& p){ return std::abs(p) == LEPTOQUARK; }
407 
408 template<class T> inline bool isPythia8Specific(const T& p){return isPythia8Specific(p->pdg_id());}
409 template<> inline bool isPythia8Specific(const DecodedPID& p){ return (p.ndigits() == 7 && p(0) == 9 && p(1) == 9);}
410 template<> inline bool isPythia8Specific(const int& p){ auto value_digits = DecodedPID(p); return isPythia8Specific(value_digits);}
411 
417 template<class T> inline bool isNeutrinoRH(const T& p){return isNeutrinoRH(p->pdg_id());}
418 template<> inline bool isNeutrinoRH(const int& p){ return (std::abs(p) == RH_NU_E || std::abs(p) == RH_NU_MU|| std::abs(p) == RH_NU_TAU);}
419 
422 template<class T> inline bool isGenSpecific(const T& p){return isGenSpecific(p->pdg_id());}
423 template<> inline bool isGenSpecific(const int& p){
424  int ap = std::abs(p);
425  if (ap >= 81 && ap <= 100) return true;
426  if (ap >= 901 && ap <= 930) return true;
427  if (ap >= 998 && ap <= 999) return true;
428  if (ap >= 1901 && ap <= 1930) return true;
429  if (ap >= 2901 && ap <= 2930) return true;
430  if (ap >= 3901 && ap <= 3930) return true;
431  return false;
432 }
433 
434 template<class T> inline bool isGeantino(const T& p){return isGeantino(p->pdg_id());}
435 template<> inline bool isGeantino(const int& p){ return (std::abs(p) == GEANTINO0 || std::abs(p) == GEANTINOPLUS);}
436 
438 template<class T> inline bool isGlueball(const T& p) { return isGlueball(p->pdg_id()); }
439 template<> inline bool isGlueball(const DecodedPID& p) {
440  if (p.ndigits() > 4) return false; // APID avoid classifying R-Glueballs as SM Glueballs
441  return
442  ( ( p.ndigits() == 3 && p(0) == COMPOSITEGLUON && p(1) == COMPOSITEGLUON && (p.last() == 1 || p.last() == 5) ) ||
443  ( p.ndigits() == 4 && p(0) == COMPOSITEGLUON && p(1) == COMPOSITEGLUON && p(2) == COMPOSITEGLUON && (p.last() == 3 || p.last() == 7) ) );
444 }
445 template<> inline bool isGlueball(const int& p) { auto value_digits = DecodedPID(p); return isGlueball(value_digits); }
446 
447 
453 
454 // APID: Super-partners of standard model quarks only
455 template<class T> inline bool isSquark(const T& p) { return isSquark(p->pdg_id()); }
456 template<> inline bool isSquark(const DecodedPID& p){
457  auto pp = p.shift(1); return (p.ndigits() == 7 && (p(0) == 1 || p(0) == 2) && isSMQuark(pp));
458 }
459 template<> inline bool isSquark(const int& p){ auto value_digits = DecodedPID(p); return isSquark(value_digits);}
460 
461 
462 // APID: Super-partners of left-handed standard model quarks only
463 template<class T> inline bool isSquarkLH(const T& p) { return isSquarkLH(p->pdg_id()); }
464 template<> inline bool isSquarkLH(const DecodedPID& p){
465  auto pp = p.shift(1); return (p.ndigits() == 7 && p(0) == 1 && isSMQuark(pp));
466 }
467 template<> inline bool isSquarkLH(const int& p){ auto value_digits = DecodedPID(p); return isSquarkLH(value_digits);}
468 
469 
470 // APID: Super-partners of right-handed standard model quarks only
471 template<class T> inline bool isSquarkRH(const T& p) { return isSquarkRH(p->pdg_id()); }
472 template<> inline bool isSquarkRH(const DecodedPID& p){
473  auto pp = p.shift(1); return (p.ndigits() == 7 && p(0) == 2 && isSMQuark(pp));
474 }
475 template<> inline bool isSquarkRH(const int& p){ auto value_digits = DecodedPID(p); return isSquarkRH(value_digits);}
476 
477 
478 // APID: Super-partners of standard model leptons only
479 template<class T> inline bool isSlepton(const T& p) { return isSlepton(p->pdg_id()); }
480 template<> inline bool isSlepton(const DecodedPID& p){ auto pp = p.shift(1); return (p.ndigits() == 7 && (p(0) == 1 || p(0) == 2) && isSMLepton(pp));}
481 template<> inline bool isSlepton(const int& p){ auto value_digits = DecodedPID(p); return isSlepton(value_digits);}
482 
483 
484 // APID: Super-partners of left-handed standard model leptons only
485 template<class T> inline bool isSleptonLH(const T& p) { return isSleptonLH(p->pdg_id()); }
486 template<> inline bool isSleptonLH(const DecodedPID& p){
487  auto pp = p.shift(1); return (p.ndigits() == 7 && p(0) == 1 && isSMLepton(pp));
488 }
489 template<> inline bool isSleptonLH(const int& p){ auto value_digits = DecodedPID(p); return isSleptonLH(value_digits);}
490 
491 
492 // APID: Super-partners of right-handed standard model leptons only
493 template<class T> inline bool isSleptonRH(const T& p) { return isSleptonRH(p->pdg_id()); }
494 template<> inline bool isSleptonRH(const DecodedPID& p){
495  auto pp = p.shift(1); return (p.ndigits() == 7 && p(0) == 2 && isSMLepton(pp));
496 }
497 template<> inline bool isSleptonRH(const int& p){ auto value_digits = DecodedPID(p); return isSleptonRH(value_digits);}
498 
499 
500 // APID: Super-partners of gauge bosons including gravitons
501 template<class T> inline bool isGaugino(const T& p) { return isGaugino(p->pdg_id()); }
502 template<> inline bool isGaugino(const DecodedPID& p){
503  auto pp = p.shift(1); return (p.ndigits() == 7 && p(0) == 1 && isBoson(pp.pid()));
504 }
505 template<> inline bool isGaugino(const int& p){ auto value_digits = DecodedPID(p); return isGaugino(value_digits);}
506 
507 
508 // APID: Super-partners of fundamental particles
509 template<class T> inline bool isSuperpartner(const T& p) { return isSuperpartner(p->pdg_id()); }
510 template<> inline bool isSuperpartner(const DecodedPID& p){
511  return isSlepton(p) || isSquark(p) || isGaugino(p);
512 }
513 template<> inline bool isSuperpartner(const int& p){ auto value_digits = DecodedPID(p); return isSuperpartner(value_digits);}
514 
515 
521 template<class T> inline bool isTechnicolor(const T& p){return isTechnicolor(p->pdg_id());}
522 template <>
523 inline bool isTechnicolor(const DecodedPID& p) {
524  const auto& pp = (p.ndigits() == 7) ? p.shift(2) : DecodedPID(0);
525  return (p.ndigits() == 7 && p(0) == 3 && (p(1) == 0 || p(0) == 1) &&
526  (isQuark(pp) || isLepton(pp) || isBoson(pp) || isGlueball(pp) ||
527  isDiquark(pp) || isHadron(pp)));
528 }
529 template<> inline bool isTechnicolor(const int& p){ auto value_digits = DecodedPID(p); return isTechnicolor(value_digits);}
530 
533 template<class T> inline bool isExcited(const T& p){return isExcited(p->pdg_id());}
534 template <>
535 inline bool isExcited(const DecodedPID& p) {
536  const auto& pp = (p.ndigits() == 7) ? p.shift(2) : DecodedPID(0);
537  return (p.ndigits() == 7 && (p(0) == 4 && p(1) == 0) &&
538  (isLepton(pp) || isQuark(pp)));
539 }
540 template<> inline bool isExcited(const int& p){ auto value_digits = DecodedPID(p); return isExcited(value_digits);}
541 
554 
559 template<class T> inline bool isRGlueball(const T& p) { return isRGlueball(p->pdg_id()); }
560 template<> inline bool isRGlueball(const DecodedPID& p) {
561  if (p.ndigits() != 7 || p(0) != 1) return false;
562  auto pp = p.shift(1);
563  return
564  ( ( pp.ndigits() == 3 && pp(0) == COMPOSITEGLUON && pp(1) == COMPOSITEGLUON && (pp.last() == 1 || pp.last() == 3) ) ||
565  ( pp.ndigits() == 4 && pp(0) == COMPOSITEGLUON && pp(1) == COMPOSITEGLUON && pp(2) == COMPOSITEGLUON && (pp.last() == 1 || pp.last() == 5) ) );
566 }
567 template<> inline bool isRGlueball(const int& p) { auto value_digits = DecodedPID(p); return isRGlueball(value_digits); }
568 
569 // APID Define R-Mesons as gluino-quark-antiquark and squark-antiquark bound states (ignore 4th generation squarks/quarks)
570 // NB Current models only allow gluino-quark-antiquark, stop-antiquark and sbottom-antiquark states
571 template<class T> inline bool isRMeson(const T& p) { return isRMeson(p->pdg_id()); }
572 template<> inline bool isRMeson(const DecodedPID& p) {
573  if (!(p.ndigits() == 7 && (p(0) == 1 || p(0) == 2))) return false;
574  auto pp = p.shift(1);
575  return (
576  // Handle ~gluino-quark-antiquark states
577  (pp.ndigits() == 4 && p(0) == 1 && pp(0) == COMPOSITEGLUON && pp.min_digit(1,3) > 0 && pp.max_digit(1,3) < QUARK_LIMIT && pp(2) <= pp(1) && (pp.last() == 1 || pp.last() == 3)) ||
578  // Handle squark-antiquark states (previously called Smeson/mesoninos)
579  (pp.ndigits() == 3 && pp.min_digit(1,3) > 0 && pp.max_digit(1,3) < QUARK_LIMIT && pp(1) <= pp(0) && pp.last() == 2)
580  );
581 }
582 template<> inline bool isRMeson(const int& p) { auto value_digits = DecodedPID(p); return isRMeson(value_digits); }
583 
584 // APID Define R-Baryons as gluino-quark-quark-quark and squark-quark-quark bound states (ignore 4th generation squarks/quarks)
585 // NB Current models only allow gluino-quark-quark-quark, stop-quark-quark and sbottom-quark-quark states
586 template<class T> inline bool isRBaryon(const T& p) { return isRBaryon(p->pdg_id()); }
587 template<> inline bool isRBaryon(const DecodedPID& p) {
588  if (!(p.ndigits() == 7 && (p(0) == 1 || p(0) == 2))) return false;
589  auto pp = p.shift(1);
590  return (
591  // Handle ~gluino-quark-quark-quark states
592  (pp.ndigits() == 5 && p(0) == 1 && pp(0) == COMPOSITEGLUON && pp.min_digit(1,4) > 0 && pp.max_digit(1,4) < QUARK_LIMIT && pp(2) <= pp(1) && pp(3) <= pp(2) && (pp.last() == 2 || pp.last() == 4)) ||
593  // Handle squark-quark-quark states (previously called Sbaryons)
594  (pp.ndigits() == 4 && pp.min_digit(1,4) > 0 && pp.max_digit(1,4) < QUARK_LIMIT && pp(1) <= pp(0) && pp(2) <= pp(1) && (pp.last() == 1 || pp.last() == 3))
595  );
596 }
597 template<> inline bool isRBaryon(const int& p) { auto value_digits = DecodedPID(p); return isRBaryon(value_digits); }
598 
599 
600 template<class T> inline bool isRHadron(const T& p) { return isRHadron(p->pdg_id()); }
601 template<> inline bool isRHadron(const DecodedPID& p) {
602  return (isRBaryon(p) || isRMeson(p) || isRGlueball(p));
603 }
604 template<> inline bool isRHadron(const int& p) { auto value_digits = DecodedPID(p); return isRHadron(value_digits); }
605 
606 
607 template<class T> inline bool hasSquark(const T& p, const int& q) { return hasSquark(p->pdg_id(), q); }
608 template<> inline bool hasSquark(const DecodedPID& p, const int& q){
609  auto pp = p.shift(1); return (
610  (isSquark(p) || isRHadron(p))
611  && pp.ndigits() != 2 // skip lepton and boson super-partners by vetoing ndigits==2
612  && pp(0) == q // After shifting, the first digit will always represent the squark in R-Hadron (and squark) PIDs
613  );
614 }
615 template<> inline bool hasSquark(const int& p, const int& q){ auto value_digits = DecodedPID(p); return hasSquark(value_digits, q);}
616 
617 
618 // APID Define isSUSY to return true for Superpartners of standard model particles and R-Hadrons
619 // TODO Should this include the MSSM extended Higgs sector??
620 template<class T> inline bool isSUSY(const T& p){return isSUSY(p->pdg_id());}
621 template<> inline bool isSUSY(const DecodedPID& p){return (isSuperpartner(p) || isRHadron(p));}
622 template<> inline bool isSUSY(const int& p){ auto value_digits = DecodedPID(p); return isSUSY(value_digits);}
623 
624 
632 template<class T> inline bool isKK(const T& p){return isKK(p->pdg_id());}
633 template<> inline bool isKK(const DecodedPID& p){return (p.ndigits() == 7 && (p(0) == 5 || p(0) == 6 ) );}
634 template<> inline bool isKK(const int& p){ auto value_digits = DecodedPID(p); return isKK(value_digits);}
635 
642 template<class T> inline bool isMonopole(const T& p){return isMonopole(p->pdg_id());}
643 template<> inline bool isMonopole(const DecodedPID& p){return (p.ndigits() == 7 && p(0) == 4 && p(1) == 1 && (p(2) == 1 || p(2) == 2 ) && p(6) == 0);}
644 template<> inline bool isMonopole(const int& p){ auto value_digits = DecodedPID(p); return isMonopole(value_digits);}
645 
656 template<class T> inline bool isDM(const T& p){return isDM(p->pdg_id());}
657 template<> inline bool isDM(const int& p){ auto sp = std::abs(p); return (sp >= 51 && sp <= 60) || sp == DARKPHOTON; }
658 
663 template<class T> inline bool isHiddenValley(const T& p){return isHiddenValley(p->pdg_id());}
664 template <>
665 inline bool isHiddenValley(const DecodedPID& p) {
666  const auto& pp = (p.ndigits() == 7) ? p.shift(2) : DecodedPID(0);
667  return (p.ndigits() == 7 && p(0) == 4 && p(1) == 9 &&
668  (isQuark(pp) || isLepton(pp) || isBoson(pp) || isGlueball(pp) ||
669  isDiquark(pp) || isHadron(pp)));
670 }
671 template<> inline bool isHiddenValley(const int& p){ auto value_digits = DecodedPID(p); return isHiddenValley(value_digits);}
672 
679 template<class T> inline bool isGenericMultichargedParticle(const T& p){return isGenericMultichargedParticle(p->pdg_id());}
680 template<> inline bool isGenericMultichargedParticle(const DecodedPID& p){return (p.ndigits() == 8 && (p(0) == 1 || p(0) == 2) && p(1) == 0 && p(2) == 0 && p(7) == 0);}
681 template<> inline bool isGenericMultichargedParticle(const int& p){ auto value_digits = DecodedPID(p); return isGenericMultichargedParticle(value_digits);}
682 
697 template<class T> inline bool isNucleus(const T& p){return isNucleus(p->pdg_id());}
698 template<> inline bool isNucleus(const DecodedPID& p){
699  if (std::abs(p.pid()) == PROTON) return true;
700  return (p.ndigits() == 10 && p(0) == 1 && p(1) == 0 );
701 }
702 template<> inline bool isNucleus(const int& p){ auto value_digits = DecodedPID(p); return isNucleus(value_digits);}
703 
704 
705 template<class T> inline bool hasQuark(const T& p, const int& q);
706 template<> inline bool hasQuark(const DecodedPID& p, const int& q){
707  if (isQuark(p.pid())) { return (std::abs(p.pid()) == q );}
708  if (isMeson(p)) { return *(p.second.rbegin() + 1) == q ||*(p.second.rbegin()+2) ==q;}
709  if (isDiquark(p)) { auto i = std::find(p.second.rbegin() + 2,p.second.rbegin()+4,q); return (i!=p.second.rbegin()+4);}
710  if (isBaryon(p)) { auto i = std::find(p.second.rbegin() + 1,p.second.rbegin()+4,q); return (i!=p.second.rbegin()+4);}
711  if (isTetraquark(p)) { auto i = std::find(p.second.rbegin() + 1,p.second.rbegin()+5,q); return (i!=p.second.rbegin()+5);}
712  if (isPentaquark(p)) { auto i = std::find(p.second.rbegin() + 1,p.second.rbegin()+6,q); return (i!=p.second.rbegin()+6);}
713  if (isNucleus(p) && std::abs(p.pid()) != PROTON) { return (q == 1 || q == 2 || (q==3 && p(2) > 0));}
714  if (isSUSY(p)) { // APID SUSY case
715  auto pp = p.shift(1);
716  if ( pp.ndigits() == 1 ) { return false; } // Handle squarks
717  if ( pp.ndigits() == 3 ) { return (pp(1) == q); } // Handle ~q qbar pairs
718  if ( pp.ndigits() == 4 ) { return (pp(1) == q || pp(2) == q); } // Ignore gluinos and squarks
719  if ( pp.ndigits() == 5 ) { return (pp(1) == q || pp(2) == q || pp(3) == q); } // Ignore gluinos and squarks
720  if ( pp.ndigits() > 5 ) { pp = pp.shift(1); } // Drop gluinos and squarks
721  return hasQuark(pp, q); }
722  return false;
723 }
724 template<> inline bool hasQuark(const int& p, const int& q){ auto value_digits = DecodedPID(p); return hasQuark(value_digits, q);}
725 
726 template<class T> inline bool hasStrange(const T& p) { return hasQuark(p,SQUARK); }
727 template<class T> inline bool hasCharm(const T& p) { return hasQuark(p,CQUARK); }
728 template<class T> inline bool hasBottom(const T& p) { return hasQuark(p,BQUARK); }
729 template<class T> inline bool hasTop(const T& p) { return hasQuark(p,TQUARK); }
730 
731 
732 // APID: The baryon number is defined as:
733 // B = (1/3)*( n_q - n_{qbar} )
734 // where n_q⁠ is the number of quarks, and ⁠n_{qbar} is the number of
735 // antiquarks. By convention, squarks have the same quantum numbers as
736 // the corresponding quarks (modulo spin and R), so have baryon number
737 // 1/3.
738 template<class T> inline int baryonNumber3(const T& p) {return baryonNumber3(p->pdg_id());}
739 template<> inline int baryonNumber3(const DecodedPID& p){
740  if (isQuark(p.pid())) { return (p.pid() > 0) ? 1 : - 1;}
741  if (isDiquark(p)) { return (p.pid() > 0) ? 2 : -2; }
742  if (isMeson(p) || isTetraquark(p)) { return 0; }
743  if (isBaryon(p) || isPentaquark(p)){ return (p.pid() > 0) ? 3 : -3; }
744  if (isNucleus(p)) {
745  const int result = 3*p(8) + 30*p(7) + 300*p(6);
746  return (p.pid() > 0) ? result : -result;
747  }
748  if (isSUSY(p)) {
749  auto pp = p.shift(1);
750  if (pp.ndigits() < 3 ) { return baryonNumber3(pp); } // super-partners of fundamental particles
751  if (pp(0) == COMPOSITEGLUON) {
752  if (pp(1) == COMPOSITEGLUON) { return 0; } // R-Glueballs
753  if ( pp.ndigits() == 4 ) { return 0; } // states with gluino-quark-antiquark
754  if ( pp.ndigits() == 5) { return (p.pid() > 0) ? 3 : -3; } // states with gluino-quark-quark-quark
755  }
756  if (pp.ndigits() == 3) { return 0; } // squark-antiquark
757  if (pp.ndigits() == 4) { return (p.pid() > 0) ? 3 : -3; } // states with squark-quark-quark
758  }
759  return 0;
760 }
761 template<> inline int baryonNumber3(const int& p){ auto value_digits = DecodedPID(p); return baryonNumber3(value_digits);}
762 
763 template<class T> inline double baryonNumber(const T& p) {return baryonNumber(p->pdg_id());}
764 template<> inline double baryonNumber(const DecodedPID& p){ return static_cast<double>(baryonNumber3(p))/3.0;}
765 template<> inline double baryonNumber(const int& p){ auto value_digits = DecodedPID(p); return static_cast<double>(baryonNumber3(value_digits))/3.0;}
766 
767 
768 // APID: The strangeness of a particle is defined as:
769 // S = − ( n_s − n_{sbar} )
770 // where n_s represents the number of strange quarks and n_{sbar}
771 // represents the number of strange antiquarks. By convention, strange
772 // squarks have the same quantum numbers as strange quarks (modulo
773 // spin and R), so have strangeness -1.
774 static const std::array<int,10> is_strange = {
775  +0, +0, +0, -1, +0, +0, +0, +0, +0, +0 };
776 template<class T> inline int strangeness(const T& p) {return strangeness(p->pdg_id());}
777 template<> inline int strangeness(const DecodedPID& p){
778  if (isNucleus(p) && p.ndigits() == 10) { return (p.pid() > 0) ? -p(2) : p(2); }
779  if (isStrange(p.pid())) { return (p.pid() > 0) ? -1 : 1; }
780  if (!hasStrange(p) && !hasSquark(p,SQUARK)) { return 0; }
781  if (std::abs(p.pid()) == K0) { return (p.pid() > 0) ? 1 : -1; }
782  size_t nq = 0;
783  int sign = 1;
784  int signmult = 1;
785  int result=0;
786  bool classified = false;
787  if (!classified && isMeson(p)) { classified = true; nq = 2; if ((*(p.second.rbegin()+2)) == 2||(*(p.second.rbegin()+2)) == 4 ) { sign=-1;} signmult =-1; }
788  if (!classified && isDiquark(p)) {return is_strange.at(p(0))+is_strange.at(p(1)); }
789  if (!classified && isBaryon(p)) { classified = true; nq = 3; }
790  if (!classified && isTetraquark(p)){ return is_strange.at(p(3)) + is_strange.at(p(4)) - is_strange.at(p(6)) - is_strange.at(p(7)); }
791  if (!classified && isPentaquark(p)){ return is_strange.at(p(3)) + is_strange.at(p(4)) + is_strange.at(p(5)) + is_strange.at(p(6)) - is_strange.at(p(7)); }
792  if (!classified && isSUSY(p)) {
793  nq = 0;
794  auto pp = p.shift(1);
795  if (pp.ndigits() < 3 ) { return strangeness(pp); } // super-partners of fundamental particles
796  if (pp(0) == COMPOSITEGLUON) {
797  if (pp(1) == COMPOSITEGLUON) { return 0; } // R-Glueballs
798  if ( pp.ndigits() == 4 || pp.ndigits() == 5) {
799  pp = pp.shift(1); // Remove gluino
800  }
801  }
802  if (pp.ndigits() == 3) { classified = true; nq = 2; if (p.last()%2==0) {sign = -1;} signmult = -1; } // states with quark-antiquark or squark-antiquark
803  if (pp.ndigits() == 4) { classified = true; nq = 3; } // states with quark-quark-quark or squark-quark-quark
804  }
805  for (auto r = p.second.rbegin() + 1; r != p.second.rbegin() + 1 + nq; ++r) {
806  result += is_strange.at(*r)*sign;
807  sign*=signmult;
808  }
809  return p.pid() > 0 ? result : -result;
810 }
811 template<> inline int strangeness(const int& p){ auto value_digits = DecodedPID(p); return strangeness(value_digits);}
812 
813 
814 template<class T> inline int numberOfLambdas(const T& p) {return numberOfLambdas(p->pdg_id());}
815 template<> inline int numberOfLambdas(const DecodedPID& p){
816  if (std::abs(p.pid()) == LAMBDA0) { return (p.pid() > 0) ? 1 : -1; }
817  if (isNucleus(p) && p.ndigits() == 10) { return (p.pid() > 0) ? p(2) : -p(2); }
818  return 0;
819 }
820 template<> inline int numberOfLambdas(const int& p){ auto value_digits = DecodedPID(p); return numberOfLambdas(value_digits);}
821 
822 
823 template<class T> inline int numberOfProtons(const T& p) {return numberOfProtons(p->pdg_id());}
824 template<> inline int numberOfProtons(const DecodedPID& p){
825  if (std::abs(p.pid()) == PROTON) { return (p.pid() > 0) ? 1 : -1; }
826  if (isNucleus(p)) {
827  const int result = p(5) + 10*p(4) + 100*p(3);
828  return (p.pid() > 0) ? result : -result;
829  }
830  return 0;
831 }
832 template<> inline int numberOfProtons(const int& p){ auto value_digits = DecodedPID(p); return numberOfProtons(value_digits);}
833 
834 
836 template<class T> inline bool isBSM(const T& p){return isBSM(p->pdg_id());}
837 template<> inline bool isBSM(const DecodedPID& p){
838  if (p.pid() == GRAVITON || std::abs(p.pid()) == MAVTOP || p.pid() == DARKPHOTON) return true;
839  if (std::abs(p.pid()) > 16 && std::abs(p.pid()) < 19) return true;
840  if (std::abs(p.pid()) > 31 && std::abs(p.pid()) < 39) return true;
841  if (std::abs(p.pid()) > 39 && std::abs(p.pid()) < 81) return true;
842  if (std::abs(p.pid()) > 6 && std::abs(p.pid()) < 9) return true;
843  if (isSUSY(p)) return true;
844  if (isNeutrinoRH(p.pid())) return true;
845  if (isGenericMultichargedParticle(p)) return true;
846  if (isTechnicolor(p)) return true;
847  if (isExcited(p)) return true;
848  if (isKK(p)) return true;
849  if (isHiddenValley(p)) return true;
850  if (isMonopole(p)) return true;
851  return false;
852 }
853 template<> inline bool isBSM(const int& p){
854  if (p == GRAVITON || std::abs(p) == MAVTOP || p == DARKPHOTON) return true;
855  if (std::abs(p) > 16 && std::abs(p) < 19) return true;
856  if (std::abs(p) > 31 && std::abs(p) < 38) return true;
857  if (std::abs(p) > 39 && std::abs(p) < 81) return true;
858  if (std::abs(p) > 6 && std::abs(p) < 9) return true;
859  auto value_digits = DecodedPID(p); return isBSM(value_digits);
860 }
861 
862 template<class T> inline bool isTransportable(const T& p){return isTransportable(p->pdg_id());}
863 template<> inline bool isTransportable(const DecodedPID& p){ return isPhoton(p.pid()) || isGeantino(p.pid()) || isHadron(p) || isLepton(p.pid()) || p.pid() == DARKPHOTON;}
864 template<> inline bool isTransportable(const int& p){ auto value_digits = DecodedPID(p); return isTransportable(value_digits);}
865 
867 template<class T> inline bool isValid(const T& p){return isValid(p->pdg_id());}
868 template<> inline bool isValid(const DecodedPID& p){
869  return p.pid() !=0 && ( isQuark(p) || isLepton(p) || isBoson(p) || isGlueball(p) ||
870  isTrajectory(p.pid()) || isGenSpecific(p.pid()) || isDiquark(p) ||
871  isBSM(p) || isHadron(p) || isNucleus(p) || isGeantino(p.pid()) ||
872  isPythia8Specific(p) ); }
873 template<> inline bool isValid(const int& p){ if (!p) return false; if (std::abs(p) < 42) return true;
874  if (isGenSpecific(p)) return true;
875  auto value_digits = DecodedPID(p); return isValid(value_digits);
876 }
877 
878 template<class T> inline int leadingQuark(const T& p) {return leadingQuark(p->pdg_id());}
879 template<> inline int leadingQuark(const DecodedPID& p){
880  if (isQuark(p.pid())) { return std::abs(p.pid());}
881  if (isMeson(p)) { return p.max_digit(1,3);}
882  if (isDiquark(p)) { return p.max_digit(2,4);}
883  if (isBaryon(p)) { return p.max_digit(1,4);}
884  if (isTetraquark(p)) { return p.max_digit(1,5);}
885  if (isPentaquark(p)) { return p.max_digit(1,6);}
886  if (isSUSY(p)) { // APID SUSY case
887  auto pp = p.shift(1);
888  if ( pp.ndigits() == 1 ) { return 0; } // Handle squarks
889  if ( pp.ndigits() == 3 ) { pp = DecodedPID(pp(1)); } // Handle ~q qbar pairs
890  if ( pp.ndigits() > 3 ) { pp = pp.shift(1); } // Drop gluinos and squarks
891  return leadingQuark(pp); }
892  return 0;
893 }
894 
895 template<> inline int leadingQuark(const int& p){ auto value_digits = DecodedPID(p); return leadingQuark(value_digits);}
896 
897 template<class T> inline bool isLightHadron(const T& p) { auto lq = leadingQuark(p); return (lq == DQUARK || lq == UQUARK||lq == SQUARK) && isHadron(p); }
898 template<class T> inline bool isHeavyHadron(const T& p) { auto lq = leadingQuark(p); return (lq == CQUARK || lq == BQUARK || lq == TQUARK ) && isHadron(p); }
899 template<class T> inline bool isStrangeHadron(const T& p) { return leadingQuark(p) == SQUARK && isHadron(p); }
900 template<class T> inline bool isCharmHadron(const T& p) { return leadingQuark(p) == CQUARK && isHadron(p); }
901 template<class T> inline bool isBottomHadron(const T& p) { return leadingQuark(p) == BQUARK && isHadron(p); }
902 template<class T> inline bool isTopHadron(const T& p) { return leadingQuark(p) == TQUARK && isHadron(p); }
903 
904 template<class T> inline bool isLightMeson(const T& p) { auto lq = leadingQuark(p); return (lq == DQUARK || lq == UQUARK||lq == SQUARK) && isMeson(p); }
905 template<class T> inline bool isHeavyMeson(const T& p) { auto lq = leadingQuark(p); return (lq == CQUARK || lq == BQUARK || lq == TQUARK) && isMeson(p); }
906 template<class T> inline bool isStrangeMeson(const T& p) { return leadingQuark(p) == SQUARK && isMeson(p); }
907 template<class T> inline bool isCharmMeson(const T& p) { return leadingQuark(p) == CQUARK && isMeson(p); }
908 template<class T> inline bool isBottomMeson(const T& p) { return leadingQuark(p) == BQUARK && isMeson(p); }
909 template<class T> inline bool isTopMeson(const T& p) { return leadingQuark(p) == TQUARK && isMeson(p); }
910 
911 template<class T> inline bool isCCbarMeson(const T& p) { return isCCbarMeson(p->pdg_id());}
912 template<> inline bool isCCbarMeson(const DecodedPID& p) { return leadingQuark(p) == CQUARK && isMeson(p) && (*(p.second.rbegin()+2)) == CQUARK && (*(p.second.rbegin()+1)) == CQUARK; }
913 template<> inline bool isCCbarMeson(const int& p) { return isCCbarMeson(DecodedPID(p)); }
914 
915 template<class T> inline bool isBBbarMeson(const T& p){ return isBBbarMeson(p->pdg_id());}
916 template<> inline bool isBBbarMeson(const DecodedPID& p) { return leadingQuark(p) == BQUARK && isMeson(p) && (*(p.second.rbegin()+2)) == BQUARK && (*(p.second.rbegin()+1)) == BQUARK; }
917 template<> inline bool isBBbarMeson(const int& p) { return isBBbarMeson(DecodedPID(p)); }
918 
919 
920 template<class T> inline bool isLightBaryon(const T& p) { auto lq = leadingQuark(p); return (lq == DQUARK || lq == UQUARK||lq == SQUARK) && isBaryon(p); }
921 template<class T> inline bool isHeavyBaryon(const T& p) { auto lq = leadingQuark(p); return (lq == CQUARK || lq == BQUARK || lq == TQUARK) && isBaryon(p); }
922 template<class T> inline bool isStrangeBaryon(const T& p) { return leadingQuark(p) == SQUARK && isBaryon(p); }
923 template<class T> inline bool isCharmBaryon(const T& p) { return leadingQuark(p) == CQUARK && isBaryon(p); }
924 template<class T> inline bool isBottomBaryon(const T& p) { return leadingQuark(p) == BQUARK && isBaryon(p); }
925 template<class T> inline bool isTopBaryon(const T& p) { return leadingQuark(p) == TQUARK && isBaryon(p); }
926 
927 
928 // APID: This function selects B-Hadrons which predominantly decay weakly. (Commonly used definition in GeneratorFilters package.)
929 // 5[1-4]1 L = J = 0, S = 0
930 // 5[1-5][1-4]2 J = 1/2, n_r = 0, n_L =0
931 template<class T> inline bool isWeaklyDecayingBHadron(const T& p) {return isWeaklyDecayingBHadron(p->pdg_id());}
932 template<> inline bool isWeaklyDecayingBHadron(const int& p) {
933  const int pid = std::abs(p);
934  return ( pid == 511 || // B0
935  pid == 521 || // B+
936  pid == 531 || // B_s0
937  pid == 541 || // B_c+
938  pid == 5122 || // Lambda_b0
939  pid == 5132 || // Xi_b-
940  pid == 5232 || // Xi_b0
941  pid == 5112 || // Sigma_b-
942  pid == 5212 || // Sigma_b0
943  pid == 5222 || // Sigma_b+
944  pid == 5332 || // Omega_b-
945  pid == 5142 || // Xi_bc0
946  pid == 5242 || // Xi_bc+
947  pid == 5412 || // Xi'_bc0
948  pid == 5422 || // Xi'_bc+
949  pid == 5342 || // Omega_bc0
950  pid == 5432 || // Omega'_bc0
951  pid == 5442 || // Omega_bcc+
952  pid == 5512 || // Xi_bb-
953  pid == 5522 || // Xi_bb0
954  pid == 5532 || // Omega_bb-
955  pid == 5542 ); // Omega_bbc0
956 }
957 template<> inline bool isWeaklyDecayingBHadron(const DecodedPID& p){ return isWeaklyDecayingBHadron(p.pid()); }
958 
959 
960 // APID: This function selects C-Hadrons which predominantly decay weakly. (Commonly used definition in GeneratorFilters package.)
961 // 4[1-3]1 L = J = 0, S = 0
962 // 4[1-4][1-3]2 J = 1/2, n_r = 0, n_L =0
963 // NB Omitting pid = 4322 (Xi'_C+) a this undergoes an EM rather than
964 // weak decay. (There was an old version of Herwig that decayed it
965 // weakly, but this was fixed in Herwig 7.)
966 template<class T> inline bool isWeaklyDecayingCHadron(const T& p) {return isWeaklyDecayingCHadron(p->pdg_id());}
967 template<> inline bool isWeaklyDecayingCHadron(const int& p) {
968  const int pid = std::abs(p);
969  return ( pid == 411 || // D+
970  pid == 421 || // D0
971  pid == 431 || // Ds+
972  pid == 4122 || // Lambda_c+
973  pid == 4132 || // Xi_c0
974  pid == 4232 || // Xi_c+
975  pid == 4212 || // Xi_c0
976  pid == 4332 || // Omega_c0
977  pid == 4412 || // Xi_cc+
978  pid == 4422 || // Xi_cc++
979  pid == 4432 ); // Omega_cc+
980 }
981 template<> inline bool isWeaklyDecayingCHadron(const DecodedPID& p){ return isWeaklyDecayingCHadron(p.pid()); }
982 
983 
984 template<class T> inline int charge3( const T& p){return charge3(p->pdg_id());}
985 template<class T> inline double fractionalCharge(const T& p){return fractionalCharge(p->pdg_id());}
986 template<class T> inline double charge( const T& p){
987  if (isGenericMultichargedParticle(p)) // BSM multi-charged particles might have a fractional charge that's not a multiple of 1/3
988  return fractionalCharge(p);
989  else
990  return 1.0*charge3(p)/3.0;
991 }
992 template<class T> inline double threeCharge( const T& p){ return charge3(p);}
993 template<class T> inline bool isCharged( const T& p){ return charge3(p) != 0;}
994 
995 
996 template<> inline int charge3(const DecodedPID& p) {
997  auto ap = std::abs(p.pid());
998  if (ap < TABLESIZE ) return p.pid() > 0 ? triple_charge.at(ap) : -triple_charge.at(ap);
999  if (ap == K0) return 0;
1000  if (ap == GEANTINO0) return 0;
1001  if (ap == GEANTINOPLUS) return p.pid() > 0 ? 3 : -3;
1002  if (ap == MAVTOP) return p.pid() > 0 ? 2 : -2;
1003  size_t nq = 0;
1004  int sign = 1;
1005  int signmult = 1;
1006  int result=0;
1007  bool classified = false;
1008  if (!classified && isMeson(p)) { classified = true; nq = 2; if ((*(p.second.rbegin()+2)) == 2||(*(p.second.rbegin()+2)) == 4 ) { sign=-1;} signmult =-1; }
1009  if (!classified && isDiquark(p)) {return triple_charge.at(p(0))+triple_charge.at(p(1)); }
1010  if (!classified && isBaryon(p)) { classified = true; nq = 3; }
1011  if (!classified && isTetraquark(p)){ return triple_charge.at(p(3)) + triple_charge.at(p(4)) - triple_charge.at(p(6)) - triple_charge.at(p(7)); }
1012  if (!classified && isPentaquark(p)){ return triple_charge.at(p(3)) + triple_charge.at(p(4)) + triple_charge.at(p(5)) + triple_charge.at(p(6)) - triple_charge.at(p(7)); }
1013  if (!classified && isNucleus(p)) { return 3*numberOfProtons(p);}
1014  if (!classified && isSUSY(p)) {
1015  nq = 0;
1016  auto pp = p.shift(1);
1017  if (pp.ndigits() < 3 ) { return charge3(pp); } // super-partners of fundamental particles
1018  if (pp(0) == COMPOSITEGLUON) {
1019  if (pp(1) == COMPOSITEGLUON) { return 0; } // R-Glueballs
1020  if ( pp.ndigits() == 4 || pp.ndigits() == 5) {
1021  pp = pp.shift(1); // Remove gluino
1022  }
1023  }
1024  if (pp.ndigits() == 3) { classified = true; nq = 2; if (p.last()%2==0) {sign = -1;} signmult = -1; } // states with squark-antiquark or quark-anti-quark
1025  if (pp.ndigits() == 4) { classified = true; nq = 3; } // states with squark-quark-quark or quark-quark-quark
1026  }
1027  if (!classified && isMonopole(p)) {
1030  result = 3*(p(3)*100 + p(4)*10 + p(5));
1031  return ( (p.pid() > 0 && p(2) == 1) || (p.pid() < 0 && p(2) == 2) ) ? result : -result;
1032  }
1033  if (!classified && isGenericMultichargedParticle(p)) {
1034  double abs_charge = 0.0;
1035  if (p(0) == 1) abs_charge = p(3)*100. + p(4)*10. + p(5)*1 + p(6)*0.1; // multi-charged particle PDG ID is +/-100XXXY0, where the charge is XXX.Y
1036  if (p(0) == 2) abs_charge = (p(3)*10. + p(4))/(p(5)*10.0 + p(6)); // multi-charged particle PDG ID is +/-200XXYY0, where the charge is XX/YY
1037  int abs_threecharge = static_cast<int>(std::round(abs_charge * 3.)); // the multi-charged particles might have a fractional charge that's not a multiple of 1/3, in that case round to the closest multiple of 1/3 for charge3 and threecharge
1038  return p.pid() > 0 ? abs_threecharge : -1 * abs_threecharge;
1039  }
1040  for (auto r = p.second.rbegin() + 1; r != p.second.rbegin() + 1 + nq; ++r) {
1041  result += triple_charge.at(*r)*sign;
1042  sign*=signmult;
1043  }
1044  return p.pid() > 0 ? result : -result;
1045 }
1046 template<> inline int charge3(const int& p){
1047  int ap = std::abs(p);
1048  if (ap < TABLESIZE) return p > 0 ? triple_charge.at(ap):-triple_charge.at(ap);
1049  auto value_digits = DecodedPID(p);
1050  return charge3(value_digits);
1051 }
1052 
1053 
1054 template<class T> inline bool isNeutral( const T& p){ return p->pdg_id() != 0 && charge3(p) == 0;}
1055 template<> inline bool isNeutral(const DecodedPID& p){ return p.pid() != 0 && charge3(p) == 0;}
1056 template<> inline bool isNeutral(const int& p){ auto value_digits = DecodedPID(p); return isNeutral(value_digits);}
1057 
1058 
1059 template<> inline double fractionalCharge(const DecodedPID& p) {
1060  if(!isGenericMultichargedParticle(p)) return 1.0*charge3(p)/3.0; // this method is written for multi-charged particles, still make sure other cases are handled properly
1061  double abs_charge = 0;
1062  if (p(0) == 1) abs_charge = p(3)*100. + p(4)*10. + p(5)*1 + p(6)*0.1; // multi-charged particle PDG ID is +/-100XXXY0, where the charge is XXX.Y
1063  if (p(0) == 2) abs_charge = (p(3)*10. + p(4))/(p(5)*10.0 + p(6)); // multi-charged particle PDG ID is +/-200XXYY0, where the charge is XX/YY
1064  return p.pid() > 0 ? abs_charge : -1 * abs_charge;
1065 }
1066 template<> inline double fractionalCharge(const int& p){auto value_digits = DecodedPID(p); return fractionalCharge(value_digits);}
1067 
1068 // APID: Including Z' and Z'' as EM interacting.
1069 template<class T> inline bool isEMInteracting(const T& p){return isEMInteracting(p->pdg_id());}
1070 template<> inline bool isEMInteracting(const int& p) {return (isPhoton(p) || isZ(p) || p == ZPRIME || p == ZDBLPRIME || std::abs(charge(p))>std::numeric_limits<double>::epsilon() || isMonopole(p));}
1071 
1072 template<class T> inline bool isParton(const T& p) { return isQuark(p)||isGluon(p);}
1073 
1074 // APID: Intended to return 2J
1075 // Useful for G4ParticleDefinition constructor
1076 template<class T> inline int spin2(const T& p) { return spin2(p->pdg_id()); }
1077 template<> inline int spin2(const DecodedPID& p) {
1078  if (isSUSY(p)) {
1079  auto pp = p.shift(1);
1080  auto ap = std::abs(pp.pid());
1081  if (ap < TABLESIZE ) { return std::abs(double_spin.at(ap)-1); } // sparticles (0->1, 1 -> 0, 2->1, 4->3)
1082  return p.last()-1; // R-Hadrons (p.last() == 2J +1)
1083  }
1084  auto ap = std::abs(p.pid());
1085  if (ap == K0S) { return 0; }
1086  if (ap == K0L) { return 0; }
1087  if (ap == MAVTOP) { return 1; } // TODO check this
1088  if (ap == DARKPHOTON) { return 2; } // TODO check this
1089  if (ap < TABLESIZE ) { return double_spin.at(ap); } // fundamental particles
1090  if (isHadron(p)) { return p.last()-1; } // Hadrons (p.last == 2J+1 - special cases handled above)
1091  if (isMonopole(p)) { return 0; } // PDG 11i - For now no spin information is provided. Also matches the definition in the G4Extensions/Monopole package.
1092  if (isGenericMultichargedParticle(p)) { return 0; } // APID Matches the definition in the G4Extensions/Monopole package.
1093  if (isNucleus(p)) { return 1; } // TODO need to explicitly deal with nuclei
1094  return p.last() > 0 ? 1 : 0; // Anything else - best guess
1095 }
1096 template<> inline int spin2(const int& p){ auto value_digits = DecodedPID(p); return spin2(value_digits);}
1097 
1098 template<class T> inline double spin(const T& p) { return spin(p->pdg_id()); }
1099 template<> inline double spin(const DecodedPID& p) { return 1.0*spin2(p)/2.0; }
1100 template<> inline double spin(const int& p){ auto value_digits = DecodedPID(p); return spin(value_digits);}
1101 
1102 
1103 // APID: Returns an unordered list of the quarks contained by the current particle
1104 template<class T> inline std::vector<int> containedQuarks(const T& p) { return containedQuarks(p->pdg_id()); }
1105 template<> inline std::vector<int> containedQuarks(const int& p) {
1106  auto pp = DecodedPID(p);
1107  std::vector<int> quarks;
1108  if (isQuark(pp.pid())) { quarks.push_back(std::abs(pp.pid())); }
1109  else if (isDiquark(pp)) { quarks.push_back(pp(0)); quarks.push_back(pp(1)); }
1110  else if (isMeson(pp)) { quarks.push_back(*(pp.second.rbegin() + 1)); quarks.push_back(*(pp.second.rbegin()+2)); }
1111  else if (isBaryon(pp)) { for (size_t digit = 1; digit < 4; ++digit) { quarks.push_back(*(pp.second.rbegin() + digit)); } }
1112  else if (isTetraquark(pp)) { for (size_t digit = 1; digit < 5; ++digit) { quarks.push_back(*(pp.second.rbegin() + digit)); } }
1113  else if (isPentaquark(pp)) { for (size_t digit = 1; digit < 6; ++digit) { quarks.push_back(*(pp.second.rbegin() + digit)); } }
1114  else if (isNucleus(pp)) { const int A = std::abs(baryonNumber3(pp)/3); const int Z = std::abs(numberOfProtons(pp)); const int L = std::abs(numberOfLambdas(pp));
1115  const int n_uquarks = A + Z; const int n_dquarks = 2*A - Z - L; const int n_squarks = L;
1116  quarks.reserve(3*A); quarks.insert(quarks.end(), n_dquarks, 1); quarks.insert(quarks.end(), n_uquarks, 2); quarks.insert(quarks.end(), n_squarks, 3); }
1117  else if (isSUSY(pp)) { // APID SUSY case
1118  pp = pp.shift(1);
1119  if ( pp.ndigits() > 1 ) { // skip squarks
1120  if ( pp.ndigits() == 3 ) { pp = DecodedPID(pp(1)); } // Handle ~q qbar pairs
1121  if ( pp.ndigits() > 3 ) { pp = pp.shift(1); } // Drop gluinos and squarks
1122  return containedQuarks(pp.pid());
1123  }
1124  }
1125  return quarks;
1126 }
1127 template<> inline std::vector<int> containedQuarks(const DecodedPID& p) { return containedQuarks(p.pid()); }
1128 
1129 template<class T> inline bool isStrongInteracting(const T& p){return isStrongInteracting(p->pdg_id());}
1130 template<> inline bool isStrongInteracting(const int& p) { return (isGluon(p) || isQuark(p) || isDiquark(p) || isGlueball(p) || isLeptoQuark(p) || isHadron(p) || isRHadron(p));} // APID: Glueballs and R-Hadrons are also strong-interacting
1131 
1132 #endif
isStrange
bool isStrange(const T &p)
Definition: AtlasPID.h:173
isStrangeMeson
bool isStrangeMeson(const T &p)
Definition: AtlasPID.h:906
beamspotman.r
def r
Definition: beamspotman.py:674
isGaugino
bool isGaugino(const T &p)
Definition: AtlasPID.h:501
isBottomMeson
bool isBottomMeson(const T &p)
Definition: AtlasPID.h:908
isStrongInteracting
bool isStrongInteracting(const T &p)
Definition: AtlasPID.h:1129
isNucleus
bool isNucleus(const T &p)
PDG rule 16 Nuclear codes are given as 10-digit numbers ±10LZZZAAAI.
Definition: AtlasPID.h:697
get_generator_info.result
result
Definition: get_generator_info.py:21
isHeavyBoson
bool isHeavyBoson(const T &p)
APID: Additional "Heavy"/"prime" versions of W and Z bosons (Used in MCTruthClassifier)
Definition: AtlasPID.h:383
isHeavyHadron
bool isHeavyHadron(const T &p)
Definition: AtlasPID.h:898
isTetraquark
bool isTetraquark(const T &p)
PDG rule 14 The 9-digit tetra-quark codes are ±1nrnLnq1nq20nq3nq4nJ.
Definition: AtlasPID.h:322
isRMeson
bool isRMeson(const T &p)
Definition: AtlasPID.h:571
isBottomBaryon
bool isBottomBaryon(const T &p)
Definition: AtlasPID.h:924
find
std::string find(const std::string &s)
return a remapped string
Definition: hcg.cxx:135
isSleptonRH
bool isSleptonRH(const T &p)
Definition: AtlasPID.h:493
isTopBaryon
bool isTopBaryon(const T &p)
Definition: AtlasPID.h:925
hasCharm
bool hasCharm(const T &p)
Definition: AtlasPID.h:727
isBSM
bool isBSM(const T &p)
APID: graviton and all Higgs extensions are BSM.
Definition: AtlasPID.h:836
baryonNumber
double baryonNumber(const T &p)
Definition: AtlasPID.h:763
isHiddenValley
bool isHiddenValley(const T &p)
PDG rule 11k Hidden Valley particles have n = 4 and n_r = 9, and trailing numbers in agreement with t...
Definition: AtlasPID.h:663
isBoson
bool isBoson(const T &p)
PDG rule 9: Two-digit numbers in the range 21–30 are provided for the Standard Model gauge and Higgs ...
Definition: AtlasPID.h:366
isGenSpecific
bool isGenSpecific(const T &p)
Main Table for MC internal use 81–100,901–930,998-999,1901–1930,2901–2930, and 3901–3930.
Definition: AtlasPID.h:422
Monitored::Z
@ Z
Definition: HistogramFillerUtils.h:24
DecodedPID::shift
DecodedPID shift(const size_t n) const
Definition: AtlasPID.h:25
isMeson
bool isMeson(const T &p)
Table 43.1 PDG rule 5a: The numbers specifying the meson’s quark content conform to the convention nq...
Definition: AtlasPID.h:241
hasQuark
bool hasQuark(const T &p, const int &q)
isKK
bool isKK(const T &p)
PDG rule 11h A black hole in models with extra dimensions has code 5000040.
Definition: AtlasPID.h:632
isSquarkRH
bool isSquarkRH(const T &p)
Definition: AtlasPID.h:471
threeCharge
double threeCharge(const T &p)
Definition: AtlasPID.h:992
isCharmMeson
bool isCharmMeson(const T &p)
Definition: AtlasPID.h:907
MuonGM::round
float round(const float toRound, const unsigned int decimals)
Definition: Mdt.cxx:27
isRGlueball
bool isRGlueball(const T &p)
PDG rule 11g: Within several scenarios of new physics, it is possible to have colored particles suffici...
Definition: AtlasPID.h:559
hasBottom
bool hasBottom(const T &p)
Definition: AtlasPID.h:728
numberOfLambdas
int numberOfLambdas(const T &p)
Definition: AtlasPID.h:814
python.SystemOfUnits.second
float second
Definition: SystemOfUnits.py:135
isNeutrino
bool isNeutrino(const T &p)
APID: the fourth generation neutrinos are neutrinos.
Definition: AtlasPID.h:209
isResonance
bool isResonance(const T &p)
Definition: AtlasPID.h:397
isParton
bool isParton(const T &p)
Definition: AtlasPID.h:1072
isSquark
bool isSquark(const T &p)
PDG rule 11d Fundamental supersymmetric particles are identified by adding a nonzero n to the particl...
Definition: AtlasPID.h:455
isSleptonLH
bool isSleptonLH(const T &p)
Definition: AtlasPID.h:485
isHeavyBaryon
bool isHeavyBaryon(const T &p)
Definition: AtlasPID.h:921
isValid
bool isValid(const T &p)
Av: we implement here an ATLAS-sepcific convention: all particles which are 99xxxxx are fine.
Definition: AtlasPID.h:867
checkRpcDigits.digit
digit
Definition: checkRpcDigits.py:186
isLightBaryon
bool isLightBaryon(const T &p)
Definition: AtlasPID.h:920
python.AtlRunQueryParser.ap
ap
Definition: AtlRunQueryParser.py:825
xAOD::Muon_v1
Class describing a Muon.
Definition: Muon_v1.h:38
isBottomHadron
bool isBottomHadron(const T &p)
Definition: AtlasPID.h:901
isLightMeson
bool isLightMeson(const T &p)
Definition: AtlasPID.h:904
isGenericMultichargedParticle
bool isGenericMultichargedParticle(const T &p)
In addition, there is a need to identify ”Q-ball” and similar very exotic (multi-charged) particles w...
Definition: AtlasPID.h:679
isSMLepton
bool isSMLepton(const T &p)
APID: the fourth generation leptons are not standard model leptons.
Definition: AtlasPID.h:191
isSMQuark
bool isSMQuark(const T &p)
Definition: AtlasPID.h:169
isGluon
bool isGluon(const T &p)
Definition: AtlasPID.h:370
DeMoUpdate.reverse
reverse
Definition: DeMoUpdate.py:563
isHiggs
bool isHiggs(const T &p)
APID: HIGGS boson is only one particle.
Definition: AtlasPID.h:387
fractionalCharge
double fractionalCharge(const T &p)
Definition: AtlasPID.h:985
isQuark
bool isQuark(const T &p)
PDG rule 2: Quarks and leptons are numbered consecutively starting from 1 and 11 respectively; to do ...
Definition: AtlasPID.h:164
A
python.utils.AtlRunQueryDQUtils.p
p
Definition: AtlRunQueryDQUtils.py:209
isLightHadron
bool isLightHadron(const T &p)
Definition: AtlasPID.h:897
isCCbarMeson
bool isCCbarMeson(const T &p)
Definition: AtlasPID.h:911
strangeness
int strangeness(const T &p)
Definition: AtlasPID.h:776
isBottom
bool isBottom(const T &p)
Definition: AtlasPID.h:179
isTopHadron
bool isTopHadron(const T &p)
Definition: AtlasPID.h:902
lumiFormat.i
int i
Definition: lumiFormat.py:85
leadingQuark
int leadingQuark(const T &p)
Definition: AtlasPID.h:878
isPythia8Specific
bool isPythia8Specific(const T &p)
Definition: AtlasPID.h:408
beamspotman.n
n
Definition: beamspotman.py:729
isMonopole
bool isMonopole(const T &p)
PDG rule 11i Magnetic monopoles and dyons are assumed to have one unit of Dirac monopole charge and a...
Definition: AtlasPID.h:642
isWeaklyDecayingCHadron
bool isWeaklyDecayingCHadron(const T &p)
Definition: AtlasPID.h:966
hasTop
bool hasTop(const T &p)
Definition: AtlasPID.h:729
isZ
bool isZ(const T &p)
Definition: AtlasPID.h:376
DecodedPID::pid
const int & pid() const
Definition: AtlasPID.h:28
isPentaquark
bool isPentaquark(const T &p)
PDG rule 15 The 9-digit penta-quark codes are ±1nrnLnq1nq2nq3nq4nq5nJ, sorted such that nq1≥nq2≥nq3≥n...
Definition: AtlasPID.h:339
isCharmBaryon
bool isCharmBaryon(const T &p)
Definition: AtlasPID.h:923
isGlueball
bool isGlueball(const T &p)
APID: Definition of Glueballs: SM glueballs 99X (X=1,5), 999Y (Y=3,7)
Definition: AtlasPID.h:438
DecodedPID::max_digit
int max_digit(const int m, const int n) const
Definition: AtlasPID.h:29
isSuperpartner
bool isSuperpartner(const T &p)
Definition: AtlasPID.h:509
ParticleGun_EoverP_Config.pid
pid
Definition: ParticleGun_EoverP_Config.py:62
sign
int sign(int a)
Definition: TRT_StrawNeighbourSvc.h:107
isBBbarMeson
bool isBBbarMeson(const T &p)
Definition: AtlasPID.h:915
isNeutral
bool isNeutral(const T &p)
Definition: AtlasPID.h:1054
hasStrange
bool hasStrange(const T &p)
Definition: AtlasPID.h:726
isChLepton
bool isChLepton(const T &p)
APID: the fourth generation leptons are leptons.
Definition: AtlasPID.h:196
isSUSY
bool isSUSY(const T &p)
Definition: AtlasPID.h:620
isTau
bool isTau(const T &p)
Definition: AtlasPID.h:205
isGraviton
bool isGraviton(const T &p)
Definition: AtlasPID.h:394
DecodedPID::DecodedPID
DecodedPID(const int &p)
Definition: AtlasPID.h:18
DecodedPID::min_digit
int min_digit(const int m, const int n) const
Definition: AtlasPID.h:30
isStrangeHadron
bool isStrangeHadron(const T &p)
Definition: AtlasPID.h:899
isDM
bool isDM(const T &p)
PDG rule 11j: The nature of Dark Matter (DM) is not known, and therefore a definitive classificationi...
Definition: AtlasPID.h:656
DecodedPID
Implementation of classification functions according to PDG2022.
Definition: AtlasPID.h:16
isExcited
bool isExcited(const T &p)
PDG rule 11f Excited (composite) quarks and leptons are identified by setting n= 4 and nr= 0.
Definition: AtlasPID.h:533
isTopMeson
bool isTopMeson(const T &p)
Definition: AtlasPID.h:909
isHadron
bool isHadron(const T &p)
Definition: AtlasPID.h:348
isNeutrinoRH
bool isNeutrinoRH(const T &p)
PDG Rule 12: APID: Helper function for right-handed neutrino states These are generator defined PDG I...
Definition: AtlasPID.h:417
charge
double charge(const T &p)
Definition: AtlasPID.h:986
DecodedPID::ndigits
size_t ndigits() const
Definition: AtlasPID.h:31
isBaryon
bool isBaryon(const T &p)
Table 43.2 APID: states with fourth generation quarks are not baryons.
Definition: AtlasPID.h:279
isGeantino
bool isGeantino(const T &p)
Definition: AtlasPID.h:434
isSlepton
bool isSlepton(const T &p)
Definition: AtlasPID.h:479
isWeaklyDecayingBHadron
bool isWeaklyDecayingBHadron(const T &p)
Definition: AtlasPID.h:931
isFourthGeneration
bool isFourthGeneration(const T &p)
Is this a 4th generation fermion? APID: 4th generation fermions are not standard model particles.
Definition: AtlasPID.h:217
spin2
int spin2(const T &p)
Definition: AtlasPID.h:1076
isW
bool isW(const T &p)
Definition: AtlasPID.h:379
isRBaryon
bool isRBaryon(const T &p)
Definition: AtlasPID.h:586
isCharged
bool isCharged(const T &p)
Definition: AtlasPID.h:993
isTop
bool isTop(const T &p)
Definition: AtlasPID.h:182
isTransportable
bool isTransportable(const T &p)
Definition: AtlasPID.h:862
numberOfProtons
int numberOfProtons(const T &p)
Definition: AtlasPID.h:823
DeMoScan.first
bool first
Definition: DeMoScan.py:534
isLepton
bool isLepton(const T &p)
APID: the fourth generation leptons are leptons.
Definition: AtlasPID.h:186
isCharm
bool isCharm(const T &p)
Definition: AtlasPID.h:176
isPhoton
bool isPhoton(const T &p)
Definition: AtlasPID.h:373
charge3
int charge3(const T &p)
Definition: AtlasPID.h:984
hasSquark
bool hasSquark(const T &p, const int &q)
Definition: AtlasPID.h:607
isSquarkLH
bool isSquarkLH(const T &p)
Definition: AtlasPID.h:463
isSMNeutrino
bool isSMNeutrino(const T &p)
Definition: AtlasPID.h:212
isCharmHadron
bool isCharmHadron(const T &p)
Definition: AtlasPID.h:900
isRHadron
bool isRHadron(const T &p)
Definition: AtlasPID.h:600
extractSporadic.q
list q
Definition: extractSporadic.py:97
isQuarkonium
bool isQuarkonium(const T &p)
Is this a heavy-flavour quarkonium meson?
Definition: AtlasPID.h:270
isTechnicolor
bool isTechnicolor(const T &p)
PDG rule 11e Technicolor states have n = 3, with technifermions treated like ordinary fermions.
Definition: AtlasPID.h:521
isTrajectory
bool isTrajectory(const T &p)
PDG rule 8: The pomeron and odderon trajectories and a generic reggeon trajectory of states in QCD ar...
Definition: AtlasPID.h:356
isHeavyMeson
bool isHeavyMeson(const T &p)
Definition: AtlasPID.h:905
isStrangeBaryon
bool isStrangeBaryon(const T &p)
Definition: AtlasPID.h:922
DecodedPID::last
const int & last() const
Definition: AtlasPID.h:27
isElectron
bool isElectron(const T &p)
Definition: AtlasPID.h:199
DecodedPID::operator()
const int & operator()(const size_t n) const
Definition: AtlasPID.h:26
isDiquark
bool isDiquark(const T &p)
PDG rule 4 Diquarks have 4-digit numbers with nq1 >= nq2 and nq3 = 0 APID: states with top quarks are...
Definition: AtlasPID.h:224
baryonNumber3
int baryonNumber3(const T &p)
Definition: AtlasPID.h:738
pow
constexpr int pow(int base, int exp) noexcept
Definition: ap_fixedTest.cxx:15
isEMInteracting
bool isEMInteracting(const T &p)
Definition: AtlasPID.h:1069
containedQuarks
std::vector< int > containedQuarks(const T &p)
Definition: AtlasPID.h:1104
spin
double spin(const T &p)
Definition: AtlasPID.h:1098
TSU::T
unsigned long long T
Definition: L1TopoDataTypes.h:35
isMSSMHiggs
bool isMSSMHiggs(const T &p)
APID: Additional Higgs bosons for MSSM (Used in MCTruthClassifier)
Definition: AtlasPID.h:391
python.SystemOfUnits.m
float m
Definition: SystemOfUnits.py:106
python.SystemOfUnits.L
float L
Definition: SystemOfUnits.py:92
isMuon
bool isMuon(const T &p)
Definition: AtlasPID.h:202
isLeptoQuark
bool isLeptoQuark(const T &p)
PDG rule 11c: “One-of-a-kind” exotic particles are assigned numbers in the range 41–80.
Definition: AtlasPID.h:405