 |
Version: SMASH-1.8
|
|
#include <particletype.h>
Particle type contains the static properties of a particle species
Before creation of Experiment, SMASH initializes the list of particles (list_all). After construction these values are immutable.
The list of particles is stored in such a way that look up of a ParticleType object (find) for a given PDG code is as efficient as possible ( \(\mathcal O(\log N)\)). This is still not efficient enough to use PdgCode as a substitute for storing information about a particle type, though. Use ParticleTypePtr instead.
Definition at line 97 of file particletype.h.
|
| ParticleType (std::string n, double m, double w, Parity p, PdgCode id) |
| Creates a fully initialized ParticleType object. More...
|
|
| ParticleType (const ParticleType &)=delete |
| Copies are not allowed as they break intended use. More...
|
|
ParticleType & | operator= (const ParticleType &)=delete |
| assignment is not allowed, see copy constructor above More...
|
|
| ParticleType (ParticleType &&)=default |
| move ctors are needed for std::sort More...
|
|
ParticleType & | operator= (ParticleType &&)=default |
| move ctors are needed for std::sort More...
|
|
const DecayModes & | decay_modes () const |
|
const std::string & | name () const |
|
double | mass () const |
|
double | mass_sqr () const |
|
double | width_at_pole () const |
|
Parity | parity () const |
|
PdgCode | pdgcode () const |
|
bool | has_antiparticle () const |
|
ParticleTypePtr | get_antiparticle () const |
|
int | antiparticle_sign () const |
|
int | isospin () const |
| Returns twice the isospin vector length \(I\). More...
|
|
int | isospin3 () const |
|
double | isospin3_rel () const |
|
IsoParticleType * | iso_multiplet () const |
|
int32_t | charge () const |
| The charge of the particle. More...
|
|
unsigned int | spin () const |
|
bool | is_hadron () const |
|
bool | is_lepton () const |
|
bool | is_baryon () const |
|
bool | is_meson () const |
|
int | baryon_number () const |
|
int | strangeness () const |
|
bool | is_nucleon () const |
|
bool | is_Delta () const |
|
bool | is_rho () const |
|
bool | is_Nstar () const |
|
bool | is_Nstar1535 () const |
|
bool | is_Deltastar () const |
|
bool | is_stable () const |
|
bool | is_nucleus () const |
|
bool | is_deuteron () const |
|
bool | is_dprime () const |
|
double | min_mass_kinematic () const |
| The minimum mass of the resonance that is kinematically allowed. More...
|
|
double | min_mass_spectral () const |
| The minimum mass of the resonance, where the spectral function is non-zero. More...
|
|
double | partial_width (const double m, const DecayBranch *mode) const |
| Get the mass-dependent partial decay width of a particle with mass m in a particular decay mode. More...
|
|
double | total_width (const double m) const |
| Get the mass-dependent total width of a particle with mass m. More...
|
|
bool | wanted_decaymode (const DecayType &t, WhichDecaymodes wh) const |
| Helper Function that containes the if-statement logic that decides if a decay mode is either a hadronic and dilepton decay mode. More...
|
|
DecayBranchList | get_partial_widths (const FourVector p, const ThreeVector x, WhichDecaymodes wh) const |
| Get all the mass-dependent partial decay widths of a particle with mass m. More...
|
|
double | get_partial_width (const double m, const ParticleType &t_a, const ParticleType &t_b) const |
| Get the mass-dependent partial width of a resonance with mass m, decaying into two given daughter particles. More...
|
|
double | get_partial_in_width (const double m, const ParticleData &p_a, const ParticleData &p_b) const |
| Get the mass-dependent partial in-width of a resonance with mass m, decaying into two given daughter particles. More...
|
|
double | spectral_function (double m) const |
| Full spectral function \( A(m) = \frac{2}{\pi} N \frac{m^2\Gamma(m)}{(m^2-m_0^2)^2+(m\Gamma(m))^2} \) of the resonance (relativistic Breit-Wigner distribution with mass-dependent width, where N is a normalization factor). More...
|
|
double | spectral_function_no_norm (double m) const |
| Full spectral function without normalization factor. More...
|
|
double | spectral_function_const_width (double m) const |
|
double | spectral_function_simple (double m) const |
| This one is the most simple form of the spectral function, using a Cauchy distribution (non-relativistic Breit-Wigner with constant width). More...
|
|
double | sample_resonance_mass (const double mass_stable, const double cms_energy, int L=0) const |
| Resonance mass sampling for 2-particle final state with one resonance (type given by 'this') and one stable particle. More...
|
|
std::pair< double, double > | sample_resonance_masses (const ParticleType &t2, const double cms_energy, int L=0) const |
| Resonance mass sampling for 2-particle final state with two resonances. More...
|
|
void | dump_width_and_spectral_function () const |
| Prints out width and spectral function versus mass to the standard output. More...
|
|
bool | operator== (const ParticleType &rhs) const |
|
bool | operator!= (const ParticleType &rhs) const |
|
bool | operator< (const ParticleType &rhs) const |
| "Less than" operator for sorting the ParticleType list (by PDG code) More...
|
|
ParticleTypePtr | operator& () const |
| Returns an object that acts like a pointer, except that it requires only 2 bytes and inhibits pointer arithmetics. More...
|
|
|
static constexpr double | width_cutoff = 1e-5 |
| Decay width cutoff for considering a particle as stable. More...
|
|
|
std::string | name_ |
| name of the particle More...
|
|
double | mass_ |
| pole mass of the particle More...
|
|
double | width_ |
| width of the particle More...
|
|
Parity | parity_ |
| Parity of the particle. More...
|
|
PdgCode | pdgcode_ |
| PDG Code of the particle. More...
|
|
double | min_mass_kinematic_ |
| minimum kinematically allowed mass of the particle Mutable, because it is initialized at first call of minimum mass function, so it's logically const, but not physically const, which is a classical case for using mutable. More...
|
|
double | min_mass_spectral_ |
| minimum mass, where the spectral function is non-zero Mutable, because it is initialized at first call of minimum mass function, so it's logically const, but not physically const, which is a classical case for using mutable. More...
|
|
double | norm_factor_ = -1. |
| This normalization factor ensures that the spectral function is normalized to unity, when integrated over its full domain. More...
|
|
int32_t | charge_ |
| Charge of the particle; filled automatically from pdgcode_. More...
|
|
int | isospin_ |
| Isospin of the particle; filled automatically from pdgcode_. More...
|
|
int | I3_ |
| Isospin projection of the particle; filled automatically from pdgcode_. More...
|
|
IsoParticleType * | iso_multiplet_ = nullptr |
| Container for the isospin multiplet information. More...
|
|
double | max_factor1_ = 1. |
| Maximum factor for single-res mass sampling, cf. sample_resonance_mass. More...
|
|
double | max_factor2_ = 1. |
| Maximum factor for double-res mass sampling, cf. sample_resonance_masses. More...
|
|
◆ ParticleType() [1/3]
smash::ParticleType::ParticleType |
( |
std::string |
n, |
|
|
double |
m, |
|
|
double |
w, |
|
|
Parity |
p, |
|
|
PdgCode |
id |
|
) |
| |
Creates a fully initialized ParticleType object.
- Parameters
-
[in] | n | The name of the particle. |
[in] | m | The mass of the particle. |
[in] | w | The width of the particle. |
[in] | p | The parity of the particle. |
[in] | id | The PDG code of the particle. |
- Note
- The remaining properties ParticleType provides are derived from the PDG code and therefore cannot be set explicitly (this avoids the chance of introducing inconsistencies).
Definition at line 128 of file particletype.cc.
◆ ParticleType() [2/3]
Copies are not allowed as they break intended use.
Instead use a const-ref or ParticleTypePtr (as returned from operator&).
◆ ParticleType() [3/3]
move ctors are needed for std::sort
◆ operator=() [1/2]
assignment is not allowed, see copy constructor above
◆ operator=() [2/2]
move ctors are needed for std::sort
◆ decay_modes()
const DecayModes & smash::ParticleType::decay_modes |
( |
| ) |
const |
- Returns
- the DecayModes object for this particle type.
Definition at line 428 of file particletype.cc.
429 const auto offset =
this - std::addressof(
list_all()[0]);
430 const auto &modes = (*DecayModes::all_decay_modes)[offset];
431 assert(
is_stable() || !modes.is_empty());
◆ name()
const std::string& smash::ParticleType::name |
( |
| ) |
const |
|
inline |
◆ mass()
double smash::ParticleType::mass |
( |
| ) |
const |
|
inline |
◆ mass_sqr()
double smash::ParticleType::mass_sqr |
( |
| ) |
const |
|
inline |
◆ width_at_pole()
double smash::ParticleType::width_at_pole |
( |
| ) |
const |
|
inline |
- Returns
- the particle width (at the mass pole).
Definition at line 150 of file particletype.h.
◆ parity()
Parity smash::ParticleType::parity |
( |
| ) |
const |
|
inline |
◆ pdgcode()
PdgCode smash::ParticleType::pdgcode |
( |
| ) |
const |
|
inline |
◆ has_antiparticle()
bool smash::ParticleType::has_antiparticle |
( |
| ) |
const |
|
inline |
- Returns
- whether a particle has a distinct antiparticle (or whether it is its own antiparticle).
Definition at line 159 of file particletype.h.
◆ get_antiparticle()
◆ antiparticle_sign()
int smash::ParticleType::antiparticle_sign |
( |
| ) |
const |
|
inline |
- Returns
- -1 for antiparticles and +1 for particles.
Definition at line 165 of file particletype.h.
◆ isospin()
int smash::ParticleType::isospin |
( |
| ) |
const |
Returns twice the isospin vector length \(I\).
This returns e.g. 1 for nucleons, 2 for pions and 3 for Deltas. It is always positive.
Definition at line 410 of file particletype.cc.
◆ isospin3()
int smash::ParticleType::isospin3 |
( |
| ) |
const |
|
inline |
- Returns
- twice the isospin-3 component \(I_3\).
This is calculated from the sum of net_quark_number of up and down.
Definition at line 176 of file particletype.h.
◆ isospin3_rel()
double smash::ParticleType::isospin3_rel |
( |
| ) |
const |
|
inline |
- Returns
- the isospin-3 component relative to the total isospin.
Definition at line 179 of file particletype.h.
181 return (I == 0) ? 0 : static_cast<double>(
isospin3()) / I;
◆ iso_multiplet()
- Returns
- a pointer to the Isospin-multiplet of this PDG Code.
Definition at line 185 of file particletype.h.
◆ charge()
int32_t smash::ParticleType::charge |
( |
| ) |
const |
|
inline |
The charge of the particle.
The charge is calculated from the quark content (for hadrons) or basically tabulated; currently leptons, neutrinos and the standard model gauge bosons are known; unknown particles return a charge of 0.
- Returns
- charge of the particle
Definition at line 188 of file particletype.h.
◆ spin()
unsigned int smash::ParticleType::spin |
( |
| ) |
const |
|
inline |
- Todo:
- (oliiny): take care of spin for nuclei
- Returns
- twice the spin of a particle.
The code is good for hadrons, leptons and spin-1-bosons. It returns 2 (meaning spin=1) for the Higgs, though.
Definition at line 191 of file particletype.h.
◆ is_hadron()
bool smash::ParticleType::is_hadron |
( |
| ) |
const |
|
inline |
- Returns
- true if this is a baryon, antibaryon or meson.
Definition at line 194 of file particletype.h.
◆ is_lepton()
bool smash::ParticleType::is_lepton |
( |
| ) |
const |
|
inline |
◆ is_baryon()
bool smash::ParticleType::is_baryon |
( |
| ) |
const |
|
inline |
- Returns
- whether this PDG code identifies a baryon.
Definition at line 200 of file particletype.h.
◆ is_meson()
bool smash::ParticleType::is_meson |
( |
| ) |
const |
|
inline |
- Returns
- whether this PDG code identifies a meson.
Definition at line 203 of file particletype.h.
◆ baryon_number()
int smash::ParticleType::baryon_number |
( |
| ) |
const |
|
inline |
- Returns
- the baryon number of the particle.
Definition at line 206 of file particletype.h.
◆ strangeness()
int smash::ParticleType::strangeness |
( |
| ) |
const |
|
inline |
- Returns
- the net number of \(\bar s\) quarks.
For particles with one strange quark, -1 is returned.
Definition at line 209 of file particletype.h.
◆ is_nucleon()
bool smash::ParticleType::is_nucleon |
( |
| ) |
const |
|
inline |
- Returns
- whether this is a nucleon/anti-nucleon (p, n, -p, -n)
Definition at line 212 of file particletype.h.
◆ is_Delta()
bool smash::ParticleType::is_Delta |
( |
| ) |
const |
|
inline |
- Returns
- whether this is a Delta(1232) (with anti-delta)
Definition at line 215 of file particletype.h.
◆ is_rho()
bool smash::ParticleType::is_rho |
( |
| ) |
const |
|
inline |
- Returns
- whether this is a rho meson (rho+/rho0/rho-)
Definition at line 218 of file particletype.h.
◆ is_Nstar()
bool smash::ParticleType::is_Nstar |
( |
| ) |
const |
|
inline |
- Returns
- Is this a nucleon resonance (N*)?
Definition at line 221 of file particletype.h.
◆ is_Nstar1535()
bool smash::ParticleType::is_Nstar1535 |
( |
| ) |
const |
|
inline |
- Returns
- whether this is a N*(1535) (+/0)
Definition at line 227 of file particletype.h.
◆ is_Deltastar()
bool smash::ParticleType::is_Deltastar |
( |
| ) |
const |
|
inline |
- Returns
- Is this a Delta resonance (Delta*)?
Definition at line 230 of file particletype.h.
◆ is_stable()
bool smash::ParticleType::is_stable |
( |
| ) |
const |
|
inline |
◆ is_nucleus()
bool smash::ParticleType::is_nucleus |
( |
| ) |
const |
|
inline |
- Returns
- whether the particle is a nucleus
Definition at line 239 of file particletype.h.
◆ is_deuteron()
bool smash::ParticleType::is_deuteron |
( |
| ) |
const |
|
inline |
- Returns
- whether the particle is an (anti-)deuteron
Definition at line 242 of file particletype.h.
◆ is_dprime()
bool smash::ParticleType::is_dprime |
( |
| ) |
const |
|
inline |
- Returns
- whether the particle is an artificial d' resonance
Definition at line 245 of file particletype.h.
◆ min_mass_kinematic()
double smash::ParticleType::min_mass_kinematic |
( |
| ) |
const |
The minimum mass of the resonance that is kinematically allowed.
Calculate the minimum rest energy the resonance must have for any decay channel to be kinematically available. (In other words, find the smallest sum of final-state particle masses.)
- Returns
- The minimum mass that a particle of this type can assume, where at least one decay is possible.
Definition at line 360 of file particletype.cc.
366 for (
const auto &mode :
decay_modes().decay_mode_list()) {
◆ min_mass_spectral()
double smash::ParticleType::min_mass_spectral |
( |
| ) |
const |
The minimum mass of the resonance, where the spectral function is non-zero.
Calculate the the smallest mass where the spectral function still has a contribution. This value can be different from min_mass_kinematic, if the spectral function becomes zero at masses higher than min_mass_kinematic, since the width is put to zero due to the width_cutoff.
The distinction between it and min_mass_kinematic() might be necessary in edge cases, where a reaction is very close to the kinematic threshold or for optimizations.
- Returns
- The minimum mass that a particle of this type can assume, where the spectral function still has a non-zero value.
Definition at line 374 of file particletype.cc.
385 const double m_step = 0.01;
386 double right_bound_bis;
387 for (
unsigned int i = 0;; i++) {
394 const double precision = 1E-6;
395 double left_bound_bis = right_bound_bis - m_step;
396 while (right_bound_bis - left_bound_bis > precision) {
397 const double mid = (left_bound_bis + right_bound_bis) / 2.0;
399 right_bound_bis = mid;
401 left_bound_bis = mid;
◆ partial_width()
double smash::ParticleType::partial_width |
( |
const double |
m, |
|
|
const DecayBranch * |
mode |
|
) |
| const |
Get the mass-dependent partial decay width of a particle with mass m in a particular decay mode.
- Parameters
-
[in] | m | Invariant mass of the decaying particle. |
[in] | mode | Decay mode to consider. |
- Returns
- the partial width of this specific mode for this mass
Definition at line 419 of file particletype.cc.
421 if (m < mode->threshold()) {
424 double partial_width_at_pole =
width_at_pole() * mode->weight();
425 return mode->type().width(
mass(), partial_width_at_pole, m);
◆ total_width()
double smash::ParticleType::total_width |
( |
const double |
m | ) |
const |
Get the mass-dependent total width of a particle with mass m.
- Parameters
-
[in] | m | Invariant mass of the decaying particle. |
- Returns
- the total width for all modes for this mass
Definition at line 435 of file particletype.cc.
442 for (
unsigned int i = 0; i < modes.size(); i++) {
◆ wanted_decaymode()
Helper Function that containes the if-statement logic that decides if a decay mode is either a hadronic and dilepton decay mode.
- Parameters
-
[in] | t | type of decay. |
[in] | wh | enum that decides which decay modes are wanted. |
- Returns
- true if a decay branch is wanted and false if not.
Definition at line 460 of file particletype.cc.
462 const auto FinalTypes = t.particle_types();
470 (t.particle_number() == 2 &&
473 (t.particle_number() == 3 &&
480 (t.particle_number() == 2 &&
482 (t.particle_number() == 3 &&
487 throw std::runtime_error(
488 "Problem in selecting decaymodes in wanted_decaymode()");
◆ get_partial_widths()
Get all the mass-dependent partial decay widths of a particle with mass m.
This function needs to know the 4-momentum and the position of the decaying particle to calculate the square root of s of the final state particles and mass
- Parameters
-
[in] | p | 4-momentum of the decaying particle. |
[in] | x | position of the decaying particle. |
[in] | wh | enum that decides which decaymodes are returned. |
- Returns
- a list of process branches, whose weights correspond to the actual partial widths.
Definition at line 492 of file particletype.cc.
496 if (decay_mode_list.size() == 0 ||
503 FourVector UB = FourVector();
504 FourVector UI3 = FourVector();
512 DecayBranchList partial;
513 partial.reserve(decay_mode_list.size());
514 for (
unsigned int i = 0; i < decay_mode_list.size(); i++) {
516 const auto FinalTypes = decay_mode_list[i]->type().particle_types();
517 double scale_B = 0.0;
518 double scale_I3 = 0.0;
522 for (
const auto finaltype : FinalTypes) {
525 finaltype->isospin3_rel();
528 double sqrt_s = (
p + UB * scale_B + UI3 * scale_I3).abs();
530 const double w =
partial_width(sqrt_s, decay_mode_list[i].get());
534 make_unique<DecayBranch>(decay_mode_list[i]->type(), w));
◆ get_partial_width()
double smash::ParticleType::get_partial_width |
( |
const double |
m, |
|
|
const ParticleType & |
t_a, |
|
|
const ParticleType & |
t_b |
|
) |
| const |
Get the mass-dependent partial width of a resonance with mass m, decaying into two given daughter particles.
- Parameters
-
[in] | m | Invariant mass of the decaying resonance. |
[in] | t_a | Type of first daughter particle. |
[in] | t_b | Type of second daughter particle. |
- Returns
- the partial width for this mass and this specific decay channel
Definition at line 541 of file particletype.cc.
548 for (
const auto &mode : decaymodes) {
549 double partial_width_at_pole =
width_at_pole() * mode->weight();
550 const ParticleTypePtrList l = {&t_a, &t_b};
551 if (mode->type().has_particles(l)) {
552 w += mode->type().width(
mass(), partial_width_at_pole, m);
◆ get_partial_in_width()
double smash::ParticleType::get_partial_in_width |
( |
const double |
m, |
|
|
const ParticleData & |
p_a, |
|
|
const ParticleData & |
p_b |
|
) |
| const |
Get the mass-dependent partial in-width of a resonance with mass m, decaying into two given daughter particles.
For stable daughter particles, the in-width equals the 'normal' partial decay width (i.e. the 'out-width').
- Parameters
-
[in] | m | Invariant mass of the decaying resonance. |
[in] | p_a | First daughter particle. |
[in] | p_b | Second daughter particle. |
- Returns
- the partial in-width for this mass and this specific decay channel
Definition at line 558 of file particletype.cc.
566 for (
const auto &mode : decaymodes) {
567 double partial_width_at_pole =
width_at_pole() * mode->weight();
568 const ParticleTypePtrList l = {&p_a.type(), &p_b.type()};
569 if (mode->type().has_particles(l)) {
570 w += mode->type().in_width(
mass(), partial_width_at_pole, m,
571 p_a.effective_mass(), p_b.effective_mass());
◆ spectral_function()
double smash::ParticleType::spectral_function |
( |
double |
m | ) |
const |
Full spectral function \( A(m) = \frac{2}{\pi} N \frac{m^2\Gamma(m)}{(m^2-m_0^2)^2+(m\Gamma(m))^2} \) of the resonance (relativistic Breit-Wigner distribution with mass-dependent width, where N is a normalization factor).
- Parameters
-
[in] | m | Actual off-shell mass of the resonance, where the spectral function is to be evaluated. |
- Returns
- the value of the spectral function for this mass
- Note
- The normalization factor N ensures that the spectral function is normalized to unity.
Definition at line 577 of file particletype.cc.
583 const double m_pole =
mass();
588 const double tanx = std::tan(x);
589 const double m_x = m_pole + width * tanx;
590 const double jacobian = width * (1.0 + tanx * tanx);
◆ spectral_function_no_norm()
double smash::ParticleType::spectral_function_no_norm |
( |
double |
m | ) |
const |
Full spectral function without normalization factor.
- See also
- spectral_function
- Parameters
-
[in] | m | Actual off-shell mass of the resonance, where the spectral function is to be evaluated. |
- Returns
- the value of the non-normalized spectral function for this mass
Definition at line 597 of file particletype.cc.
◆ spectral_function_const_width()
double smash::ParticleType::spectral_function_const_width |
( |
double |
m | ) |
const |
- Todo:
- unused The spectral function with a constant width (= width at pole).
It is guaranteed to be normalized to 1, when integrated from 0 to inf.
Definition at line 607 of file particletype.cc.
◆ spectral_function_simple()
double smash::ParticleType::spectral_function_simple |
( |
double |
m | ) |
const |
This one is the most simple form of the spectral function, using a Cauchy distribution (non-relativistic Breit-Wigner with constant width).
It can be integrated analytically, and is normalized to 1 when integrated from -inf to inf.
- Parameters
-
[in] | m | Actual off-shell mass of the resonance, where the spectral function is to be evaluated. |
- Returns
- the Cauchy spectral function at mass m
Definition at line 617 of file particletype.cc.
◆ sample_resonance_mass()
double smash::ParticleType::sample_resonance_mass |
( |
const double |
mass_stable, |
|
|
const double |
cms_energy, |
|
|
int |
L = 0 |
|
) |
| const |
Resonance mass sampling for 2-particle final state with one resonance (type given by 'this') and one stable particle.
- Parameters
-
[in] | mass_stable | Mass of the stable particle. |
[in] | cms_energy | center-of-mass energy of the 2-particle final state. |
[in] | L | relative angular momentum of the final-state particles |
- Returns
- The mass of the resonance particle.
Definition at line 622 of file particletype.cc.
627 const double max_mass = std::nextafter(cms_energy - mass_stable, 0.);
633 const double pcm_max =
pCM(cms_energy, mass_stable, min_mass);
640 const double sf_ratio_max =
644 double mass_res, val;
648 const double max = blw_max * q_max;
655 const double pcm =
pCM(cms_energy, mass_stable, mass_res);
666 "maximum is being increased in sample_resonance_mass: ",
667 this->max_factor1_,
" ", val / max,
" ", this->
pdgcode(),
" ",
668 mass_stable,
" ", cms_energy,
" ", mass_res);
669 this->max_factor1_ *= val / max;
◆ sample_resonance_masses()
std::pair< double, double > smash::ParticleType::sample_resonance_masses |
( |
const ParticleType & |
t2, |
|
|
const double |
cms_energy, |
|
|
int |
L = 0 |
|
) |
| const |
Resonance mass sampling for 2-particle final state with two resonances.
- Parameters
-
[in] | t2 | Type of the second resonance (the first resonance is given by 'this'). |
[in] | cms_energy | center-of-mass energy of the 2-particle final state. |
[in] | L | relative angular momentum of the final-state particles |
- Returns
- The masses of the resonance particles.
Definition at line 679 of file particletype.cc.
684 const double max_mass_1 =
685 std::nextafter(cms_energy - t2.min_mass_spectral(), 0.);
686 const double max_mass_2 =
687 std::nextafter(cms_energy - t1.min_mass_spectral(), 0.);
689 const double pcm_max =
690 pCM(cms_energy, t1.min_mass_spectral(), t2.min_mass_spectral());
693 double mass_1, mass_2, val;
697 const double max = blw_max * t1.max_factor2_;
702 t1.min_mass_spectral(), max_mass_1);
704 t2.min_mass_spectral(), max_mass_2);
706 const double pcm =
pCM(cms_energy, mass_1, mass_2);
710 t1.spectral_function(mass_1) / t1.spectral_function_simple(mass_1);
712 t2.spectral_function(mass_2) / t2.spectral_function_simple(mass_2);
718 "maximum is being increased in sample_resonance_masses: ",
719 t1.max_factor2_,
" ", val / max,
" ", t1.pdgcode(),
" ", t2.pdgcode(),
720 " ", cms_energy,
" ", mass_1,
" ", mass_2);
721 t1.max_factor2_ *= val / max;
727 return {mass_1, mass_2};
◆ dump_width_and_spectral_function()
void smash::ParticleType::dump_width_and_spectral_function |
( |
| ) |
const |
Prints out width and spectral function versus mass to the standard output.
This is useful for debugging and analysis.
- Exceptions
-
if | the particle type is stable |
Definition at line 730 of file particletype.cc.
732 std::stringstream err;
733 err <<
"Particle " << *
this <<
" is stable, so it makes no"
734 <<
" sense to print its spectral function, etc.";
735 throw std::runtime_error(err.str());
738 double rightmost_pole = 0.0;
740 for (
const auto &mode : decaymodes) {
741 double pole_mass_sum = 0.0;
742 for (
const ParticleTypePtr
p : mode->type().particle_types()) {
743 pole_mass_sum +=
p->mass();
745 if (pole_mass_sum > rightmost_pole) {
746 rightmost_pole = pole_mass_sum;
750 std::cout <<
"# mass m[GeV], width w(m) [GeV],"
751 <<
" spectral function(m^2)*m [GeV^-1] of " << *
this << std::endl;
752 constexpr
double m_step = 0.02;
756 constexpr
double spectral_function_threshold = 8.e-3;
757 std::cout << std::fixed << std::setprecision(5);
758 for (
unsigned int i = 0;; i++) {
759 const double m = m_min + m_step * i;
761 if (m > rightmost_pole * 2 && sf < spectral_function_threshold) {
764 std::cout << m <<
" " << w <<
" " << sf << std::endl;
◆ list_all()
const ParticleTypeList & smash::ParticleType::list_all |
( |
| ) |
|
|
static |
- Returns
- a list of all ParticleType objects.
- Note
- This list is currently sorted, but do not rely on it.
-
It might make sense to inline this function to optimize runtime performance.
Definition at line 57 of file particletype.cc.
◆ list_nucleons()
ParticleTypePtrList & smash::ParticleType::list_nucleons |
( |
| ) |
|
|
static |
- Returns
- a list of all nucleons (i.e. proton and neutron).
Definition at line 75 of file particletype.cc.
◆ list_anti_nucleons()
ParticleTypePtrList & smash::ParticleType::list_anti_nucleons |
( |
| ) |
|
|
static |
- Returns
- a list of all anti-nucleons (i.e. anti-proton and anti-neutron).
Definition at line 77 of file particletype.cc.
◆ list_Deltas()
ParticleTypePtrList & smash::ParticleType::list_Deltas |
( |
| ) |
|
|
static |
- Returns
- a list of the Delta(1232) baryons (i.e. all four charge states).
Definition at line 81 of file particletype.cc.
◆ list_anti_Deltas()
ParticleTypePtrList & smash::ParticleType::list_anti_Deltas |
( |
| ) |
|
|
static |
- Returns
- a list of the anti-Delta(1232) baryons (i.e. all four charge states).
Definition at line 83 of file particletype.cc.
◆ list_baryon_resonances()
ParticleTypePtrList & smash::ParticleType::list_baryon_resonances |
( |
| ) |
|
|
static |
- Returns
- a list of all baryon resonances, i.e. unstable baryons (not including antibaryons).
Definition at line 87 of file particletype.cc.
◆ list_light_nuclei()
ParticleTypePtrList & smash::ParticleType::list_light_nuclei |
( |
| ) |
|
|
static |
- Returns
- a list of all light nuclei from SMASH particle list. Nucleons are not included into light nuclei by convention.
Definition at line 91 of file particletype.cc.
◆ try_find()
Returns the ParticleTypePtr for the given pdgcode
.
If the particle type is not found, an invalid ParticleTypePtr is returned. You can convert a ParticleTypePtr to a bool to check whether it is valid.
- Note
- The complexity of the search is \(\mathcal O(\log N)\). Therefore, do not use this function except for user input that selects a particle type. All other internal references for a particle type should use ParticleTypePtr instead.
- Parameters
-
[in] | pdgcode | the unique pdg code to try to find |
- Returns
- the ParticleTypePtr that corresponds to this pdg code, or an invalid pointer
Definition at line 95 of file particletype.cc.
96 const auto found = std::lower_bound(
98 [](
const ParticleType &l,
const PdgCode &r) {
return l.pdgcode() < r; });
◆ find()
Returns the ParticleType object for the given pdgcode
.
If the particle is not found, a PdgNotFoundFailure is thrown.
- Note
- The complexity of the search is \(\mathcal O(\log N)\). Therefore, do not use this function except for user input that selects a particle type. All other internal references for a particle type should use ParticleTypePtr instead.
- Parameters
-
[in] | pdgcode | the unique pdg code to try to find |
- Returns
- the ParticleTypePtr that corresponds to this pdg code
- Exceptions
-
Definition at line 105 of file particletype.cc.
108 throw PdgNotFoundFailure(
"PDG code " +
pdgcode.
string() +
" not found!");
◆ exists() [1/2]
bool smash::ParticleType::exists |
( |
PdgCode |
pdgcode | ) |
|
|
static |
- Parameters
-
- Returns
- whether the ParticleType with the given
pdgcode
exists.
- Note
- The complexity of the search is \(\mathcal O(\log N)\).
Definition at line 113 of file particletype.cc.
◆ exists() [2/2]
bool smash::ParticleType::exists |
( |
const std::string & |
name | ) |
|
|
static |
- Parameters
-
[in] | name | the name to look for |
- Returns
- whether the ParticleType with the given
name
exists.
- Note
- The complexity of the search is \(\mathcal O(N)\).
Definition at line 118 of file particletype.cc.
◆ create_type_list()
void smash::ParticleType::create_type_list |
( |
const std::string & |
particles | ) |
|
|
static |
Initialize the global ParticleType list (list_all) from the given input data.
This function must only be called once (will fail on second invocation).
- Parameters
-
[in] | particles | A string that contains the definition of ParticleTypes to be created. |
- Exceptions
-
LoadFailure | if a line in the particle file could not be read, or if there are duplicates in it |
runtime_error | if the mass of of nucleons, kaons and deltas are different from the hardcoded masses, or if this function is called more than once |
Definition at line 205 of file particletype.cc.
206 static ParticleTypeList type_list;
210 std::istringstream lineinput(line.text);
213 std::string parity_string;
214 std::vector<std::string> pdgcode_strings;
216 pdgcode_strings.reserve(4);
217 lineinput >>
name >>
mass >> width >> parity_string;
220 if (parity_string ==
"+") {
222 }
else if (parity_string ==
"-") {
227 if (lineinput.fail() || fail) {
229 "While loading the ParticleType data:\nFailed to convert the input "
230 "string to the expected data types.",
234 while (!lineinput.eof()) {
235 pdgcode_strings.push_back(
"");
236 lineinput >> pdgcode_strings.back();
237 if (lineinput.fail()) {
239 "While loading the ParticleType data:\nFailed to convert the input "
240 "string to the expected data types.",
244 if (pdgcode_strings.size() < 1) {
246 "While loading the ParticleType data:\nFailed to convert the input "
247 "string due to missing PDG code.",
251 pdgcode.resize(pdgcode_strings.size());
252 std::transform(pdgcode_strings.begin(), pdgcode_strings.end(),
254 [](
const std::string &s) {
return PdgCode(s); });
260 throw std::runtime_error(
261 "Nucleon mass in input file different from 0.938");
264 throw std::runtime_error(
"Kaon mass in input file different from 0.494");
267 throw std::runtime_error(
"Delta mass in input file different from 1.232");
270 throw std::runtime_error(
"d mass in input file different from 1.8756");
274 for (
size_t i = 0; i <
pdgcode.size(); i++) {
275 std::string full_name =
name;
282 <<
"Setting particle type: " << type_list.back();
287 const auto anti_parity = (anti.spin() % 2 == 0) ?
parity : -
parity;
289 type_list.emplace_back(full_name,
mass, width, anti_parity, anti);
291 <<
"Setting antiparticle type: " << type_list.back();
295 type_list.shrink_to_fit();
298 std::sort(type_list.begin(), type_list.end());
301 PdgCode prev_pdg = 0;
302 for (
const auto &t : type_list) {
303 if (t.pdgcode() == prev_pdg) {
304 throw ParticleType::LoadFailure(
"Duplicate PdgCode in particles.txt: " +
305 t.pdgcode().string());
307 prev_pdg = t.pdgcode();
311 throw std::runtime_error(
"Error: Type list was already built!");
318 for (
const auto &t : type_list) {
322 for (
auto &t : type_list) {
345 if (type_resonance.is_stable() ||
346 type_resonance.pdgcode().baryon_number() != 1) {
354 if (type.is_nucleus()) {
◆ operator==()
bool smash::ParticleType::operator== |
( |
const ParticleType & |
rhs | ) |
const |
|
inline |
◆ operator!=()
bool smash::ParticleType::operator!= |
( |
const ParticleType & |
rhs | ) |
const |
|
inline |
◆ operator<()
bool smash::ParticleType::operator< |
( |
const ParticleType & |
rhs | ) |
const |
|
inline |
"Less than" operator for sorting the ParticleType list (by PDG code)
- Parameters
-
- Returns
- whether the PDG code of rhs is larger than the one from *this
Definition at line 544 of file particletype.h.
545 return pdgcode() < rhs.pdgcode();
◆ check_consistency()
void smash::ParticleType::check_consistency |
( |
| ) |
|
|
static |
- Exceptions
-
runtime_error | if unstable particles have no decay modes |
Note that the particles and decay modes have to be initialized, otherwise calling this is undefined behavior.
Definition at line 451 of file particletype.cc.
453 if (!ptype.is_stable() && ptype.decay_modes().is_empty()) {
454 throw std::runtime_error(
"Unstable particle " + ptype.name() +
455 " has no decay chanels!");
◆ operator&()
Returns an object that acts like a pointer, except that it requires only 2 bytes and inhibits pointer arithmetics.
This is an optimization for creating references to ParticleType objects. With a normal pointer you would require the original object to stay in its place in memory, as otherwise the pointer would dangle. The ParticleTypePtr does not have this problem as it only stores the index of the ParticleType in the global vector of ParticleType objects.
In addition to returning a more efficient reference type, the overload of operator& effectively inhibits passing ParticleType objects by pointer. This is an intended restriction since ParticleType objects should only be passed by const-ref. (You can now pass by ParticleTypePtr instead, but please prefer not to. It might be useful when you want to store a reference inside the function anyway, but that can just as well be done with a const-ref parameter. A ParticleTypePtr can be invalid, a const-ref is always a valid reference semantically.)
- Pre-condition:
- The operator expects that the ParticleType object is stored in the vector returned by ParticleType::list_all. Therefore, never create new ParticleType objects (that includes copies and moves)!
- Note on distributed execution:
- At some point we might want to have several copies of the ParticleType vector - on different machines or NUMA nodes. In that case ParticleType::list_all will return the local vector. This operator will continue to work as expected as long as the ParticleType object is an entry of this local vector. The ParticleTypePtr can then be used to communicate a ParticleType over node / NUMA boundaries (an actual pointer would not work, though).
- Note
- It might make sense to inline this function to optimize runtime performance.
- Returns
- A pointer-like object referencing this ParticleType object.
- See also
- ParticleTypePtr
Definition at line 63 of file particletype.cc.
65 const auto offset =
this - std::addressof(
list_all()[0]);
70 assert(offset >= 0 && offset < 0xffff);
72 return ParticleTypePtr(static_cast<uint16_t>(offset));
◆ width_cutoff
constexpr double smash::ParticleType::width_cutoff = 1e-5 |
|
staticconstexpr |
Decay width cutoff for considering a particle as stable.
We currently regard a particle type as stable if its on-shell width is less than 10 keV. The cutoff is chosen such that the η and the η' are stable.
Definition at line 107 of file particletype.h.
◆ name_
std::string smash::ParticleType::name_ |
|
private |
◆ mass_
double smash::ParticleType::mass_ |
|
private |
◆ width_
double smash::ParticleType::width_ |
|
private |
◆ parity_
Parity smash::ParticleType::parity_ |
|
private |
◆ pdgcode_
PdgCode smash::ParticleType::pdgcode_ |
|
private |
◆ min_mass_kinematic_
double smash::ParticleType::min_mass_kinematic_ |
|
mutableprivate |
minimum kinematically allowed mass of the particle Mutable, because it is initialized at first call of minimum mass function, so it's logically const, but not physically const, which is a classical case for using mutable.
Definition at line 621 of file particletype.h.
◆ min_mass_spectral_
double smash::ParticleType::min_mass_spectral_ |
|
mutableprivate |
minimum mass, where the spectral function is non-zero Mutable, because it is initialized at first call of minimum mass function, so it's logically const, but not physically const, which is a classical case for using mutable.
Definition at line 628 of file particletype.h.
◆ norm_factor_
double smash::ParticleType::norm_factor_ = -1. |
|
mutableprivate |
This normalization factor ensures that the spectral function is normalized to unity, when integrated over its full domain.
Definition at line 631 of file particletype.h.
◆ charge_
int32_t smash::ParticleType::charge_ |
|
private |
Charge of the particle; filled automatically from pdgcode_.
Definition at line 633 of file particletype.h.
◆ isospin_
int smash::ParticleType::isospin_ |
|
mutableprivate |
Isospin of the particle; filled automatically from pdgcode_.
Definition at line 635 of file particletype.h.
◆ I3_
int smash::ParticleType::I3_ |
|
private |
Isospin projection of the particle; filled automatically from pdgcode_.
Definition at line 637 of file particletype.h.
◆ iso_multiplet_
Container for the isospin multiplet information.
Definition at line 640 of file particletype.h.
◆ max_factor1_
double smash::ParticleType::max_factor1_ = 1. |
|
mutableprivate |
Maximum factor for single-res mass sampling, cf. sample_resonance_mass.
Definition at line 643 of file particletype.h.
◆ max_factor2_
double smash::ParticleType::max_factor2_ = 1. |
|
mutableprivate |
Maximum factor for double-res mass sampling, cf. sample_resonance_masses.
Definition at line 645 of file particletype.h.
The documentation for this class was generated from the following files:
double breit_wigner(double m, double pole, double width)
Returns a relativistic Breit-Wigner distribution.
build_vector_< Line > line_parser(const std::string &input)
Helper function for parsing particles.txt and decaymodes.txt.
PdgCode get_antiparticle() const
Construct the antiparticle to a given PDG code.
double spectral_function_no_norm(double m) const
Full spectral function without normalization factor.
T cauchy(T pole, T width, T min, T max)
Draws a random number from a Cauchy distribution (sometimes also called Lorentz or non-relativistic B...
int isospin() const
Returns twice the total isospin of the multiplet.
double width_
width of the particle
static const IsoParticleType & find(const std::string &name)
Returns the IsoParticleType object for the given name.
double width_at_pole() const
double min_mass_spectral() const
The minimum mass of the resonance, where the spectral function is non-zero.
bool wanted_decaymode(const DecayType &t, WhichDecaymodes wh) const
Helper Function that containes the if-statement logic that decides if a decay mode is either a hadron...
std::string build_error_string(std::string message, const Line &line)
Builds a meaningful error message.
ParticleTypePtrList anti_deltas_list
Global pointer to the Particle Type list of anti-deltas.
int baryon_number() const
const DecayBranchList & decay_mode_list() const
void ensure_all_read(std::istream &input, const Line &line)
Makes sure that nothing is left to read from this line.
static std::pair< double, int > force_scale(const ParticleType &data)
Evaluates the scaling factor of the forces acting on the particles.
int32_t charge_
Charge of the particle; filled automatically from pdgcode_.
unsigned int spin() const
#define unlikely(x)
Tell the branch predictor that this expression is likely false.
ParticleTypePtrList deltas_list
Global pointer to the Particle Type list of deltas.
const ParticleTypeList * all_particle_types
Global pointer to the Particle Type list.
int antiparticle_sign() const
double blatt_weisskopf_sqr(const double p_ab, const int L)
double spectral_function(double m) const
Full spectral function of the resonance (relativistic Breit-Wigner distribution with mass-dependent ...
int isospin_
Isospin of the particle; filled automatically from pdgcode_.
constexpr double nucleon_mass
Nucleon mass in GeV.
std::string name_
name of the particle
double min_mass_kinematic_
minimum kinematically allowed mass of the particle Mutable, because it is initialized at first call o...
std::array< einhard::Logger<>, std::tuple_size< LogArea::AreaTuple >::value > logg
An array that stores all pre-configured Logger objects.
static constexpr int LResonances
T pCM(const T sqrts, const T mass_a, const T mass_b) noexcept
double total_width(const double m) const
Get the mass-dependent total width of a particle with mass m.
constexpr double really_small
Numerical error tolerance.
static const ParticleType & find(PdgCode pdgcode)
Returns the ParticleType object for the given pdgcode.
int32_t get_decimal() const
bool has_antiparticle() const
const DecayModes & decay_modes() const
RectangularLattice< FourVector > * UI3_lat_pointer
Pointer to the symmmetry potential on the lattice.
std::string string() const
double norm_factor_
This normalization factor ensures that the spectral function is normalized to unity,...
IsoParticleType * iso_multiplet_
Container for the isospin multiplet information.
static Integrator integrate
static constexpr int LParticleType
int isospin() const
Returns twice the isospin vector length .
int I3_
Isospin projection of the particle; filled automatically from pdgcode_.
Ignore dilepton decay modes widths.
int charge() const
The charge of the particle.
RectangularLattice< FourVector > * UB_lat_pointer
Pointer to the skyrme potential on the lattice.
static bool exists(const std::string &name)
Returns whether the ParticleType with the given pdgcode exists.
double breit_wigner_nonrel(double m, double pole, double width)
Returns a non-relativistic Breit-Wigner distribution, which is essentially a Cauchy distribution with...
ParticleTypePtrList nucleons_list
Global pointer to the Particle Type list of nucleons.
double isospin3_rel() const
double spectral_function_simple(double m) const
This one is the most simple form of the spectral function, using a Cauchy distribution (non-relativis...
const std::string & name() const
bool almost_equal(const N x, const N y)
Checks two numbers for relative approximate equality.
ParticleTypePtrList baryon_resonances_list
Global pointer to the Particle Type list of baryon resonances.
ParticleTypePtrList light_nuclei_list
Global pointer to the Particle Type list of light nuclei.
PdgCode pdgcode_
PDG Code of the particle.
static const ParticleTypePtr try_find(PdgCode pdgcode)
Returns the ParticleTypePtr for the given pdgcode.
ParticleType(std::string n, double m, double w, Parity p, PdgCode id)
Creates a fully initialized ParticleType object.
Parity parity_
Parity of the particle.
static std::string chargestr(int charge)
Construct a charge string, given the charge as integer.
static std::string antiname(const std::string &name, PdgCode code)
Construct an antiparticle name-string from the given name-string for the particle and its PDG code.
double partial_width(const double m, const DecayBranch *mode) const
Get the mass-dependent partial decay width of a particle with mass m in a particular decay mode.
bool has_lepton_pair(const PdgCode pdg1, const PdgCode pdg2, const PdgCode pdg3)
double min_mass_spectral_
minimum mass, where the spectral function is non-zero Mutable, because it is initialized at first cal...
constexpr double delta_mass
Delta mass in GeV.
double max_factor1_
Maximum factor for single-res mass sampling, cf. sample_resonance_mass.
Parity
Represent the parity of a particle type.
bool has_antiparticle() const
double min_mass_kinematic() const
The minimum mass of the resonance that is kinematically allowed.
Potentials * pot_pointer
Pointer to a Potential class.
double mass_
pole mass of the particle
constexpr double deuteron_mass
Deuteron mass in GeV.
Only return dilepton decays widths.
constexpr double kaon_mass
Kaon mass in GeV.
ParticleTypePtrList anti_nucs_list
Global pointer to the Particle Type list of anti-nucleons.
bool is_dilepton(const PdgCode pdg1, const PdgCode pdg2)
bool is_Nstar1535() const
ParticleTypePtrList get_states() const
Returns list of states that form part of the multiplet.
static constexpr double width_cutoff
Decay width cutoff for considering a particle as stable.
static const ParticleTypeList & list_all()
static void create_multiplet(const ParticleType &type)
Add a new multiplet to the global list of IsoParticleTypes, which contains type.
int32_t charge() const
The charge of the particle.