Version: SMASH-3.3
smash::ParticleData Class Reference

#include <particledata.h>

ParticleData contains the dynamic information of a certain particle.

Each particle has its momentum, position and other relevant physical data entry.

Definition at line 59 of file particledata.h.

Public Member Functions

 ParticleData (const ParticleType &particle_type, int unique_id=-1)
 Create a new particle with the given particle_type and optionally a specific unique_id. More...
 
int32_t id () const
 Get the id of the particle. More...
 
void set_id (int i)
 Set id of the particle. More...
 
PdgCode pdgcode () const
 Get the pdgcode of the particle. More...
 
bool is_hadron () const
 
bool is_baryon () const
 
bool is_nucleus () const
 
bool is_rho () const
 
bool is_proton () const
 
bool is_neutron () const
 
bool is_pion () const
 
bool is_sigmastar () const
 
double pole_mass () const
 Get the particle's pole mass ("on-shell"). More...
 
double effective_mass () const
 Get the particle's effective mass. More...
 
const ParticleTypetype () const
 Get the type of the particle. More...
 
uint32_t id_process () const
 Get the id of the last action. More...
 
HistoryData get_history () const
 Get history information. More...
 
void set_history (HistoryData &&history)
 Set history_ from rvalue reference. More...
 
void set_history (int ncoll, uint32_t pid, ProcessType pt, double time_last_coll, const ParticleList &plist)
 Store history information. More...
 
const FourVectormomentum () const
 Get the particle's 4-momentum. More...
 
void set_4momentum (const FourVector &momentum_vector)
 Set the particle's 4-momentum directly. More...
 
void set_4momentum (double mass, const ThreeVector &mom)
 Set the momentum of the particle given its mass and momentum three-vector. More...
 
void set_4momentum (double mass, double px, double py, double pz)
 Set the momentum of the particle. More...
 
void set_3momentum (const ThreeVector &mom)
 Set the momentum of the particle without modifying the energy. More...
 
const FourVectorposition () const
 Get the particle's position in Minkowski space. More...
 
void set_4position (const FourVector &pos)
 Set the particle's 4-position directly. More...
 
void set_3position (const ThreeVector &pos)
 Set particle's 3-position. More...
 
ParticleData translated (const ThreeVector &delta) const
 Translate the particle position. More...
 
double formation_time () const
 Get the absolute formation time of the particle. More...
 
double begin_formation_time () const
 Get the absolute time, where the cross section scaling factor slowly starts increasing from the given scaling factor to 1. More...
 
void set_formation_time (double form_time)
 Set the absolute formation time. More...
 
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 cross section scaling factor to 1. More...
 
const double & initial_xsec_scaling_factor () const
 Get the initially assigned cross section scaling factor. More...
 
void set_cross_section_scaling_factor (const double &xsec_scal)
 Set the particle's initial cross_section_scaling_factor. More...
 
ThreeVector velocity () const
 Get the velocity 3-vector. More...
 
double inverse_gamma () const
 Get the inverse of the gamma factor from the current velocity of the particle. More...
 
void boost (const ThreeVector &v)
 Apply a full Lorentz boost of momentum and position. More...
 
void boost_momentum (const ThreeVector &v)
 Apply a Lorentz-boost to only the momentum. More...
 
int spin () const
 Get the (maximum positive) spin s of a particle in multiples of 1/2. More...
 
const FourVectorspin_vector () const
 Get the mean spin 4-vector (Pauli–Lubanski vector) of the particle (const reference, no copy). More...
 
FourVectorspin_vector ()
 Get the mean spin 4-vector (Pauli–Lubanski vector) of the particle (non const reference). More...
 
void set_spin_vector (const FourVector &s)
 Set the mean spin 4-vector (Pauli–Lubanski vector) of the particle. More...
 
void set_spin_vector_component (int index, double value)
 Set a single component of the mean spin 4-vector (Pauli-Lubanski vector). More...
 
void set_unpolarized_spin_vector ()
 Set the 4 components of the spin vector such that the particle is unpolarized. More...
 
void set_belongs_to (BelongsTo label)
 Setter for belongs_to label. More...
 
BelongsTo belongs_to () const
 Getter for belongs_to label. More...
 
void fluidize ()
 Fluidize the particle. More...
 
bool is_core () const
 Check whether the particle is core. More...
 
double hyperbolic_time () const
 Particle \(tau\) (hyperbolic time) More...
 
double spatial_rapidity () const
 Particle spacetime rapidity \(\eta_s\). More...
 
double transverse_mass () const
 Particle \(m_T\). More...
 
double rapidity () const
 Particle momentum rapidity \(y_\mathrm{rap}\). More...
 
void set_perturbative_weight (const double weight)
 Set the perturbative weight. More...
 
double perturbative_weight () const
 Get the perturbative weight. More...
 
