Version: SMASH-3.3
pdgcode.h
Go to the documentation of this file.
1 /*
2  *
3  * Copyright (c) 2014-2025
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 <cstdint>
17 #include <cstdlib>
18 #include <functional>
19 #include <iomanip>
20 #include <iosfwd>
21 #include <iostream>
22 #include <sstream>
23 #include <stdexcept>
24 #include <string>
25 #include <type_traits>
26 
27 #include "pdgcode_constants.h"
28 
29 namespace smash {
30 
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 
216  template <typename T>
217  PdgCode(T codenumber, // NOLINT(runtime/explicit)
218  typename std::enable_if_t<std::is_integral_v<T> && 4 < sizeof(T),
219  bool> = true)
220  : PdgCode{std::invoke([&codenumber]() {
221  std::stringstream stream;
222  char sign = '+';
223  if (codenumber < 0) {
224  sign = '-';
225  codenumber = -codenumber;
226  }
227  stream << sign << std::hex << codenumber;
228  return stream.str();
229  })} {}
230 
231  /****************************************************************************
232  * *
233  * test function and export functions *
234  * *
235  ****************************************************************************/
236 
250  inline int test_code() const {
251  int fail = 0;
252  if (digits_.n_ > 9) {
253  fail |= 1 << 6;
254  }
255  if (digits_.n_R_ > 9) {
256  fail |= 1 << 5;
257  }
258  if (digits_.n_L_ > 9) {
259  fail |= 1 << 4;
260  }
261  if (digits_.n_q1_ > 6) {
262  fail |= 1 << 3;
263  }
264  if (digits_.n_q2_ > 6) {
265  fail |= 1 << 2;
266  }
267  if (digits_.n_q3_ > 6) {
268  fail |= 1 << 1;
269  }
270  if (digits_.n_J_ > 15) {
271  fail |= 1;
272  }
273  return fail;
274  }
275 
284  void check() const {
285  // n_J must be odd for mesons and even for baryons (and cannot be zero)
286  if (is_hadron()) {
287  if (baryon_number() == 0) {
288  // mesons: special cases K0_L=0x130 and K0_S=0x310
289  if ((digits_.n_J_ % 2 == 0) && dump() != 0x130 && dump() != 0x310) {
290  throw InvalidPdgCode("Invalid PDG code " + string() +
291  " (meson with even n_J)");
292  }
293  } else {
294  if ((digits_.n_J_ % 2 != 0) || digits_.n_J_ == 0) {
295  throw InvalidPdgCode("Invalid PDG code " + string() +
296  " (baryon with odd n_J)");
297  }
298  }
299  } else {
300  if (digits_.n_J_ == 0 && dump() != 0x0) {
301  throw InvalidPdgCode("Invalid PDG code " + string() + " (n_J==0)");
302  }
303  }
304  /* The antiparticle flag only makes sense for particle types
305  * that have an antiparticle. */
306  if (digits_.antiparticle_ && !has_antiparticle()) {
307  throw InvalidPdgCode("Invalid PDG code " + string() +
308  " (cannot be negative)");
309  }
310  }
311 
313  inline std::uint32_t dump() const {
314  // this cuts the three unused bits.
315  return (dump_ & 0x8fffffff);
316  }
317 
319  inline std::int32_t code() const { return antiparticle_sign() * ucode(); }
320 
322  inline std::string string() const {
323  std::stringstream ss;
324  ss << get_decimal();
325  return ss.str();
326  }
327 
330  PdgCode result = *this;
331  result.digits_.antiparticle_ = !digits_.antiparticle_;
332  return result;
333  }
334 
339  static PdgCode from_decimal(const int pdgcode_decimal) {
340  // Nucleus and special codes with 2J+1 > 9
341  if (std::abs(pdgcode_decimal) > 1E7) {
342  return PdgCode(std::to_string(pdgcode_decimal));
343  }
344  int a = pdgcode_decimal;
345  int hex_pdg = 0, tmp = 1;
346  while (a) {
347  hex_pdg += (a % 10) * tmp;
348  tmp *= 16;
349  a = a / 10;
350  }
351  return PdgCode(hex_pdg);
352  }
353 
354  /****************************************************************************
355  * *
356  * accessors of various properties *
357  * *
358  ****************************************************************************/
359 
361  inline bool is_nucleus() const {
362  assert(digits_.is_nucleus_ == nucleus_.is_nucleus_);
363  return nucleus_.is_nucleus_;
364  }
365 
367  inline bool is_hadron() const {
368  return (digits_.n_q3_ != 0 && digits_.n_q2_ != 0 && !is_nucleus());
369  }
370 
372  inline bool is_lepton() const {
373  return (digits_.n_q1_ == 0 && digits_.n_q2_ == 0 && digits_.n_q3_ == 1 &&
374  !is_nucleus());
375  }
376 
378  inline bool is_neutrino() const {
379  return (is_lepton() && digits_.n_J_ % 2 == 0);
380  }
381 
383  inline bool is_charmonia() const {
384  return is_meson() && digits_.n_q2_ == 4 && digits_.n_q3_ == 4;
385  }
386 
388  inline int baryon_number() const {
389  if (is_nucleus()) {
390  return static_cast<int>(nucleus_.A_) * antiparticle_sign();
391  }
392  if (!is_hadron() || digits_.n_q1_ == 0) {
393  return 0;
394  }
395  return antiparticle_sign();
396  }
398  inline bool is_baryon() const { return is_hadron() && digits_.n_q1_ != 0; }
399 
401  inline bool is_meson() const { return is_hadron() && digits_.n_q1_ == 0; }
402 
404  inline bool is_nucleon() const {
405  const auto abs_code = std::abs(code());
406  return (abs_code == pdg::p || abs_code == pdg::n);
407  }
408 
410  inline bool is_proton() const {
411  const auto abs_code = std::abs(code());
412  return (abs_code == pdg::p);
413  }
414 
416  inline bool is_neutron() const {
417  const auto abs_code = std::abs(code());
418  return (abs_code == pdg::n);
419  }
420 
422  inline bool is_Nstar1535() const {
423  const auto abs_code = std::abs(code());
424  return (abs_code == pdg::N1535_p || abs_code == pdg::N1535_z);
425  }
426 
428  inline bool is_Delta() const {
429  const auto abs_code = std::abs(code());
430  return (abs_code == pdg::Delta_pp || abs_code == pdg::Delta_p ||
431  abs_code == pdg::Delta_z || abs_code == pdg::Delta_m);
432  }
433 
435  inline bool is_hyperon() const { return is_hadron() && digits_.n_q1_ == 3; }
436 
438  inline bool is_Omega() const {
439  return is_hyperon() && digits_.n_q2_ == 3 && digits_.n_q3_ == 3;
440  }
441 
443  inline bool is_Xi() const {
444  return is_hyperon() && digits_.n_q2_ == 3 && digits_.n_q3_ != 3;
445  }
446 
448  inline bool is_Lambda() const {
449  return is_hyperon() && digits_.n_q2_ == 1 && digits_.n_q3_ == 2;
450  }
451 
453  inline bool is_Sigma() const {
454  return is_hyperon() && digits_.n_q2_ != 3 && !is_Lambda();
455  }
456 
458  inline bool is_Sigmastar() const {
459  const auto abs_code = std::abs(code());
460  return (abs_code == pdg::Sigma_star_m || abs_code == pdg::Sigma_star_z ||
461  abs_code == pdg::Sigma_star_p);
462  }
463 
465  inline bool is_kaon() const {
466  const auto abs_code = std::abs(code());
467  return (abs_code == pdg::K_p) || (abs_code == pdg::K_z);
468  }
469 
471  inline bool is_pion() const {
472  const auto c = code();
473  return (c == pdg::pi_z) || (c == pdg::pi_p) || (c == pdg::pi_m);
474  }
475 
477  inline bool is_omega() const {
478  const auto c = code();
479  return c == pdg::omega;
480  }
481 
483  inline bool is_rho() const {
484  const auto c = code();
485  return (c == pdg::rho_z) || (c == pdg::rho_p) || (c == pdg::rho_m);
486  }
487 
489  inline bool is_deuteron() const {
490  return is_nucleus() && nucleus_.A_ == 2 && nucleus_.Z_ == 1 &&
491  nucleus_.n_Lambda_ == 0 && nucleus_.I_ == 0;
492  }
493 
495  inline bool is_triton() const {
496  return is_nucleus() && nucleus_.A_ == 3 && nucleus_.Z_ == 1 &&
497  nucleus_.n_Lambda_ == 0 && nucleus_.I_ == 0;
498  }
499 
504  bool has_antiparticle() const {
505  if (is_nucleus()) {
506  return true;
507  }
508  if (is_hadron()) {
509  return (baryon_number() != 0) || (digits_.n_q2_ != digits_.n_q3_);
510  } else {
511  return digits_.n_q3_ == 1; // leptons!
512  }
513  }
514 
520  inline int isospin3() const {
521  /* net_quark_number(2) is the number of u quarks,
522  * net_quark_number(1) is the number of d quarks. */
523  return net_quark_number(2) - net_quark_number(1);
524  }
525 
533  inline double frac_strange() const {
534  if (is_baryon()) {
535  return std::abs(strangeness()) / 3.;
536  } else if (is_meson()) {
537  /* The quarkonium state has 0 net strangeness
538  * but there are actually 2 strange quarks out of 2 total */
539  if (digits_.n_q3_ == 3 && digits_.n_q2_ == 3) {
540  return 1.;
541  } else {
542  return std::abs(strangeness()) / 2.;
543  }
544  } else {
545  /* If not baryon or meson, this should be 0, as AQM does not
546  * extend to non-hadrons */
547  return 0.;
548  }
549  }
550 
558  inline double frac_charm() const {
559  if (is_baryon()) {
560  return std::abs(charmness()) / 3.;
561  } else if (is_meson()) {
562  /* The charmonium state has 0 net charmness
563  * but there are actually 2 charm quarks out of 2 total */
564  if (digits_.n_q3_ == 4 && digits_.n_q2_ == 4) {
565  return 1.;
566  } else {
567  return std::abs(charmness()) / 2.;
568  }
569  } else {
570  /* If not baryon or meson, this should be 0, as AQM does not
571  * extend to non-hadrons */
572  return 0.;
573  }
574  }
575 
583  inline double frac_bottom() const {
584  if (is_baryon()) {
585  return std::abs(bottomness()) / 3.;
586  } else if (is_meson()) {
587  /* The bottomonium state has 0 net bottomness
588  * but there are actually 2 bottom quarks out of 2 total */
589  if (digits_.n_q3_ == 5 && digits_.n_q2_ == 5) {
590  return 1.;
591  } else {
592  return std::abs(bottomness()) / 2.;
593  }
594  } else {
595  /* If not baryon or meson, this should be 0, as AQM does not
596  * extend to non-hadrons */
597  return 0.;
598  }
599  }
600 
602  inline bool is_heavy_flavor() const {
603  return (frac_charm() != 0) || (frac_bottom() != 0);
604  }
605 
611  inline int strangeness() const { return -net_quark_number(3); }
612 
618  inline int charmness() const { return +net_quark_number(4); }
619 
625  inline int bottomness() const { return -net_quark_number(5); }
626 
635  int charge() const {
636  if (is_hadron() || is_nucleus()) {
637  // Q will accumulate 3*charge (please excuse the upper case. I
638  // want to distinguish this from q which might be interpreted as
639  // shorthand for "quark".)
640  int Q = 0;
641  /* This loops over d,u,s,c,b,t quarks (the latter can be safely ignored,
642  * but I don't think this will be a bottle neck. */
643  for (int i = 1; i < 7; i++) {
644  /* u,c,t quarks have charge = 2/3 e, while d,s,b quarks have -1/3 e.
645  * The antiparticle sign is already in net_quark_number. */
646  Q += (i % 2 == 0 ? 2 : -1) * net_quark_number(i);
647  }
648  return Q / 3;
649  }
650  /* non-hadron:
651  * Leptons: 11, 13, 15 are e, μ, τ and have a charge -1, while
652  * 12, 14, 16 are the neutrinos that have no charge. */
653  if (digits_.n_q3_ == 1) {
654  return -1 * (digits_.n_J_ % 2) * antiparticle_sign();
655  }
656  /* Bosons: 24 is the W+, all else is uncharged.
657  * we ignore the first digits so that this also finds strange gauge
658  * boson "resonances" (in particular, \f$\tilde \chi_1^+\f$ with PDG
659  * Code 1000024). */
660  if ((dump_ & 0x0000ffff) == 0x24) {
661  return antiparticle_sign();
662  }
663  // default (this includes all other Bosons) is 0.
664  return 0;
665  }
666 
676  inline unsigned int spin() const {
677  if (is_nucleus()) {
678  // Generally spin of a nucleus cannot be simply guessed, it should be
679  // provided from some table. However, here we only care about a
680  // limited set of light nuclei with A <= 4.
681  if (nucleus_.A_ == 2) {
682  // Deuteron spin is 1
683  return 2;
684  } else if (nucleus_.A_ == 3) {
685  // Tritium and He-3 spin are 1/2
686  // Hypertriton spin is not firmly determined, but decay branching ratios
687  // indicate spin 1/2
688  return 1;
689  } else if (nucleus_.A_ == 4) {
690  // He-4 spin is 0
691  return 0;
692  }
693  throw std::runtime_error("Unknown spin of nucleus.");
694  // Alternative possibility is to guess 1/2 for fermions and 0 for bosons
695  // as 2 * (nucleus_.A_ % 2).
696  }
697 
698  if (is_hadron()) {
699  if (digits_.n_J_ == 0) {
700  return 0; // special cases: K0_L=0x130 & K0_S=0x310
701  } else {
702  return digits_.n_J_ - 1;
703  }
704  }
705  /* this assumes that we only have white particles (no single
706  * quarks): Electroweak fermions have 11-17, so the
707  * second-to-last-digit is the spin. The same for the Bosons: they
708  * have 21-29 and 2spin = 2 (this fails for the Higgs). */
709  return digits_.n_q3_;
710  }
712  inline unsigned int spin_degeneracy() const {
713  if (is_hadron() && digits_.n_J_ > 0) {
714  return digits_.n_J_;
715  }
716  return spin() + 1;
717  }
719  inline int antiparticle_sign() const {
720  return (digits_.antiparticle_ ? -1 : +1);
721  }
723  inline std::int32_t quarks() const {
724  if (!is_hadron() || is_nucleus()) {
725  return 0;
726  }
727  return chunks_.quarks_;
728  }
729 
739  std::array<int, 3> quark_content() const {
740  std::array<int, 3> result = {static_cast<int>(digits_.n_q1_),
741  static_cast<int>(digits_.n_q2_),
742  static_cast<int>(digits_.n_q3_)};
743  if (is_hadron()) {
744  // Antibaryons
745  if (digits_.n_q1_ != 0 && digits_.antiparticle_) {
746  for (size_t i = 0; i < 3; i++) {
747  result[i] = -result[i];
748  }
749  }
750  // Mesons
751  if (digits_.n_q1_ == 0) {
752  // Own antiparticle
753  if (digits_.n_q2_ == digits_.n_q3_) {
754  result[2] = -result[2];
755  } else {
756  // Like pi-
757  if (digits_.antiparticle_) {
758  result[1] = -result[1];
759  // Like pi+
760  } else {
761  result[2] = -result[2];
762  }
763  }
764  // add extra minus sign according to the pdg convention
765  if (digits_.n_q2_ != digits_.n_q3_ && digits_.n_q2_ % 2 == 1) {
766  for (int i = 1; i <= 2; i++) {
767  result[i] = -result[i];
768  }
769  }
770  }
771  } else {
772  result = {0, 0, 0};
773  }
774  return result;
775  }
776 
787  bool contains_enough_valence_quarks(int valence_quarks_required) const;
788 
789  /****************************************************************************
790  * *
791  * operations with more than one PDG Code *
792  * *
793  ****************************************************************************/
794 
799  inline bool operator<(const PdgCode rhs) const {
800  return dump_ < rhs.dump_;
801  /* the complex thing to do here is to calculate:
802  * code() < rhs.code()
803  * but for getting a total order that's overkill. The uint32_t value in
804  * dump_ works just fine. */
805  }
806 
808  inline bool operator==(const PdgCode rhs) const { return dump_ == rhs.dump_; }
809 
811  inline bool operator!=(const PdgCode rhs) const { return !(*this == rhs); }
812 
814  inline bool is_antiparticle_of(const PdgCode rhs) const {
815  return code() == -rhs.code();
816  }
817 
819  friend std::istream& operator>>(std::istream& is, PdgCode& code);
820 
826  static PdgCode invalid() { return PdgCode(0x0); }
827 
837  int32_t get_decimal() const {
838  if (is_nucleus()) {
839  // ±10LZZZAAAI
840  return antiparticle_sign() *
841  (nucleus_.I_ + 10 * nucleus_.A_ + 10000 * nucleus_.Z_ +
842  10000000 * nucleus_.n_Lambda_ + 1000000000);
843  }
844  int n_J_1 = 0;
845  int n_J_2 = digits_.n_J_;
846  if (n_J_2 > 9) {
847  n_J_1 = n_J_2 - 9;
848  n_J_2 = 9;
849  }
850  return antiparticle_sign() *
851  (n_J_2 + digits_.n_q3_ * 10 + digits_.n_q2_ * 100 +
852  digits_.n_q1_ * 1000 + digits_.n_L_ * 10000 +
853  digits_.n_R_ * 100000 + digits_.n_ * 1000000 + n_J_1 * 10000000);
854  }
855 
857  void deexcite() {
858  if (!is_nucleus()) {
859  chunks_.excitation_ = 0;
860  } else {
861  nucleus_.I_ = 0;
862  }
863  }
864 
875  int net_quark_number(const int quark) const;
876 
878  int nucleus_p() const {
879  return (is_nucleus() && !nucleus_.antiparticle_) ? nucleus_.Z_ : 0;
880  }
882  int nucleus_n() const {
883  return (is_nucleus() && !nucleus_.antiparticle_)
884  ? nucleus_.A_ - nucleus_.Z_ - nucleus_.n_Lambda_
885  : 0;
886  }
888  int nucleus_La() const {
889  return (is_nucleus() && !nucleus_.antiparticle_) ? nucleus_.n_Lambda_ : 0;
890  }
892  int nucleus_ap() const {
893  return (is_nucleus() && nucleus_.antiparticle_) ? nucleus_.Z_ : 0;
894  }
896  int nucleus_an() const {
897  return (is_nucleus() && nucleus_.antiparticle_)
898  ? nucleus_.A_ - nucleus_.Z_ - nucleus_.n_Lambda_
899  : 0;
900  }
902  int nucleus_aLa() const {
903  return (is_nucleus() && nucleus_.antiparticle_) ? nucleus_.n_Lambda_ : 0;
904  }
906  int nucleus_A() const { return is_nucleus() ? nucleus_.A_ : 0; }
907 
908  private:
914  union {
919  struct {
920 #if defined(LITTLE_ENDIAN_ARCHITECTURE) || defined(DOXYGEN)
922  std::uint32_t n_J_ : 4;
924  std::uint32_t n_q3_ : 4;
926  std::uint32_t n_q2_ : 4;
928  std::uint32_t n_q1_ : 4;
930  std::uint32_t n_L_ : 4;
932  std::uint32_t n_R_ : 4;
934  std::uint32_t n_ : 4, : 2;
936  bool is_nucleus_ : 1;
938  bool antiparticle_ : 1;
939 #elif defined(BIG_ENDIAN_ARCHITECTURE) // reverse ordering
940  bool antiparticle_ : 1;
941  bool is_nucleus_ : 1, : 2;
942  std::uint32_t n_ : 4;
943  std::uint32_t n_R_ : 4;
944  std::uint32_t n_L_ : 4;
945  std::uint32_t n_q1_ : 4;
946  std::uint32_t n_q2_ : 4;
947  std::uint32_t n_q3_ : 4;
948  std::uint32_t n_J_ : 4;
949 #else
950 #error Endianness macro of the machine not defined.
951 #endif
952  } digits_;
957  std::uint32_t dump_;
962  struct {
963 #if defined(LITTLE_ENDIAN_ARCHITECTURE) || defined(DOXYGEN)
964  std::uint32_t : 4;
966  std::uint32_t quarks_ : 12;
968  std::uint32_t excitation_ : 12, : 4;
969 #elif defined(BIG_ENDIAN_ARCHITECTURE) // reverse ordering
970  std::uint32_t : 4, excitation_ : 12;
971  std::uint32_t quarks_ : 12, : 4;
972 #else
973 #error Endianness macro of the machine not defined.
974 #endif
975  } chunks_;
977  struct {
978 #if defined(LITTLE_ENDIAN_ARCHITECTURE) || defined(DOXYGEN)
979  std::uint32_t n_Lambda_ : 6;
980  std::uint32_t Z_ : 10;
981  std::uint32_t A_ : 10;
982  std::uint32_t I_ : 4;
983  bool is_nucleus_ : 1;
984  bool antiparticle_ : 1;
985 #elif defined(BIG_ENDIAN_ARCHITECTURE) // reverse ordering
986  bool antiparticle_ : 1;
987  bool is_nucleus_ : 1;
988  std::uint32_t I_ : 4;
989  std::uint32_t A_ : 10;
990  std::uint32_t Z_ : 10;
991  std::uint32_t n_Lambda_ : 6;
992 #else
993 #error Endianness macro of the machine not defined.
994 #endif
995  } nucleus_;
996  };
997 
1002  inline std::uint32_t ucode() const { return (dump_ & 0x0fffffff); }
1003 
1010  inline std::uint32_t get_digit_from_char(const char inp) const {
1011  // Decimal digit
1012  if (48 <= inp && inp <= 57) {
1013  return inp - 48;
1014  }
1015  // Hexdecimal digit, uppercase
1016  if (65 <= inp && inp <= 70) {
1017  return inp - 65 + 10;
1018  }
1019  // Hexdecimal digit, lowercase
1020  if (97 <= inp && inp <= 102) {
1021  return inp - 97 + 10;
1022  }
1023  throw InvalidPdgCode("PdgCode: Invalid character " + std::string(&inp, 1) +
1024  " found.\n");
1025  }
1026 
1048  inline void set_from_string(const std::string& codestring) {
1049  dump_ = 0;
1050  // Implicit with the above: digits_.antiparticle_ = false;
1051  digits_.n_ = digits_.n_R_ = digits_.n_L_ = digits_.n_q1_ = digits_.n_q2_ =
1052  digits_.n_q3_ = digits_.n_J_ = digits_.is_nucleus_ = 0;
1053  size_t length = codestring.size();
1054  if (length < 1) {
1055  throw InvalidPdgCode("Empty string does not contain PDG Code\n");
1056  }
1057  int c = 0;
1058  /* Look at current character; if it is a + or minus sign, read it
1059  * and advance to next char. */
1060  if (codestring[c] == '-') {
1061  digits_.antiparticle_ = true;
1062  ++c;
1063  } else if (codestring[c] == '+') {
1064  digits_.antiparticle_ = false;
1065  ++c;
1066  }
1067  // Save if the first character was a sign:
1068  unsigned int sign = c;
1069 
1070  // Nucleus
1071  if (length == 10 + sign) {
1072  nucleus_.is_nucleus_ = true;
1073  if (codestring[c] != '1' || codestring[c + 1] != '0') {
1074  throw InvalidPdgCode("Pdg code of nucleus \"" + codestring +
1075  "\" should start with 10\n");
1076  }
1077  c += 2;
1078  // ±10LZZZAAAI is the standard for nuclei
1079  std::array<int, 8> digits;
1080  for (int i = 0; i < 8; i++) {
1081  digits[i] = get_digit_from_char(codestring[c + i]);
1082  }
1083  nucleus_.n_Lambda_ = digits[0];
1084  nucleus_.Z_ = 100 * digits[1] + 10 * digits[2] + digits[3];
1085  nucleus_.A_ = 100 * digits[4] + 10 * digits[5] + digits[6];
1086  nucleus_.I_ = digits[7];
1087  return;
1088  }
1089 
1090  // Codestring shouldn't be longer than 8 + sign, except for nuclei
1091  if (length > 8 + sign) {
1092  throw InvalidPdgCode("String \"" + codestring +
1093  "\" too long for PDG Code\n");
1094  }
1095  /* Please note that in what follows, we actually need c++, not ++c.
1096  * first digit is used for n_J if the last digit is not enough. */
1097  if (length > 7 + sign) {
1098  digits_.n_J_ += get_digit_from_char(codestring[c++]);
1099  }
1100  // Codestring has 7 digits? 7th from last goes in n_.
1101  if (length > 6 + sign) {
1102  digits_.n_ = get_digit_from_char(codestring[c++]);
1103  }
1104  // It has 6 or 7 digits? 6th from last is n_R_.
1105  if (length > 5 + sign) {
1106  digits_.n_R_ = get_digit_from_char(codestring[c++]);
1107  }
1108  // 5th from last is n_L_.
1109  if (length > 4 + sign) {
1110  digits_.n_L_ = get_digit_from_char(codestring[c++]);
1111  }
1112  // 4th from last is n_q1_.
1113  if (length > 3 + sign) {
1114  digits_.n_q1_ = get_digit_from_char(codestring[c++]);
1115  if (digits_.n_q1_ > 6) {
1116  throw InvalidPdgCode("Invalid PDG code " + codestring + " (n_q1>6)");
1117  }
1118  }
1119  // 3rd from last is n_q2_.
1120  if (length > 2 + sign) {
1121  digits_.n_q2_ = get_digit_from_char(codestring[c++]);
1122  if (digits_.n_q2_ > 6) {
1123  throw InvalidPdgCode("Invalid PDG code " + codestring + " (n_q2>6)");
1124  }
1125  }
1126  // Next to last is n_q3_.
1127  if (length > 1 + sign) {
1128  digits_.n_q3_ = get_digit_from_char(codestring[c++]);
1129  if (digits_.n_q3_ > 6) {
1130  throw InvalidPdgCode("Invalid PDG code " + codestring + " (n_q3>6)");
1131  }
1132  }
1133  // Last digit is the spin degeneracy.
1134  if (length > sign) {
1135  digits_.n_J_ += get_digit_from_char(codestring[c++]);
1136  } else {
1137  throw InvalidPdgCode(
1138  "String \"" + codestring +
1139  "\" only consists of a sign, that is no valid PDG Code\n");
1140  }
1141  check();
1142  }
1143 
1153  inline void set_fields(std::uint32_t abscode) {
1154  /* "dump_ =" overwrites antiparticle_, but this needs to have been set
1155  * already, so we carry it around the assignment. */
1156  bool ap = digits_.antiparticle_;
1157  dump_ = abscode & 0x0fffffff;
1158  digits_.antiparticle_ = ap;
1159  int test = test_code();
1160  if (test > 0) {
1161  throw InvalidPdgCode("Invalid digits " + std::to_string(test) +
1162  " in PDG Code " + string());
1163  }
1164  check();
1165  }
1166 };
1167 
1168 static_assert(sizeof(PdgCode) == 4, "should fit into 32 bit integer");
1169 
1176 std::istream& operator>>(std::istream& is, PdgCode& code);
1182 std::ostream& operator<<(std::ostream& is, const PdgCode& code);
1183 
1185 inline bool is_dilepton(const PdgCode pdg1, const PdgCode pdg2) {
1186  const auto c1 = pdg1.code();
1187  const auto c2 = pdg2.code();
1188  const auto min = std::min(c1, c2);
1189  const auto max = std::max(c1, c2);
1190  return (max == 0x11 && min == -0x11) || (max == 0x13 && min == -0x13) ||
1191  (max == 0x12 && min == -0x11) || (max == 0x11 && min == -0x12);
1192 }
1193 
1198 inline bool has_lepton_pair(const PdgCode pdg1, const PdgCode pdg2,
1199  const PdgCode pdg3) {
1200  return is_dilepton(pdg1, pdg2) || is_dilepton(pdg1, pdg3) ||
1201  is_dilepton(pdg2, pdg3);
1202 }
1203 
1204 } // namespace smash
1205 
1206 #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:966
bool is_Nstar1535() const
Definition: pdgcode.h:422
std::int32_t code() const
Definition: pdgcode.h:319
bool is_omega() const
Definition: pdgcode.h:477
bool is_rho() const
Definition: pdgcode.h:483
int antiparticle_sign() const
Definition: pdgcode.h:719
PdgCode(T codenumber, typename std::enable_if_t< std::is_integral_v< T > &&4< sizeof(T), bool >=true)
The creation of PdgCode instances for nuclei that have a 10-digits code cannot be done using the inte...
Definition: pdgcode.h:217
bool operator!=(const PdgCode rhs) const
Definition: pdgcode.h:811
int nucleus_n() const
Number of neutrons in nucleus.
Definition: pdgcode.h:882
int baryon_number() const
Definition: pdgcode.h:388
bool is_meson() const
Definition: pdgcode.h:401
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:928
bool is_Sigma() const
Definition: pdgcode.h:453
bool is_charmonia() const
Definition: pdgcode.h:383
PdgCode(const std::string &codestring)
Initialize using a string The string is interpreted as a hexadecimal number, i.e.,...
Definition: pdgcode.h:131
int nucleus_ap() const
Number of antiprotons in nucleus.
Definition: pdgcode.h:892
int nucleus_an() const
Number of antineutrons in nucleus.
Definition: pdgcode.h:896
int bottomness() const
Definition: pdgcode.h:625
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:739
std::uint32_t n_q3_
third quark field
Definition: pdgcode.h:924
unsigned int spin() const
Definition: pdgcode.h:676
bool is_pion() const
Definition: pdgcode.h:471
bool is_kaon() const
Definition: pdgcode.h:465
bool antiparticle_
first bit: stores the sign.
Definition: pdgcode.h:938
bool is_Sigmastar() const
Definition: pdgcode.h:458
std::uint32_t dump() const
Dumps the bitfield into an unsigned integer.
Definition: pdgcode.h:313
std::uint32_t Z_
Definition: pdgcode.h:980
std::uint32_t n_
first field: "counter"
Definition: pdgcode.h:934
bool is_lepton() const
Definition: pdgcode.h:372
double frac_charm() const
Definition: pdgcode.h:558
bool is_hyperon() const
Definition: pdgcode.h:435
std::uint32_t n_J_
spin quantum number .
Definition: pdgcode.h:922
bool is_nucleus() const
Definition: pdgcode.h:361
bool is_deuteron() const
Definition: pdgcode.h:489
bool is_proton() const
Definition: pdgcode.h:410
static PdgCode from_decimal(const int pdgcode_decimal)
Construct PDG code from decimal number.
Definition: pdgcode.h:339
bool is_nucleus_
Definition: pdgcode.h:983
bool is_antiparticle_of(const PdgCode rhs) const
Definition: pdgcode.h:814
std::uint32_t bool is_nucleus_
1 for nuclei, 0 for the rest
Definition: pdgcode.h:934
int charmness() const
Definition: pdgcode.h:618
int strangeness() const
Definition: pdgcode.h:611
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:878
int32_t get_decimal() const
Definition: pdgcode.h:837
bool is_baryon() const
Definition: pdgcode.h:398
void set_fields(std::uint32_t abscode)
Sets the bitfield from an unsigned integer.
Definition: pdgcode.h:1153
int isospin3() const
Definition: pdgcode.h:520
void check() const
Do all sorts of validity checks.
Definition: pdgcode.h:284
bool is_nucleon() const
Definition: pdgcode.h:404
std::string string() const
Definition: pdgcode.h:322
int test_code() const
Checks the integer for invalid hex digits.
Definition: pdgcode.h:250
PdgCode get_antiparticle() const
Construct the antiparticle to a given PDG code.
Definition: pdgcode.h:329
int nucleus_La() const
Number of Lambdas in nucleus.
Definition: pdgcode.h:888
double frac_bottom() const
Definition: pdgcode.h:583
bool is_hadron() const
Definition: pdgcode.h:367
bool operator==(const PdgCode rhs) const
Definition: pdgcode.h:808
PdgCode()
Standard initializer.
Definition: pdgcode.h:125
int nucleus_A() const
Nucleus mass number.
Definition: pdgcode.h:906
static PdgCode invalid()
PdgCode 0x0 is guaranteed not to be valid by the PDG standard, but it passes all tests here,...
Definition: pdgcode.h:826
int nucleus_aLa() const
Number of anti-Lambdas in nucleus.
Definition: pdgcode.h:902
bool operator<(const PdgCode rhs) const
Sorts PDG Codes according to their numeric value.
Definition: pdgcode.h:799
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:438
std::uint32_t A_
Definition: pdgcode.h:981
bool is_neutron() const
Definition: pdgcode.h:416
void deexcite()
Remove all excitation, except spin. Sign and quark content remains.
Definition: pdgcode.h:857
std::uint32_t n_L_
"angular momentum"
Definition: pdgcode.h:930
void set_from_string(const std::string &codestring)
Set the PDG code from the given string.
Definition: pdgcode.h:1048
bool is_triton() const
Definition: pdgcode.h:495
std::uint32_t ucode() const
Definition: pdgcode.h:1002
std::uint32_t n_R_
"radial excitation"
Definition: pdgcode.h:932
std::uint32_t dump_
The bitfield dumped into a single integer.
Definition: pdgcode.h:957
bool is_neutrino() const
Definition: pdgcode.h:378
bool is_Delta() const
Definition: pdgcode.h:428
std::uint32_t n_Lambda_
Definition: pdgcode.h:979
double frac_strange() const
Definition: pdgcode.h:533
bool has_antiparticle() const
Definition: pdgcode.h:504
bool is_Lambda() const
Definition: pdgcode.h:448
std::uint32_t get_digit_from_char(const char inp) const
Definition: pdgcode.h:1010
bool is_Xi() const
Definition: pdgcode.h:443
bool is_heavy_flavor() const
Definition: pdgcode.h:602
std::int32_t quarks() const
Definition: pdgcode.h:723
std::uint32_t I_
Definition: pdgcode.h:982
std::uint32_t n_q2_
second quark field
Definition: pdgcode.h:926
unsigned int spin_degeneracy() const
Definition: pdgcode.h:712
std::uint32_t excitation_
The excitation digits n_, n_R_, n_L_.
Definition: pdgcode.h:968
int charge() const
The charge of the particle.
Definition: pdgcode.h:635
std::ostream & operator<<(std::ostream &is, const PdgCode &code)
Writes the textual representation of the PDG code to the output stream.
Definition: pdgcode.cc:88
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⁰.
constexpr int Sigma_star_z
Σ*⁰
constexpr int p
Proton.
constexpr int omega
ω.
constexpr int N1535_z
N(1535)⁰.
constexpr int Sigma_star_p
Σ*⁺
constexpr int pi_z
π⁰.
constexpr int n
Neutron.
constexpr int Delta_m
Δ⁻.
constexpr int Delta_z
Δ⁰.
constexpr int rho_z
ρ⁰.
constexpr int Sigma_star_m
Σ*⁻
constexpr int pi_m
π⁻.
constexpr int N1535_p
N(1535)⁺.
Definition: action.h:24
bool is_dilepton(const PdgCode pdg1, const PdgCode pdg2)
Definition: pdgcode.h:1185
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:1198
std::string to_string(ThermodynamicQuantity quantity)
Convert a ThermodynamicQuantity enum value to its corresponding string.
Definition: stringify.cc:26
thrown for invalid inputs
Definition: pdgcode.h:114