Version: SMASH-2.2
pdgcode.h
Go to the documentation of this file.
1 /*
2  *
3  * Copyright (c) 2014-2022
4  * SMASH Team
5  *
6  * GNU General Public License (GPLv3 or later)
7  *
8  */
9 
10 #ifndef SRC_INCLUDE_SMASH_PDGCODE_H_
11 #define SRC_INCLUDE_SMASH_PDGCODE_H_
12 
13 #include <algorithm>
14 #include <array>
15 #include <cassert>
16 #include <cstdlib>
17 #include <iosfwd>
18 #include <sstream>
19 #include <stdexcept>
20 #include <string>
21 
22 #include <iostream>
23 
24 #include "pdgcode_constants.h"
25 
26 namespace smash {
27 
108 class PdgCode {
109  public:
114  struct InvalidPdgCode : public std::invalid_argument {
115  using std::invalid_argument::invalid_argument;
116  };
117 
118  /****************************************************************************
119  * *
120  * First, the constructors *
121  * *
122  ****************************************************************************/
123 
125  PdgCode() : dump_(0x0) {}
131  explicit PdgCode(const std::string& codestring) {
132  set_from_string(codestring);
133  }
134 
143  PdgCode(std::int32_t codenumber) : dump_(0x0) { // NOLINT(runtime/explicit)
144  digits_.antiparticle_ = false;
145  if (codenumber < 0) {
146  digits_.antiparticle_ = true;
147  codenumber = -codenumber;
148  }
149  set_fields(codenumber);
150  }
155  explicit PdgCode(const std::uint32_t abscode) : dump_(0x0) {
156  // use the first bit for the antiparticle_ boolean.
157  digits_.antiparticle_ = ((abscode & 0x80000000u) != 0);
158  set_fields(abscode);
159  }
160 
161  /****************************************************************************
162  * *
163  * test function and export functions *
164  * *
165  ****************************************************************************/
166 
180  inline int test_code() const {
181  int fail = 0;
182  if (digits_.n_ > 9) {
183  fail |= 1 << 6;
184  }
185  if (digits_.n_R_ > 9) {
186  fail |= 1 << 5;
187  }
188  if (digits_.n_L_ > 9) {
189  fail |= 1 << 4;
190  }
191  if (digits_.n_q1_ > 6) {
192  fail |= 1 << 3;
193  }
194  if (digits_.n_q2_ > 6) {
195  fail |= 1 << 2;
196  }
197  if (digits_.n_q3_ > 6) {
198  fail |= 1 << 1;
199  }
200  if (digits_.n_J_ > 15) {
201  fail |= 1;
202  }
203  return fail;
204  }
205 
214  void check() const {
215  // n_J must be odd for mesons and even for baryons (and cannot be zero)
216  if (is_hadron()) {
217  if (baryon_number() == 0) {
218  // mesons: special cases K0_L=0x130 and K0_S=0x310
219  if ((digits_.n_J_ % 2 == 0) && dump() != 0x130 && dump() != 0x310) {
220  throw InvalidPdgCode("Invalid PDG code " + string() +
221  " (meson with even n_J)");
222  }
223  } else {
224  if ((digits_.n_J_ % 2 != 0) || digits_.n_J_ == 0) {
225  throw InvalidPdgCode("Invalid PDG code " + string() +
226  " (baryon with odd n_J)");
227  }
228  }
229  } else {
230  if (digits_.n_J_ == 0 && dump() != 0x0) {
231  throw InvalidPdgCode("Invalid PDG code " + string() + " (n_J==0)");
232  }
233  }
234  /* The antiparticle flag only makes sense for particle types
235  * that have an antiparticle. */
236  if (digits_.antiparticle_ && !has_antiparticle()) {
237  throw InvalidPdgCode("Invalid PDG code " + string() +
238  " (cannot be negative)");
239  }
240  }
241 
243  inline std::uint32_t dump() const {
244  // this cuts the three unused bits.
245  return (dump_ & 0x8fffffff);
246  }
247 
249  inline std::int32_t code() const { return antiparticle_sign() * ucode(); }
250 
252  inline std::string string() const {
253  std::stringstream ss;
254  ss << get_decimal();
255  return ss.str();
256  }
257 
260  PdgCode result = *this;
261  result.digits_.antiparticle_ = !digits_.antiparticle_;
262  return result;
263  }
264 
269  static PdgCode from_decimal(const int pdgcode_decimal) {
270  // Nucleus and special codes with 2J+1 > 9
271  if (std::abs(pdgcode_decimal) > 1E7) {
272  return PdgCode(std::to_string(pdgcode_decimal));
273  }
274  int a = pdgcode_decimal;
275  int hex_pdg = 0, tmp = 1;
276  while (a) {
277  hex_pdg += (a % 10) * tmp;
278  tmp *= 16;
279  a = a / 10;
280  }
281  return PdgCode(hex_pdg);
282  }
283 
284  /****************************************************************************
285  * *
286  * accessors of various properties *
287  * *
288  ****************************************************************************/
289 
291  inline bool is_nucleus() const {
292  assert(digits_.is_nucleus_ == nucleus_.is_nucleus_);
293  return nucleus_.is_nucleus_;
294  }
295 
297  inline bool is_hadron() const {
298  return (digits_.n_q3_ != 0 && digits_.n_q2_ != 0 && !is_nucleus());
299  }
300 
302  inline bool is_lepton() const {
303  return (digits_.n_q1_ == 0 && digits_.n_q2_ == 0 && digits_.n_q3_ == 1 &&
304  !is_nucleus());
305  }
306 
308  inline int baryon_number() const {
309  if (is_nucleus()) {
310  return static_cast<int>(nucleus_.A_) * antiparticle_sign();
311  }
312  if (!is_hadron() || digits_.n_q1_ == 0) {
313  return 0;
314  }
315  return antiparticle_sign();
316  }
318  inline bool is_baryon() const { return is_hadron() && digits_.n_q1_ != 0; }
319 
321  inline bool is_meson() const { return is_hadron() && digits_.n_q1_ == 0; }
322 
324  inline bool is_nucleon() const {
325  const auto abs_code = std::abs(code());
326  return (abs_code == pdg::p || abs_code == pdg::n);
327  }
328 
330  inline bool is_proton() const {
331  const auto abs_code = std::abs(code());
332  return (abs_code == pdg::p);
333  }
334 
336  inline bool is_neutron() const {
337  const auto abs_code = std::abs(code());
338  return (abs_code == pdg::n);
339  }
340 
342  inline bool is_Nstar1535() const {
343  const auto abs_code = std::abs(code());
344  return (abs_code == pdg::N1535_p || abs_code == pdg::N1535_z);
345  }
346 
348  inline bool is_Delta() const {
349  const auto abs_code = std::abs(code());
350  return (abs_code == pdg::Delta_pp || abs_code == pdg::Delta_p ||
351  abs_code == pdg::Delta_z || abs_code == pdg::Delta_m);
352  }
353 
355  inline bool is_hyperon() const { return is_hadron() && digits_.n_q1_ == 3; }
356 
358  inline bool is_Omega() const {
359  return is_hyperon() && digits_.n_q2_ == 3 && digits_.n_q3_ == 3;
360  }
361 
363  inline bool is_Xi() const {
364  return is_hyperon() && digits_.n_q2_ == 3 && digits_.n_q3_ != 3;
365  }
366 
368  inline bool is_Lambda() const {
369  return is_hyperon() && digits_.n_q2_ == 1 && digits_.n_q3_ == 2;
370  }
371 
373  inline bool is_Sigma() const {
374  return is_hyperon() && digits_.n_q2_ != 3 && !is_Lambda();
375  }
376 
378  inline bool is_kaon() const {
379  const auto abs_code = std::abs(code());
380  return (abs_code == pdg::K_p) || (abs_code == pdg::K_z);
381  }
382 
384  inline bool is_pion() const {
385  const auto c = code();
386  return (c == pdg::pi_z) || (c == pdg::pi_p) || (c == pdg::pi_m);
387  }
388 
390  inline bool is_rho() const {
391  const auto c = code();
392  return (c == pdg::rho_z) || (c == pdg::rho_p) || (c == pdg::rho_m);
393  }
394 
396  inline bool is_deuteron() const {
397  return is_nucleus() && nucleus_.A_ == 2 && nucleus_.Z_ == 1 &&
398  nucleus_.n_Lambda_ == 0 && nucleus_.I_ == 0;
399  }
400 
402  inline bool is_triton() const {
403  return is_nucleus() && nucleus_.A_ == 3 && nucleus_.Z_ == 1 &&
404  nucleus_.n_Lambda_ == 0 && nucleus_.I_ == 0;
405  }
406 
411  bool has_antiparticle() const {
412  if (is_nucleus()) {
413  return true;
414  }
415  if (is_hadron()) {
416  return (baryon_number() != 0) || (digits_.n_q2_ != digits_.n_q3_);
417  } else {
418  return digits_.n_q3_ == 1; // leptons!
419  }
420  }
421 
427  inline int isospin3() const {
428  /* net_quark_number(2) is the number of u quarks,
429  * net_quark_number(1) is the number of d quarks. */
430  return net_quark_number(2) - net_quark_number(1);
431  }
432 
440  inline double frac_strange() const {
441  /* The quarkonium state has 0 net strangeness
442  * but there are actually 2 strange quarks out of 2 total */
443  if (is_hadron() && digits_.n_q3_ == 3 && digits_.n_q2_ == 3) {
444  return 1.;
445  } else {
446  // For all other cases, there isn't both a strange and anti-strange
447  if (is_baryon()) {
448  return std::abs(strangeness()) / 3.;
449  } else if (is_meson()) {
450  return std::abs(strangeness()) / 2.;
451  } else {
452  /* If not baryon or meson, this should be 0, as AQM does not
453  * extend to non-hadrons */
454  return 0.;
455  }
456  }
457  }
458 
464  inline int strangeness() const { return -net_quark_number(3); }
465 
471  inline int charmness() const { return +net_quark_number(4); }
472 
478  inline int bottomness() const { return -net_quark_number(5); }
479 
488  int charge() const {
489  if (is_hadron() || is_nucleus()) {
490  // Q will accumulate 3*charge (please excuse the upper case. I
491  // want to distinguish this from q which might be interpreted as
492  // shorthand for "quark".)
493  int Q = 0;
494  /* This loops over d,u,s,c,b,t quarks (the latter can be safely ignored,
495  * but I don't think this will be a bottle neck. */
496  for (int i = 1; i < 7; i++) {
497  /* u,c,t quarks have charge = 2/3 e, while d,s,b quarks have -1/3 e.
498  * The antiparticle sign is already in net_quark_number. */
499  Q += (i % 2 == 0 ? 2 : -1) * net_quark_number(i);
500  }
501  return Q / 3;
502  }
503  /* non-hadron:
504  * Leptons: 11, 13, 15 are e, μ, τ and have a charge -1, while
505  * 12, 14, 16 are the neutrinos that have no charge. */
506  if (digits_.n_q3_ == 1) {
507  return -1 * (digits_.n_J_ % 2) * antiparticle_sign();
508  }
509  /* Bosons: 24 is the W+, all else is uncharged.
510  * we ignore the first digits so that this also finds strange gauge
511  * boson "resonances" (in particular, \f$\tilde \chi_1^+\f$ with PDG
512  * Code 1000024). */
513  if ((dump_ & 0x0000ffff) == 0x24) {
514  return antiparticle_sign();
515  }
516  // default (this includes all other Bosons) is 0.
517  return 0;
518  }
519 
529  inline unsigned int spin() const {
530  if (is_nucleus()) {
531  // Generally spin of a nucleus cannot be simply guessed, it should be
532  // provided from some table. However, here we only care about a
533  // limited set of light nuclei with A <= 4.
534  if (nucleus_.A_ == 2) {
535  // Deuteron spin is 1
536  return 2;
537  } else if (nucleus_.A_ == 3) {
538  // Tritium and He-3 spin are 1/2
539  // Hypertriton spin is not firmly determined, but decay branching ratios
540  // indicate spin 1/2
541  return 1;
542  } else if (nucleus_.A_ == 4) {
543  // He-4 spin is 0
544  return 0;
545  }
546  throw std::runtime_error("Unknown spin of nucleus.");
547  // Alternative possibility is to guess 1/2 for fermions and 0 for bosons
548  // as 2 * (nucleus_.A_ % 2).
549  }
550 
551  if (is_hadron()) {
552  if (digits_.n_J_ == 0) {
553  return 0; // special cases: K0_L=0x130 & K0_S=0x310
554  } else {
555  return digits_.n_J_ - 1;
556  }
557  }
558  /* this assumes that we only have white particles (no single
559  * quarks): Electroweak fermions have 11-17, so the
560  * second-to-last-digit is the spin. The same for the Bosons: they
561  * have 21-29 and 2spin = 2 (this fails for the Higgs). */
562  return digits_.n_q3_;
563  }
565  inline unsigned int spin_degeneracy() const {
566  if (is_hadron() && digits_.n_J_ > 0) {
567  return digits_.n_J_;
568  }
569  return spin() + 1;
570  }
572  inline int antiparticle_sign() const {
573  return (digits_.antiparticle_ ? -1 : +1);
574  }
576  inline std::int32_t quarks() const {
577  if (!is_hadron() || is_nucleus()) {
578  return 0;
579  }
580  return chunks_.quarks_;
581  }
582 
592  std::array<int, 3> quark_content() const {
593  std::array<int, 3> result = {static_cast<int>(digits_.n_q1_),
594  static_cast<int>(digits_.n_q2_),
595  static_cast<int>(digits_.n_q3_)};
596  if (is_hadron()) {
597  // Antibaryons
598  if (digits_.n_q1_ != 0 && digits_.antiparticle_) {
599  for (size_t i = 0; i < 3; i++) {
600  result[i] = -result[i];
601  }
602  }
603  // Mesons
604  if (digits_.n_q1_ == 0) {
605  // Own antiparticle
606  if (digits_.n_q2_ == digits_.n_q3_) {
607  result[2] = -result[2];
608  } else {
609  // Like pi-
610  if (digits_.antiparticle_) {
611  result[1] = -result[1];
612  // Like pi+
613  } else {
614  result[2] = -result[2];
615  }
616  }
617  // add extra minus sign according to the pdg convention
618  if (digits_.n_q2_ != digits_.n_q3_ && digits_.n_q2_ % 2 == 1) {
619  for (int i = 1; i <= 2; i++) {
620  result[i] = -result[i];
621  }
622  }
623  }
624  } else {
625  result = {0, 0, 0};
626  }
627  return result;
628  }
629 
640  bool contains_enough_valence_quarks(int valence_quarks_required) const;
641 
642  /****************************************************************************
643  * *
644  * operations with more than one PDG Code *
645  * *
646  ****************************************************************************/
647 
652  inline bool operator<(const PdgCode rhs) const {
653  return dump_ < rhs.dump_;
654  /* the complex thing to do here is to calculate:
655  * code() < rhs.code()
656  * but for getting a total order that's overkill. The uint32_t value in
657  * dump_ works just fine. */
658  }
659 
661  inline bool operator==(const PdgCode rhs) const { return dump_ == rhs.dump_; }
662 
664  inline bool operator!=(const PdgCode rhs) const { return !(*this == rhs); }
665 
667  inline bool is_antiparticle_of(const PdgCode rhs) const {
668  return code() == -rhs.code();
669  }
670 
672  friend std::istream& operator>>(std::istream& is, PdgCode& code);
673 
679  static PdgCode invalid() { return PdgCode(0x0); }
680 
690  int32_t get_decimal() const {
691  if (is_nucleus()) {
692  // ±10LZZZAAAI
693  return antiparticle_sign() *
694  (nucleus_.I_ + 10 * nucleus_.A_ + 10000 * nucleus_.Z_ +
695  10000000 * nucleus_.n_Lambda_ + 1000000000);
696  }
697  int n_J_1 = 0;
698  int n_J_2 = digits_.n_J_;
699  if (n_J_2 > 9) {
700  n_J_1 = n_J_2 - 9;
701  n_J_2 = 9;
702  }
703  return antiparticle_sign() *
704  (n_J_2 + digits_.n_q3_ * 10 + digits_.n_q2_ * 100 +
705  digits_.n_q1_ * 1000 + digits_.n_L_ * 10000 +
706  digits_.n_R_ * 100000 + digits_.n_ * 1000000 + n_J_1 * 10000000);
707  }
708 
710  void deexcite() {
711  if (!is_nucleus()) {
712  chunks_.excitation_ = 0;
713  } else {
714  nucleus_.I_ = 0;
715  }
716  }
717 
728  int net_quark_number(const int quark) const;
729 
731  int nucleus_p() const {
732  return (is_nucleus() && !nucleus_.antiparticle_) ? nucleus_.Z_ : 0;
733  }
735  int nucleus_n() const {
736  return (is_nucleus() && !nucleus_.antiparticle_)
737  ? nucleus_.A_ - nucleus_.Z_ - nucleus_.n_Lambda_
738  : 0;
739  }
741  int nucleus_La() const {
742  return (is_nucleus() && !nucleus_.antiparticle_) ? nucleus_.n_Lambda_ : 0;
743  }
745  int nucleus_ap() const {
746  return (is_nucleus() && nucleus_.antiparticle_) ? nucleus_.Z_ : 0;
747  }
749  int nucleus_an() const {
750  return (is_nucleus() && nucleus_.antiparticle_)
751  ? nucleus_.A_ - nucleus_.Z_ - nucleus_.n_Lambda_
752  : 0;
753  }
755  int nucleus_aLa() const {
756  return (is_nucleus() && nucleus_.antiparticle_) ? nucleus_.n_Lambda_ : 0;
757  }
759  int nucleus_A() const { return is_nucleus() ? nucleus_.A_ : 0; }
760 
761  private:
767  union {
772  struct {
773 #if defined(LITTLE_ENDIAN_ARCHITECTURE) || defined(DOXYGEN)
775  std::uint32_t n_J_ : 4;
777  std::uint32_t n_q3_ : 4;
779  std::uint32_t n_q2_ : 4;
781  std::uint32_t n_q1_ : 4;
783  std::uint32_t n_L_ : 4;
785  std::uint32_t n_R_ : 4;
787  std::uint32_t n_ : 4, : 2;
789  bool is_nucleus_ : 1;
791  bool antiparticle_ : 1;
792 #elif defined(BIG_ENDIAN_ARCHITECTURE) // reverse ordering
793  bool antiparticle_ : 1;
794  bool is_nucleus_ : 1, : 2;
795  std::uint32_t n_ : 4;
796  std::uint32_t n_R_ : 4;
797  std::uint32_t n_L_ : 4;
798  std::uint32_t n_q1_ : 4;
799  std::uint32_t n_q2_ : 4;
800  std::uint32_t n_q3_ : 4;
801  std::uint32_t n_J_ : 4;
802 #else
803 #error Endianness macro of the machine not defined.
804 #endif
810  std::uint32_t dump_;
815  struct {
816 #if defined(LITTLE_ENDIAN_ARCHITECTURE) || defined(DOXYGEN)
817  std::uint32_t : 4;
819  std::uint32_t quarks_ : 12;
821  std::uint32_t excitation_ : 12, : 4;
822 #elif defined(BIG_ENDIAN_ARCHITECTURE) // reverse ordering
823  std::uint32_t : 4, excitation_ : 12;
824  std::uint32_t quarks_ : 12, : 4;
825 #else
826 #error Endianness macro of the machine not defined.
827 #endif
830  struct {
831 #if defined(LITTLE_ENDIAN_ARCHITECTURE) || defined(DOXYGEN)
832  std::uint32_t n_Lambda_ : 6;
833  std::uint32_t Z_ : 10;
834  std::uint32_t A_ : 10;
835  std::uint32_t I_ : 4;
836  bool is_nucleus_ : 1;
837  bool antiparticle_ : 1;
838 #elif defined(BIG_ENDIAN_ARCHITECTURE) // reverse ordering
839  bool antiparticle_ : 1;
840  bool is_nucleus_ : 1;
841  std::uint32_t I_ : 4;
842  std::uint32_t A_ : 10;
843  std::uint32_t Z_ : 10;
844  std::uint32_t n_Lambda_ : 6;
845 #else
846 #error Endianness macro of the machine not defined.
847 #endif
849  };
850 
855  inline std::uint32_t ucode() const { return (dump_ & 0x0fffffff); }
856 
863  inline std::uint32_t get_digit_from_char(const char inp) const {
864  // Decimal digit
865  if (48 <= inp && inp <= 57) {
866  return inp - 48;
867  }
868  // Hexdecimal digit, uppercase
869  if (65 <= inp && inp <= 70) {
870  return inp - 65 + 10;
871  }
872  // Hexdecimal digit, lowercase
873  if (97 <= inp && inp <= 102) {
874  return inp - 97 + 10;
875  }
876  throw InvalidPdgCode("PdgCode: Invalid character " + std::string(&inp, 1) +
877  " found.\n");
878  }
879 
901  inline void set_from_string(const std::string& codestring) {
902  dump_ = 0;
903  // Implicit with the above: digits_.antiparticle_ = false;
904  digits_.n_ = digits_.n_R_ = digits_.n_L_ = digits_.n_q1_ = digits_.n_q2_ =
905  digits_.n_q3_ = digits_.n_J_ = digits_.is_nucleus_ = 0;
906  size_t length = codestring.size();
907  if (length < 1) {
908  throw InvalidPdgCode("Empty string does not contain PDG Code\n");
909  }
910  int c = 0;
911  /* Look at current character; if it is a + or minus sign, read it
912  * and advance to next char. */
913  if (codestring[c] == '-') {
914  digits_.antiparticle_ = true;
915  ++c;
916  } else if (codestring[c] == '+') {
917  digits_.antiparticle_ = false;
918  ++c;
919  }
920  // Save if the first character was a sign:
921  unsigned int sign = c;
922 
923  // Nucleus
924  if (length == 10 + sign) {
925  nucleus_.is_nucleus_ = true;
926  if (codestring[c] != '1' || codestring[c + 1] != '0') {
927  throw InvalidPdgCode("Pdg code of nucleus \"" + codestring +
928  "\" should start with 10\n");
929  }
930  c += 2;
931  // ±10LZZZAAAI is the standard for nuclei
932  std::array<int, 8> digits;
933  for (int i = 0; i < 8; i++) {
934  digits[i] = get_digit_from_char(codestring[c + i]);
935  }
936  nucleus_.n_Lambda_ = digits[0];
937  nucleus_.Z_ = 100 * digits[1] + 10 * digits[2] + digits[3];
938  nucleus_.A_ = 100 * digits[4] + 10 * digits[5] + digits[6];
939  nucleus_.I_ = digits[7];
940  return;
941  }
942 
943  // Codestring shouldn't be longer than 8 + sign, except for nuclei
944  if (length > 8 + sign) {
945  throw InvalidPdgCode("String \"" + codestring +
946  "\" too long for PDG Code\n");
947  }
948  /* Please note that in what follows, we actually need c++, not ++c.
949  * first digit is used for n_J if the last digit is not enough. */
950  if (length > 7 + sign) {
951  digits_.n_J_ += get_digit_from_char(codestring[c++]);
952  }
953  // Codestring has 7 digits? 7th from last goes in n_.
954  if (length > 6 + sign) {
955  digits_.n_ = get_digit_from_char(codestring[c++]);
956  }
957  // It has 6 or 7 digits? 6th from last is n_R_.
958  if (length > 5 + sign) {
959  digits_.n_R_ = get_digit_from_char(codestring[c++]);
960  }
961  // 5th from last is n_L_.
962  if (length > 4 + sign) {
963  digits_.n_L_ = get_digit_from_char(codestring[c++]);
964  }
965  // 4th from last is n_q1_.
966  if (length > 3 + sign) {
967  digits_.n_q1_ = get_digit_from_char(codestring[c++]);
968  if (digits_.n_q1_ > 6) {
969  throw InvalidPdgCode("Invalid PDG code " + codestring + " (n_q1>6)");
970  }
971  }
972  // 3rd from last is n_q2_.
973  if (length > 2 + sign) {
974  digits_.n_q2_ = get_digit_from_char(codestring[c++]);
975  if (digits_.n_q2_ > 6) {
976  throw InvalidPdgCode("Invalid PDG code " + codestring + " (n_q2>6)");
977  }
978  }
979  // Next to last is n_q3_.
980  if (length > 1 + sign) {
981  digits_.n_q3_ = get_digit_from_char(codestring[c++]);
982  if (digits_.n_q3_ > 6) {
983  throw InvalidPdgCode("Invalid PDG code " + codestring + " (n_q3>6)");
984  }
985  }
986  // Last digit is the spin degeneracy.
987  if (length > sign) {
988  digits_.n_J_ += get_digit_from_char(codestring[c++]);
989  } else {
990  throw InvalidPdgCode(
991  "String \"" + codestring +
992  "\" only consists of a sign, that is no valid PDG Code\n");
993  }
994  check();
995  }
996 
1006  inline void set_fields(std::uint32_t abscode) {
1007  /* "dump_ =" overwrites antiparticle_, but this needs to have been set
1008  * already, so we carry it around the assignment. */
1009  bool ap = digits_.antiparticle_;
1010  dump_ = abscode & 0x0fffffff;
1011  digits_.antiparticle_ = ap;
1012  int test = test_code();
1013  if (test > 0) {
1014  throw InvalidPdgCode("Invalid digits " + std::to_string(test) +
1015  " in PDG Code " + string());
1016  }
1017  check();
1018  }
1019 };
1020 
1021 static_assert(sizeof(PdgCode) == 4, "should fit into 32 bit integer");
1022 
1029 std::istream& operator>>(std::istream& is, PdgCode& code);
1035 std::ostream& operator<<(std::ostream& is, const PdgCode& code);
1036 
1038 inline bool is_dilepton(const PdgCode pdg1, const PdgCode pdg2) {
1039  const auto c1 = pdg1.code();
1040  const auto c2 = pdg2.code();
1041  const auto min = std::min(c1, c2);
1042  const auto max = std::max(c1, c2);
1043  return (max == 0x11 && min == -0x11) || (max == 0x13 && min == -0x13);
1044 }
1045 
1050 inline bool has_lepton_pair(const PdgCode pdg1, const PdgCode pdg2,
1051  const PdgCode pdg3) {
1052  return is_dilepton(pdg1, pdg2) || is_dilepton(pdg1, pdg3) ||
1053  is_dilepton(pdg2, pdg3);
1054 }
1055 
1059 namespace pdg {
1061 const PdgCode d(PdgCode::from_decimal(1000010020));
1071 const PdgCode he3(PdgCode::from_decimal(1000020030));
1078 
1079 } // namespace pdg
1080 
1081 } // namespace smash
1082 
1083 #endif // SRC_INCLUDE_SMASH_PDGCODE_H_
PdgCode stores a Particle Data Group Particle Numbering Scheme particle type number.
Definition: pdgcode.h:108
std::uint32_t quarks_
The quark digits n_q{1,2,3}_.
Definition: pdgcode.h:819
bool is_Nstar1535() const
Definition: pdgcode.h:342
std::int32_t code() const
Definition: pdgcode.h:249
bool is_rho() const
Definition: pdgcode.h:390
int antiparticle_sign() const
Definition: pdgcode.h:572
bool operator!=(const PdgCode rhs) const
Definition: pdgcode.h:664
int nucleus_n() const
Number of neutrons in nucleus.
Definition: pdgcode.h:735
int baryon_number() const
Definition: pdgcode.h:308
bool is_meson() const
Definition: pdgcode.h:321
PdgCode(const std::uint32_t abscode)
receive an unsigned integer and process it into a PDG Code.
Definition: pdgcode.h:155
std::uint32_t n_q1_
first quark field. 0 for mesons.
Definition: pdgcode.h:781
bool is_Sigma() const
Definition: pdgcode.h:373
PdgCode(const std::string &codestring)
Initialize using a string The string is interpreted as a hexadecimal number, i.e.,...
Definition: pdgcode.h:131
int net_quark_number(const int quark) const
Returns the net number of quarks with given flavour number For public use, see strangeness(),...
Definition: pdgcode.cc:31
int nucleus_ap() const
Number of antiprotons in nucleus.
Definition: pdgcode.h:745
int nucleus_an() const
Number of antineutrons in nucleus.
Definition: pdgcode.h:749
friend std::istream & operator>>(std::istream &is, PdgCode &code)
istream >> PdgCode assigns the PDG Code from an istream.
Definition: pdgcode.cc:14
int bottomness() const
Definition: pdgcode.h:478
std::array< int, 3 > quark_content() const
The return is always an array of three numbers, which are pdgcodes of quarks: 1 - d,...
Definition: pdgcode.h:592
std::uint32_t n_q3_
third quark field
Definition: pdgcode.h:777
unsigned int spin() const
Definition: pdgcode.h:529
bool is_pion() const
Definition: pdgcode.h:384
bool is_kaon() const
Definition: pdgcode.h:378
bool antiparticle_
first bit: stores the sign.
Definition: pdgcode.h:791
std::uint32_t dump() const
Dumps the bitfield into an unsigned integer.
Definition: pdgcode.h:243
std::uint32_t Z_
Definition: pdgcode.h:833
std::uint32_t n_
first field: "counter"
Definition: pdgcode.h:787
bool is_lepton() const
Definition: pdgcode.h:302
bool is_hyperon() const
Definition: pdgcode.h:355
std::uint32_t n_J_
spin quantum number .
Definition: pdgcode.h:775
bool is_nucleus() const
Definition: pdgcode.h:291
bool is_deuteron() const
Definition: pdgcode.h:396
bool is_proton() const
Definition: pdgcode.h:330
static PdgCode from_decimal(const int pdgcode_decimal)
Construct PDG code from decimal number.
Definition: pdgcode.h:269
bool is_nucleus_
Definition: pdgcode.h:836
bool is_antiparticle_of(const PdgCode rhs) const
Definition: pdgcode.h:667
std::uint32_t bool is_nucleus_
1 for nuclei, 0 for the rest
Definition: pdgcode.h:787
int charmness() const
Definition: pdgcode.h:471
int strangeness() const
Definition: pdgcode.h:464
struct smash::PdgCode::@0::@2 digits_
The single digits collection of the code.
int nucleus_p() const
Number of protons in nucleus.
Definition: pdgcode.h:731
int32_t get_decimal() const
Definition: pdgcode.h:690
bool is_baryon() const
Definition: pdgcode.h:318
void set_fields(std::uint32_t abscode)
Sets the bitfield from an unsigned integer.
Definition: pdgcode.h:1006
int isospin3() const
Definition: pdgcode.h:427
void check() const
Do all sorts of validity checks.
Definition: pdgcode.h:214
bool is_nucleon() const
Definition: pdgcode.h:324
std::string string() const
Definition: pdgcode.h:252
int test_code() const
Checks the integer for invalid hex digits.
Definition: pdgcode.h:180
PdgCode get_antiparticle() const
Construct the antiparticle to a given PDG code.
Definition: pdgcode.h:259
int nucleus_La() const
Number of Lambdas in nucleus.
Definition: pdgcode.h:741
bool is_hadron() const
Definition: pdgcode.h:297
bool operator==(const PdgCode rhs) const
Definition: pdgcode.h:661
PdgCode()
Standard initializer.
Definition: pdgcode.h:125
int nucleus_A() const
Nucleus mass number.
Definition: pdgcode.h:759
static PdgCode invalid()
PdgCode 0x0 is guaranteed not to be valid by the PDG standard, but it passes all tests here,...
Definition: pdgcode.h:679
int nucleus_aLa() const
Number of anti-Lambdas in nucleus.
Definition: pdgcode.h:755
bool contains_enough_valence_quarks(int valence_quarks_required) const
Definition: pdgcode.cc:92
bool operator<(const PdgCode rhs) const
Sorts PDG Codes according to their numeric value.
Definition: pdgcode.h:652
PdgCode(std::int32_t codenumber)
Receive a signed integer and process it into a PDG Code.
Definition: pdgcode.h:143
bool is_Omega() const
Definition: pdgcode.h:358
std::uint32_t A_
Definition: pdgcode.h:834
bool is_neutron() const
Definition: pdgcode.h:336
void deexcite()
Remove all excitation, except spin. Sign and quark content remains.
Definition: pdgcode.h:710
std::uint32_t n_L_
"angular momentum"
Definition: pdgcode.h:783
void set_from_string(const std::string &codestring)
Set the PDG code from the given string.
Definition: pdgcode.h:901
bool is_triton() const
Definition: pdgcode.h:402
std::uint32_t ucode() const
Definition: pdgcode.h:855
std::uint32_t n_R_
"radial excitation"
Definition: pdgcode.h:785
std::uint32_t dump_
The bitfield dumped into a single integer.
Definition: pdgcode.h:810
struct smash::PdgCode::@0::@4 nucleus_
Structure for the nuclei.
struct smash::PdgCode::@0::@3 chunks_
Chunk collection: here, the chunks with and are directly accessible.
bool is_Delta() const
Definition: pdgcode.h:348
std::uint32_t n_Lambda_
Definition: pdgcode.h:832
double frac_strange() const
Definition: pdgcode.h:440
bool has_antiparticle() const
Definition: pdgcode.h:411
bool is_Lambda() const
Definition: pdgcode.h:368
std::uint32_t get_digit_from_char(const char inp) const
Definition: pdgcode.h:863
bool is_Xi() const
Definition: pdgcode.h:363
std::int32_t quarks() const
Definition: pdgcode.h:576
std::uint32_t I_
Definition: pdgcode.h:835
std::uint32_t n_q2_
second quark field
Definition: pdgcode.h:779
unsigned int spin_degeneracy() const
Definition: pdgcode.h:565
std::uint32_t excitation_
The excitation digits n_, n_R_, n_L_.
Definition: pdgcode.h:821
int charge() const
The charge of the particle.
Definition: pdgcode.h:488
std::ostream & operator<<(std::ostream &out, const ActionPtr &action)
Convenience: dereferences the ActionPtr to Action.
Definition: action.h:532
constexpr int pi_p
π⁺.
constexpr int rho_p
ρ⁺.
constexpr int Delta_p
Δ⁺.
constexpr int rho_m
ρ⁻.
constexpr int Delta_pp
Δ⁺⁺.
constexpr int K_p
K⁺.
constexpr int K_z
K⁰.
const PdgCode triton(PdgCode::from_decimal(1000010030))
Triton.
constexpr int p
Proton.
constexpr int N1535_z
N(1535)⁰.
const PdgCode d(PdgCode::from_decimal(1000010020))
Deuteron.
const PdgCode he3(PdgCode::from_decimal(1000020030))
He-3.
const PdgCode antihe3(PdgCode::from_decimal(-1000020030))
Anti-He-3.
constexpr int pi_z
π⁰.
constexpr int n
Neutron.
constexpr int Delta_m
Δ⁻.
constexpr int Delta_z
Δ⁰.
constexpr int rho_z
ρ⁰.
const PdgCode antitriton(PdgCode::from_decimal(-1000010030))
Anti-triton.
constexpr int pi_m
π⁻.
constexpr int N1535_p
N(1535)⁺.
const PdgCode hypertriton(PdgCode::from_decimal(1010010030))
Hypertriton.
const PdgCode antid(PdgCode::from_decimal(-1000010020))
Anti-deuteron in decimal digits.
const PdgCode antihypertriton(PdgCode::from_decimal(-1010010030))
Anti-Hypertriton.
const PdgCode dprime(PdgCode::from_decimal(1000010021))
Deuteron-prime resonance.
Definition: action.h:24
bool is_dilepton(const PdgCode pdg1, const PdgCode pdg2)
Definition: pdgcode.h:1038
std::istream & operator>>(std::istream &is, PdgCode &code)
Sets the PDG code from the textual representation in the input stream.
Definition: pdgcode.cc:14
bool has_lepton_pair(const PdgCode pdg1, const PdgCode pdg2, const PdgCode pdg3)
Definition: pdgcode.h:1050
thrown for invalid inputs
Definition: pdgcode.h:114