bool operator== (const ParticleData &a) const
 Check whether two particles have the same id. More...
 
bool operator< (const ParticleData &a) const
 Check if this particle has a smaller id than another particle. More...
 
bool operator== (int id_a) const
 Check if the particle has a given id. More...
 
bool operator< (int id_a) const
 Check whether the particle's id is smaller than the given id. More...
 
 ParticleData (const ParticleType &ptype, int uid, int index)
 Construct a particle with the given type, id and index in Particles. More...
 
double xsec_scaling_factor (double delta_time=0.) const
 Return the cross section scaling factor at a given time. More...
 

Static Public Attributes

static double formation_power_ = 0.0
 Power with which the cross section scaling factor grows in time. More...
 

Private Member Functions

 ParticleData ()=default
 Default constructor. More...
 
void copy_to (ParticleData &dst) const
 Copies some information of the particle to the given particle dst. More...
 

Private Attributes

int32_t id_ = -1
 Each particle has a unique identifier. More...
 
unsigned index_ = std::numeric_limits<unsigned>::max()
 Internal index in the Particles list. More...
 
ParticleTypePtr type_
 A reference to the ParticleType object for this particle (this contains all the static information). More...
 
bool hole_ = false
 If true, the object is an entry in Particles::data_ and does not hold valid particle data. More...
 
bool core_ = false
 If the particle is part of a pseudofluid. More...
 
FourVector momentum_
 momenta of the particle: x0, x1, x2, x3 as E, px, py, pz More...
 
FourVector position_
 position in space: x0, x1, x2, x3 as t, x, y, z More...
 
FourVector spin_vector_
 Pauli-Lubanski vector (mean spin 4-vector) of the particle. More...
 
double formation_time_ = smash_NaN<double>
 Formation time at which the particle is fully formed given as an absolute value in the computational frame. More...
 
double begin_formation_time_ = smash_NaN<double>
 time when the cross section scaling factor starts to increase to 1 More...
 
double initial_xsec_scaling_factor_ = 1.0
 Initial cross section scaling factor. More...
 
double perturbative_weight_ = 1.0
 Perturbative weight attributed to heavy flavor particles. More...
 
HistoryData history_
 history information More...
 
BelongsTo belongs_to_ = BelongsTo::Nothing
 is it part of projectile or target nuclei? More...
 

Friends

class Particles
 

Constructor & Destructor Documentation

◆ ParticleData() [1/3]

smash::ParticleData::ParticleData ( const ParticleType particle_type,
int  unique_id = -1 
)
inlineexplicit

Create a new particle with the given particle_type and optionally a specific unique_id.

All other values are initialized to unphysical values.

Parameters
[in]particle_typeType of particle to be created
[in]unique_idid of particle to be created

Definition at line 70 of file particledata.h.

71  : id_(unique_id), type_(&particle_type) {}
int32_t id_
Each particle has a unique identifier.
Definition: particledata.h:505
ParticleTypePtr type_
A reference to the ParticleType object for this particle (this contains all the static information).
Definition: particledata.h:521

◆ ParticleData() [2/3]

smash::ParticleData::ParticleData ( const ParticleType ptype,
int  uid,
int  index 
)
inline

Construct a particle with the given type, id and index in Particles.

This constructor may only be called (directly or indirectly) from Particles. This constructor should be private, but can't be in order to support vector::emplace_back.

Parameters
[in]ptypeType of the particle to be constructed
[in]uidid of the particle to be constructed
[in]indexindex of the particle to be constructed

Definition at line 458 of file particledata.h.

459  : id_(uid), index_(index), type_(&ptype) {}
unsigned index_
Internal index in the Particles list.
Definition: particledata.h:515

◆ ParticleData() [3/3]

smash::ParticleData::ParticleData ( )
privatedefault

Default constructor.

Member Function Documentation

◆ id()

int32_t smash::ParticleData::id ( ) const
inline

Get the id of the particle.

Returns
particle id

Definition at line 77 of file particledata.h.

77 { return id_; }

◆ set_id()

void smash::ParticleData::set_id ( int  i)
inline

Set id of the particle.

Parameters
[in]iid to be assigned to the particle

Definition at line 82 of file particledata.h.

82 { id_ = i; }

◆ pdgcode()

PdgCode smash::ParticleData::pdgcode ( ) const
inline

Get the pdgcode of the particle.

Returns
pdgcode of the particle

Definition at line 88 of file particledata.h.

88 { return type_->pdgcode(); }
PdgCode pdgcode() const
Definition: particletype.h:157

◆ is_hadron()

bool smash::ParticleData::is_hadron ( ) const
inline

Returns
true if this is a baryon, antibaryon or meson.

Definition at line 92 of file particledata.h.

92 { return type_->is_hadron(); }
bool is_hadron() const
Definition: particletype.h:198

◆ is_baryon()

bool smash::ParticleData::is_baryon ( ) const
inline

