Version: SMASH-3.1
smash::ThermLatticeNode Class Reference

#include <grandcan_thermalizer.h>

The ThermLatticeNode class is intended to compute thermodynamical quantities in a cell given a set of particles.

It accumulates the upper row of the energy-momentum tensor \( T^{\mu 0}\), net baryon density nb and net strangeness densities in the computational frame. From these quantities it allows to compute the local rest frame quantites: temperature T, chemical potentials mub, mus and muq, the velocity of the local rest frame with respect to the computational frame.

An example of the intended use is:

ThermLatticeNode node; for (const ParticleData &particle: some_particles_list) { const double some_smearing_factor = ...; node.add_particle(particle, some_smearing_factor); } HadronGasEos eos = HadronGasEos(false); node.compute_rest_frame_quantities(eos);

Note that before calling the compute_rest_frame_quantities T, mu, p and e are set to zero.

Definition at line 48 of file grandcan_thermalizer.h.

Public Member Functions

 ThermLatticeNode ()
 Default constructor of thermal quantities on the lattice returning thermodynamic quantities in computational frame. More...
 
void add_particle (const ParticleData &p, double factor)
 Add particle contribution to Tmu0, nb, ns and nq May look like unused at first glance, but it is actually used by update_lattice, where the node type of the lattice is templated. More...
 
void add_particle_for_derivatives (const ParticleData &, double, ThreeVector)
 dummy function for update_lattice More...
 
void compute_rest_frame_quantities (HadronGasEos &eos)
 Temperature, chemical potentials and rest frame velocity are calculated given the hadron gas equation of state object. More...
 
void set_rest_frame_quantities (double T0, double mub0, double mus0, double muq0, const ThreeVector v0)
 Set all the rest frame quantities to some values, this is useful for testing. More...
 
FourVector Tmu0 () const
 Get Four-momentum flow of the cell. More...
 
double nb () const
 Get net baryon density of the cell in the computational frame. More...
 
double ns () const
 Get net strangeness density of the cell in the computational frame. More...
 
double nq () const
 Get net charge density of the cell in the computational frame. More...
 
double e () const
 Get energy density in the rest frame. More...
 
double p () const
 Get pressure in the rest frame. More...
 
ThreeVector v () const
 Get 3-velocity of the rest frame. More...
 
double T () const
 Get the temperature. More...
 
double mub () const
 Get the net baryon chemical potential. More...
 
double mus () const
 Get the net strangeness chemical potential. More...
 
double muq () const
 Get the net charge chemical potential. More...
 

Private Attributes

FourVector Tmu0_
 Four-momentum flow of the cell. More...
 
double nb_
 Net baryon density of the cell in the computational frame. More...
 
double ns_
 Net strangeness density of the cell in the computational frame. More...
 
double nq_
 Net charge density of the cell in the computational frame. More...
 
double e_
 Energy density in the rest frame. More...
 
double p_
 Pressure in the rest frame. More...
 
ThreeVector v_
 Velocity of the rest frame. More...
 
double T_
 Temperature. More...
 
double mub_
 Net baryon chemical potential. More...
 
double mus_
 Net strangeness chemical potential. More...
 
double muq_
 Net charge chemical potential. More...
 

Constructor & Destructor Documentation

◆ ThermLatticeNode()

smash::ThermLatticeNode::ThermLatticeNode ( )

Default constructor of thermal quantities on the lattice returning thermodynamic quantities in computational frame.

Returns
Tmu0_ Four vector \(T^{\mu 0}\), nb_ Net baryon density at this location, ns_ Net strangeness density, nq_ Net charge density, e_ Energy density, v_ 3 vector velocity of local rest frame, T_ Temperature mub_ Net baryon chemical potential, mus_ Net strangeness chemical potential, muq_ Net charge chemical potential

Definition at line 22 of file grandcan_thermalizer.cc.

23  : Tmu0_(FourVector()),
24  nb_(0.0),
25  ns_(0.0),
26  nq_(0.0),
27  e_(0.0),
28  p_(0.0),
29  v_(ThreeVector()),
30  T_(0.0),
31  mub_(0.0),
32  mus_(0.0),
33  muq_(0.0) {}
double ns_
Net strangeness density of the cell in the computational frame.
double p_
Pressure in the rest frame.
double mus_
Net strangeness chemical potential.
double e_
Energy density in the rest frame.
double nb_
Net baryon density of the cell in the computational frame.
double muq_
Net charge chemical potential.
double mub_
Net baryon chemical potential.
FourVector Tmu0_
Four-momentum flow of the cell.
ThreeVector v_
Velocity of the rest frame.
double nq_
Net charge density of the cell in the computational frame.

