Version: SMASH-3.2
pdgcode.h
Go to the documentation of this file.
1 /*
2  *
3  * Copyright (c) 2014-2024
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 
111 class PdgCode {
112  public:
117  struct InvalidPdgCode : public std::invalid_argument {
118  using std::invalid_argument::invalid_argument;
119  };
120 
121  /****************************************************************************
122  * *
123  * First, the constructors *
124  * *
125  ****************************************************************************/
126 
128  PdgCode() : dump_(0x0) {}
134  explicit PdgCode(const std::string& codestring) {
135  set_from_string(codestring);
136  }
137 
146  PdgCode(std::int32_t codenumber) : dump_(0x0) { // NOLINT(runtime/explicit)
147  digits_.antiparticle_ = false;
148  if (codenumber < 0) {
149  digits_.antiparticle_ = true;
150  codenumber = -codenumber;
151  }
152  set_fields(codenumber);
153  }
158  explicit PdgCode(const std::uint32_t abscode) : dump_(0x0) {
159  // use the first bit for the antiparticle_ boolean.
160  digits_.antiparticle_ = ((abscode & 0x80000000u) != 0);
161  set_fields(abscode);
162  }
163 
219  template <typename T>
220  PdgCode(T codenumber, // NOLINT(runtime/explicit)
221  typename std::enable_if_t<std::is_integral_v<T> && 4 < sizeof(T),
222  bool> = true)
223  : PdgCode{std::invoke([&codenumber]() {
224  std::stringstream stream;
225  char sign = '+';
226  if (codenumber < 0) {
227  sign = '-';
228  codenumber = -codenumber;
229  }
230  stream << sign << std::hex << codenumber;
231  return stream.str();
232  })} {}
233 
234  /****************************************************************************
235  * *
236  * test function and export functions *
237  * *
238  ****************************************************************************/
239 
253  inline int test_code() const {
254  int fail = 0;
255  if (digits_.n_ > 9) {
256  fail |= 1 << 6;
257  }
258  if (digits_.n_R_ > 9) {
259  fail |= 1 << 5;
260  }
261  if (digits_.n_L_ > 9) {
262  fail |= 1 << 4;
263  }
264  if (digits_.n_q1_ > 6) {
265  fail |= 1 << 3;
266  }
267  if (digits_.n_q2_ > 6) {
268  fail |= 1 << 2;
269  }
270  if (digits_.n_q3_ > 6) {
271  fail |= 1 << 1;
272  }
273  if (digits_.n_J_ > 15) {
274  fail |= 1;
275  }
276  return fail;
277  }
278 
287  void check() const {
288  // n_J must be odd for mesons and even for baryons (and cannot be zero)
289  if (is_hadron()) {
290  if (baryon_number() == 0) {
291  // mesons: special cases K0_L=0x130 and K0_S=0x310
292  if ((digits_.n_J_ % 2 == 0) && dump() != 0x130 && dump() != 0x310) {
293  throw InvalidPdgCode("Invalid PDG code " + string() +
294  " (meson with even n_J)");
295  }
296  } else {
297  if ((digits_.n_J_ % 2 != 0) || digits_.n_J_ == 0) {
298  throw InvalidPdgCode("Invalid PDG code " + string() +
299  " (baryon with odd n_J)");
300  }
301  }
302  } else {
303  if (digits_.n_J_ == 0 && dump() != 0x0) {
304  throw InvalidPdgCode("Invalid PDG code " + string() + " (n_J==0)");
305  }
306  }
307  /* The antiparticle flag only makes sense for particle types
308  * that have an antiparticle. */
309  if (digits_.antiparticle_ && !has_antiparticle()) {
310  throw InvalidPdgCode("Invalid PDG code " + string() +
311  " (cannot be negative)");
312  }
313  }
314 
316  inline std::uint32_t dump() const {
317  // this cuts the three unused bits.
318  return (dump_ & 0x8fffffff);
319  }
320 
322  inline std::int32_t code() const { return antiparticle_sign() * ucode(); }
323 
325  inline std::string string() const {
326  std::stringstream ss;
327  ss << get_decimal();
328  return ss.str();
329  }
330 
333  PdgCode result = *this;
334  result.digits_.antiparticle_ = !digits_.antiparticle_;
335  return result;
336  }
337 
342  static PdgCode from_decimal(const int pdgcode_decimal) {
343  // Nucleus and special codes with 2J+1 > 9
344  if (std::abs(pdgcode_decimal) > 1E7) {
345  return PdgCode(std::to_string(pdgcode_decimal));
346  }
347  int a = pdgcode_decimal;
348  int hex_pdg = 0, tmp = 1;
349  while (a) {
350  hex_pdg += (a % 10) * tmp;
351  tmp *= 16;
352  a = a / 10;
353  }
354  return PdgCode(hex_pdg);
355  }
356 
357  /****************************************************************************
358  * *
359  * accessors of various properties *
360  * *
361  ****************************************************************************/
362 
364  inline bool is_nucleus() const {
365  assert(digits_.is_nucleus_ == nucleus_.is_nucleus_);
366  return nucleus_.is_nucleus_;
367  }
368 
370  inline bool is_hadron() const {
371  return (digits_.n_q3_ != 0 && digits_.n_q2_ != 0 && !is_nucleus());
372  }
373 
375  inline bool is_lepton() const {
376  return (digits_.n_q1_ == 0 && digits_.n_q2_ == 0 && digits_.n_q3_ == 1 &&
377  !is_nucleus());
378  }
379 
381  inline bool is_neutrino() const {
382  return (is_lepton() && digits_.n_J_ % 2 == 0);
383  }
384 
386  inline bool is_charmonia() const {
387  return is_meson() && digits_.n_q2_ == 4 && digits_.n_q3_ == 4;
388  }
389 
391  inline int baryon_number() const {
392  if (is_nucleus()) {
393  return static_cast<int>(nucleus_.A_) * antiparticle_sign();
394  }
395  if (!is_hadron() || digits_.n_q1_ == 0) {
396  return 0;
397  }
398  return antiparticle_sign();
399  }
401  inline bool is_baryon() const { return is_hadron() && digits_.n_q1_ != 0; }
402 
404  inline bool is_meson() const { return is_hadron() && digits_.n_q1_ == 0; }
405 
407  inline bool is_nucleon() const {
408  const auto abs_code = std::abs(code());
409  return (abs_code == pdg::p || abs_code == pdg::n);
410  }
411 
413  inline bool is_proton() const {
414  const auto abs_code = std::abs(code());
415  return (abs_code == pdg::p);
416  }
417 
419  inline bool is_neutron() const {
420  const auto abs_code = std::abs(code());
421  return (abs_code == pdg::n);
422  }
423 
425  inline bool is_Nstar1535() const {
426  const auto abs_code = std::abs(code());
427  return (abs_code == pdg::N1535_p || abs_code == pdg::N1535_z);
428  }
429 
431  inline bool is_Delta() const {
432  const auto abs_code = std::abs(code());
433  return (abs_code == pdg::Delta_pp || abs_code == pdg::Delta_p ||
434  abs_code == pdg::Delta_z || abs_code == pdg::Delta_m);
435  }
436 
438  inline bool is_hyperon() const { return is_hadron() && digits_.n_q1_ == 3; }
439 
441  inline bool is_Omega() const {
442  return is_hyperon() && digits_.n_q2_ == 3 && digits_.n_q3_ == 3;
443  }
444 
446  inline bool is_Xi() const {
447  return is_hyperon() && digits_.n_q2_ == 3 && digits_.n_q3_ != 3;
448  }
449 
451  inline bool is_Lambda() const {
452  return is_hyperon() && digits_.n_q2_ == 1 && digits_.n_q3_ == 2;
453  }
454 
456  inline bool is_Sigma() const {
457  return is_hyperon() && digits_.n_q2_ != 3 && !is_Lambda();
458  }
459 
461  inline bool is_kaon() const {
462  const auto abs_code = std::abs(code());
463  return (abs_code == pdg::K_p) || (abs_code == pdg::K_z);
464  }
465 
467  inline bool is_pion() const {
468  const auto c = code();
469  return (c == pdg::pi_z) || (c == pdg::pi_p) || (c == pdg::pi_m);
470  }
471 
473  inline bool is_omega() const {
474  const auto c = code();
475  return c == pdg::omega;
476  }
477 
479  inline bool is_rho() const {
480  const auto c = code();
481  return (c == pdg::rho_z) || (c == pdg::rho_p) || (c == pdg::rho_m);
482  }
483 
485  inline bool is_deuteron() const {
486  return is_nucleus() && nucleus_.A_ == 2 && nucleus_.Z_ == 1 &&
487  nucleus_.n_Lambda_ == 0 && nucleus_.I_ == 0;
488  }
489 
491  inline bool is_triton() const {
492  return is_nucleus() && nucleus_.A_ == 3 && nucleus_.Z_ == 1 &&
493  nucleus_.n_Lambda_ == 0 && nucleus_.I_ == 0;
494  }
495 
500  bool has_antiparticle() const {
501  if (is_nucleus()) {
502  return true;
503  }
504  if (is_hadron()) {
505  return (baryon_number() != 0) || (digits_.n_q2_ != digits_.n_q3_);
506  } else {
507  return digits_.n_q3_ == 1; // leptons!
508  }
509  }
510 
516  inline int isospin3() const {
517  /* net_quark_number(2) is the number of u quarks,
518  * net_quark_number(1) is the number of d quarks. */
519  return net_quark_number(2) - net_quark_number(1);
520  }
521 
529  inline double frac_strange() const {
530  /* The quarkonium state has 0 net strangeness
531  * but there are actually 2 strange quarks out of 2 total */
532  if (is_hadron() && digits_.n_q3_ == 3 && digits_.n_q2_ == 3) {
533  return 1.;
534  } else {
535  // For all other cases, there isn't both a strange and anti-strange
536  if (is_baryon()) {
537  return std::abs(strangeness()) / 3.;
538  } else if (is_meson()) {
539  return std::abs(strangeness()) / 2.;
540  } else {
541  /* If not baryon or meson, this should be 0, as AQM does not
542  * extend to non-hadrons */
543  return 0.;
544  }
545  }
546  }
547 
555  inline double frac_charm() const {
556  /* The charmonium state has 0 net charmness
557  * but there are actually 2 charm quarks out of 2 total */
558  if (is_hadron() && digits_.n_q3_ == 4 && digits_.n_q2_ == 4) {
559  return 1.;
560  } else {
561  // For all other cases, there isn't both a charm and anti-charm
562  if (is_baryon()) {
563  return std::abs(charmness()) / 3.;
564  } else if (is_meson()) {
565  return std::abs(charmness()) / 2.;
566  } else {
567  /* If not baryon or meson, this should be 0, as AQM does not
568  * extend to non-hadrons */
569  return 0.;
570  }
571  }
572  }
573 
581  inline double frac_bottom() const {
582  /* The bottomonium state has 0 net bottomness
583  * but there are actually 2 bottom quarks out of 2 total */
584  if (is_hadron() && digits_.n_q3_ == 3 && digits_.n_q2_ == 3) {
585  return 1.;
586  } else {
587  // For all other cases, there isn't both a bottom and anti-bottom
588  if (is_baryon()) {
589  return std::abs(bottomness()) / 3.;
590  } else if (is_meson()) {
591  return std::abs(bottomness()) / 2.;
592  } else {
593  /* If not baryon or meson, this should be 0, as AQM does not
594  * extend to non-hadrons */
595  return 0.;
596  }
597  }
598  }
599 
605  inline int strangeness() const { return -net_quark_number(3); }
606 
612  inline int charmness() const { return +net_quark_number(4); }
613 
619  inline int bottomness() const { return -net_quark_number(5); }
620 
629  int charge() const {
630  if (is_hadron() || is_nucleus()) {
631  // Q will accumulate 3*charge (please excuse the upper case. I
632  // want to distinguish this from q which might be interpreted as
633  // shorthand for "quark".)
634  int Q = 0;
635  /* This loops over d,u,s,c,b,t quarks (the latter can be safely ignored,
636  * but I don't think this will be a bottle neck. */
637  for (int i = 1; i < 7; i++) {
638  /* u,c,t quarks have charge = 2/3 e, while d,s,b quarks have -1/3 e.
639  * The antiparticle sign is already in net_quark_number. */
640  Q += (i % 2 == 0 ? 2 : -1) * net_quark_number(i);
641  }
642  return Q / 3;
643  }
644  /* non-hadron:
645  * Leptons: 11, 13, 15 are e, μ, τ and have a charge -1, while
646  * 12, 14, 16 are the neutrinos that have no charge. */
647  if (digits_.n_q3_ == 1) {
648  return -1 * (digits_.n_J_ % 2) * antiparticle_sign();
649  }
650  /* Bosons: 24 is the W+, all else is uncharged.
651  * we ignore the first digits so that this also finds strange gauge
652  * boson "resonances" (in particular, \f$\tilde \chi_1^+\f$ with PDG
653  * Code 1000024). */
654  if ((dump_ & 0x0000ffff) == 0x24) {
655  return antiparticle_sign();
656  }
657  // default (this includes all other Bosons) is 0.
658  return 0;
659  }
660 
670  inline unsigned int spin() const {
671  if (is_nucleus()) {
672  // Generally spin of a nucleus cannot be simply guessed, it should be
673  // provided from some table. However, here we only care about a
674  // limited set of light nuclei with A <= 4.
675  if (nucleus_.A_ == 2) {
676  // Deuteron spin is 1
677  return 2;
678  } else if (nucleus_.A_ == 3) {
679  // Tritium and He-3 spin are 1/2
680  // Hypertriton spin is not firmly determined, but decay branching ratios
681  // indicate spin 1/2
682  return 1;
683  } else if (nucleus_.A_ == 4) {
684  // He-4 spin is 0
685  return 0;
686  }
687  throw std::runtime_error("Unknown spin of nucleus.");
688  // Alternative possibility is to guess 1/2 for fermions and 0 for bosons
689  // as 2 * (nucleus_.A_ % 2).
690  }
691 
692  if (is_hadron()) {
693  if (digits_.n_J_ == 0) {
694  return 0; // special cases: K0_L=0x130 & K0_S=0x310
695  } else {
696  return digits_.n_J_ - 1;
697  }
698  }
699  /* this assumes that we only have white particles (no single
700  * quarks): Electroweak fermions have 11-17, so the
701  * second-to-last-digit is the spin. The same for the Bosons: they
702  * have 21-29 and 2spin = 2 (this fails for the Higgs). */
703  return digits_.n_q3_;
704  }
706  inline unsigned int spin_degeneracy() const {
707  if (is_hadron() && digits_.n_J_ > 0) {
708  return digits_.n_J_;
709  }
710  return spin() + 1;
711  }
713  inline int antiparticle_sign() const {
714  return (digits_.antiparticle_ ? -1 : +1);
715  }
717  inline std::int32_t quarks() const {
718  if (!is_hadron() || is_nucleus()) {
719  return 0;
720  }
721  return chunks_.quarks_;
722  }
723 
733  std::array<int, 3> quark_content() const {
734  std::array<int, 3> result = {static_cast<int>(digits_.n_q1_),
735  static_cast<int>(digits_.n_q2_),
736  static_cast<int>(digits_.n_q3_)};
737  if (is_hadron()) {
738  // Antibaryons
739  if (digits_.n_q1_ != 0 && digits_.antiparticle_) {
740  for (size_t i = 0; i < 3; i++) {
741  result[i] = -result[i];
742  }
743  }
744  // Mesons
745  if (digits_.n_q1_ == 0) {
746  // Own antiparticle
747  if (digits_.n_q2_ == digits_.n_q3_) {
748  result[2] = -result[2];
749  } else {
750  // Like pi-
751  if (digits_.antiparticle_) {
752  result[1] = -result[1];
753  // Like pi+
754  } else {
755  result[2] = -result[2];
756  }
757  }
758  // add extra minus sign according to the pdg convention
759  if (digits_.n_q2_ != digits_.n_q3_ && digits_.n_q2_ % 2 == 1) {
760  for (int i = 1; i <= 2; i++) {
761  result[i] = -result[i];
762  }
763  }
764  }
765  } else {
766  result = {0, 0, 0};
767  }
768  return result;
769  }
770 
781  bool contains_enough_valence_quarks(int valence_quarks_required) const;
782 
783  /****************************************************************************
784  * *
785  * operations with more than one PDG Code *
786  * *
787  ****************************************************************************/
788 
793  inline bool operator<(const PdgCode rhs) const {
794  return dump_ < rhs.dump_;
795  /* the complex thing to do here is to calculate:
796  * code() < rhs.code()
797  * but for getting a total order that's overkill. The uint32_t value in
798  * dump_ works just fine. */
799  }
800 
802  inline bool operator==(const PdgCode rhs) const { return dump_ == rhs.dump_; }
803 
805  inline bool operator!=(const PdgCode rhs) const { return !(*this == rhs); }
806 
808  inline bool is_antiparticle_of(const PdgCode rhs) const {
809  return code() == -rhs.code();
810  }
811 
813  friend std::istream& operator>>(std::istream& is, PdgCode& code);
814 
820  static PdgCode invalid() { return PdgCode(0x0); }
821 
831  int32_t get_decimal() const {
832  if (is_nucleus()) {
833  // ±10LZZZAAAI
834  return antiparticle_sign() *
835  (nucleus_.I_ + 10 * nucleus_.A_ + 10000 * nucleus_.Z_ +
836  10000000 * nucleus_.n_Lambda_ + 1000000000);
837  }
838  int n_J_1 = 0;
839  int n_J_2 = digits_.n_J_;
840  if (n_J_2 > 9) {
841  n_J_1 = n_J_2 - 9;
842  n_J_2 = 9;
843  }
844  return antiparticle_sign() *
845  (n_J_2 + digits_.n_q3_ * 10 + digits_.n_q2_ * 100 +
846  digits_.n_q1_ * 1000 + digits_.n_L_ * 10000 +
847  digits_.n_R_ * 100000 + digits_.n_ * 1000000 + n_J_1 * 10000000);
848  }
849 
851  void deexcite() {
852  if (!is_nucleus()) {
853  chunks_.excitation_ = 0;
854  } else {
855  nucleus_.I_ = 0;
856  }
857  }
858 
869  int net_quark_number(const int quark) const;
870 
872  int nucleus_p() const {
873  return (is_nucleus() && !nucleus_.antiparticle_) ? nucleus_.Z_ : 0;
874  }
876  int nucleus_n() const {
877  return (is_nucleus() && !nucleus_.antiparticle_)
878  ? nucleus_.A_ - nucleus_.Z_ - nucleus_.n_Lambda_
879  : 0;
880  }
882  int nucleus_La() const {
883  return (is_nucleus() && !nucleus_.antiparticle_) ? nucleus_.n_Lambda_ : 0;
884  }
886  int nucleus_ap() const {
887  return (is_nucleus() && nucleus_.antiparticle_) ? nucleus_.Z_ : 0;
888  }
890  int nucleus_an() const {
891  return (is_nucleus() && nucleus_.antiparticle_)
892  ? nucleus_.A_ - nucleus_.Z_ - nucleus_.n_Lambda_
893  : 0;
894  }
896  int nucleus_aLa() const {
897  return (is_nucleus() && nucleus_.antiparticle_) ? nucleus_.n_Lambda_ : 0;
898  }
900  int nucleus_A() const { return is_nucleus() ? nucleus_.A_ : 0; }
901 
902  private:
908  union {
913  struct {
914 #if defined(LITTLE_ENDIAN_ARCHITECTURE) || defined(DOXYGEN)
916  std::uint32_t n_J_ : 4;
918  std::uint32_t n_q3_ : 4;
920  std::uint32_t n_q2_ : 4;
922  std::uint32_t n_q1_ : 4;
924  std::uint32_t n_L_ : 4;
926  std::uint32_t n_R_ : 4;
928  std::uint32_t n_ : 4, : 2;
930  bool is_nucleus_ : 1;
932  bool antiparticle_ : 1;
933 #elif defined(BIG_ENDIAN_ARCHITECTURE) // reverse ordering
934  bool antiparticle_ : 1;
935  bool is_nucleus_ : 1, : 2;
936  std::uint32_t n_ : 4;
937  std::uint32_t n_R_ : 4;
938  std::uint32_t n_L_ : 4;
939  std::uint32_t n_q1_ : 4;
940  std::uint32_t n_q2_ : 4;
941  std::uint32_t n_q3_ : 4;
942  std::uint32_t n_J_ : 4;
943 #else
944 #error Endianness macro of the machine not defined.
945 #endif
946  } digits_;
951  std::uint32_t dump_;
956  struct {
957 #if defined(LITTLE_ENDIAN_ARCHITECTURE) || defined(DOXYGEN)
958  std::uint32_t : 4;
960  std::uint32_t quarks_ : 12;
962  std::uint32_t excitation_ : 12, : 4;
963 #elif defined(BIG_ENDIAN_ARCHITECTURE) // reverse ordering
964  std::uint32_t : 4, excitation_ : 12;
965  std::uint32_t quarks_ : 12, : 4;
966 #else
967 #error Endianness macro of the machine not defined.
968 #endif
969  } chunks_;
971  struct {
972 #if defined(LITTLE_ENDIAN_ARCHITECTURE) || defined(DOXYGEN)
973  std::uint32_t n_Lambda_ : 6;
974  std::uint32_t Z_ : 10;
975  std::uint32_t A_ : 10;
976  std::uint32_t I_ : 4;
977  bool is_nucleus_ : 1;
978  bool antiparticle_ : 1;
979 #elif defined(BIG_ENDIAN_ARCHITECTURE) // reverse ordering
980  bool antiparticle_ : 1;
981  bool is_nucleus_ : 1;
982  std::uint32_t I_ : 4;
983  std::uint32_t A_ : 10;
984  std::uint32_t Z_ : 10;
985  std::uint32_t n_Lambda_ : 6;
986 #else
987 #error Endianness macro of the machine not defined.
988 #endif
989  } nucleus_;
990  };
991 
996  inline std::uint32_t ucode() const { return (dump_ & 0x0fffffff); }
997 
1004  inline std::uint32_t get_digit_from_char(const char inp) const {
1005  // Decimal digit
1006  if (48 <= inp && inp <= 57) {
1007  return inp - 48;
1008  }
1009  // Hexdecimal digit, uppercase
1010  if (65 <= inp && inp <= 70) {
1011  return inp - 65 + 10;
1012  }
1013  // Hexdecimal digit, lowercase
1014  if (97 <= inp && inp <= 102) {
1015  return inp - 97 + 10;
1016  }
1017  throw InvalidPdgCode("PdgCode: Invalid character " + std::string(&inp, 1) +
1018  " found.\n");
1019  }
1020 
1042  inline void set_from_string(const std::string& codestring) {
1043  dump_ = 0;
1044  // Implicit with the above: digits_.antiparticle_ = false;
1045  digits_.n_ = digits_.n_R_ = digits_.n_L_ = digits_.n_q1_ = digits_.n_q2_ =
1046  digits_.n_q3_ = digits_.n_J_ = digits_.is_nucleus_ = 0;
1047  size_t length = codestring.size();
1048  if (length < 1) {
1049  throw InvalidPdgCode("Empty string does not contain PDG Code\n");
1050  }
1051  int c = 0;
1052  /* Look at current character; if it is a + or minus sign, read it
1053  * and advance to next char. */
1054  if (codestring[c] == '-') {
1055  digits_.antiparticle_ = true;
1056  ++c;
1057  } else if (codestring[c] == '+') {
1058  digits_.antiparticle_ = false;
1059  ++c;
1060  }
1061  // Save if the first character was a sign:
1062  unsigned int sign = c;
1063 
1064  // Nucleus
1065  if (length == 10 + sign) {
1066  nucleus_.is_nucleus_ = true;
1067  if (codestring[c] != '1' || codestring[c + 1] != '0') {
1068  throw InvalidPdgCode("Pdg code of nucleus \"" + codestring +
1069  "\" should start with 10\n");
1070  }
1071  c += 2;
1072  // ±10LZZZAAAI is the standard for nuclei
1073  std::array<int, 8> digits;
1074  for (int i = 0; i < 8; i++) {
1075  digits[i] = get_digit_from_char(codestring[c + i]);
1076  }
1077  nucleus_.n_Lambda_ = digits[0];
1078  nucleus_.Z_ = 100 * digits[1] + 10 * digits[2] + digits[3];
1079  nucleus_.A_ = 100 * digits[4] + 10 * digits[5] + digits[6];
1080  nucleus_.I_ = digits[7];
1081  return;
1082  }
1083 
1084  // Codestring shouldn't be longer than 8 + sign, except for nuclei
1085  if (length > 8 + sign) {
1086  throw InvalidPdgCode("String \"" + codestring +
1087  "\" too long for PDG Code\n");
1088  }
1089  /* Please note that in what follows, we actually need c++, not ++c.
1090  * first digit is used for n_J if the last digit is not enough. */
1091  if (length > 7 + sign) {
1092  digits_.n_J_ += get_digit_from_char(codestring[c++]);
1093  }
1094  // Codestring has 7 digits? 7th from last goes in n_.
1095  if (length > 6 + sign) {
1096  digits_.n_ = get_digit_from_char(codestring[c++]);
1097  }
1098  // It has 6 or 7 digits? 6th from last is n_R_.
1099  if (length > 5 + sign) {
1100  digits_.n_R_ = get_digit_from_char(codestring[c++]);
1101  }
1102  // 5th from last is n_L_.
1103  if (length > 4 + sign) {
1104  digits_.n_L_ = get_digit_from_char(codestring[c++]);
1105  }
1106  // 4th from last is n_q1_.
1107  if (length > 3 + sign) {
1108  digits_.n_q1_ = get_digit_from_char(codestring[c++]);
1109  if (digits_.n_q1_ > 6) {
1110  throw InvalidPdgCode("Invalid PDG code " + codestring + " (n_q1>6)");
1111  }
1112  }
1113  // 3rd from last is n_q2_.
1114  if (length > 2 + sign) {
1115  digits_.n_q2_ = get_digit_from_char(codestring[c++]);
1116  if (digits_.n_q2_ > 6) {
1117  throw InvalidPdgCode("Invalid PDG code " + codestring + " (n_q2>6)");
1118  }
1119  }
1120  // Next to last is n_q3_.
1121  if (length > 1 + sign) {
1122  digits_.n_q3_ = get_digit_from_char(codestring[c++]);
1123  if (digits_.n_q3_ > 6) {
1124  throw InvalidPdgCode("Invalid PDG code " + codestring + " (n_q3>6)");
1125  }
1126  }
1127  // Last digit is the spin degeneracy.
1128  if (length > sign) {
1129  digits_.n_J_ += get_digit_from_char(codestring[c++]);
1130  } else {
1131  throw InvalidPdgCode(
1132  "String \"" + codestring +
1133  "\" only consists of a sign, that is no valid PDG Code\n");
1134  }
1135  check();
1136  }
1137 
1147  inline void set_fields(std::uint32_t abscode) {
1148  /* "dump_ =" overwrites antiparticle_, but this needs to have been set
1149  * already, so we carry it around the assignment. */
1150  bool ap = digits_.antiparticle_;
1151  dump_ = abscode & 0x0fffffff;
1152  digits_.antiparticle_ = ap;
1153  int test = test_code();
1154  if (test > 0) {
1155  throw InvalidPdgCode("Invalid digits " + std::to_string(test) +
1156  " in PDG Code " + string());
1157  }
1158  check();
1159  }
1160 };
1161 
1162 static_assert(sizeof(PdgCode) == 4, "should fit into 32 bit integer");
1163 
1170 std::istream& operator>>(std::istream& is, PdgCode& code);
1176 std::ostream& operator<<(std::ostream& is, const PdgCode& code);
1177 
1179 inline bool is_dilepton(const PdgCode pdg1, const PdgCode pdg2) {
1180  const auto c1 = pdg1.code();
1181  const auto c2 = pdg2.code();
1182  const auto min = std::min(c1, c2);
1183  const auto max = std::max(c1, c2);
1184  return (max == 0x11 && min == -0x11) || (max == 0x13 && min == -0x13) ||
1185  (max == 0x12 && min == -0x11) || (max == 0x11 && min == -0x12);
1186 }
1187 
1192 inline bool has_lepton_pair(const PdgCode pdg1, const PdgCode pdg2,
1193  const PdgCode pdg3) {
1194  return is_dilepton(pdg1, pdg2) || is_dilepton(pdg1, pdg3) ||
1195  is_dilepton(pdg2, pdg3);
1196 }
1197 
1198 } // namespace smash
1199 
1200 #endif // SRC_INCLUDE_SMASH_PDGCODE_H_
PdgCode stores a Particle Data Group Particle Numbering Scheme particle type number.
Definition: pdgcode.h:111
std::uint32_t quarks_
The quark digits n_q{1,2,3}_.
Definition: pdgcode.h:960
bool is_Nstar1535() const
Definition: pdgcode.h:425
std::int32_t code() const
Definition: pdgcode.h:322
bool is_omega() const
Definition: pdgcode.h:473
bool is_rho() const
Definition: pdgcode.h:479
int antiparticle_sign() const
Definition: pdgcode.h:713
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:220
bool operator!=(const PdgCode rhs) const
Definition: pdgcode.h:805
int nucleus_n() const
Number of neutrons in nucleus.
Definition: pdgcode.h:876
int baryon_number() const
Definition: pdgcode.h:391
bool is_meson() const
Definition: pdgcode.h:404
PdgCode(const std::uint32_t abscode)
Receive an unsigned integer and process it into a PDG Code.
Definition: pdgcode.h:158
std::uint32_t n_q1_
first quark field. 0 for mesons.
Definition: pdgcode.h:922
bool is_Sigma() const
Definition: pdgcode.h:456
bool is_charmonia() const
Definition: pdgcode.h:386
PdgCode(const std::string &codestring)
Initialize using a string The string is interpreted as a hexadecimal number, i.e.,...
Definition: pdgcode.h:134
int nucleus_ap() const
Number of antiprotons in nucleus.
Definition: pdgcode.h:886
int nucleus_an() const
Number of antineutrons in nucleus.
Definition: pdgcode.h:890
int bottomness() const
Definition: pdgcode.h:619
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:733
std::uint32_t n_q3_
third quark field
Definition: pdgcode.h:918
unsigned int spin() const
Definition: pdgcode.h:670
bool is_pion() const
Definition: pdgcode.h:467
bool is_kaon() const
Definition: pdgcode.h:461
bool antiparticle_
first bit: stores the sign.
Definition: pdgcode.h:932
std::uint32_t dump() const
Dumps the bitfield into an unsigned integer.
Definition: pdgcode.h:316
std::uint32_t Z_
Definition: pdgcode.h:974
std::uint32_t n_
first field: "counter"
Definition: pdgcode.h:928
bool is_lepton() const
Definition: pdgcode.h:375
double frac_charm() const
Definition: pdgcode.h:555
bool is_hyperon() const
Definition: pdgcode.h:438
std::uint32_t n_J_
spin quantum number .
Definition: pdgcode.h:916
bool is_nucleus() const
Definition: pdgcode.h:364
bool is_deuteron() const
Definition: pdgcode.h:485
bool is_proton() const
Definition: pdgcode.h:413
static PdgCode from_decimal(const int pdgcode_decimal)
Construct PDG code from decimal number.
Definition: pdgcode.h:342
bool is_nucleus_
Definition: pdgcode.h:977
bool is_antiparticle_of(const PdgCode rhs) const
Definition: pdgcode.h:808
std::uint32_t bool is_nucleus_
1 for nuclei, 0 for the rest
Definition: pdgcode.h:928
int charmness() const
Definition: pdgcode.h:612
int strangeness() const
Definition: pdgcode.h:605
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:872
int32_t get_decimal() const
Definition: pdgcode.h:831
bool is_baryon() const
Definition: pdgcode.h:401
void set_fields(std::uint32_t abscode)
Sets the bitfield from an unsigned integer.
Definition: pdgcode.h:1147
int isospin3() const
Definition: pdgcode.h:516
void check() const
Do all sorts of validity checks.
Definition: pdgcode.h:287
bool is_nucleon() const
Definition: pdgcode.h:407
std::string string() const
Definition: pdgcode.h:325
int test_code() const
Checks the integer for invalid hex digits.
Definition: pdgcode.h:253
PdgCode get_antiparticle() const
Construct the antiparticle to a given PDG code.
Definition: pdgcode.h:332
int nucleus_La() const
Number of Lambdas in nucleus.
Definition: pdgcode.h:882
double frac_bottom() const
Definition: pdgcode.h:581
bool is_hadron() const
Definition: pdgcode.h:370
bool operator==(const PdgCode rhs) const
Definition: pdgcode.h:802
PdgCode()
Standard initializer.
Definition: pdgcode.h:128
int nucleus_A() const
Nucleus mass number.
Definition: pdgcode.h:900
static PdgCode invalid()
PdgCode 0x0 is guaranteed not to be valid by the PDG standard, but it passes all tests here,...
Definition: pdgcode.h:820
int nucleus_aLa() const
Number of anti-Lambdas in nucleus.
Definition: pdgcode.h:896
bool operator<(const PdgCode rhs) const
Sorts PDG Codes according to their numeric value.
Definition: pdgcode.h:793
PdgCode(std::int32_t codenumber)
Receive a signed integer and process it into a PDG Code.
Definition: pdgcode.h:146
bool is_Omega() const
Definition: pdgcode.h:441
std::uint32_t A_
Definition: pdgcode.h:975
bool is_neutron() const
Definition: pdgcode.h:419
void deexcite()
Remove all excitation, except spin. Sign and quark content remains.
Definition: pdgcode.h:851
std::uint32_t n_L_
"angular momentum"
Definition: pdgcode.h:924
void set_from_string(const std::string &codestring)
Set the PDG code from the given string.
Definition: pdgcode.h:1042
bool is_triton() const
Definition: pdgcode.h:491
std::uint32_t ucode() const
Definition: pdgcode.h:996
std::uint32_t n_R_
"radial excitation"
Definition: pdgcode.h:926
std::uint32_t dump_
The bitfield dumped into a single integer.
Definition: pdgcode.h:951
bool is_neutrino() const
Definition: pdgcode.h:381
bool is_Delta() const
Definition: pdgcode.h:431
std::uint32_t n_Lambda_
Definition: pdgcode.h:973
double frac_strange() const
Definition: pdgcode.h:529
bool has_antiparticle() const
Definition: pdgcode.h:500
bool is_Lambda() const
Definition: pdgcode.h:451
std::uint32_t get_digit_from_char(const char inp) const
Definition: pdgcode.h:1004
bool is_Xi() const
Definition: pdgcode.h:446
std::int32_t quarks() const
Definition: pdgcode.h:717
std::uint32_t I_
Definition: pdgcode.h:976
std::uint32_t n_q2_
second quark field
Definition: pdgcode.h:920
unsigned int spin_degeneracy() const
Definition: pdgcode.h:706
std::uint32_t excitation_
The excitation digits n_, n_R_, n_L_.
Definition: pdgcode.h:962
int charge() const
The charge of the particle.
Definition: pdgcode.h:629
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 p
Proton.
constexpr int omega
ω.
constexpr int N1535_z
N(1535)⁰.
constexpr int pi_z
π⁰.
constexpr int n
Neutron.
constexpr int Delta_m
Δ⁻.
constexpr int Delta_z
Δ⁰.
constexpr int rho_z
ρ⁰.
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:1179
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:1192
thrown for invalid inputs
Definition: pdgcode.h:117