Returns
whether this PDG code identifies a baryon.

Definition at line 95 of file particledata.h.

95 { return pdgcode().is_baryon(); }
PdgCode pdgcode() const
Get the pdgcode of the particle.
Definition: particledata.h:88
bool is_baryon() const
Definition: pdgcode.h:398

◆ is_nucleus()

bool smash::ParticleData::is_nucleus ( ) const
inline

Returns
true if this is a nucleus, false otherwise

Definition at line 98 of file particledata.h.

98 { return pdgcode().is_nucleus(); }
bool is_nucleus() const
Definition: pdgcode.h:361

◆ is_rho()

bool smash::ParticleData::is_rho ( ) const
inline

Returns
whether this is a rho meson (rho+/rho0/rho-)

Definition at line 101 of file particledata.h.

101 { return type_->is_rho(); }
bool is_rho() const
Definition: particletype.h:228

◆ is_proton()

bool smash::ParticleData::is_proton ( ) const
inline

Returns
whether this is a proton/anti-proton

Definition at line 104 of file particledata.h.

104 { return pdgcode().is_proton(); }
bool is_proton() const
Definition: pdgcode.h:410

◆ is_neutron()

bool smash::ParticleData::is_neutron ( ) const
inline

Returns
whether this is a neutron/anti-neutron

Definition at line 107 of file particledata.h.

107 { return pdgcode().is_neutron(); }
bool is_neutron() const
Definition: pdgcode.h:416

◆ is_pion()

bool smash::ParticleData::is_pion ( ) const
inline

Returns
whether this is a pion (pi+/pi0/pi-)

Definition at line 110 of file particledata.h.

110 { return pdgcode().is_pion(); }
bool is_pion() const
Definition: pdgcode.h:471

◆ is_sigmastar()

bool smash::ParticleData::is_sigmastar ( ) const
inline

Returns
whether this is a Sigma* with pdgcodes 3224, 3214, 3114

Definition at line 113 of file particledata.h.

113 { return pdgcode().is_Sigmastar(); }
bool is_Sigmastar() const
Definition: pdgcode.h:458

◆ pole_mass()

double smash::ParticleData::pole_mass ( ) const
inline

Get the particle's pole mass ("on-shell").

Returns
pole mass of the particle [GeV]

Definition at line 119 of file particledata.h.

119 { return type_->mass(); }
double mass() const
Definition: particletype.h:145

◆ effective_mass()

double smash::ParticleData::effective_mass ( ) const

Get the particle's effective mass.

Determined from the 4-momentum \(m=\sqrt{p_\mu p^\mu}\). Possibly "off-shell".

Returns
Effective mass [GeV]

Definition at line 25 of file particledata.cc.

25  {
26  const double m_pole = pole_mass();
27  if (m_pole < really_small) {
28  // prevent numerical problems with massless or very light particles
29  return m_pole;
30  } else {
31  return momentum().abs();
32  }
33 }
double abs() const
calculate the lorentz invariant absolute value
Definition: fourvector.h:464
const FourVector & momentum() const
Get the particle's 4-momentum.
Definition: particledata.h:171
double pole_mass() const
Get the particle's pole mass ("on-shell").
Definition: particledata.h:119
constexpr double really_small
Numerical error tolerance.
Definition: constants.h:41

◆ type()

const ParticleType& smash::ParticleData::type ( ) const
inline

Get the type of the particle.

Returns
ParticleType object associated to this particle.

Definition at line 132 of file particledata.h.

132 { return *type_; }

◆ id_process()

uint32_t smash::ParticleData::id_process ( ) const
inline

Get the id of the last action.

Returns
id of particle's latest collision

Definition at line 138 of file particledata.h.

138 { return history_.id_process; }
HistoryData history_
history information
Definition: particledata.h:567
int32_t id_process
id of the last action
Definition: particledata.h:35

◆ get_history()

HistoryData smash::ParticleData::get_history ( ) const
inline

Get history information.

Returns
particle history struct

Definition at line 143 of file particledata.h.

143 { return history_; }

◆ set_history() [1/2]

void smash::ParticleData::set_history ( HistoryData &&  history)
inline

Set history_ from rvalue reference.

Meant to be used only in special situations e.g. in the ListModus, where a temporary HistoryData is constructed from the user input.

Parameters
[in]historyobject to be moved from.

Definition at line 151 of file particledata.h.

151 { history_ = std::move(history); }

◆ set_history() [2/2]

void smash::ParticleData::set_history ( int  ncoll,
uint32_t  pid,
ProcessType  pt,
double  time_last_coll,
const ParticleList &  plist 
)

Store history information.

The history contains the type of process and possibly the PdgCodes of the parent particles (plist). Note that history is not set for dileptons and photons.

