Version: SMASH-3.3
particledata.h
Go to the documentation of this file.
1 /*
2  * Copyright (c) 2012-2020,2022-2025
3  * SMASH Team
4  *
5  * GNU General Public License (GPLv3 or later)
6  */
7 #ifndef SRC_INCLUDE_SMASH_PARTICLEDATA_H_
8 #define SRC_INCLUDE_SMASH_PARTICLEDATA_H_
9 
10 #include <limits>
11 #include <utility>
12 
13 #include "forwarddeclarations.h"
14 #include "fourvector.h"
15 #include "particletype.h"
16 #include "pdgcode.h"
17 #include "processbranch.h"
18 
19 namespace smash {
20 
21 enum class BelongsTo : uint8_t {
22  Nothing = 0,
23  Projectile = 1,
24  Target = 2,
25 };
26 
31 struct HistoryData {
35  int32_t id_process = 0;
44  double time_last_collision = smash_NaN<double>;
46  PdgCode p1 = 0x0;
48  PdgCode p2 = 0x0;
49 };
50 
59 class ParticleData {
60  public:
70  explicit ParticleData(const ParticleType &particle_type, int unique_id = -1)
71  : id_(unique_id), type_(&particle_type) {}
72 
77  int32_t id() const { return id_; }
82  void set_id(int i) { id_ = i; }
83 
88  PdgCode pdgcode() const { return type_->pdgcode(); }
89 
90  // Convenience accessors to PdgCode:
92  bool is_hadron() const { return type_->is_hadron(); }
93 
95  bool is_baryon() const { return pdgcode().is_baryon(); }
96 
98  bool is_nucleus() const { return pdgcode().is_nucleus(); }
99 
101  bool is_rho() const { return type_->is_rho(); }
102 
104  bool is_proton() const { return pdgcode().is_proton(); }
105 
107  bool is_neutron() const { return pdgcode().is_neutron(); }
108 
110  bool is_pion() const { return pdgcode().is_pion(); }
111 
113  bool is_sigmastar() const { return pdgcode().is_Sigmastar(); }
114 
119  double pole_mass() const { return type_->mass(); }
127  double effective_mass() const;
132  const ParticleType &type() const { return *type_; }
133 
138  uint32_t id_process() const { return history_.id_process; }
143  HistoryData get_history() const { return history_; }
151  void set_history(HistoryData &&history) { history_ = std::move(history); }
152 
164  void set_history(int ncoll, uint32_t pid, ProcessType pt,
165  double time_last_coll, const ParticleList &plist);
166 
171  const FourVector &momentum() const { return momentum_; }
172 
177  void set_4momentum(const FourVector &momentum_vector) {
178  momentum_ = momentum_vector;
179  }
180 
187  void set_4momentum(double mass, const ThreeVector &mom) {
188  momentum_ = FourVector(std::sqrt(mass * mass + mom * mom), mom);
189  }
190 
199  void set_4momentum(double mass, double px, double py, double pz) {
200  momentum_ = FourVector(std::sqrt(mass * mass + px * px + py * py + pz * pz),
201  px, py, pz);
202  }
209  void set_3momentum(const ThreeVector &mom) {
210  momentum_ = FourVector(momentum_.x0(), mom);
211  }
212 
217  const FourVector &position() const { return position_; }
222  void set_4position(const FourVector &pos) { position_ = pos; }
229  void set_3position(const ThreeVector &pos) {
230  position_ = FourVector(position_.x0(), pos);
231  }
232 
237  ParticleData translated(const ThreeVector &delta) const {
238  ParticleData p = *this;
239  p.position_[1] += delta[0];
240  p.position_[2] += delta[1];
241  p.position_[3] += delta[2];
242  return p;
243  }
244 
249  double formation_time() const { return formation_time_; }
255  double begin_formation_time() const { return begin_formation_time_; }
256 
264  void set_formation_time(double form_time) {
265  formation_time_ = form_time;
266  // cross section scaling factor will be a step function in time
267  begin_formation_time_ = form_time;
268  // if time of the last collision is NaN set it to the formation time
269  if (std::isnan(history_.time_last_collision)) {
270  history_.time_last_collision = form_time;
271  }
272  }
282  void set_slow_formation_times(double begin_form_time, double form_time) {
283  begin_formation_time_ = begin_form_time;
284  formation_time_ = form_time;
285  if (std::isnan(history_.time_last_collision)) {
286  history_.time_last_collision = form_time;
287  }
288  }
289 
298  const double &initial_xsec_scaling_factor() const {
300  }
313  void set_cross_section_scaling_factor(const double &xsec_scal) {
314  initial_xsec_scaling_factor_ = xsec_scal;
315  }
316 
321  ThreeVector velocity() const { return momentum_.velocity(); }
322 
335  double inverse_gamma() const {
336  return std::sqrt(1. - momentum_.sqr3() / (momentum_.x0() * momentum_.x0()));
337  }
338 
343  void boost(const ThreeVector &v) {
346  }
347 
352  void boost_momentum(const ThreeVector &v) {
354  }
360  int spin() const { return pdgcode().spin(); }
365  const FourVector &spin_vector() const { return spin_vector_; }
375  void set_spin_vector(const FourVector &s) { spin_vector_ = s; }
381  void set_spin_vector_component(int index, double value) {
382  if (index < 0 || index > 3) {
383  throw std::out_of_range("Invalid spin vector component index");
384  }
385  spin_vector_[index] = value;
386  }
387 
395  void set_belongs_to(BelongsTo label) { belongs_to_ = label; }
397  BelongsTo belongs_to() const { return belongs_to_; }
398 
400  void fluidize() { core_ = true; }
402  bool is_core() const { return core_; }
404  double hyperbolic_time() const { return position_.tau(); }
406  double spatial_rapidity() const { return position_.eta(); }
408  double transverse_mass() const { return momentum_.tau(); }
410  double rapidity() const { return momentum_.eta(); }
411 
415  void set_perturbative_weight(const double weight) {
416  perturbative_weight_ = weight;
417  }
421  double perturbative_weight() const { return perturbative_weight_; }
422 
428  bool operator==(const ParticleData &a) const { return this->id_ == a.id_; }
434  bool operator<(const ParticleData &a) const { return this->id_ < a.id_; }
435 
441  bool operator==(int id_a) const { return this->id_ == id_a; }
446  bool operator<(int id_a) const { return this->id_ < id_a; }
447 
458  ParticleData(const ParticleType &ptype, int uid, int index)
459  : id_(uid), index_(index), type_(&ptype) {}
460 
468  double xsec_scaling_factor(double delta_time = 0.) const;
469 
471  static double formation_power_;
472 
473  private:
474  friend class Particles;
476  ParticleData() = default;
477 
484  void copy_to(ParticleData &dst) const {
485  dst.history_ = history_;
486  dst.momentum_ = momentum_;
487  dst.position_ = position_;
492  dst.belongs_to_ = belongs_to_;
493  dst.core_ = core_;
495  }
496 
505  int32_t id_ = -1;
506 
515  unsigned index_ = std::numeric_limits<unsigned>::max();
516 
522 
523  // this leaves us two Bytes padding to use for "free"
524  static_assert(sizeof(ParticleTypePtr) == 2, "");
525  // make sure we don't exceed that space
526  static_assert(sizeof(bool) <= 2, "");
535  bool hole_ = false;
536 
538  // A particle cannot be un-core, and any children it produces inherits
539  // this trait.
540  bool core_ = false;
541 
551  FourVector spin_vector_ = FourVector(smash_NaN<double>, smash_NaN<double>,
552  smash_NaN<double>, smash_NaN<double>);
556  double formation_time_ = smash_NaN<double>;
558  double begin_formation_time_ = smash_NaN<double>;
565  double perturbative_weight_ = 1.0;
570 };
571 
576 std::ostream &operator<<(std::ostream &s, const ParticleData &p);
577 
583 std::ostream &operator<<(std::ostream &out, const ParticleList &particle_list);
584 
592  const ParticleList &list;
593 };
599 inline PrintParticleListDetailed detailed(const ParticleList &list) {
600  return {list};
601 }
602 
609 std::ostream &operator<<(std::ostream &out,
610  const PrintParticleListDetailed &particle_list);
611 
642  PdgCode pdgcode, double mass, const FourVector &four_position,
643  const FourVector &four_momentum, int log_area, bool &mass_warning,
644  bool &on_shell_warning);
645 
660 bool are_particles_identical_at_given_time(const ParticleData &p1,
661  const ParticleData &p2, double time);
662 
663 } // namespace smash
664 
665 #endif // SRC_INCLUDE_SMASH_PARTICLEDATA_H_
The FourVector class holds relevant values in Minkowski spacetime with (+, −, −, −) metric signature.
Definition: fourvector.h:33
FourVector lorentz_boost(const ThreeVector &v) const
Returns the FourVector boosted with velocity v.
Definition: fourvector.cc:17
double sqr3() const
calculate the square of the spatial three-vector
Definition: fourvector.h:474
double x0() const
Definition: fourvector.h:313
ThreeVector velocity() const
Get the velocity (3-vector divided by zero component).
Definition: fourvector.h:333
double eta() const
calculate the space-time rapidity from the given four vector
Definition: fourvector.h:482
double tau() const
calculate the proper time from the given four vector
Definition: fourvector.h:478
ParticleData contains the dynamic information of a certain particle.
Definition: particledata.h:59
double formation_time_
Formation time at which the particle is fully formed given as an absolute value in the computational ...
Definition: particledata.h:556
PdgCode pdgcode() const
Get the pdgcode of the particle.
Definition: particledata.h:88
void set_history(HistoryData &&history)
Set history_ from rvalue reference.
Definition: particledata.h:151
unsigned index_
Internal index in the Particles list.
Definition: particledata.h:515
void set_4momentum(const FourVector &momentum_vector)
Set the particle's 4-momentum directly.
Definition: particledata.h:177
bool core_
If the particle is part of a pseudofluid.
Definition: particledata.h:540
bool is_proton() const
Definition: particledata.h:104
const ParticleType & type() const
Get the type of the particle.
Definition: particledata.h:132
BelongsTo belongs_to_
is it part of projectile or target nuclei?
Definition: particledata.h:569
void set_3position(const ThreeVector &pos)
Set particle's 3-position.
Definition: particledata.h:229
double inverse_gamma() const
Get the inverse of the gamma factor from the current velocity of the particle.
Definition: particledata.h:335
bool is_pion() const
Definition: particledata.h:110
ParticleData translated(const ThreeVector &delta) const
Translate the particle position.
Definition: particledata.h:237
int32_t id_
Each particle has a unique identifier.
Definition: particledata.h:505
double begin_formation_time_
time when the cross section scaling factor starts to increase to 1
Definition: particledata.h:558
void set_4position(const FourVector &pos)
Set the particle's 4-position directly.
Definition: particledata.h:222
bool is_nucleus() const
Definition: particledata.h:98
double xsec_scaling_factor(double delta_time=0.) const
Return the cross section scaling factor at a given time.
double hyperbolic_time() const
Particle (hyperbolic time)
Definition: particledata.h:404
bool operator<(int id_a) const
Check whether the particle's id is smaller than the given id.
Definition: particledata.h:446
bool is_hadron() const
Definition: particledata.h:92
static double formation_power_
Power with which the cross section scaling factor grows in time.
Definition: particledata.h:471
bool hole_
If true, the object is an entry in Particles::data_ and does not hold valid particle data.
Definition: particledata.h:535
uint32_t id_process() const
Get the id of the last action.
Definition: particledata.h:138
ParticleData(const ParticleType &particle_type, int unique_id=-1)
Create a new particle with the given particle_type and optionally a specific unique_id.
Definition: particledata.h:70
double begin_formation_time() const
Get the absolute time, where the cross section scaling factor slowly starts increasing from the given...
Definition: particledata.h:255
bool is_sigmastar() const
Definition: particledata.h:113
void set_id(int i)
Set id of the particle.
Definition: particledata.h:82
bool is_core() const
Check whether the particle is core.
Definition: particledata.h:402
const FourVector & momentum() const
Get the particle's 4-momentum.
Definition: particledata.h:171
BelongsTo belongs_to() const
Getter for belongs_to label.
Definition: particledata.h:397
void copy_to(ParticleData &dst) const
Copies some information of the particle to the given particle dst.
Definition: particledata.h:484
bool operator==(int id_a) const
Check if the particle has a given id.
Definition: particledata.h:441
FourVector momentum_
momenta of the particle: x0, x1, x2, x3 as E, px, py, pz
Definition: particledata.h:543
double formation_time() const
Get the absolute formation time of the particle.
Definition: particledata.h:249
bool is_neutron() const
Definition: particledata.h:107
ThreeVector velocity() const
Get the velocity 3-vector.
Definition: particledata.h:321
double initial_xsec_scaling_factor_
Initial cross section scaling factor.
Definition: particledata.h:563
bool is_baryon() const
Definition: particledata.h:95
FourVector & spin_vector()
Get the mean spin 4-vector (Pauli–Lubanski vector) of the particle (non const reference).
Definition: particledata.h:370
FourVector spin_vector_
Pauli-Lubanski vector (mean spin 4-vector) of the particle.
Definition: particledata.h:551
bool is_rho() const
Definition: particledata.h:101
bool operator<(const ParticleData &a) const
Check if this particle has a smaller id than another particle.
Definition: particledata.h:434
double effective_mass() const
Get the particle's effective mass.
Definition: particledata.cc:25
void set_formation_time(double form_time)
Set the absolute formation time.
Definition: particledata.h:264
int32_t id() const
Get the id of the particle.
Definition: particledata.h:77
ParticleData()=default
Default constructor.
double perturbative_weight_
Perturbative weight attributed to heavy flavor particles.
Definition: particledata.h:565
const FourVector & spin_vector() const
Get the mean spin 4-vector (Pauli–Lubanski vector) of the particle (const reference,...
Definition: particledata.h:365
double perturbative_weight() const
Get the perturbative weight.
Definition: particledata.h:421
HistoryData get_history() const
Get history information.
Definition: particledata.h:143
void fluidize()
Fluidize the particle.
Definition: particledata.h:400
const double & initial_xsec_scaling_factor() const
Get the initially assigned cross section scaling factor.
Definition: particledata.h:298
void set_4momentum(double mass, double px, double py, double pz)
Set the momentum of the particle.
Definition: particledata.h:199
void set_cross_section_scaling_factor(const double &xsec_scal)
Set the particle's initial cross_section_scaling_factor.
Definition: particledata.h:313
void set_spin_vector_component(int index, double value)
Set a single component of the mean spin 4-vector (Pauli-Lubanski vector).
Definition: particledata.h:381
FourVector position_
position in space: x0, x1, x2, x3 as t, x, y, z
Definition: particledata.h:545
double pole_mass() const
Get the particle's pole mass ("on-shell").
Definition: particledata.h:119
bool operator==(const ParticleData &a) const
Check whether two particles have the same id.
Definition: particledata.h:428
double rapidity() const
Particle momentum rapidity .
Definition: particledata.h:410
int spin() const
Get the (maximum positive) spin s of a particle in multiples of 1/2.
Definition: particledata.h:360
ParticleTypePtr type_
A reference to the ParticleType object for this particle (this contains all the static information).
Definition: particledata.h:521
void set_3momentum(const ThreeVector &mom)
Set the momentum of the particle without modifying the energy.
Definition: particledata.h:209
void set_slow_formation_times(double begin_form_time, double form_time)
Set the time, when the cross section scaling factor begins, and finishes to increase from the given c...
Definition: particledata.h:282
void set_spin_vector(const FourVector &s)
Set the mean spin 4-vector (Pauli–Lubanski vector) of the particle.
Definition: particledata.h:375
double transverse_mass() const
Particle .
Definition: particledata.h:408
void boost_momentum(const ThreeVector &v)
Apply a Lorentz-boost to only the momentum.
Definition: particledata.h:352
ParticleData(const ParticleType &ptype, int uid, int index)
Construct a particle with the given type, id and index in Particles.
Definition: particledata.h:458
double spatial_rapidity() const
Particle spacetime rapidity .
Definition: particledata.h:406
void set_belongs_to(BelongsTo label)
Setter for belongs_to label.
Definition: particledata.h:395
HistoryData history_
history information
Definition: particledata.h:567
void set_unpolarized_spin_vector()
Set the 4 components of the spin vector such that the particle is unpolarized.
Definition: particledata.cc:88
void boost(const ThreeVector &v)
Apply a full Lorentz boost of momentum and position.
Definition: particledata.h:343
const FourVector & position() const
Get the particle's position in Minkowski space.
Definition: particledata.h:217
void set_perturbative_weight(const double weight)
Set the perturbative weight.
Definition: particledata.h:415
void set_4momentum(double mass, const ThreeVector &mom)
Set the momentum of the particle given its mass and momentum three-vector.
Definition: particledata.h:187
A pointer-like interface to global references to ParticleType objects.
Definition: particletype.h:682
Particle type contains the static properties of a particle species.
Definition: particletype.h:98
PdgCode pdgcode() const
Definition: particletype.h:157
bool is_rho() const
Definition: particletype.h:228
bool is_hadron() const
Definition: particletype.h:198
double mass() const
Definition: particletype.h:145
The Particles class abstracts the storage and manipulation of particles.
Definition: particles.h:33
PdgCode stores a Particle Data Group Particle Numbering Scheme particle type number.
Definition: pdgcode.h:108
unsigned int spin() const
Definition: pdgcode.h:676
bool is_pion() const
Definition: pdgcode.h:471
bool is_Sigmastar() const
Definition: pdgcode.h:458
bool is_nucleus() const
Definition: pdgcode.h:361
bool is_proton() const
Definition: pdgcode.h:410
bool is_baryon() const
Definition: pdgcode.h:398
bool is_neutron() const
Definition: pdgcode.h:416
The ThreeVector class represents a physical three-vector with the components .
Definition: threevector.h:31
std::ostream & operator<<(std::ostream &out, const ActionPtr &action)
Convenience: dereferences the ActionPtr to Action.
Definition: action.h:555
PrintParticleListDetailed detailed(const ParticleList &list)
Request the ParticleList to be printed in full detail (i.e.
Definition: particledata.h:599
constexpr int p
Proton.
Definition: action.h:24
ProcessType
ProcessTypes are used to identify the type of the process.
Definition: processbranch.h:39
@ None
See here for a short description.
bool are_particles_identical_at_given_time(const ParticleData &p1, const ParticleData &p2, double time)
Utility function to compare two ParticleData instances with respect to their PDG code,...
ParticleData create_valid_smash_particle_matching_provided_quantities(PdgCode pdgcode, double mass, const FourVector &four_position, const FourVector &four_momentum, int log_area, bool &mass_warning, bool &on_shell_warning)
This function creates a SMASH particle validating the provided information.
A structure to hold information about the history of the particle, e.g.
Definition: particledata.h:31
double time_last_collision
Time of the last action (excluding walls), time of kinetic freeze_out for HBT analysis this time shou...
Definition: particledata.h:44
int32_t id_process
id of the last action
Definition: particledata.h:35
PdgCode p2
PdgCode of the second parent particles.
Definition: particledata.h:48
PdgCode p1
PdgCode of the first parent particles.
Definition: particledata.h:46
int32_t collisions_per_particle
Collision counter per particle, zero only for initially present particles.
Definition: particledata.h:33
ProcessType process_type
type of the last action
Definition: particledata.h:37
const ParticleList & list
Particle list.
Definition: particledata.h:592