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 
48 //Note: The PDG rules assign the range 51-60 for generic DM, with explicit spins defined for 51-55.
49 // To keep the rest of this range consistent within ATLAS, 56, 57 and 58 are assigned spins 0, +1/2, +1 (respectively) as seen in arXiv:2504.10597v2
50 static const std::array<int,TABLESIZE> double_spin = {
51  +0, +1, +1, +1, +1, +1, +1, +1, +1, +0,
52  +0, +1, +1, +1, +1, +1, +1, +1, +1, +0,
53  +2, +2, +2, +2, +2, +0, +0, +0, +0, +0,
54  +0, +0, +2, +2, +2, +0, +0, +0, +0, +4,
55  +0, +0, +0, +0, +0, +0, +0, +0, +0, +0,
56  +0, +0, +1, +2, +0, +2, +0, +1, +2, +0,
57  +0, +0, +0, +0, +0, +0, +0, +0, +0, +0,
58  +0, +0, +0, +0, +0, +0, +0, +0, +0, +0,
59  +0, +0, +0, +0, +0, +0, +0, +0, +0, +0,
60  +0, +0, +0, +0, +0, +0, +0, +0, +0, +0
61 };
62 
63 
64 static const int DQUARK = 1;
65 static const int UQUARK = 2;
66 static const int SQUARK = 3;
67 static const int CQUARK = 4;
68 static const int BQUARK = 5;
69 static const int TQUARK = 6;
70 static const int BPRIME = 7; // 4th Generation quark
71 static const int TPRIME = 8; // 4th Generation quark
72 static const int QUARK_LIMIT = BPRIME; // Quark pdg_ids less than this are considered in (R-)Hadrons and Diquarks
73 
74 static const int ELECTRON = 11;
75 static const int POSITRON = -ELECTRON;
76 static const int NU_E = 12;
77 static const int MUON = 13;
78 static const int NU_MU = 14;
79 static const int TAU = 15;
80 static const int NU_TAU = 16;
81 static const int LPRIME = 17; // 4th Generation charged lepton
82 static const int NUPRIME = 18; // 4th Generation neutrino
83 
84 static const int GLUON = 21;
85 // APID: 9 rather than 21 is used to denote a gluon/gluino in composite states. (From PDG 11g)
86 static const int COMPOSITEGLUON = 9;
87 static const int PHOTON = 22;
88 static const int Z0BOSON = 23;
89 static const int WPLUSBOSON = 24;
90 static const int HIGGSBOSON = 25;
91 static const int ZPRIME = 32; // Z′/Z^0_2
92 static const int ZDBLPRIME = 33; // Z′′/Z^0_3
93 static const int WPLUSPRIME = 34; // W ′/W^+_2
94 static const int HIGGS2 = 35; // H^0/H^0_2 FIXME Any better ideas?
95 static const int HIGGS3 = 36; // A^0/H^0_3 FIXME Any better ideas?
96 static const int HIGGSPLUS = 37; // H^+
97 static const int HIGGSPLUSPLUS = 38; // H^++
98 static const int GRAVITON = 39;
99 static const int HIGGS4 = 40; // a^0/H^0_4 FIXME Any better ideas?
100 static const int LEPTOQUARK = 42;
101 
106 static const int DARKPHOTON = 60000;
107 static const int MAVTOP = 60001;
108 
109 static const int PIPLUS = 211;
110 static const int PIMINUS = -PIPLUS;
111 static const int PI0 = 111;
112 static const int K0L = 130;
113 
114 static const int K0S = 310;
115 static const int K0 = 311;
116 static const int KPLUS = 321;
117 static const int DPLUS = 411;
118 static const int DSTAR = 413;
119 static const int D0 = 421;
120 static const int DSPLUS = 431;
121 static const int JPSI = 443;
122 static const int B0 = 511;
123 static const int BCPLUS = 541;
124 static const int PROTON = 2212;
125 static const int NEUTRON = 2112;
126 static const int LAMBDA0 = 3122;
127 static const int LAMBDACPLUS = 4122;
128 static const int LAMBDAB0 = 5122;
129 static const int PSI2S = 20443;
130 
136 static const int RH_NU_E = 9900012;
137 static const int RH_NU_MU = 9900014;
138 static const int RH_NU_TAU = 9900016;
139 static const int WBOSON_LRSM = 9900024;
140 
141 static const int LEAD = 1000822080;
142 static const int OXYGEN = 1000080160;
143 static const int NEON = 1000100200;
144 static const int HELIUM = 1000020040;
145 
149 static const int POMERON = 990;
150 static const int ODDERON = 9990;
151 static const int REGGEON = 110;
152 
158 static const int GEANTINOPLUS = 998;
159 static const int GEANTINO0 = 999;
160 
161 
167 template<class T> inline bool isQuark(const T& p) {return isQuark(p->pdg_id());}
168 template<> inline bool isQuark(const int& p) { return p != 0 && (std::abs(p) <= TPRIME || std::abs(p) == MAVTOP);}
169 template<> inline bool isQuark(const DecodedPID& p){ return isQuark(p.pid()); }
170 
171 // APID: the fourth generation quarks are not standard model quarks
172 template<class T> inline bool isSMQuark(const T& p) {return isSMQuark(p->pdg_id());}
173 template<> inline bool isSMQuark(const int& p) { return p != 0 && std::abs(p) <= TQUARK;}
174 template<> inline bool isSMQuark(const DecodedPID& p){ return isSMQuark(p.pid()); }
175 
176 template<class T> inline bool isStrange(const T& p) {return isStrange(p->pdg_id());}
177 template<> inline bool isStrange(const int& p){ return std::abs(p) == SQUARK;}
178 
179 template<class T> inline bool isCharm(const T& p){return isCharm(p->pdg_id());}
180 template<> inline bool isCharm(const int& p){ return std::abs(p) == CQUARK;}
181 
182 template<class T> inline bool isBottom(const T& p){return isBottom(p->pdg_id());}
183 template<> inline bool isBottom(const int& p){ return std::abs(p) == BQUARK;}
184 
185 template<class T> inline bool isTop(const T& p){return isTop(p->pdg_id());}
186 template<> inline bool isTop(const int& p){ return std::abs(p) == TQUARK;}
187 
189 template<class T> inline bool isLepton(const T& p){return isLepton(p->pdg_id());}
190 template<> inline bool isLepton(const int& p){ auto sp = std::abs(p); return sp >= ELECTRON && sp <= NUPRIME; }
191 template<> inline bool isLepton(const DecodedPID& p){ return isLepton(p.pid()); }
192 
194 template<class T> inline bool isSMLepton(const T& p){return isSMLepton(p->pdg_id());}
195 template<> inline bool isSMLepton(const int& p){ auto sp = std::abs(p); return sp >= ELECTRON && sp <= NU_TAU; }
196 template<> inline bool isSMLepton(const DecodedPID& p){ return isSMLepton(p.pid()); }
197 
199 template<class T> inline bool isChLepton(const T& p){return isChLepton(p->pdg_id());}
200 template<> inline bool isChLepton(const int& p){ auto sp = std::abs(p); return sp >= ELECTRON && sp <= LPRIME && sp%2 == 1; }
201 
202 template<class T> inline bool isElectron(const T& p){return isElectron(p->pdg_id());}
203 template<> inline bool isElectron(const int& p){ return std::abs(p) == ELECTRON;}
204 
205 template<class T> inline bool isMuon(const T& p){return isMuon(p->pdg_id());}
206 template<> inline bool isMuon(const int& p){ return std::abs(p) == MUON;}
207 
208 template<class T> inline bool isTau(const T& p){return isTau(p->pdg_id());}
209 template<> inline bool isTau(const int& p){ return std::abs(p) == TAU;}
210 
212 template<class T> inline bool isNeutrino(const T& p){return isNeutrino(p->pdg_id());}
213 template<> inline bool isNeutrino(const int& p){ auto sp = std::abs(p); return sp == NU_E || sp == NU_MU || sp == NU_TAU || sp == NUPRIME; }
214 
215 template<class T> inline bool isSMNeutrino(const T& p){return isSMNeutrino(p->pdg_id());}
216 template<> inline bool isSMNeutrino(const int& p){ auto sp = std::abs(p); return sp == NU_E || sp == NU_MU || sp == NU_TAU; }
217 
220 template<class T> inline bool isFourthGeneration(const T& p){return isFourthGeneration(p->pdg_id());}
221 template<> inline bool isFourthGeneration(const int& p) {return std::abs(p) == BPRIME || std::abs(p) == TPRIME || std::abs(p) == LPRIME || std::abs(p) == NUPRIME;}
222 
227 template<class T> inline bool isDiquark(const T& p){return isDiquark(p->pdg_id());}
228 template<> inline bool isDiquark(const DecodedPID& p){
229  if ( p.ndigits() == 4 && p(0) >= p(1) && p(1) !=0 && p(2) == 0 && (p.last() == 1 || p.last() == 3)
230  && p.max_digit(2,4) < QUARK_LIMIT
231  ) return true;
232  return false;
233 }
234 template<> inline bool isDiquark(const int& p){ auto value_digits = DecodedPID(p); return isDiquark(value_digits);}
235 
244 template<class T> inline bool isMeson(const T& p){return isMeson(p->pdg_id());}
245 template<> inline bool isMeson(const DecodedPID& p){
246  if (p.ndigits() < 3 ) return false;
247  if (p.ndigits() == 7 && (p(0) == 1 || p(0) == 2)) return false; // APID don't match SUSY particles
248  if (std::abs(p.pid()) == K0S) return true;
249  if (std::abs(p.pid()) == K0L) return true;
250  if (std::abs(p.pid()) == K0) return true;
251  if (p.last() % 2 != 1 ) return false;
252  if (p.max_digit(1,3) >= QUARK_LIMIT ) return false; // Ignore pdg_ids which would describe states including fourth generation quarks
253  if (p.min_digit(1,3) == 0 ) return false;
254  if (*(p.second.rbegin() + 2) < *(p.second.rbegin() + 1) ) return false; // Quark ordering (nq2 >= nq3)
255  if (*(p.second.rbegin() + 2) == *(p.second.rbegin() + 1) && p.pid() < 0 ) return false; // Illegal antiparticle check (nq2 == nq3)
256  if (p.ndigits() == 3 ) return true;
257  if (*(p.second.rbegin() + 3) != 0 ) return false; // Only two quarks! (nq1 == 0)
258  if (p.ndigits() == 5 && p(0) == 1 ) return true;
259  if (p.ndigits() == 5 && p(0) == 2 && p.last() > 1 ) return true;
260  if (p.ndigits() == 5 && p(0) == 3 && p.last() > 1 ) return true;
261  if (p.ndigits() == 6 && p.last() % 2 == 1 ) return true;
262  if (p.ndigits() == 7 && p(0) == 9 && p(1) == 0 ) return true;
263 
264  return false;
265 }
266 template<> inline bool isMeson(const int& p){ auto value_digits = DecodedPID(p); return isMeson(value_digits);}
267 
273 template<class T> inline bool isQuarkonium(const T& p){return isQuarkonium(p->pdg_id());}
274 template<> inline bool isQuarkonium(const DecodedPID& p) {
275  if (!isMeson(p)) return false; //< all quarkonia are mesons
276  return (*(p.second.rbegin() + 2) > SQUARK && p.last() > 0 && *(p.second.rbegin() + 1) == *(p.second.rbegin() + 2));
277 }
278 template<> inline bool isQuarkonium(const int& p){ auto value_digits = DecodedPID(p); return isQuarkonium(value_digits);}
279 
282 template<class T> inline bool isBaryon(const T& p){return isBaryon(p->pdg_id());}
283 template<> inline bool isBaryon(const DecodedPID& p){
284  if (p.ndigits() < 4 ) return false;
285  if (p.max_digit(1,4) >= QUARK_LIMIT ) return false; // Ignore pdg_ids which would describe states including fourth generation quarks
286  if (p.min_digit(1,4) == 0) return false; // Ignore pdg_ids with zero for nq1, nq2, nq3
287  if (p.ndigits() == 4 && (p.last() == 2 || p.last() == 4|| p.last() == 6|| p.last() == 8) ) return true;
288 
289  if (p.ndigits() == 5 && p(0) == 1 && (p.last() == 2 || p.last() == 4) ) return true;
290  if (p.ndigits() == 5 && p(0) == 3 && (p.last() == 2 || p.last() == 4) ) return true;
291 
292  if (p.ndigits() == 6 ) {
293  if (p(0) == 1 && p(1) == 0 && p.last() == 2 ) return true;
294  if (p(0) == 1 && p(1) == 0 && p.last() == 4 ) return true;
295  if (p(0) == 1 && p(1) == 0 && p.last() == 6 ) return true;
296  if (p(0) == 1 && p(1) == 1 && p.last() == 2 ) return true;
297  if (p(0) == 1 && p(1) == 2 && p.last() == 4 ) return true;
298 
299  if (p(0) == 2 && p(1) == 0 && p.last() == 2 ) return true;
300  if (p(0) == 2 && p(1) == 0 && p.last() == 4 ) return true;
301  if (p(0) == 2 && p(1) == 0 && p.last() == 6 ) return true;
302  if (p(0) == 2 && p(1) == 0 && p.last() == 8 ) return true;
303  if (p(0) == 2 && p(1) == 1 && p.last() == 2 ) return true;
304  }
305 
306  if (p.ndigits() == 5 ) {
307  if (p(0) == 2 && p.last() == 2 ) return true;
308  if (p(0) == 2 && p.last() == 4 ) return true;
309  if (p(0) == 2 && p.last() == 6 ) return true;
310  if (p(0) == 5 && p.last() == 2 ) return true;
311  if (p(0) == 1 && p.last() == 6 ) return true;
312  if (p(0) == 4 && p.last() == 2 ) return true;
313  }
314  return false;
315 }
316 template<> inline bool isBaryon(const int& p){ auto value_digits = DecodedPID(p); return isBaryon(value_digits);}
317 
325 template<class T> inline bool isTetraquark(const T& p){return isTetraquark(p->pdg_id());}
326 template<> inline bool isTetraquark(const DecodedPID& p){
327  return (p.ndigits() == 9 && p(0) == 1 && p(5) == 0 &&
328  p.max_digit(1,3) < QUARK_LIMIT && p.min_digit(1,3) > 0 && // ignore 4th generation quarks for nq3 and nq4
329  p.max_digit(4,6) < QUARK_LIMIT && p.min_digit(4,6) > 0 && // ignore 4th generation quarks for nq1 and nq2
330  ( p(3) >= p(4) && p(6) >= p(7) ) && ( ( p(3) > p(6) ) || ( p(3) == p(6) && (p(4) >= p(7))))
331  );
332 }
333 template<> inline bool isTetraquark(const int& p){ auto value_digits = DecodedPID(p); return isTetraquark(value_digits);}
334 
341 // APID: states with fourth generation quarks are not pentaquarks
342 template<class T> inline bool isPentaquark(const T& p){return isPentaquark(p->pdg_id());}
343 template<> inline bool isPentaquark(const DecodedPID& p){
344  return (p.ndigits() == 9 && p(0) == 1 &&
345  p.max_digit(1,6) < QUARK_LIMIT && p.min_digit(1,6) > 0 && // ignore 4th generation (anti-)quarks
346  ( p(3) >= p(4) && p(4) >= p(5) && p(5) >= p(6)) );
347 }
348 template<> inline bool isPentaquark(const int& p){ auto value_digits = DecodedPID(p); return isPentaquark(value_digits);}
349 
350 // APID Mesons, Baryons, Tetraquarks and Pentaquarks are Hadrons
351 template<class T> inline bool isHadron(const T& p){return isHadron(p->pdg_id());}
352 template<> inline bool isHadron(const DecodedPID& p){ return isMeson(p) || isBaryon(p) || isTetraquark(p) || isPentaquark(p); }
353 template<> inline bool isHadron(const int& p){ auto value_digits = DecodedPID(p); return isHadron(value_digits);}
354 
355 
359 template<class T> inline bool isTrajectory(const T& p){return isTrajectory(p->pdg_id());}
360 template<> inline bool isTrajectory(const int& p){ return std::abs(p) == POMERON || std::abs(p) == ODDERON || std::abs(p) == REGGEON; }
361 
362 
369 template<class T> inline bool isBoson(const T& p){return isBoson(p->pdg_id());}
370 template<> inline bool isBoson(const int& p){ auto sp = std::abs(p); return sp > 20 && sp < 41; }
371 template<> inline bool isBoson(const DecodedPID& p){ return isBoson(p.pid()); }
372 
373 template<class T> inline bool isGluon(const T& p){return isGluon(p->pdg_id());}
374 template<> inline bool isGluon(const int& p){ return p == GLUON; }
375 
376 template<class T> inline bool isPhoton(const T& p){return isPhoton(p->pdg_id());}
377 template<> inline bool isPhoton(const int& p){ return p == PHOTON; }
378 
379 template<class T> inline bool isZ(const T& p){return isZ(p->pdg_id());}
380 template<> inline bool isZ(const int& p){ return p == Z0BOSON; }
381 
382 template<class T> inline bool isW(const T& p){return isW(p->pdg_id());}
383 template<> inline bool isW(const int& p){ return std::abs(p) == WPLUSBOSON; }
384 
386 template<class T> inline bool isHeavyBoson(const T& p){return isHeavyBoson(p->pdg_id());}
387 template<> inline bool isHeavyBoson(const int& p){ return p == ZPRIME || p == ZDBLPRIME || std::abs(p) == WPLUSPRIME; }
388 
390 template<class T> inline bool isHiggs(const T& p){return isHiggs(p->pdg_id());}
391 template<> inline bool isHiggs(const int& p){ return p == HIGGSBOSON; }
392 
394 template<class T> inline bool isMSSMHiggs(const T& p){return isMSSMHiggs(p->pdg_id());}
395 template<> inline bool isMSSMHiggs(const int& p){ return p == HIGGS2 || p == HIGGS3 || std::abs(p) == HIGGSPLUS; }
396 
397 template<class T> inline bool isGraviton(const T& p) {return isGraviton(p->pdg_id());}
398 template<> inline bool isGraviton(const int& p){ return p == GRAVITON; }
399 
400 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
401 
408 template<class T> inline bool isLeptoQuark(const T& p){return isLeptoQuark(p->pdg_id());}
409 template<> inline bool isLeptoQuark(const int& p){ return std::abs(p) == LEPTOQUARK; }
410 
411 template<class T> inline bool isPythia8Specific(const T& p){return isPythia8Specific(p->pdg_id());}
412 template<> inline bool isPythia8Specific(const DecodedPID& p){ return (p.ndigits() == 7 && p(0) == 9 && p(1) == 9);}
413 template<> inline bool isPythia8Specific(const int& p){ auto value_digits = DecodedPID(p); return isPythia8Specific(value_digits);}
414 
420 template<class T> inline bool isNeutrinoRH(const T& p){return isNeutrinoRH(p->pdg_id());}
421 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);}
422 
425 template<class T> inline bool isGenSpecific(const T& p){return isGenSpecific(p->pdg_id());}
426 template<> inline bool isGenSpecific(const int& p){
427  int ap = std::abs(p);
428  if (ap >= 81 && ap <= 100) return true;
429  if (ap >= 901 && ap <= 930) return true;
430  if (ap >= 998 && ap <= 999) return true;
431  if (ap >= 1901 && ap <= 1930) return true;
432  if (ap >= 2901 && ap <= 2930) return true;
433  if (ap >= 3901 && ap <= 3930) return true;
434  return false;
435 }
436 
437 template<class T> inline bool isGeantino(const T& p){return isGeantino(p->pdg_id());}
438 template<> inline bool isGeantino(const int& p){ return (std::abs(p) == GEANTINO0 || std::abs(p) == GEANTINOPLUS);}
439 
441 template<class T> inline bool isGlueball(const T& p) { return isGlueball(p->pdg_id()); }
442 template<> inline bool isGlueball(const DecodedPID& p) {
443  if (p.ndigits() > 4) return false; // APID avoid classifying R-Glueballs as SM Glueballs
444  return
445  ( ( p.ndigits() == 3 && p(0) == COMPOSITEGLUON && p(1) == COMPOSITEGLUON && (p.last() == 1 || p.last() == 5) ) ||
446  ( p.ndigits() == 4 && p(0) == COMPOSITEGLUON && p(1) == COMPOSITEGLUON && p(2) == COMPOSITEGLUON && (p.last() == 3 || p.last() == 7) ) );
447 }
448 template<> inline bool isGlueball(const int& p) { auto value_digits = DecodedPID(p); return isGlueball(value_digits); }
449 
450 
456 
457 // APID: Super-partners of standard model quarks only
458 template<class T> inline bool isSquark(const T& p) { return isSquark(p->pdg_id()); }
459 template<> inline bool isSquark(const DecodedPID& p){
460  auto pp = p.shift(1); return (p.ndigits() == 7 && (p(0) == 1 || p(0) == 2) && isSMQuark(pp));
461 }
462 template<> inline bool isSquark(const int& p){ auto value_digits = DecodedPID(p); return isSquark(value_digits);}
463 
464 
465 // APID: Super-partners of left-handed standard model quarks only
466 template<class T> inline bool isSquarkLH(const T& p) { return isSquarkLH(p->pdg_id()); }
467 template<> inline bool isSquarkLH(const DecodedPID& p){
468  auto pp = p.shift(1); return (p.ndigits() == 7 && p(0) == 1 && isSMQuark(pp));
469 }
470 template<> inline bool isSquarkLH(const int& p){ auto value_digits = DecodedPID(p); return isSquarkLH(value_digits);}
471 
472 
473 // APID: Super-partners of right-handed standard model quarks only
474 template<class T> inline bool isSquarkRH(const T& p) { return isSquarkRH(p->pdg_id()); }
475 template<> inline bool isSquarkRH(const DecodedPID& p){
476  auto pp = p.shift(1); return (p.ndigits() == 7 && p(0) == 2 && isSMQuark(pp));
477 }
478 template<> inline bool isSquarkRH(const int& p){ auto value_digits = DecodedPID(p); return isSquarkRH(value_digits);}
479 
480 
481 // APID: Super-partners of standard model leptons only
482 template<class T> inline bool isSlepton(const T& p) { return isSlepton(p->pdg_id()); }
483 template<> inline bool isSlepton(const DecodedPID& p){ auto pp = p.shift(1); return (p.ndigits() == 7 && (p(0) == 1 || p(0) == 2) && isSMLepton(pp));}
484 template<> inline bool isSlepton(const int& p){ auto value_digits = DecodedPID(p); return isSlepton(value_digits);}
485 
486 
487 // APID: Super-partners of left-handed standard model leptons only
488 template<class T> inline bool isSleptonLH(const T& p) { return isSleptonLH(p->pdg_id()); }
489 template<> inline bool isSleptonLH(const DecodedPID& p){
490  auto pp = p.shift(1); return (p.ndigits() == 7 && p(0) == 1 && isSMLepton(pp));
491 }
492 template<> inline bool isSleptonLH(const int& p){ auto value_digits = DecodedPID(p); return isSleptonLH(value_digits);}
493 
494 
495 // APID: Super-partners of right-handed standard model leptons only
496 template<class T> inline bool isSleptonRH(const T& p) { return isSleptonRH(p->pdg_id()); }
497 template<> inline bool isSleptonRH(const DecodedPID& p){
498  auto pp = p.shift(1); return (p.ndigits() == 7 && p(0) == 2 && isSMLepton(pp));
499 }
500 template<> inline bool isSleptonRH(const int& p){ auto value_digits = DecodedPID(p); return isSleptonRH(value_digits);}
501 
502 
503 // APID: Super-partners of gauge bosons including gravitons
504 template<class T> inline bool isGaugino(const T& p) { return isGaugino(p->pdg_id()); }
505 template<> inline bool isGaugino(const DecodedPID& p){
506  auto pp = p.shift(1); return (p.ndigits() == 7 && p(0) == 1 && isBoson(pp.pid()));
507 }
508 template<> inline bool isGaugino(const int& p){ auto value_digits = DecodedPID(p); return isGaugino(value_digits);}
509 
510 
511 // APID: Super-partners of fundamental particles
512 template<class T> inline bool isSuperpartner(const T& p) { return isSuperpartner(p->pdg_id()); }
513 template<> inline bool isSuperpartner(const DecodedPID& p){
514  return isSlepton(p) || isSquark(p) || isGaugino(p);
515 }
516 template<> inline bool isSuperpartner(const int& p){ auto value_digits = DecodedPID(p); return isSuperpartner(value_digits);}
517 
518 
524 template<class T> inline bool isTechnicolor(const T& p){return isTechnicolor(p->pdg_id());}
525 template <>
526 inline bool isTechnicolor(const DecodedPID& p) {
527  const auto& pp = (p.ndigits() == 7) ? p.shift(2) : DecodedPID(0);
528  return (p.ndigits() == 7 && p(0) == 3 && (p(1) == 0 || p(1) == 1) &&
529  (isQuark(pp) || isLepton(pp) || isBoson(pp) || isGlueball(pp) ||
530  isDiquark(pp) || isHadron(pp)));
531 }
532 template<> inline bool isTechnicolor(const int& p){ auto value_digits = DecodedPID(p); return isTechnicolor(value_digits);}
533 
536 template<class T> inline bool isExcited(const T& p){return isExcited(p->pdg_id());}
537 template <>
538 inline bool isExcited(const DecodedPID& p) {
539  const auto& pp = (p.ndigits() == 7) ? p.shift(2) : DecodedPID(0);
540  return (p.ndigits() == 7 && (p(0) == 4 && p(1) == 0) &&
541  (isLepton(pp) || isQuark(pp)));
542 }
543 template<> inline bool isExcited(const int& p){ auto value_digits = DecodedPID(p); return isExcited(value_digits);}
544 
557 
562 template<class T> inline bool isRGlueball(const T& p) { return isRGlueball(p->pdg_id()); }
563 template<> inline bool isRGlueball(const DecodedPID& p) {
564  if (p.ndigits() != 7 || p(0) != 1) return false;
565  auto pp = p.shift(1);
566  return
567  ( ( pp.ndigits() == 3 && pp(0) == COMPOSITEGLUON && pp(1) == COMPOSITEGLUON && (pp.last() == 1 || pp.last() == 3) ) ||
568  ( pp.ndigits() == 4 && pp(0) == COMPOSITEGLUON && pp(1) == COMPOSITEGLUON && pp(2) == COMPOSITEGLUON && (pp.last() == 1 || pp.last() == 5) ) );
569 }
570 template<> inline bool isRGlueball(const int& p) { auto value_digits = DecodedPID(p); return isRGlueball(value_digits); }
571 
572 // APID Define R-Mesons as gluino-quark-antiquark and squark-antiquark bound states (ignore 4th generation squarks/quarks)
573 // NB Current models only allow gluino-quark-antiquark, stop-antiquark and sbottom-antiquark states
574 template<class T> inline bool isRMeson(const T& p) { return isRMeson(p->pdg_id()); }
575 template<> inline bool isRMeson(const DecodedPID& p) {
576  if (!(p.ndigits() == 7 && (p(0) == 1 || p(0) == 2))) return false;
577  auto pp = p.shift(1);
578  return (
579  // Handle ~gluino-quark-antiquark states
580  (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)) ||
581  // Handle squark-antiquark states (previously called Smeson/mesoninos)
582  (pp.ndigits() == 3 && pp.min_digit(1,3) > 0 && pp.max_digit(1,3) < QUARK_LIMIT && pp(1) <= pp(0) && pp.last() == 2)
583  );
584 }
585 template<> inline bool isRMeson(const int& p) { auto value_digits = DecodedPID(p); return isRMeson(value_digits); }
586 
587 // APID Define R-Baryons as gluino-quark-quark-quark and squark-quark-quark bound states (ignore 4th generation squarks/quarks)
588 // NB Current models only allow gluino-quark-quark-quark, stop-quark-quark and sbottom-quark-quark states
589 template<class T> inline bool isRBaryon(const T& p) { return isRBaryon(p->pdg_id()); }
590 template<> inline bool isRBaryon(const DecodedPID& p) {
591  if (!(p.ndigits() == 7 && (p(0) == 1 || p(0) == 2))) return false;
592  auto pp = p.shift(1);
593  return (
594  // Handle ~gluino-quark-quark-quark states
595  (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)) ||
596  // Handle squark-quark-quark states (previously called Sbaryons)
597  (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))
598  );
599 }
600 template<> inline bool isRBaryon(const int& p) { auto value_digits = DecodedPID(p); return isRBaryon(value_digits); }
601 
602 
603 template<class T> inline bool isRHadron(const T& p) { return isRHadron(p->pdg_id()); }
604 template<> inline bool isRHadron(const DecodedPID& p) {
605  return (isRBaryon(p) || isRMeson(p) || isRGlueball(p));
606 }
607 template<> inline bool isRHadron(const int& p) { auto value_digits = DecodedPID(p); return isRHadron(value_digits); }
608 
609 
610 template<class T> inline bool hasSquark(const T& p, const int& q) { return hasSquark(p->pdg_id(), q); }
611 template<> inline bool hasSquark(const DecodedPID& p, const int& q){
612  auto pp = p.shift(1); return (
613  (isSquark(p) || isRHadron(p))
614  && pp.ndigits() != 2 // skip lepton and boson super-partners by vetoing ndigits==2
615  && pp(0) == q // After shifting, the first digit will always represent the squark in R-Hadron (and squark) PIDs
616  );
617 }
618 template<> inline bool hasSquark(const int& p, const int& q){ auto value_digits = DecodedPID(p); return hasSquark(value_digits, q);}
619 
620 
621 // APID Define isSUSY to return true for Superpartners of standard model particles and R-Hadrons
622 // TODO Should this include the MSSM extended Higgs sector??
623 template<class T> inline bool isSUSY(const T& p){return isSUSY(p->pdg_id());}
624 template<> inline bool isSUSY(const DecodedPID& p){return (isSuperpartner(p) || isRHadron(p));}
625 template<> inline bool isSUSY(const int& p){ auto value_digits = DecodedPID(p); return isSUSY(value_digits);}
626 
627 
635 template<class T> inline bool isKK(const T& p){return isKK(p->pdg_id());}
636 template<> inline bool isKK(const DecodedPID& p){return (p.ndigits() == 7 && (p(0) == 5 || p(0) == 6 ) && (p(1) != 9) );}
637 template<> inline bool isKK(const int& p){ auto value_digits = DecodedPID(p); return isKK(value_digits);}
638 
645 template<class T> inline bool isMonopole(const T& p){return isMonopole(p->pdg_id());}
646 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);}
647 template<> inline bool isMonopole(const int& p){ auto value_digits = DecodedPID(p); return isMonopole(value_digits);}
648 
658 template<class T> inline bool isDM(const T& p){return isDM(p->pdg_id());}
659 template<> inline bool isDM(const int& p){
660  auto sp = std::abs(p);
661  auto value_digits = DecodedPID(p);
662  return (sp >= 51 && sp <= 60) || (value_digits.ndigits() == 7 && value_digits(0) == 5 && value_digits(1) == 9) || sp == DARKPHOTON; }
663 
668 template<class T> inline bool isHiddenValley(const T& p){return isHiddenValley(p->pdg_id());}
669 template <>
670 inline bool isHiddenValley(const DecodedPID& p) {
671  const auto& pp = (p.ndigits() == 7) ? p.shift(2) : DecodedPID(0);
672  return (p.ndigits() == 7 && p(0) == 4 && p(1) == 9 &&
673  (isQuark(pp) || isLepton(pp) || isBoson(pp) || isGlueball(pp) ||
674  isDiquark(pp) || isHadron(pp)));
675 }
676 template<> inline bool isHiddenValley(const int& p){ auto value_digits = DecodedPID(p); return isHiddenValley(value_digits);}
677 
684 template<class T> inline bool isGenericMultichargedParticle(const T& p){return isGenericMultichargedParticle(p->pdg_id());}
685 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);}
686 template<> inline bool isGenericMultichargedParticle(const int& p){ auto value_digits = DecodedPID(p); return isGenericMultichargedParticle(value_digits);}
687 
702 template<class T> inline bool isNucleus(const T& p){return isNucleus(p->pdg_id());}
703 template<> inline bool isNucleus(const DecodedPID& p){
704  if (std::abs(p.pid()) == PROTON) return true;
705  if (p.ndigits() != 10) return false;
706  // charge should always be less than or equal to baryon number
707  // the following line is A >= Z
708  const int A = p(8) + 10*p(7) + 100*p(6);
709  const int Z = p(5) + 10*p(4) + 100*p(3);
710  return ( A >= Z && p(0) == 1 && p(1) == 0 );
711 }
712 template<> inline bool isNucleus(const int& p){ auto value_digits = DecodedPID(p); return isNucleus(value_digits);}
713 
714 
715 template<class T> inline bool hasQuark(const T& p, const int& q);
716 template<> inline bool hasQuark(const DecodedPID& p, const int& q){
717  if (isQuark(p.pid())) { return (std::abs(p.pid()) == q );}
718  if (isMeson(p)) { return *(p.second.rbegin() + 1) == q ||*(p.second.rbegin()+2) ==q;}
719  if (isDiquark(p)) { auto i = std::find(p.second.rbegin() + 2,p.second.rbegin()+4,q); return (i!=p.second.rbegin()+4);}
720  if (isBaryon(p)) { auto i = std::find(p.second.rbegin() + 1,p.second.rbegin()+4,q); return (i!=p.second.rbegin()+4);}
721  if (isTetraquark(p)) { auto i = std::find(p.second.rbegin() + 1,p.second.rbegin()+5,q); return (i!=p.second.rbegin()+5);}
722  if (isPentaquark(p)) { auto i = std::find(p.second.rbegin() + 1,p.second.rbegin()+6,q); return (i!=p.second.rbegin()+6);}
723  if (isNucleus(p) && std::abs(p.pid()) != PROTON) { return (q == 1 || q == 2 || (q==3 && p(2) > 0));}
724  if (isSUSY(p)) { // APID SUSY case
725  auto pp = p.shift(1);
726  if ( pp.ndigits() == 1 ) { return false; } // Handle squarks
727  if ( pp.ndigits() == 3 ) { return (pp(1) == q); } // Handle ~q qbar pairs
728  if ( pp.ndigits() == 4 ) { return (pp(1) == q || pp(2) == q); } // Ignore gluinos and squarks
729  if ( pp.ndigits() == 5 ) { return (pp(1) == q || pp(2) == q || pp(3) == q); } // Ignore gluinos and squarks
730  if ( pp.ndigits() > 5 ) { pp = pp.shift(1); } // Drop gluinos and squarks
731  return hasQuark(pp, q); }
732  return false;
733 }
734 template<> inline bool hasQuark(const int& p, const int& q){ auto value_digits = DecodedPID(p); return hasQuark(value_digits, q);}
735 
736 template<class T> inline bool hasStrange(const T& p) { return hasQuark(p,SQUARK); }
737 template<class T> inline bool hasCharm(const T& p) { return hasQuark(p,CQUARK); }
738 template<class T> inline bool hasBottom(const T& p) { return hasQuark(p,BQUARK); }
739 template<class T> inline bool hasTop(const T& p) { return hasQuark(p,TQUARK); }
740 
741 
742 // APID: The baryon number is defined as:
743 // B = (1/3)*( n_q - n_{qbar} )
744 // where n_q⁠ is the number of quarks, and ⁠n_{qbar} is the number of
745 // antiquarks. By convention, squarks have the same quantum numbers as
746 // the corresponding quarks (modulo spin and R), so have baryon number
747 // 1/3.
748 template<class T> inline int baryonNumber3(const T& p) {return baryonNumber3(p->pdg_id());}
749 template<> inline int baryonNumber3(const DecodedPID& p){
750  if (isQuark(p.pid())) { return (p.pid() > 0) ? 1 : - 1;}
751  if (isDiquark(p)) { return (p.pid() > 0) ? 2 : -2; }
752  if (isMeson(p) || isTetraquark(p)) { return 0; }
753  if (isBaryon(p) || isPentaquark(p)){ return (p.pid() > 0) ? 3 : -3; }
754  if (isNucleus(p)) {
755  const int result = 3*p(8) + 30*p(7) + 300*p(6);
756  return (p.pid() > 0) ? result : -result;
757  }
758  if (isSUSY(p)) {
759  auto pp = p.shift(1);
760  if (pp.ndigits() < 3 ) { return baryonNumber3(pp); } // super-partners of fundamental particles
761  if (pp(0) == COMPOSITEGLUON) {
762  if (pp(1) == COMPOSITEGLUON) { return 0; } // R-Glueballs
763  if ( pp.ndigits() == 4 ) { return 0; } // states with gluino-quark-antiquark
764  if ( pp.ndigits() == 5) { return (p.pid() > 0) ? 3 : -3; } // states with gluino-quark-quark-quark
765  }
766  if (pp.ndigits() == 3) { return 0; } // squark-antiquark
767  if (pp.ndigits() == 4) { return (p.pid() > 0) ? 3 : -3; } // states with squark-quark-quark
768  }
769  return 0;
770 }
771 template<> inline int baryonNumber3(const int& p){ auto value_digits = DecodedPID(p); return baryonNumber3(value_digits);}
772 
773 template<class T> inline double baryonNumber(const T& p) {return baryonNumber(p->pdg_id());}
774 template<> inline double baryonNumber(const DecodedPID& p){ return static_cast<double>(baryonNumber3(p))/3.0;}
775 template<> inline double baryonNumber(const int& p){ auto value_digits = DecodedPID(p); return static_cast<double>(baryonNumber3(value_digits))/3.0;}
776 
777 
778 // APID: The strangeness of a particle is defined as:
779 // S = − ( n_s − n_{sbar} )
780 // where n_s represents the number of strange quarks and n_{sbar}
781 // represents the number of strange antiquarks. By convention, strange
782 // squarks have the same quantum numbers as strange quarks (modulo
783 // spin and R), so have strangeness -1.
784 static const std::array<int,10> is_strange = {
785  +0, +0, +0, -1, +0, +0, +0, +0, +0, +0 };
786 template<class T> inline int strangeness(const T& p) {return strangeness(p->pdg_id());}
787 template<> inline int strangeness(const DecodedPID& p){
788  if (isNucleus(p) && p.ndigits() == 10) { return (p.pid() > 0) ? -p(2) : p(2); }
789  if (isStrange(p.pid())) { return (p.pid() > 0) ? -1 : 1; }
790  if (!hasStrange(p) && !hasSquark(p,SQUARK)) { return 0; }
791  if (std::abs(p.pid()) == K0) { return (p.pid() > 0) ? 1 : -1; }
792  size_t nq = 0;
793  int sign = 1;
794  int signmult = 1;
795  int result=0;
796  bool classified = false;
797  if (!classified && isMeson(p)) { classified = true; nq = 2; if ((*(p.second.rbegin()+2)) == 2||(*(p.second.rbegin()+2)) == 4 ) { sign=-1;} signmult =-1; }
798  if (!classified && isDiquark(p)) {return is_strange.at(p(0))+is_strange.at(p(1)); }
799  if (!classified && isBaryon(p)) { classified = true; nq = 3; }
800  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)); }
801  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)); }
802  if (!classified && isSUSY(p)) {
803  nq = 0;
804  auto pp = p.shift(1);
805  if (pp.ndigits() < 3 ) { return strangeness(pp); } // super-partners of fundamental particles
806  if (pp(0) == COMPOSITEGLUON) {
807  if (pp(1) == COMPOSITEGLUON) { return 0; } // R-Glueballs
808  if ( pp.ndigits() == 4 || pp.ndigits() == 5) {
809  pp = pp.shift(1); // Remove gluino
810  }
811  }
812  if (pp.ndigits() == 3) { classified = true; nq = 2; if (p.last()%2==0) {sign = -1;} signmult = -1; } // states with quark-antiquark or squark-antiquark
813  if (pp.ndigits() == 4) { classified = true; nq = 3; } // states with quark-quark-quark or squark-quark-quark
814  }
815  for (auto r = p.second.rbegin() + 1; r != p.second.rbegin() + 1 + nq; ++r) {
816  result += is_strange.at(*r)*sign;
817  sign*=signmult;
818  }
819  return p.pid() > 0 ? result : -result;
820 }
821 template<> inline int strangeness(const int& p){ auto value_digits = DecodedPID(p); return strangeness(value_digits);}
822 
823 
824 template<class T> inline int numberOfLambdas(const T& p) {return numberOfLambdas(p->pdg_id());}
825 template<> inline int numberOfLambdas(const DecodedPID& p){
826  if (std::abs(p.pid()) == LAMBDA0) { return (p.pid() > 0) ? 1 : -1; }
827  if (isNucleus(p) && p.ndigits() == 10) { return (p.pid() > 0) ? p(2) : -p(2); }
828  return 0;
829 }
830 template<> inline int numberOfLambdas(const int& p){ auto value_digits = DecodedPID(p); return numberOfLambdas(value_digits);}
831 
832 
833 template<class T> inline int numberOfProtons(const T& p) {return numberOfProtons(p->pdg_id());}
834 template<> inline int numberOfProtons(const DecodedPID& p){
835  if (std::abs(p.pid()) == PROTON) { return (p.pid() > 0) ? 1 : -1; }
836  if (isNucleus(p)) {
837  const int result = p(5) + 10*p(4) + 100*p(3);
838  return (p.pid() > 0) ? result : -result;
839  }
840  return 0;
841 }
842 template<> inline int numberOfProtons(const int& p){ auto value_digits = DecodedPID(p); return numberOfProtons(value_digits);}
843 
844 
846 template<class T> inline bool isBSM(const T& p){return isBSM(p->pdg_id());}
847 template<> inline bool isBSM(const DecodedPID& p){
848  if (p.pid() == GRAVITON || std::abs(p.pid()) == MAVTOP || p.pid() == DARKPHOTON) return true;
849  if (std::abs(p.pid()) > 16 && std::abs(p.pid()) < 19) return true;
850  if (std::abs(p.pid()) > 31 && std::abs(p.pid()) < 39) return true;
851  if (std::abs(p.pid()) > 39 && std::abs(p.pid()) < 81) return true;
852  if (std::abs(p.pid()) > 6 && std::abs(p.pid()) < 9) return true;
853  if (isSUSY(p)) return true;
854  if (isNeutrinoRH(p.pid())) return true;
855  if (isGenericMultichargedParticle(p)) return true;
856  if (isTechnicolor(p)) return true;
857  if (isExcited(p)) return true;
858  if (isKK(p)) return true;
859  if (isHiddenValley(p)) return true;
860  if (isMonopole(p)) return true;
861  if (isDM(p.pid())) return true;
862  return false;
863 }
864 template<> inline bool isBSM(const int& p){
865  if (p == GRAVITON || std::abs(p) == MAVTOP || p == DARKPHOTON) return true;
866  if (std::abs(p) > 16 && std::abs(p) < 19) return true;
867  if (std::abs(p) > 31 && std::abs(p) < 38) return true;
868  if (std::abs(p) > 39 && std::abs(p) < 81) return true;
869  if (std::abs(p) > 6 && std::abs(p) < 9) return true;
870  auto value_digits = DecodedPID(p); return isBSM(value_digits);
871 }
872 
873 template<class T> inline bool isTransportable(const T& p){return isTransportable(p->pdg_id());}
874 template<> inline bool isTransportable(const DecodedPID& p){ return isPhoton(p.pid()) || isGeantino(p.pid()) || isHadron(p) || isLepton(p.pid()) || p.pid() == DARKPHOTON;}
875 template<> inline bool isTransportable(const int& p){ auto value_digits = DecodedPID(p); return isTransportable(value_digits);}
876 
878 template<class T> inline bool isValid(const T& p){return isValid(p->pdg_id());}
879 template<> inline bool isValid(const DecodedPID& p){
880  return p.pid() !=0 && ( isQuark(p) || isLepton(p) || isBoson(p) || isGlueball(p) ||
881  isTrajectory(p.pid()) || isGenSpecific(p.pid()) || isDiquark(p) ||
882  isBSM(p) || isHadron(p) || isNucleus(p) || isGeantino(p.pid()) ||
883  isPythia8Specific(p) ); }
884 template<> inline bool isValid(const int& p){ if (!p) return false; if (std::abs(p) < 42) return true;
885  if (isGenSpecific(p)) return true;
886  auto value_digits = DecodedPID(p); return isValid(value_digits);
887 }
888 
889 template<class T> inline int leadingQuark(const T& p) {return leadingQuark(p->pdg_id());}
890 template<> inline int leadingQuark(const DecodedPID& p){
891  if (isQuark(p.pid())) { return std::abs(p.pid());}
892  if (isMeson(p)) { return p.max_digit(1,3);}
893  if (isDiquark(p)) { return p.max_digit(2,4);}
894  if (isBaryon(p)) { return p.max_digit(1,4);}
895  if (isTetraquark(p)) { return p.max_digit(1,5);}
896  if (isPentaquark(p)) { return p.max_digit(1,6);}
897  if (isSUSY(p)) { // APID SUSY case
898  auto pp = p.shift(1);
899  if ( pp.ndigits() == 1 ) { return 0; } // Handle squarks
900  if ( pp.ndigits() == 3 ) { pp = DecodedPID(pp(1)); } // Handle ~q qbar pairs
901  if ( pp.ndigits() > 3 ) { pp = pp.shift(1); } // Drop gluinos and squarks
902  return leadingQuark(pp); }
903  return 0;
904 }
905 
906 template<> inline int leadingQuark(const int& p){ auto value_digits = DecodedPID(p); return leadingQuark(value_digits);}
907 
908 template<class T> inline bool isLightHadron(const T& p) { auto lq = leadingQuark(p); return (lq == DQUARK || lq == UQUARK||lq == SQUARK) && isHadron(p); }
909 template<class T> inline bool isHeavyHadron(const T& p) { auto lq = leadingQuark(p); return (lq == CQUARK || lq == BQUARK || lq == TQUARK ) && isHadron(p); }
910 template<class T> inline bool isStrangeHadron(const T& p) { return leadingQuark(p) == SQUARK && isHadron(p); }
911 template<class T> inline bool isCharmHadron(const T& p) { return leadingQuark(p) == CQUARK && isHadron(p); }
912 template<class T> inline bool isBottomHadron(const T& p) { return leadingQuark(p) == BQUARK && isHadron(p); }
913 template<class T> inline bool isTopHadron(const T& p) { return leadingQuark(p) == TQUARK && isHadron(p); }
914 
915 template<class T> inline bool isLightMeson(const T& p) { auto lq = leadingQuark(p); return (lq == DQUARK || lq == UQUARK||lq == SQUARK) && isMeson(p); }
916 template<class T> inline bool isHeavyMeson(const T& p) { auto lq = leadingQuark(p); return (lq == CQUARK || lq == BQUARK || lq == TQUARK) && isMeson(p); }
917 template<class T> inline bool isStrangeMeson(const T& p) { return leadingQuark(p) == SQUARK && isMeson(p); }
918 template<class T> inline bool isCharmMeson(const T& p) { return leadingQuark(p) == CQUARK && isMeson(p); }
919 template<class T> inline bool isBottomMeson(const T& p) { return leadingQuark(p) == BQUARK && isMeson(p); }
920 template<class T> inline bool isTopMeson(const T& p) { return leadingQuark(p) == TQUARK && isMeson(p); }
921 
922 template<class T> inline bool isCCbarMeson(const T& p) { return isCCbarMeson(p->pdg_id());}
923 template<> inline bool isCCbarMeson(const DecodedPID& p) { return leadingQuark(p) == CQUARK && isMeson(p) && (*(p.second.rbegin()+2)) == CQUARK && (*(p.second.rbegin()+1)) == CQUARK; }
924 template<> inline bool isCCbarMeson(const int& p) { return isCCbarMeson(DecodedPID(p)); }
925 
926 template<class T> inline bool isBBbarMeson(const T& p){ return isBBbarMeson(p->pdg_id());}
927 template<> inline bool isBBbarMeson(const DecodedPID& p) { return leadingQuark(p) == BQUARK && isMeson(p) && (*(p.second.rbegin()+2)) == BQUARK && (*(p.second.rbegin()+1)) == BQUARK; }
928 template<> inline bool isBBbarMeson(const int& p) { return isBBbarMeson(DecodedPID(p)); }
929 
930 
931 template<class T> inline bool isLightBaryon(const T& p) { auto lq = leadingQuark(p); return (lq == DQUARK || lq == UQUARK||lq == SQUARK) && isBaryon(p); }
932 template<class T> inline bool isHeavyBaryon(const T& p) { auto lq = leadingQuark(p); return (lq == CQUARK || lq == BQUARK || lq == TQUARK) && isBaryon(p); }
933 template<class T> inline bool isStrangeBaryon(const T& p) { return leadingQuark(p) == SQUARK && isBaryon(p); }
934 template<class T> inline bool isCharmBaryon(const T& p) { return leadingQuark(p) == CQUARK && isBaryon(p); }
935 template<class T> inline bool isBottomBaryon(const T& p) { return leadingQuark(p) == BQUARK && isBaryon(p); }
936 template<class T> inline bool isTopBaryon(const T& p) { return leadingQuark(p) == TQUARK && isBaryon(p); }
937 
938 
939 // APID: This function selects B-Hadrons which predominantly decay weakly. (Commonly used definition in GeneratorFilters package.)
940 // 5[1-4]1 L = J = 0, S = 0
941 // 5[1-5][1-4]2 J = 1/2, n_r = 0, n_L =0
942 template<class T> inline bool isWeaklyDecayingBHadron(const T& p) {return isWeaklyDecayingBHadron(p->pdg_id());}
943 template<> inline bool isWeaklyDecayingBHadron(const int& p) {
944  const int pid = std::abs(p);
945  return ( pid == 511 || // B0
946  pid == 521 || // B+
947  pid == 531 || // B_s0
948  pid == 541 || // B_c+
949  pid == 5122 || // Lambda_b0
950  pid == 5132 || // Xi_b-
951  pid == 5232 || // Xi_b0
952  pid == 5112 || // Sigma_b-
953  pid == 5212 || // Sigma_b0
954  pid == 5222 || // Sigma_b+
955  pid == 5332 || // Omega_b-
956  pid == 5142 || // Xi_bc0
957  pid == 5242 || // Xi_bc+
958  pid == 5412 || // Xi'_bc0
959  pid == 5422 || // Xi'_bc+
960  pid == 5342 || // Omega_bc0
961  pid == 5432 || // Omega'_bc0
962  pid == 5442 || // Omega_bcc+
963  pid == 5512 || // Xi_bb-
964  pid == 5522 || // Xi_bb0
965  pid == 5532 || // Omega_bb-
966  pid == 5542 ); // Omega_bbc0
967 }
968 template<> inline bool isWeaklyDecayingBHadron(const DecodedPID& p){ return isWeaklyDecayingBHadron(p.pid()); }
969 
970 
971 // APID: This function selects C-Hadrons which predominantly decay weakly. (Commonly used definition in GeneratorFilters package.)
972 // 4[1-3]1 L = J = 0, S = 0
973 // 4[1-4][1-3]2 J = 1/2, n_r = 0, n_L =0
974 // NB Omitting pid = 4322 (Xi'_C+) a this undergoes an EM rather than
975 // weak decay. (There was an old version of Herwig that decayed it
976 // weakly, but this was fixed in Herwig 7.)
977 template<class T> inline bool isWeaklyDecayingCHadron(const T& p) {return isWeaklyDecayingCHadron(p->pdg_id());}
978 template<> inline bool isWeaklyDecayingCHadron(const int& p) {
979  const int pid = std::abs(p);
980  return ( pid == 411 || // D+
981  pid == 421 || // D0
982  pid == 431 || // Ds+
983  pid == 4122 || // Lambda_c+
984  pid == 4132 || // Xi_c0
985  pid == 4232 || // Xi_c+
986  pid == 4212 || // Xi_c0
987  pid == 4332 || // Omega_c0
988  pid == 4412 || // Xi_cc+
989  pid == 4422 || // Xi_cc++
990  pid == 4432 ); // Omega_cc+
991 }
992 template<> inline bool isWeaklyDecayingCHadron(const DecodedPID& p){ return isWeaklyDecayingCHadron(p.pid()); }
993 
994 
995 template<class T> inline int charge3( const T& p){return charge3(p->pdg_id());}
996 template<class T> inline double fractionalCharge(const T& p){return fractionalCharge(p->pdg_id());}
997 template<class T> inline double charge( const T& p){
998  if (isGenericMultichargedParticle(p)) // BSM multi-charged particles might have a fractional charge that's not a multiple of 1/3
999  return fractionalCharge(p);
1000  else
1001  return 1.0*charge3(p)/3.0;
1002 }
1003 template<class T> inline double threeCharge( const T& p){ return charge3(p);}
1004 template<class T> inline bool isCharged( const T& p){ return charge3(p) != 0;}
1005 
1006 
1007 template<> inline int charge3(const DecodedPID& p) {
1008  auto ap = std::abs(p.pid());
1009  if (ap < TABLESIZE ) return p.pid() > 0 ? triple_charge.at(ap) : -triple_charge.at(ap);
1010  if (ap == K0) return 0;
1011  if (ap == GEANTINO0) return 0;
1012  if (ap == GEANTINOPLUS) return p.pid() > 0 ? 3 : -3;
1013  if (ap == MAVTOP) return p.pid() > 0 ? 2 : -2;
1014  size_t nq = 0;
1015  int sign = 1;
1016  int signmult = 1;
1017  int result=0;
1018  bool classified = false;
1019  if (!classified && isMeson(p)) { classified = true; nq = 2; if ((*(p.second.rbegin()+2)) == 2||(*(p.second.rbegin()+2)) == 4 ) { sign=-1;} signmult =-1; }
1020  if (!classified && isDiquark(p)) {return triple_charge.at(p(0))+triple_charge.at(p(1)); }
1021  if (!classified && isBaryon(p)) { classified = true; nq = 3; }
1022  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)); }
1023  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)); }
1024  if (!classified && isNucleus(p)) { return 3*numberOfProtons(p);}
1025  if (!classified && isSUSY(p)) {
1026  nq = 0;
1027  auto pp = p.shift(1);
1028  if (pp.ndigits() < 3 ) { return charge3(pp); } // super-partners of fundamental particles
1029  if (pp(0) == COMPOSITEGLUON) {
1030  if (pp(1) == COMPOSITEGLUON) { return 0; } // R-Glueballs
1031  if ( pp.ndigits() == 4 || pp.ndigits() == 5) {
1032  pp = pp.shift(1); // Remove gluino
1033  }
1034  }
1035  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
1036  if (pp.ndigits() == 4) { classified = true; nq = 3; } // states with squark-quark-quark or quark-quark-quark
1037  }
1038  if (!classified && isHiddenValley(p)) { // Hidden Valley particles
1039  auto pp = p.shift(2);
1040  if (!classified && isMeson(pp)) { classified = true; nq = 2; if ((*(pp.second.rbegin()+2)) == 2||(*(pp.second.rbegin()+2)) == 4 ) { sign=-1;} signmult =-1; }
1041  if (!classified && isDiquark(pp)) {return triple_charge.at(pp(0))+triple_charge.at(pp(1)); }
1042  if (!classified && isBaryon(pp)) { classified = true; nq = 3; }
1043 
1044  }
1045  if (!classified && isKK(p)) { // Kaluza-Klein particles
1046  auto pp = p.shift(2);
1047  auto ap = std::abs(pp.pid());
1048  if (ap < TABLESIZE ) return pp.pid() > 0 ? triple_charge.at(ap) : -triple_charge.at(ap);
1049 
1050  }
1051  if (!classified && isDM(p.pid())) { //Dark Matter Particles
1052  if (p.ndigits() == 7){ // Determining the charges for the more elaborate, 7-digit DM codes
1053  auto pp = p.shift(3); // The first two digits indicate the particle is DM, the third indicates left/right-handedness (see 11(j))
1054  auto ap = std::abs(pp.pid());
1055  if (ap < TABLESIZE ) return pp.pid() > 0 ? triple_charge.at(ap) : -triple_charge.at(ap);
1056  }else if (std::abs(p.pid()) < TABLESIZE) return p.pid() > 0 ? triple_charge.at(ap) : -triple_charge.at(ap); // Just to make sure the correct charge is returned for DM 51-60
1057  }
1058  if (!classified && isMonopole(p)) {
1061  result = 3*(p(3)*100 + p(4)*10 + p(5));
1062  return ( (p.pid() > 0 && p(2) == 1) || (p.pid() < 0 && p(2) == 2) ) ? result : -result;
1063  }
1064  if (!classified && isGenericMultichargedParticle(p)) {
1065  double abs_charge = 0.0;
1066  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
1067  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
1068  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
1069  return p.pid() > 0 ? abs_threecharge : -1 * abs_threecharge;
1070  }
1071  for (auto r = p.second.rbegin() + 1; r != p.second.rbegin() + 1 + nq; ++r) {
1072  result += triple_charge.at(*r)*sign;
1073  sign*=signmult;
1074  }
1075  return p.pid() > 0 ? result : -result;
1076 }
1077 template<> inline int charge3(const int& p){
1078  int ap = std::abs(p);
1079  if (ap < TABLESIZE) return p > 0 ? triple_charge.at(ap):-triple_charge.at(ap);
1080  auto value_digits = DecodedPID(p);
1081  return charge3(value_digits);
1082 }
1083 
1084 
1085 template<class T> inline bool isNeutral( const T& p){ return p->pdg_id() != 0 && charge3(p) == 0;}
1086 template<> inline bool isNeutral(const DecodedPID& p){ return p.pid() != 0 && charge3(p) == 0;}
1087 template<> inline bool isNeutral(const int& p){ auto value_digits = DecodedPID(p); return isNeutral(value_digits);}
1088 
1089 
1090 template<> inline double fractionalCharge(const DecodedPID& p) {
1091  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
1092  double abs_charge = 0;
1093  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
1094  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
1095  return p.pid() > 0 ? abs_charge : -1 * abs_charge;
1096 }
1097 template<> inline double fractionalCharge(const int& p){auto value_digits = DecodedPID(p); return fractionalCharge(value_digits);}
1098 
1099 // APID: Including Z' and Z'' as EM interacting.
1100 template<class T> inline bool isEMInteracting(const T& p){return isEMInteracting(p->pdg_id());}
1101 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));}
1102 
1103 template<class T> inline bool isParton(const T& p) { return isQuark(p)||isGluon(p);}
1104 
1105 // APID: Intended to return 2J
1106 // Useful for G4ParticleDefinition constructor
1107 template<class T> inline int spin2(const T& p) { return spin2(p->pdg_id()); }
1108 template<> inline int spin2(const DecodedPID& p) {
1109  if (isSUSY(p)) {
1110  auto pp = p.shift(1);
1111  auto ap = std::abs(pp.pid());
1112  if (ap < TABLESIZE ) { return std::abs(double_spin.at(ap)-1); } // sparticles (0->1, 1 -> 0, 2->1, 4->3)
1113  return p.last()-1; // R-Hadrons (p.last() == 2J +1)
1114  }
1115  if (isHiddenValley(p)) { //Hidden Valley spins
1116  auto pp = p.shift(2);
1117  if (isHadron(pp)) { return pp.last()-1; } // Hadrons (p.last == 2J+1 - special cases handled above)
1118  }
1119  if (isKK(p)) { // Kaluza-Klein spins
1120  auto pp = p.shift(2);
1121  auto ap = std::abs(pp.pid());
1122  if (ap < TABLESIZE ) { return double_spin.at(ap); } // fundamental particles
1123  }
1124  if (isDM(std::abs(p.pid()))) { //DM spins
1125  if (p.ndigits() == 7) { // Determining the spins for the more elaborate, 7-digit DM codes
1126  auto pp = p.shift(3); // The first two digits indicate the particle is DM, the third indicates left/right-handedness (see 11(j))
1127  auto ap = std::abs(pp.pid());
1128  if (ap < TABLESIZE) { return double_spin.at(ap); } // fundamental particles
1129  }else if (std::abs(p.pid()) < TABLESIZE) { // Just to make sure the correct spin is returned for DM 51-60
1130  return std::abs(double_spin.at(std::abs(p.pid())));
1131  }
1132  }
1133  auto ap = std::abs(p.pid());
1134  if (ap == K0S) { return 0; }
1135  if (ap == K0L) { return 0; }
1136  if (ap == MAVTOP) { return 1; } // TODO check this
1137  if (ap == DARKPHOTON) { return 2; } // TODO check this
1138  if (ap < TABLESIZE ) { return double_spin.at(ap); } // fundamental particles
1139  if (isHadron(p)) { return p.last()-1; } // Hadrons (p.last == 2J+1 - special cases handled above)
1140  if (isMonopole(p)) { return 0; } // PDG 11i - For now no spin information is provided. Also matches the definition in the G4Extensions/Monopole package.
1141  if (isGenericMultichargedParticle(p)) { return 0; } // APID Matches the definition in the G4Extensions/Monopole package.
1142  if (isNucleus(p)) { return 1; } // TODO need to explicitly deal with nuclei
1143  return p.last() > 0 ? 1 : 0; // Anything else - best guess
1144 }
1145 template<> inline int spin2(const int& p){ auto value_digits = DecodedPID(p); return spin2(value_digits);}
1146 
1147 template<class T> inline double spin(const T& p) { return spin(p->pdg_id()); }
1148 template<> inline double spin(const DecodedPID& p) { return 1.0*spin2(p)/2.0; }
1149 template<> inline double spin(const int& p){ auto value_digits = DecodedPID(p); return spin(value_digits);}
1150 
1151 
1152 // APID: Returns an unordered list of the quarks contained by the current particle
1153 template<class T> inline std::vector<int> containedQuarks(const T& p) { return containedQuarks(p->pdg_id()); }
1154 template<> inline std::vector<int> containedQuarks(const int& p) {
1155  auto pp = DecodedPID(p);
1156  std::vector<int> quarks;
1157  if (isQuark(pp.pid())) { quarks.push_back(std::abs(pp.pid())); }
1158  else if (isDiquark(pp)) { quarks.push_back(pp(0)); quarks.push_back(pp(1)); }
1159  else if (isMeson(pp)) { quarks.push_back(*(pp.second.rbegin() + 1)); quarks.push_back(*(pp.second.rbegin()+2)); }
1160  else if (isBaryon(pp)) { for (size_t digit = 1; digit < 4; ++digit) { quarks.push_back(*(pp.second.rbegin() + digit)); } }
1161  else if (isTetraquark(pp)) { for (size_t digit = 1; digit < 5; ++digit) { quarks.push_back(*(pp.second.rbegin() + digit)); } }
1162  else if (isPentaquark(pp)) { for (size_t digit = 1; digit < 6; ++digit) { quarks.push_back(*(pp.second.rbegin() + digit)); } }
1163  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));
1164  const int n_uquarks = A + Z; const int n_dquarks = 2*A - Z - L; const int n_squarks = L;
1165  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); }
1166  else if (isSUSY(pp)) { // APID SUSY case
1167  pp = pp.shift(1);
1168  if ( pp.ndigits() > 1 ) { // skip squarks
1169  if ( pp.ndigits() == 3 ) { pp = DecodedPID(pp(1)); } // Handle ~q qbar pairs
1170  if ( pp.ndigits() > 3 ) { pp = pp.shift(1); } // Drop gluinos and squarks
1171  return containedQuarks(pp.pid());
1172  }
1173  }
1174  return quarks;
1175 }
1176 template<> inline std::vector<int> containedQuarks(const DecodedPID& p) { return containedQuarks(p.pid()); }
1177 
1178 template<class T> inline bool isStrongInteracting(const T& p){return isStrongInteracting(p->pdg_id());}
1179 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
1180 
1181 #endif
isStrange
bool isStrange(const T &p)
Definition: AtlasPID.h:176
isStrangeMeson
bool isStrangeMeson(const T &p)
Definition: AtlasPID.h:917
beamspotman.r
def r
Definition: beamspotman.py:672
isGaugino
bool isGaugino(const T &p)
Definition: AtlasPID.h:504
isBottomMeson
bool isBottomMeson(const T &p)
Definition: AtlasPID.h:919
isStrongInteracting
bool isStrongInteracting(const T &p)
Definition: AtlasPID.h:1178
isNucleus
bool isNucleus(const T &p)
PDG rule 16 Nuclear codes are given as 10-digit numbers ±10LZZZAAAI.
Definition: AtlasPID.h:702
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:386
isHeavyHadron
bool isHeavyHadron(const T &p)
Definition: AtlasPID.h:909
isTetraquark
bool isTetraquark(const T &p)
PDG rule 14 The 9-digit tetra-quark codes are ±1nrnLnq1nq20nq3nq4nJ.
Definition: AtlasPID.h:325
isRMeson
bool isRMeson(const T &p)
Definition: AtlasPID.h:574
isBottomBaryon
bool isBottomBaryon(const T &p)
Definition: AtlasPID.h:935
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:496
isTopBaryon
bool isTopBaryon(const T &p)
Definition: AtlasPID.h:936
hasCharm
bool hasCharm(const T &p)
Definition: AtlasPID.h:737
isBSM
bool isBSM(const T &p)
APID: graviton and all Higgs extensions are BSM.
Definition: AtlasPID.h:846
baryonNumber
double baryonNumber(const T &p)
Definition: AtlasPID.h:773
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:668
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:369
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:425
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:244
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:635
isSquarkRH
bool isSquarkRH(const T &p)
Definition: AtlasPID.h:474
threeCharge
double threeCharge(const T &p)
Definition: AtlasPID.h:1003
isCharmMeson
bool isCharmMeson(const T &p)
Definition: AtlasPID.h:918
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:562
CaloClusterMLCalib::epsilon
constexpr float epsilon
Definition: CaloClusterMLGaussianMixture.h:16
hasBottom
bool hasBottom(const T &p)
Definition: AtlasPID.h:738
numberOfLambdas
int numberOfLambdas(const T &p)
Definition: AtlasPID.h:824
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:212
isResonance
bool isResonance(const T &p)
Definition: AtlasPID.h:400
isParton
bool isParton(const T &p)
Definition: AtlasPID.h:1103
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:458
isSleptonLH
bool isSleptonLH(const T &p)
Definition: AtlasPID.h:488
isHeavyBaryon
bool isHeavyBaryon(const T &p)
Definition: AtlasPID.h:932
isValid
bool isValid(const T &p)
Av: we implement here an ATLAS-sepcific convention: all particles which are 99xxxxx are fine.
Definition: AtlasPID.h:878
checkRpcDigits.digit
digit
Definition: checkRpcDigits.py:186
isLightBaryon
bool isLightBaryon(const T &p)
Definition: AtlasPID.h:931
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:912
isLightMeson
bool isLightMeson(const T &p)
Definition: AtlasPID.h:915
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:684
isSMLepton
bool isSMLepton(const T &p)
APID: the fourth generation leptons are not standard model leptons.
Definition: AtlasPID.h:194
isSMQuark
bool isSMQuark(const T &p)
Definition: AtlasPID.h:172
isGluon
bool isGluon(const T &p)
Definition: AtlasPID.h:373
DeMoUpdate.reverse
reverse
Definition: DeMoUpdate.py:563
isHiggs
bool isHiggs(const T &p)
APID: HIGGS boson is only one particle.
Definition: AtlasPID.h:390
fractionalCharge
double fractionalCharge(const T &p)
Definition: AtlasPID.h:996
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:167
A
python.utils.AtlRunQueryDQUtils.p
p
Definition: AtlRunQueryDQUtils.py:209
isLightHadron
bool isLightHadron(const T &p)
Definition: AtlasPID.h:908
isCCbarMeson
bool isCCbarMeson(const T &p)
Definition: AtlasPID.h:922
strangeness
int strangeness(const T &p)
Definition: AtlasPID.h:786
isBottom
bool isBottom(const T &p)
Definition: AtlasPID.h:182
isTopHadron
bool isTopHadron(const T &p)
Definition: AtlasPID.h:913
lumiFormat.i
int i
Definition: lumiFormat.py:85
leadingQuark
int leadingQuark(const T &p)
Definition: AtlasPID.h:889
isPythia8Specific
bool isPythia8Specific(const T &p)
Definition: AtlasPID.h:411
beamspotman.n
n
Definition: beamspotman.py:727
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:645
isWeaklyDecayingCHadron
bool isWeaklyDecayingCHadron(const T &p)
Definition: AtlasPID.h:977
hasTop
bool hasTop(const T &p)
Definition: AtlasPID.h:739
isZ
bool isZ(const T &p)
Definition: AtlasPID.h:379
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:342
isCharmBaryon
bool isCharmBaryon(const T &p)
Definition: AtlasPID.h:934
isGlueball
bool isGlueball(const T &p)
APID: Definition of Glueballs: SM glueballs 99X (X=1,5), 999Y (Y=3,7)
Definition: AtlasPID.h:441
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:512
ParticleGun_EoverP_Config.pid
pid
Definition: ParticleGun_EoverP_Config.py:62
sign
int sign(int a)
Definition: TRT_StrawNeighbourSvc.h:108
isBBbarMeson
bool isBBbarMeson(const T &p)
Definition: AtlasPID.h:926
isNeutral
bool isNeutral(const T &p)
Definition: AtlasPID.h:1085
hasStrange
bool hasStrange(const T &p)
Definition: AtlasPID.h:736
isChLepton
bool isChLepton(const T &p)
APID: the fourth generation leptons are leptons.
Definition: AtlasPID.h:199
isSUSY
bool isSUSY(const T &p)
Definition: AtlasPID.h:623
isTau
bool isTau(const T &p)
Definition: AtlasPID.h:208
isGraviton
bool isGraviton(const T &p)
Definition: AtlasPID.h:397
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:910
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:658
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:536
isTopMeson
bool isTopMeson(const T &p)
Definition: AtlasPID.h:920
isHadron
bool isHadron(const T &p)
Definition: AtlasPID.h:351
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:420
charge
double charge(const T &p)
Definition: AtlasPID.h:997
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:282
isGeantino
bool isGeantino(const T &p)
Definition: AtlasPID.h:437
isSlepton
bool isSlepton(const T &p)
Definition: AtlasPID.h:482
isWeaklyDecayingBHadron
bool isWeaklyDecayingBHadron(const T &p)
Definition: AtlasPID.h:942
isFourthGeneration
bool isFourthGeneration(const T &p)
Is this a 4th generation fermion? APID: 4th generation fermions are not standard model particles.
Definition: AtlasPID.h:220
spin2
int spin2(const T &p)
Definition: AtlasPID.h:1107
isW
bool isW(const T &p)
Definition: AtlasPID.h:382
isRBaryon
bool isRBaryon(const T &p)
Definition: AtlasPID.h:589
isCharged
bool isCharged(const T &p)
Definition: AtlasPID.h:1004
isTop
bool isTop(const T &p)
Definition: AtlasPID.h:185
isTransportable
bool isTransportable(const T &p)
Definition: AtlasPID.h:873
numberOfProtons
int numberOfProtons(const T &p)
Definition: AtlasPID.h:833
DeMoScan.first
bool first
Definition: DeMoScan.py:534
isLepton
bool isLepton(const T &p)
APID: the fourth generation leptons are leptons.
Definition: AtlasPID.h:189
isCharm
bool isCharm(const T &p)
Definition: AtlasPID.h:179
isPhoton
bool isPhoton(const T &p)
Definition: AtlasPID.h:376
charge3
int charge3(const T &p)
Definition: AtlasPID.h:995
hasSquark
bool hasSquark(const T &p, const int &q)
Definition: AtlasPID.h:610
isSquarkLH
bool isSquarkLH(const T &p)
Definition: AtlasPID.h:466
isSMNeutrino
bool isSMNeutrino(const T &p)
Definition: AtlasPID.h:215
isCharmHadron
bool isCharmHadron(const T &p)
Definition: AtlasPID.h:911
isRHadron
bool isRHadron(const T &p)
Definition: AtlasPID.h:603
extractSporadic.q
list q
Definition: extractSporadic.py:97
isQuarkonium
bool isQuarkonium(const T &p)
Is this a heavy-flavour quarkonium meson?
Definition: AtlasPID.h:273
isTechnicolor
bool isTechnicolor(const T &p)
PDG rule 11e Technicolor states have n = 3, with technifermions treated like ordinary fermions.
Definition: AtlasPID.h:524
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:359
isHeavyMeson
bool isHeavyMeson(const T &p)
Definition: AtlasPID.h:916
isStrangeBaryon
bool isStrangeBaryon(const T &p)
Definition: AtlasPID.h:933
DecodedPID::last
const int & last() const
Definition: AtlasPID.h:27
isElectron
bool isElectron(const T &p)
Definition: AtlasPID.h:202
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:227
baryonNumber3
int baryonNumber3(const T &p)
Definition: AtlasPID.h:748
pow
constexpr int pow(int base, int exp) noexcept
Definition: ap_fixedTest.cxx:15
isEMInteracting
bool isEMInteracting(const T &p)
Definition: AtlasPID.h:1100
containedQuarks
std::vector< int > containedQuarks(const T &p)
Definition: AtlasPID.h:1153
spin
double spin(const T &p)
Definition: AtlasPID.h:1147
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:394
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:205
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:408