Parameters
[in]ncollparticle's number of collisions
[in]pidid of the particle's latest process
[in]ptprocess type of the particle's latest process
[in]time_last_colltime of latest collision [fm]
[in]plistlist of parent particles

Definition at line 35 of file particledata.cc.

37  {
40  history_.time_last_collision = time_last_coll;
41  }
42  history_.id_process = pid;
44  switch (pt) {
45  case ProcessType::Decay:
46  case ProcessType::Wall:
47  // only store one parent
48  history_.p1 = plist[0].pdgcode();
49  history_.p2 = 0x0;
50  break;
55  // Parent particles are not updated by the elastic scatterings,
56  // failed string processes, or fluidizations
57  break;
70  // store two parent particles
71  history_.p1 = plist[0].pdgcode();
72  history_.p2 = plist[1].pdgcode();
73  break;
80  case ProcessType::None:
81  // nullify parents
82  history_.p1 = 0x0;
83  history_.p2 = 0x0;
84  break;
85  }
86 }
@ FluidizationNoRemoval
See here for a short description.
@ FailedString
See here for a short description.
@ TwoToOne
See here for a short description.
@ MultiParticleThreeToTwo
See here for a short description.
@ StringSoftDoubleDiffractive
See here for a short description.
@ Fluidization
See here for a short description.
@ Bremsstrahlung
See here for a short description.
@ Thermalization
See here for a short description.
@ Freeforall
See here for a short description.
@ Decay
See here for a short description.
@ TwoToFive
See here for a short description.
@ None
See here for a short description.
@ StringSoftSingleDiffractiveXB
See here for a short description.
@ TwoToTwo
See here for a short description.
@ Wall
See here for a short description.
@ Elastic
See here for a short description.
@ TwoToFour
See here for a short description.
@ StringSoftAnnihilation
See here for a short description.
@ MultiParticleThreeMesonsToOne
See here for a short description.
@ StringSoftNonDiffractive
See here for a short description.
@ MultiParticleFourToTwo
See here for a short description.
@ StringSoftSingleDiffractiveAX
See here for a short description.
@ StringHard
See here for a short description.
@ TwoToThree
See here for a short description.
@ MultiParticleFiveToTwo
See here for a short description.
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
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

◆ momentum()

const FourVector& smash::ParticleData::momentum ( ) const
inline

Get the particle's 4-momentum.

Returns
particle's 4-momentum [GeV]

Definition at line 171 of file particledata.h.

171 { return momentum_; }
FourVector momentum_
momenta of the particle: x0, x1, x2, x3 as E, px, py, pz
Definition: particledata.h:543

◆ set_4momentum() [1/3]

void smash::ParticleData::set_4momentum ( const FourVector momentum_vector)
inline

Set the particle's 4-momentum directly.

Parameters
[in]momentum_vector4-vector \(p^\mu = (E,\mathbf{p})^T\)

Definition at line 177 of file particledata.h.

177  {
178  momentum_ = momentum_vector;
179  }

◆ set_4momentum() [2/3]

void smash::ParticleData::set_4momentum ( double  mass,
const ThreeVector mom 
)
inline

Set the momentum of the particle given its mass and momentum three-vector.

Parameters
[in]massthe mass of the particle (without E_kin contribution) [GeV]
[in]momthe three-momentum of the particle [GeV]

Definition at line 187 of file particledata.h.

187  {
188  momentum_ = FourVector(std::sqrt(mass * mass + mom * mom), mom);
189  }

◆ set_4momentum() [3/3]

void smash::ParticleData::set_4momentum ( double  mass,
double  px,
double  py,
double  pz 
)
inline

Set the momentum of the particle.

Parameters
[in]massthe mass of the particle (without E_kin contribution) [GeV]
[in]pxx-component of the momentum [GeV]
[in]pyy-component of the momentum [GeV]
[in]pzz-component of the momentum [GeV]

Definition at line 199 of file particledata.h.

199  {
200  momentum_ = FourVector(std::sqrt(mass * mass + px * px + py * py + pz * pz),
201  px, py, pz);
202  }

◆ set_3momentum()

void smash::ParticleData::set_3momentum ( const ThreeVector mom)
inline

Set the momentum of the particle without modifying the energy.

WARNING: Mass gets modified.

Parameters
[in]mommomentum 3-vector [GeV]

Definition at line 209 of file particledata.h.

209  {
210  momentum_ = FourVector(momentum_.x0(), mom);
211  }
double x0() const
Definition: fourvector.h:313

◆ position()

const FourVector& smash::ParticleData::position ( ) const
inline

Get the particle's position in Minkowski space.

Returns
particle's position 4-vector

Definition at line 217 of file particledata.h.

217 { return position_; }
FourVector position_
position in space: x0, x1, x2, x3 as t, x, y, z
Definition: particledata.h:545

◆ set_4position()

void smash::ParticleData::set_4position ( const FourVector pos)
inline