Member Function Documentation

◆ add_particle()

void smash::ThermLatticeNode::add_particle ( const ParticleData p,
double  factor 
)

Add particle contribution to Tmu0, nb, ns and nq May look like unused at first glance, but it is actually used by update_lattice, where the node type of the lattice is templated.

Definition at line 35 of file grandcan_thermalizer.cc.

35  {
36  Tmu0_ += part.momentum() * factor;
37  nb_ += static_cast<double>(part.type().baryon_number()) * factor;
38  ns_ += static_cast<double>(part.type().strangeness()) * factor;
39  nq_ += static_cast<double>(part.type().charge()) * factor;
40 }

◆ add_particle_for_derivatives()

void smash::ThermLatticeNode::add_particle_for_derivatives ( const ParticleData ,
double  ,
ThreeVector   
)
inline

dummy function for update_lattice

Definition at line 68 of file grandcan_thermalizer.h.

68 {}

◆ compute_rest_frame_quantities()

void smash::ThermLatticeNode::compute_rest_frame_quantities ( HadronGasEos eos)

Temperature, chemical potentials and rest frame velocity are calculated given the hadron gas equation of state object.

Parameters
[in]eos
See also
HadronGasEos based on Tmu0, nb, ns and nq
Returns
Temperature T, net baryon chemical potential mub, net strangeness chemical potential mus, net charge chemical potential muq and the velocity of the Landau rest frame, under assumption that the energy-momentum tensor has an ideal-fluid form. For more details and discussion see Oliinychenko:2015lva [40]. The advantage of this rest frame transformation is that it conserves energy and momentum, even though the dissipative part of the energy-momentum tensor is neglected.
Todo:
(oliiny): use Newton's method instead of these iterations

Definition at line 42 of file grandcan_thermalizer.cc.

42  {
44  const int max_iter = 50;
45  v_ = ThreeVector(0.0, 0.0, 0.0);
46  double e_previous_step = 0.0;
47  const double tolerance = 5.e-4;
48  int iter;
49  for (iter = 0; iter < max_iter; iter++) {
50  e_previous_step = e_;
51  e_ = Tmu0_.x0() - Tmu0_.threevec() * v_;
52  if (std::abs(e_ - e_previous_step) < tolerance) {
53  break;
54  }
55  const double gamma_inv = std::sqrt(1.0 - v_.sqr());
56  EosTable::table_element tabulated;
57  eos.from_table(tabulated, e_, gamma_inv * nb_, nq_);
58  if (!eos.is_tabulated() || tabulated.p < 0.0) {
59  auto T_mub_mus_muq =
60  eos.solve_eos(e_, gamma_inv * nb_, gamma_inv * ns_, nq_);
61  T_ = T_mub_mus_muq[0];
62  mub_ = T_mub_mus_muq[1];
63  mus_ = T_mub_mus_muq[2];
64  muq_ = T_mub_mus_muq[3];
66  } else {
67  p_ = tabulated.p;
68  T_ = tabulated.T;
69  mub_ = tabulated.mub;
70  mus_ = tabulated.mus;
71  muq_ = tabulated.muq;
72  }
73  v_ = Tmu0_.threevec() / (Tmu0_.x0() + p_);
74  }
75  if (iter == max_iter) {
76  std::cout << "Warning from solver: max iterations exceeded."
77  << " Accuracy: " << std::abs(e_ - e_previous_step)
78  << " is less than tolerance " << tolerance << std::endl;
79  }
80 }
ThreeVector threevec() const
Definition: fourvector.h:329
double x0() const
Definition: fourvector.h:313
static double pressure(double T, double mub, double mus, double muq, bool account_for_resonance_widths=false)
Compute pressure .
Definition: hadgas_eos.h:195
double sqr() const
Definition: threevector.h:266

◆ set_rest_frame_quantities()

void smash::ThermLatticeNode::set_rest_frame_quantities ( double  T0,
double  mub0,
double  mus0,
double  muq0,
const ThreeVector  v0 
)

Set all the rest frame quantities to some values, this is useful for testing.

Parameters
[out]T0Rest frame temperature
[out]mub0Rest frame net baryon chemical potential
[out]mus0Rest frame net strangeness chemical potential
[out]muq0Rest frame net charge chemical potential
[out]v0Velocity of the rest frame

Definition at line 82 of file grandcan_thermalizer.cc.