Set the particle's 4-position directly.

Parameters
[in]posposition 4-vector

Definition at line 222 of file particledata.h.

222 { position_ = pos; }

◆ set_3position()

void smash::ParticleData::set_3position ( const ThreeVector pos)
inline

Set particle's 3-position.

The time component is not changed

Parameters
[in]posposition 3-vector

Definition at line 229 of file particledata.h.

229  {
230  position_ = FourVector(position_.x0(), pos);
231  }

◆ translated()

ParticleData smash::ParticleData::translated ( const ThreeVector delta) const
inline

Translate the particle position.

Parameters
[in]delta3-vector by which the particle is translated [fm]

Definition at line 237 of file particledata.h.

237  {
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  }
ParticleData()=default
Default constructor.
constexpr int p
Proton.

◆ formation_time()

double smash::ParticleData::formation_time ( ) const
inline

Get the absolute formation time of the particle.

Returns
particle's formation time

Definition at line 249 of file particledata.h.

249 { return formation_time_; }
double formation_time_
Formation time at which the particle is fully formed given as an absolute value in the computational ...
Definition: particledata.h:556

◆ begin_formation_time()

double smash::ParticleData::begin_formation_time ( ) const
inline

Get the absolute time, where the cross section scaling factor slowly starts increasing from the given scaling factor to 1.

Returns
time, when scaling factor starts increasing

Definition at line 255 of file particledata.h.

255 { return begin_formation_time_; }
double begin_formation_time_
time when the cross section scaling factor starts to increase to 1
Definition: particledata.h:558

◆ set_formation_time()

void smash::ParticleData::set_formation_time ( double  form_time)
inline

Set the absolute formation time.

The particle's cross section scaling factor will be a Heavyside fuction of time.

Parameters
[in]form_timeabsolute formation time

Definition at line 264 of file particledata.h.

264  {
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  }

◆ set_slow_formation_times()

void smash::ParticleData::set_slow_formation_times ( double  begin_form_time,
double  form_time 
)
inline

Set the time, when the cross section scaling factor begins, and finishes to increase from the given cross section scaling factor to 1.

The cross section will only grow slowly, if the option is used.

Parameters
[in]begin_form_timetime when the cross section starts to increase
[in]form_timetime when the cross section reaches 1

Definition at line 282 of file particledata.h.

282  {
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  }

◆ initial_xsec_scaling_factor()

const double& smash::ParticleData::initial_xsec_scaling_factor ( ) const
inline

Get the initially assigned cross section scaling factor.

Depending on the config, the cross section scaling factor might change with time, while this value will not be updated.

Returns
particle's initially assigned cross section scaling factor

Definition at line 298 of file particledata.h.

298  {
300  }
double initial_xsec_scaling_factor_
Initial cross section scaling factor.
Definition: particledata.h:563

◆ set_cross_section_scaling_factor()

void smash::ParticleData::set_cross_section_scaling_factor ( const double &  xsec_scal)
inline

Set the particle's initial cross_section_scaling_factor.

All cross sections of this particle are scaled down by this factor until the formation time is over.

If the particle formation power is set to be positive, this will only be the initial scaling factor, while the actual scaling factor grows with time.

Parameters
[in]xsec_scalcross section scaling factor

Definition at line 313 of file particledata.h.

313  {
314  initial_xsec_scaling_factor_ = xsec_scal;
315  }

◆ velocity()

ThreeVector smash::ParticleData::velocity ( ) const
inline

Get the velocity 3-vector.

Returns
3-velocity of the particle

Definition at line 321 of file particledata.h.

321 { return momentum_.velocity(); }
ThreeVector velocity() const
Get the velocity (3-vector divided by zero component).
Definition: fourvector.h:333

◆ inverse_gamma()

double smash::ParticleData::inverse_gamma ( ) const
inline

Get the inverse of the gamma factor from the current velocity of the particle.

\[\frac{1}{\gamma}=\sqrt{1-v^2}\]

This functions is more efficient than calculating the gamma factor from velocity, since the velocity function must execute three divisions (for every space component of the momentum vector).

Returns
inverse gamma factor

Definition at line 335 of file particledata.h.

335  {
336  return std::sqrt(1. - momentum_.sqr3() / (momentum_.x0() * momentum_.x0()));
337  }
double sqr3() const
calculate the square of the spatial three-vector
Definition: fourvector.h:474

◆ boost()

void smash::ParticleData::boost ( const ThreeVector v)
inline

Apply a full Lorentz boost of momentum and position.

Parameters
[in]vboost 3-velocity

Definition at line 343 of file particledata.h.

343  {
346  }
FourVector lorentz_boost(const ThreeVector &v) const
Returns the FourVector boosted with velocity v.
Definition: fourvector.cc:17
void set_4momentum(const FourVector &momentum_vector)
Set the particle's 4-momentum directly.
Definition: particledata.h:177
void set_4position(const FourVector &pos)
Set the particle's 4-position directly.
Definition: particledata.h:222

◆ boost_momentum()

void smash::ParticleData::boost_momentum ( const ThreeVector v)
inline

Apply a Lorentz-boost to only the momentum.

Parameters
[in]vboost 3-veloctity

Definition at line 352 of file particledata.h.

352  {
354  }

◆ spin()

int smash::ParticleData::spin ( ) const
inline

Get the (maximum positive) spin s of a particle in multiples of 1/2.

E.g. for a spin-1 particle s=2.

Returns
particle's spin in multiples of 1/2

Definition at line 360 of file particledata.h.

360 { return pdgcode().spin(); }
unsigned int spin() const
Definition: pdgcode.h:676

◆ spin_vector() [1/2]

const FourVector& smash::ParticleData::spin_vector ( ) const
inline

Get the mean spin 4-vector (Pauli–Lubanski vector) of the particle (const reference, no copy).

Returns
particle's mean spin 4-vector

Definition at line 365 of file particledata.h.

365 { return spin_vector_; }
FourVector spin_vector_
Pauli-Lubanski vector (mean spin 4-vector) of the particle.
Definition: particledata.h:551

◆ spin_vector() [2/2]

FourVector& smash::ParticleData::spin_vector ( )
inline

Get the mean spin 4-vector (Pauli–Lubanski vector) of the particle (non const reference).

Returns
particle's mean spin 4-vector

Definition at line 370 of file particledata.h.

370 { return spin_vector_; }

◆ set_spin_vector()

void smash::ParticleData::set_spin_vector ( const FourVector s)
inline

Set the mean spin 4-vector (Pauli–Lubanski vector) of the particle.

Parameters
[in]sparticle's mean spin 4-vector

Definition at line 375 of file particledata.h.

375 { spin_vector_ = s; }

◆ set_spin_vector_component()

void smash::ParticleData::set_spin_vector_component ( int  index,
double  value 
)
inline

Set a single component of the mean spin 4-vector (Pauli-Lubanski vector).

Parameters
[in]indexcomponent index (0-3)
[in]valuecomponent value

Definition at line 381 of file particledata.h.

381  {
382  if (index < 0 || index > 3) {
383  throw std::out_of_range("Invalid spin vector component index");
384  }
385  spin_vector_[index] = value;
386  }

◆ set_unpolarized_spin_vector()

void smash::ParticleData::set_unpolarized_spin_vector ( )

Set the 4 components of the spin vector such that the particle is unpolarized.

This function is used only to initialize the spin vector of particles at creation.

For finite-spin particles, we assign unpolarized spin vectors by sampling the spatial components from a normal distribution with mean 0 in the particle rest frame, ensuring ⟨S⟩ = 0 on average. The time component S⁰ is set to 0. The resulting spin vector is then Lorentz-boosted to the lab frame.

This initialization is not physical spin quantization, but a statistically unpolarized setup. The standard deviation is arbitrary and chosen small to avoid unphysical artifacts.

Definition at line 88 of file particledata.cc.

88  {
89  // For massless particles, a rest frame does not exist,
90  // so the spin 4-vector cannot be defined via a rest-frame boost.
91  // In such cases, we assign a vanishing spin vector to ensure
92  // numerical stability and well-defined behavior.
93  // For spin-0 particles, the spin 4-vector is physically zero
94  if (pole_mass() == 0.0 || spin() == 0) {
95  spin_vector_ = FourVector(0., 0., 0., 0.);
96  return;
97  }
98 
99  // Check whether the velocity of a particle is set and not nan
100  const ThreeVector v = velocity();
101  assert(!(std::isnan(v.x1()) || std::isnan(v.x2()) || std::isnan(v.x3())));
102 
114  constexpr double mean = 0.0;
115  constexpr double sigma = 0.75;
116 
117  const FourVector rest_frame_spin(0., random::normal(mean, sigma),
118  random::normal(mean, sigma),
119  random::normal(mean, sigma));
120 
121  // Boost the spin vector from rest frame to lab frame
122  spin_vector_ = rest_frame_spin.lorentz_boost(v);
123 }
ThreeVector velocity() const
Get the velocity 3-vector.
Definition: particledata.h:321
int spin() const
Get the (maximum positive) spin s of a particle in multiples of 1/2.
Definition: particledata.h:360
double normal(const T &mean, const T &sigma)
Returns a random number drawn from a normal distribution.
Definition: random.h:250

◆ set_belongs_to()

void smash::ParticleData::set_belongs_to ( BelongsTo  label)
inline

Setter for belongs_to label.

Definition at line 395 of file particledata.h.

395 { belongs_to_ = label; }
BelongsTo belongs_to_
is it part of projectile or target nuclei?
Definition: particledata.h:569

◆ belongs_to()