84  {
85  T_ = T0;
86  mub_ = mub0;
87  mus_ = mus0;
88  muq_ = muq0;
89  v_ = v0;
95 }
static double net_charge_density(double T, double mub, double mus, double muq, bool account_for_resonance_widths=false)
Compute net charge density.
Definition: hadgas_eos.cc:366
static double net_baryon_density(double T, double mub, double mus, double muq, bool account_for_resonance_widths=false)
Compute net baryon density.
Definition: hadgas_eos.cc:328
static double energy_density(double T, double mub, double mus, double muq)
Compute energy density.
Definition: hadgas_eos.cc:281
static double net_strange_density(double T, double mub, double mus, double muq, bool account_for_resonance_widths=false)
Compute net strangeness density.
Definition: hadgas_eos.cc:347

◆ Tmu0()

FourVector smash::ThermLatticeNode::Tmu0 ( ) const
inline

Get Four-momentum flow of the cell.

Definition at line 94 of file grandcan_thermalizer.h.

94 { return Tmu0_; }

◆ nb()

double smash::ThermLatticeNode::nb ( ) const
inline

Get net baryon density of the cell in the computational frame.

Definition at line 96 of file grandcan_thermalizer.h.

96 { return nb_; }

◆ ns()

double smash::ThermLatticeNode::ns ( ) const
inline

Get net strangeness density of the cell in the computational frame.

Definition at line 98 of file grandcan_thermalizer.h.

98 { return ns_; }

◆ nq()

double smash::ThermLatticeNode::nq ( ) const
inline

Get net charge density of the cell in the computational frame.

Definition at line 100 of file grandcan_thermalizer.h.

100 { return nq_; }

◆ e()

double smash::ThermLatticeNode::e ( ) const
inline

Get energy density in the rest frame.

Definition at line 102 of file grandcan_thermalizer.h.

102 { return e_; }

◆ p()

double smash::ThermLatticeNode::p ( ) const
inline

Get pressure in the rest frame.

Definition at line 104 of file grandcan_thermalizer.h.

104 { return p_; }

◆ v()

ThreeVector smash::ThermLatticeNode::v ( ) const
inline

Get 3-velocity of the rest frame.

Definition at line 106 of file grandcan_thermalizer.h.

106 { return v_; }

◆ T()

double smash::ThermLatticeNode::T ( ) const
inline

Get the temperature.

Definition at line 108 of file grandcan_thermalizer.h.

108 { return T_; }

◆ mub()

double smash::ThermLatticeNode::mub ( ) const
inline

Get the net baryon chemical potential.

Definition at line 110 of file grandcan_thermalizer.h.

110 { return mub_; }

◆ mus()

double smash::ThermLatticeNode::mus ( ) const
inline

Get the net strangeness chemical potential.

Definition at line 112 of file grandcan_thermalizer.h.

112 { return mus_; }

◆ muq()

double smash::ThermLatticeNode::muq ( ) const
inline

Get the net charge chemical potential.

Definition at line 114 of file grandcan_thermalizer.h.

114 { return muq_; }

Member Data Documentation

◆ Tmu0_

FourVector smash::ThermLatticeNode::Tmu0_
private

Four-momentum flow of the cell.

Definition at line 118 of file grandcan_thermalizer.h.

◆ nb_

double smash::ThermLatticeNode::nb_
private

Net baryon density of the cell in the computational frame.

Definition at line 120 of file grandcan_thermalizer.h.

◆ ns_

double smash::ThermLatticeNode::ns_
private

Net strangeness density of the cell in the computational frame.

Definition at line 122 of file grandcan_thermalizer.h.

◆ nq_

double smash::ThermLatticeNode::nq_
private

Net charge density of the cell in the computational frame.

Definition at line 124 of file grandcan_thermalizer.h.

◆ e_

double smash::ThermLatticeNode::e_
private

Energy density in the rest frame.

Definition at line 126 of file grandcan_thermalizer.h.

◆ p_

double smash::ThermLatticeNode::p_
private

Pressure in the rest frame.

Definition at line 128 of file grandcan_thermalizer.h.

◆ v_

ThreeVector smash::ThermLatticeNode::v_
private

Velocity of the rest frame.

Definition at line 130 of file grandcan_thermalizer.h.

◆ T_

double smash::ThermLatticeNode::T_
private

Temperature.

Definition at line 132 of file grandcan_thermalizer.h.

◆ mub_

double smash::ThermLatticeNode::mub_
private

Net baryon chemical potential.

Definition at line 134 of file grandcan_thermalizer.h.

◆ mus_

double smash::ThermLatticeNode::mus_
private

Net strangeness chemical potential.

Definition at line 136 of file grandcan_thermalizer.h.

◆ muq_

double smash::ThermLatticeNode::muq_
private

Net charge chemical potential.

Definition at line 138 of file grandcan_thermalizer.h.


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