BelongsTo smash::ParticleData::belongs_to ( ) const
inline

Getter for belongs_to label.

Definition at line 397 of file particledata.h.

397 { return belongs_to_; }

◆ fluidize()

void smash::ParticleData::fluidize ( )
inline

Fluidize the particle.

Definition at line 400 of file particledata.h.

400 { core_ = true; }
bool core_
If the particle is part of a pseudofluid.
Definition: particledata.h:540

◆ is_core()

bool smash::ParticleData::is_core ( ) const
inline

Check whether the particle is core.

Definition at line 402 of file particledata.h.

402 { return core_; }

◆ hyperbolic_time()

double smash::ParticleData::hyperbolic_time ( ) const
inline

Particle \(tau\) (hyperbolic time)

Definition at line 404 of file particledata.h.

404 { return position_.tau(); }
double tau() const
calculate the proper time from the given four vector
Definition: fourvector.h:478

◆ spatial_rapidity()

double smash::ParticleData::spatial_rapidity ( ) const
inline

Particle spacetime rapidity \(\eta_s\).

Definition at line 406 of file particledata.h.

406 { return position_.eta(); }
double eta() const
calculate the space-time rapidity from the given four vector
Definition: fourvector.h:482

◆ transverse_mass()

double smash::ParticleData::transverse_mass ( ) const
inline

Particle \(m_T\).

Definition at line 408 of file particledata.h.

408 { return momentum_.tau(); }

◆ rapidity()

double smash::ParticleData::rapidity ( ) const
inline

Particle momentum rapidity \(y_\mathrm{rap}\).

Definition at line 410 of file particledata.h.

410 { return momentum_.eta(); }

◆ set_perturbative_weight()

void smash::ParticleData::set_perturbative_weight ( const double  weight)
inline

Set the perturbative weight.

Definition at line 415 of file particledata.h.

415  {
416  perturbative_weight_ = weight;
417  }
double perturbative_weight_
Perturbative weight attributed to heavy flavor particles.
Definition: particledata.h:565

◆ perturbative_weight()

double smash::ParticleData::perturbative_weight ( ) const
inline

Get the perturbative weight.

Definition at line 421 of file particledata.h.

421 { return perturbative_weight_; }

◆ operator==() [1/2]

bool smash::ParticleData::operator== ( const ParticleData a) const
inline

Check whether two particles have the same id.

Parameters
[in]aparticle to compare to
Returns
whether the particles have the same id

Definition at line 428 of file particledata.h.

428 { return this->id_ == a.id_; }

◆ operator<() [1/2]

bool smash::ParticleData::operator< ( const ParticleData a) const
inline

Check if this particle has a smaller id than another particle.

Parameters
[in]aparticle to compare to
Returns
whether this particle has a smaller id than other particle

Definition at line 434 of file particledata.h.

434 { return this->id_ < a.id_; }

◆ operator==() [2/2]

bool smash::ParticleData::operator== ( int  id_a) const
inline

Check if the particle has a given id.

Parameters
[in]id_aid to compare to
Returns
whether the particle has the given id

Definition at line 441 of file particledata.h.

441 { return this->id_ == id_a; }

◆ operator<() [2/2]

bool smash::ParticleData::operator< ( int  id_a) const
inline

Check whether the particle's id is smaller than the given id.

Parameters
[in]id_anumber to compare particle's id to

Definition at line 446 of file particledata.h.

446 { return this->id_ < id_a; }

◆ xsec_scaling_factor()

double smash::ParticleData::xsec_scaling_factor ( double  delta_time = 0.) const

Return the cross section scaling factor at a given time.

Parameters
[in]delta_timescaling factor at current time plus this time will be returned.
Returns
the cross section scaling factor at a specified time.

Definition at line 125 of file particledata.cc.

125  {
126  // if formation times are NaNs simply return a NaN to avoid floating-point
127  // FE_INVALID errors (that would be raised by unordered comparisons with NaNs)
128  if (std::isnan(formation_time_) || std::isnan(begin_formation_time_)) {
129  return smash_NaN<double>;
130  }
131 
132  double time_of_interest = position_.x0() + delta_time;
133  // cross section scaling factor at the time_of_interest
134  double scaling_factor;
135 
136  if (formation_power_ <= 0.) {
137  // use a step function to form particles
138  if (time_of_interest < formation_time_) {
139  // particles will not be fully formed at time of interest
140  scaling_factor = initial_xsec_scaling_factor_;
141  } else {
142  // particles are fully formed at time of interest
143  scaling_factor = 1.;
144  }
145  } else {
146  // use smooth function to scale cross section (unless particles are already
147  // fully formed at desired time or will start to form later)
148  if (formation_time_ <= time_of_interest) {
149  // particles are fully formed when colliding
150  scaling_factor = 1.;
151  } else if (begin_formation_time_ >= time_of_interest) {
152  // particles will start formimg later
153  scaling_factor = initial_xsec_scaling_factor_;
154  } else {
155  // particles are in the process of formation at the given time
156  scaling_factor =
159  std::pow((time_of_interest - begin_formation_time_) /
162  }
163  }
164  return scaling_factor;
165 }
static double formation_power_
Power with which the cross section scaling factor grows in time.
Definition: particledata.h:471

◆ copy_to()

void smash::ParticleData::copy_to ( ParticleData dst) const
inlineprivate

Copies some information of the particle to the given particle dst.

Specifically it avoids to copy id_, index_, and type_.

Parameters
[in]dstparticle values are copied to

Definition at line 484 of file particledata.h.

484  {
485  dst.history_ = history_;
486  dst.momentum_ = momentum_;
487  dst.position_ = position_;
488  dst.spin_vector_ = spin_vector_;
489  dst.formation_time_ = formation_time_;
490  dst.initial_xsec_scaling_factor_ = initial_xsec_scaling_factor_;
491  dst.begin_formation_time_ = begin_formation_time_;
492  dst.belongs_to_ = belongs_to_;
493  dst.core_ = core_;
494  dst.perturbative_weight_ = perturbative_weight_;
495  }

Friends And Related Function Documentation

◆ Particles

friend class Particles
friend

Definition at line 474 of file particledata.h.

Member Data Documentation

◆ formation_power_

double smash::ParticleData::formation_power_ = 0.0
static

Power with which the cross section scaling factor grows in time.

Definition at line 471 of file particledata.h.

◆ id_

int32_t smash::ParticleData::id_ = -1
private

Each particle has a unique identifier.

This identifier is used for identifying the particle in the output files. It is specifically not used for searching for ParticleData objects in lists of particles, though it may be used to identify two ParticleData objects as referencing the same particle. This is why the comparison operators depend only on the id_ member.

Definition at line 505 of file particledata.h.

◆ index_

unsigned smash::ParticleData::index_ = std::numeric_limits<unsigned>::max()
private

Internal index in the Particles list.

This number is used to find the Experiment-wide original of this copy.

The value is read and written from the Particles class.

See also
Particles::data_

Definition at line 515 of file particledata.h.

◆ type_

ParticleTypePtr smash::ParticleData::type_
private

A reference to the ParticleType object for this particle (this contains all the static information).

Default-initialized with an invalid index.

Definition at line 521 of file particledata.h.

◆ hole_

bool smash::ParticleData::hole_ = false
private

If true, the object is an entry in Particles::data_ and does not hold valid particle data.

Specifically iterations over Particles must skip objects with hole_ == true. All other ParticleData instances should set this member to false.

See also
Particles::data_

Definition at line 535 of file particledata.h.

◆ core_

bool smash::ParticleData::core_ = false
private

If the particle is part of a pseudofluid.

Definition at line 540 of file particledata.h.

◆ momentum_

FourVector smash::ParticleData::momentum_
private

momenta of the particle: x0, x1, x2, x3 as E, px, py, pz

Definition at line 543 of file particledata.h.

◆ position_

FourVector smash::ParticleData::position_
private

position in space: x0, x1, x2, x3 as t, x, y, z

Definition at line 545 of file particledata.h.

◆ spin_vector_

FourVector smash::ParticleData::spin_vector_
private
Initial value:
= FourVector(smash_NaN<double>, smash_NaN<double>,
smash_NaN<double>, smash_NaN<double>)

Pauli-Lubanski vector (mean spin 4-vector) of the particle.

Each component is initialized with NaN (double) to indicate that the spin vector has not been set.

Definition at line 551 of file particledata.h.

◆ formation_time_

double smash::ParticleData::formation_time_ = smash_NaN<double>
private

Formation time at which the particle is fully formed given as an absolute value in the computational frame.

Definition at line 556 of file particledata.h.

◆ begin_formation_time_

double smash::ParticleData::begin_formation_time_ = smash_NaN<double>
private

time when the cross section scaling factor starts to increase to 1

Definition at line 558 of file particledata.h.

◆ initial_xsec_scaling_factor_

double smash::ParticleData::initial_xsec_scaling_factor_ = 1.0
private

Initial cross section scaling factor.

1 by default, since a particle is fully formed in this case.

Definition at line 563 of file particledata.h.

◆ perturbative_weight_

double smash::ParticleData::perturbative_weight_ = 1.0
private

Perturbative weight attributed to heavy flavor particles.

Definition at line 565 of file particledata.h.

◆ history_

HistoryData smash::ParticleData::history_
private

history information

Definition at line 567 of file particledata.h.

◆ belongs_to_

BelongsTo smash::ParticleData::belongs_to_ = BelongsTo::Nothing
private

is it part of projectile or target nuclei?

Definition at line 569 of file particledata.h.


The documentation for this class was generated from the following files: