Version: SMASH-2.2
smash::Potentials Class Reference

#include <potentials.h>

A class that stores parameters of potentials, calculates potentials and their gradients.

Potentials are responsible for long-range interactions and stand in the left part of Boltzmann equation. Short-range interactions are taken into account in the right part of it - in the collision term.

Definition at line 32 of file potentials.h.

Collaboration diagram for smash::Potentials:
[legend]

Public Member Functions

 Potentials (Configuration conf, const DensityParameters &parameters)
 Potentials constructor. More...
 
virtual ~Potentials ()
 Standard destructor. More...
 
double skyrme_pot (const double baryon_density) const
 Evaluates skyrme potential given a baryon density. More...
 
double symmetry_pot (const double baryon_isospin_density, const double baryon_density) const
 Evaluates symmetry potential given baryon isospin density. More...
 
double symmetry_S (const double baryon_density) const
 Calculate the factor \(S(\rho)\) in the symmetry potential. More...
 
FourVector vdf_pot (double rhoB, const FourVector jmuB_net) const
 Evaluates the FourVector potential in the VDF model given the rest frame density and the computational frame baryon current. More...
 
double potential (const ThreeVector &r, const ParticleList &plist, const ParticleType &acts_on) const
 Evaluates potential (Skyrme with optional Symmetry or VDF) at point r. More...
 
std::pair< ThreeVector, ThreeVectorskyrme_force (const double rhoB, const ThreeVector grad_j0B, const ThreeVector dvecjB_dt, const ThreeVector curl_vecjB) const
 Evaluates the electric and magnetic components of the skyrme force. More...
 
std::pair< ThreeVector, ThreeVectorsymmetry_force (const double rhoI3, const ThreeVector grad_j0I3, const ThreeVector dvecjI3_dt, const ThreeVector curl_vecjI3, const double rhoB, const ThreeVector grad_j0B, const ThreeVector dvecjB_dt, const ThreeVector curl_vecjB) const
 Evaluates the electric and magnetic components of the symmetry force. More...
 
std::pair< ThreeVector, ThreeVectorvdf_force (double rhoB, const double drhoB_dt, const ThreeVector grad_rhoB, const ThreeVector gradrhoB_cross_vecjB, const double j0B, const ThreeVector grad_j0B, const ThreeVector vecjB, const ThreeVector dvecjB_dt, const ThreeVector curl_vecjB) const
 Evaluates the electric and magnetic components of force in the VDF model given the derivatives of the baryon current \(j^{\mu}\). More...
 
std::pair< ThreeVector, ThreeVectorvdf_force (const ThreeVector grad_A_0, const ThreeVector dA_dt, const ThreeVector curl_vecA) const
 Evaluates the electric and magnetic components of force in the VDF force given the derivatives of the VDF mean-field \(A^\mu\). More...
 
virtual std::tuple< ThreeVector, ThreeVector, ThreeVector, ThreeVectorall_forces (const ThreeVector &r, const ParticleList &plist) const
 Evaluates the electric and magnetic components of the forces at point r. More...
 
virtual bool use_skyrme () const
 
virtual bool use_symmetry () const
 
virtual bool use_coulomb () const
 
double skyrme_a () const
 
double skyrme_b () const
 
double skyrme_tau () const
 
double symmetry_S_pot () const
 
virtual bool use_vdf () const
 
double saturation_density () const
 
const std::vector< double > & coeffs () const
 
const std::vector< double > & powers () const
 
int number_of_terms () const
 
double coulomb_r_cut () const
 

Static Public Member Functions

static std::pair< double, int > force_scale (const ParticleType &data)
 Evaluates the scaling factor of the forces acting on the particles. More...
 
static ThreeVector E_field_integrand (ThreeVector pos, DensityOnLattice &charge_density, ThreeVector point)
 Integrand for calculating the electric field. More...
 
static ThreeVector B_field_integrand (ThreeVector pos, DensityOnLattice &charge_density, ThreeVector point)
 Integrand for calculating the magnetic field using the Biot-Savart formula. More...
 

Private Member Functions

double dVsym_drhoI3 (const double rhoB, const double rhoI3) const
 Calculate the derivative of the symmetry potential with respect to the isospin density in GeV * fm^3. More...
 
double dVsym_drhoB (const double rhoB, const double rhoI3) const
 Calculate the derivative of the symmetry potential with respect to the net baryon density in GeV * fm^3. More...
 

Private Attributes

const DensityParameters param_
 Struct that contains the gaussian smearing width \(\sigma\), the distance cutoff \(r_{\rm cut}\) and the testparticle number needed for the density calculation. More...
 
bool use_skyrme_
 Skyrme potential on/off. More...
 
bool use_symmetry_
 Symmetry potential on/off. More...
 
bool use_coulomb_
 Coulomb potential on/Off. More...
 
bool use_vdf_
 VDF potential on/off. More...
 
double skyrme_a_
 Parameter of skyrme potentials: the coefficient in front of \(\frac{\rho}{\rho_0}\) in GeV. More...
 
double skyrme_b_
 Parameters of skyrme potentials: the coefficient in front of \((\frac{\rho}{\rho_0})^\tau\) in GeV. More...
 
double skyrme_tau_
 Parameters of skyrme potentials: the power index. More...
 
double symmetry_S_Pot_
 Parameter S_Pot in the symmetry potential in MeV. More...
 
bool symmetry_is_rhoB_dependent_ = false
 Whether the baryon density dependence of the symmetry potential is included. More...
 
double symmetry_gamma_
 Power \( \gamma \) in formula for \( S(\rho) \): More...
 
double coulomb_r_cut_
 Cutoff in integration for coulomb potential. More...
 
double saturation_density_
 Saturation density of nuclear matter used in the VDF potential; it may vary between different parameterizations. More...
 
std::vector< double > coeffs_
 Parameters of the VDF potential: coefficients \(C_i\), in GeV. More...
 
std::vector< double > powers_
 Parameters of the VDF potential: exponents \(b_i\). More...
 

Constructor & Destructor Documentation

◆ Potentials()

smash::Potentials::Potentials ( Configuration  conf,
const DensityParameters parameters 
)

Potentials constructor.

Parameters
[in]confConfiguration which contains the switches determining whether to turn on the Skyrme or the symmetry potentials, and the coefficents controlling how strong the potentials are.
[in]parametersStruct that contains the gaussian smearing factor \(\sigma\), the distance cutoff \(r_{\rm cut}\) and the testparticle number needed for the density calculation.

◆ ~Potentials()

smash::Potentials::~Potentials ( )
virtual

Standard destructor.

Definition at line 232 of file potentials.cc.

232 {}

Member Function Documentation

◆ skyrme_pot()

double smash::Potentials::skyrme_pot ( const double  baryon_density) const

Evaluates skyrme potential given a baryon density.

Parameters
[in]baryon_densityBaryon density \(\rho\) evaluated in the local rest frame in fm \(^{-3}\).
Returns
Skyrme potential

\[U_B=10^{-3}\times\frac{\rho}{|\rho|} (A\frac{\rho}{\rho_0}+B(\frac{\rho}{\rho_0})^\tau)\]

in GeV

Definition at line 234 of file potentials.cc.

234  {
235  const double tmp = baryon_density / nuclear_density;
236  /* U = U(|rho|) * sgn , because the sign of the potential changes
237  * under a charge reversal transformation. */
238  const int sgn = tmp > 0 ? 1 : -1;
239  // Return in GeV
240  return mev_to_gev * sgn *
241  (skyrme_a_ * std::abs(tmp) +
242  skyrme_b_ * std::pow(std::abs(tmp), skyrme_tau_));
243 }
double skyrme_a_
Parameter of skyrme potentials: the coefficient in front of in GeV.
Definition: potentials.h:358
double skyrme_tau_
Parameters of skyrme potentials: the power index.
Definition: potentials.h:370
double skyrme_b_
Parameters of skyrme potentials: the coefficient in front of in GeV.
Definition: potentials.h:364
int sgn(T val)
Signum function.
Definition: random.h:190
constexpr double mev_to_gev
MeV to GeV conversion factor.
Definition: constants.h:34
constexpr double nuclear_density
Ground state density of symmetric nuclear matter [fm^-3].
Definition: constants.h:48
Here is the call graph for this function:

◆ symmetry_pot()

double smash::Potentials::symmetry_pot ( const double  baryon_isospin_density,
const double  baryon_density 
) const

Evaluates symmetry potential given baryon isospin density.

Note
The second term is neglected if \(\gamma\) is not specified in the config.
Parameters
[in]baryon_isospin_densityThe difference between the proton and the neutron density in the local rest frame in fm \(^{-3}\).
[in]baryon_density
Returns
Symmetry potential

\[U_I=2\times 10^{-3}S_{\rm sym} \frac{\rho_{I_3}}{\rho_0} + \left[12.3\left(\frac{\rho_B}{\rho_0}\right)^{2/3} + 20\left(\frac{\rho_B}{\rho_0}\right)^\gamma\right] \left(\frac{\rho_{I_3}}{\rho_B}\right)^2\]

in GeV

Definition at line 252 of file potentials.cc.

253  {
254  double pot = mev_to_gev * 2. * symmetry_S_Pot_ * baryon_isospin_density /
257  pot += mev_to_gev * symmetry_S(baryon_density) * baryon_isospin_density *
258  baryon_isospin_density / (baryon_density * baryon_density);
259  }
260  return pot;
261 }
double symmetry_S(const double baryon_density) const
Calculate the factor in the symmetry potential.
Definition: potentials.cc:244
double symmetry_S_Pot_
Parameter S_Pot in the symmetry potential in MeV.
Definition: potentials.h:373
bool symmetry_is_rhoB_dependent_
Whether the baryon density dependence of the symmetry potential is included.
Definition: potentials.h:379

◆ symmetry_S()

double smash::Potentials::symmetry_S ( const double  baryon_density) const

Calculate the factor \(S(\rho)\) in the symmetry potential.

Parameters
[in]baryon_densitybaryon density
Returns
factor S in symmetry potenial

Definition at line 244 of file potentials.cc.

244  {
246  return 12.3 * std::pow(baryon_density / nuclear_density, 2. / 3.) +
247  20.0 * std::pow(baryon_density / nuclear_density, symmetry_gamma_);
248  } else {
249  return 0.;
250  }
251 }
double symmetry_gamma_
Power in formula for :
Definition: potentials.h:386

◆ vdf_pot()

FourVector smash::Potentials::vdf_pot ( double  rhoB,
const FourVector  jmuB_net 
) const

Evaluates the FourVector potential in the VDF model given the rest frame density and the computational frame baryon current.

Parameters
[in]rhoBrest frame baryon density, in fm \(^{-3}\)
[in]jmuB_netnet baryon current in the computational frame, in fm \(^{-3}\)
Returns
VDF potential

\[A^{\mu} = 10^{-3}\times \sum_i C_i \left(\frac{\rho}{\rho_0}\right)^{b_i - 2} \frac{j^{\mu}}{\rho_0}\]

in GeV

Definition at line 263 of file potentials.cc.

263  {
264  // this needs to be used in order to prevent trying to calculate something
265  // like
266  // (-rho_B)^{3.4}
267  const int sgn = rhoB > 0 ? 1 : -1;
268  double abs_rhoB = std::abs(rhoB);
269  // to prevent NAN expressions
270  if (abs_rhoB < very_small_double) {
271  abs_rhoB = very_small_double;
272  }
273  // F_2 is a multiplicative factor in front of the baryon current
274  // in the VDF potential
275  double F_2 = 0.0;
276  for (int i = 0; i < number_of_terms(); i++) {
277  F_2 += coeffs_[i] * std::pow(abs_rhoB, powers_[i] - 2.0) /
278  std::pow(saturation_density_, powers_[i] - 1.0);
279  }
280  F_2 = F_2 * sgn;
281  // Return in GeV
282  return F_2 * jmuB_net;
283 }
std::vector< double > powers_
Parameters of the VDF potential: exponents .
Definition: potentials.h:399
std::vector< double > coeffs_
Parameters of the VDF potential: coefficients , in GeV.
Definition: potentials.h:397
int number_of_terms() const
Definition: potentials.h:329
double saturation_density_
Saturation density of nuclear matter used in the VDF potential; it may vary between different paramet...
Definition: potentials.h:395
constexpr double very_small_double
A very small double, used to avoid division by zero.
Definition: constants.h:40
Here is the call graph for this function:

◆ potential()

double smash::Potentials::potential ( const ThreeVector r,
const ParticleList &  plist,
const ParticleType acts_on 
) const

Evaluates potential (Skyrme with optional Symmetry or VDF) at point r.

For Skyrme and Symmetry options, potential is always taken in the local Eckart rest frame, but point r is in the computational frame.

Parameters
[in]rArbitrary space point where potential is calculated
[in]plistList of all particles to be used in \(j^{\mu}\) calculation. If the distance between particle and calculation point r, \( |r-r_i| > r_{cut} \) then particle input to density will be ignored.
[in]acts_onType of particle on which potential is going to act. It gives the charges (or more precisely, the scaling factors) of the particle moving in the potential field.
Returns
Total potential energy acting on the particle: for Skyrme and Symmetry potentials,

\[U_{\rm tot} =Q_BU_B+2I_3U_I\]

in GeV, while for the VDF potential

\[U_{\rm tot} =Q_B A^0\]

in GeV, where \(Q_B\) is the baryon charge scaled by the ratio of the light (u, d) quark to the total quark number and \(I_3\) is the third compnent of the isospin.

Definition at line 285 of file potentials.cc.

286  {
287  double total_potential = 0.0;
288  const bool compute_gradient = false;
289  const bool smearing = true;
290  const auto scale = force_scale(acts_on);
291 
292  if (!(acts_on.is_baryon() || acts_on.is_nucleus())) {
293  return total_potential;
294  }
295  const auto baryon_density_and_gradient = current_eckart(
296  r, plist, param_, DensityType::Baryon, compute_gradient, smearing);
297  const double rhoB = std::get<0>(baryon_density_and_gradient);
298  if (use_skyrme_) {
299  total_potential += scale.first * skyrme_pot(rhoB);
300  }
301  if (use_symmetry_) {
302  const double rho_iso = std::get<0>(
304  compute_gradient, smearing));
305  const double sym_pot = symmetry_pot(rho_iso, rhoB) * acts_on.isospin3_rel();
306  total_potential += scale.second * sym_pot;
307  }
308 
309  if (use_vdf_) {
310  const FourVector jmuB = std::get<1>(baryon_density_and_gradient);
311  const FourVector VDF_potential = vdf_pot(rhoB, jmuB);
312  total_potential += scale.first * VDF_potential.x0();
313  }
314 
315  return total_potential;
316 }
bool use_skyrme_
Skyrme potential on/off.
Definition: potentials.h:343
bool use_symmetry_
Symmetry potential on/off.
Definition: potentials.h:346
double skyrme_pot(const double baryon_density) const
Evaluates skyrme potential given a baryon density.
Definition: potentials.cc:234
double symmetry_pot(const double baryon_isospin_density, const double baryon_density) const
Evaluates symmetry potential given baryon isospin density.
Definition: potentials.cc:252
FourVector vdf_pot(double rhoB, const FourVector jmuB_net) const
Evaluates the FourVector potential in the VDF model given the rest frame density and the computationa...
Definition: potentials.cc:263
static std::pair< double, int > force_scale(const ParticleType &data)
Evaluates the scaling factor of the forces acting on the particles.
Definition: potentials.cc:318
const DensityParameters param_
Struct that contains the gaussian smearing width , the distance cutoff and the testparticle number n...
Definition: potentials.h:340
bool use_vdf_
VDF potential on/off.
Definition: potentials.h:352
std::tuple< double, FourVector, ThreeVector, ThreeVector, FourVector, FourVector, FourVector, FourVector > current_eckart(const ThreeVector &r, const ParticleList &plist, const DensityParameters &par, DensityType dens_type, bool compute_gradient, bool smearing)
Calculates Eckart rest frame density and 4-current of a given density type and optionally the gradien...
Definition: density.cc:167
Here is the call graph for this function:

◆ force_scale()

std::pair< double, int > smash::Potentials::force_scale ( const ParticleType data)
static

Evaluates the scaling factor of the forces acting on the particles.

The forces are equal to the product of the scaling factor and the gradient of the potential. We need these scaling factors to describe the motions of the hyperons as well as the anti-particles in the potentials. For Lambda and Sigma, since they carry 2 light (u or d) quarks, they are affected by 2/3 of the Skyrme force. Xi carries 1 light quark, it is affected by 1/3 of the Skyrme force. Omega carries no light quark, so it's not affected by the Skyrme force. Anti-baryons are affected by the force as large as the force acting on baryons but with an opposite direction.

Parameters
[in]dataType of particle on which potential is going to act.
Returns
( \(Q_B(1-\frac{|Q_S|}{3}), Q_B\)) where \(Q_B\) is the baryon charge and \(Q_S\) is the strangeness.

Definition at line 318 of file potentials.cc.

318  {
319  const auto &pdg = data.pdgcode();
320  const double skyrme_or_VDF_scale =
321  (3 - std::abs(pdg.strangeness())) / 3. * pdg.baryon_number();
322  const int symmetry_scale = pdg.baryon_number();
323  return std::make_pair(skyrme_or_VDF_scale, symmetry_scale);
324 }
Here is the call graph for this function:
Here is the caller graph for this function:

◆ skyrme_force()

std::pair< ThreeVector, ThreeVector > smash::Potentials::skyrme_force ( const double  rhoB,
const ThreeVector  grad_j0B,
const ThreeVector  dvecjB_dt,
const ThreeVector  curl_vecjB 
) const

Evaluates the electric and magnetic components of the skyrme force.

Parameters
[in]rhoBEckart baryon density [fm \(^{-3}\)].
[in]grad_j0BGradient of baryon density [fm \(^{-4}\)]. This density is evaluated in the computational frame.
[in]dvecjB_dtTime derivative of the vector baryon current density [fm \(^{-4}\)
[in]curl_vecjBCurl of the baryon vector current density [fm \(^{-4}\)
Returns
( \(E_B, B_B\)), where

\[E_B = -V_B^\prime(\rho^\ast)(\nabla\rho_B + \partial_t \vec j_B)\]

is the electro component of Skyrme force and

\[B_B = V_B^\prime(\rho^\ast) \nabla\times\vec j_B\]

is the magnetic component of the Skyrme force with \(\rho^\ast\) being the Eckart baryon density.

Definition at line 326 of file potentials.cc.

328  {
329  ThreeVector E_component(0.0, 0.0, 0.0), B_component(0.0, 0.0, 0.0);
330  if (use_skyrme_) {
331  const int sgn = rhoB > 0 ? 1 : -1;
332  const double abs_rhoB = std::abs(rhoB);
333  const double dV_drho = sgn *
335  std::pow(abs_rhoB / nuclear_density,
336  skyrme_tau_ - 1)) *
338  E_component -= dV_drho * (grad_j0B + dvecjB_dt);
339  B_component += dV_drho * curl_vecjB;
340  }
341  return std::make_pair(E_component, B_component);
342 }
Here is the call graph for this function:

◆ symmetry_force()

std::pair< ThreeVector, ThreeVector > smash::Potentials::symmetry_force ( const double  rhoI3,
const ThreeVector  grad_j0I3,
const ThreeVector  dvecjI3_dt,
const ThreeVector  curl_vecjI3,
const double  rhoB,
const ThreeVector  grad_j0B,
const ThreeVector  dvecjB_dt,
const ThreeVector  curl_vecjB 
) const

Evaluates the electric and magnetic components of the symmetry force.

Parameters
[in]rhoI3Relative isospin 3 density.
[in]grad_j0I3Gradient of I3/I density [fm \(^{-4}\)]. This density is evaluated in the computational frame.
[in]dvecjI3_dtTime derivative of the I3/I vector current density [fm \(^{-4}\)]
[in]curl_vecjI3Curl of the I3/I vector current density [fm \(^{-4}\)]
[in]rhoBNet-baryon density in the rest frame
[in]grad_j0BGradient of the net-baryon density in the computational frame
[in]dvecjB_dtTime derivative of the net-baryon vector current density
[in]curl_vecjBCurl of the net-baryon vector current density
Returns
( \(E_I3, B_I3\)) [GeV/fm], where

\[\vec{E} = - \frac{\partial V^\ast}{\partial\rho_{I_3}^\ast} (\nabla\rho_{I_3} + \partial_t \vec j_{I_3}) - \frac{\partial V^\ast}{\partial\rho_B^\ast}(\nabla\rho_B + \partial_t \vec j_B)\]

is the electrical component of symmetry force and

\[\vec{B} = \frac{\partial V^\ast}{\rho_{I_3}^\ast} \nabla\times\vec j_{I_3} + \frac{\partial V^\ast}{\rho_B^\ast} \nabla\times\vec j_B \]

is the magnetic component of the symmetry force with \(\rho^\ast\) being the respective Eckart density.

Definition at line 344 of file potentials.cc.

348  {
349  ThreeVector E_component(0.0, 0.0, 0.0), B_component(0.0, 0.0, 0.0);
350  if (use_symmetry_) {
351  E_component -= dVsym_drhoI3(rhoB, rhoI3) * (grad_j0I3 + dvecjI3_dt) +
352  dVsym_drhoB(rhoB, rhoI3) * (grad_j0B + dvecjB_dt);
353  B_component += dVsym_drhoI3(rhoB, rhoI3) * curl_vecjI3 +
354  dVsym_drhoB(rhoB, rhoI3) * curl_vecjB;
355  }
356  return std::make_pair(E_component, B_component);
357 }
double dVsym_drhoB(const double rhoB, const double rhoI3) const
Calculate the derivative of the symmetry potential with respect to the net baryon density in GeV * fm...
Definition: potentials.cc:420
double dVsym_drhoI3(const double rhoB, const double rhoI3) const
Calculate the derivative of the symmetry potential with respect to the isospin density in GeV * fm^3.
Definition: potentials.cc:410

◆ E_field_integrand()

static ThreeVector smash::Potentials::E_field_integrand ( ThreeVector  pos,
DensityOnLattice charge_density,
ThreeVector  point 
)
inlinestatic

Integrand for calculating the electric field.

The field is calculaed via

\[\vec{E}(\vec{r})= \int\frac{(\vec{r}-\vec{r}')\rho(\vec{r}')}{|\vec{r}-\vec{r}'|^3}d^3r' \]

Parameters
[in]posposition vector to be integrated over
[in]charge_densityelectric charge density at position pos
[in]pointposition where to calculate the field

Definition at line 204 of file potentials.h.

206  {
207  ThreeVector dr = point - pos;
208  if (dr.abs() < really_small) {
209  return {0., 0., 0.};
210  }
211  return elementary_charge * charge_density.rho() * dr /
212  std::pow(dr.abs(), 3);
213  }
constexpr double really_small
Numerical error tolerance.
Definition: constants.h:37
const double elementary_charge
Elementary electric charge in natural units, approximately 0.3.
Definition: constants.h:104
Here is the call graph for this function:
Here is the caller graph for this function:

◆ B_field_integrand()

static ThreeVector smash::Potentials::B_field_integrand ( ThreeVector  pos,
DensityOnLattice charge_density,
ThreeVector  point 
)
inlinestatic

Integrand for calculating the magnetic field using the Biot-Savart formula.

Parameters
[in]posposition vector to be integrated over
[in]charge_densityelectric charge density and current
[in]pointposition where the magnetic field will be calculated

Definition at line 222 of file potentials.h.

224  {
225  ThreeVector dr = point - pos;
226  if (dr.abs() < really_small) {
227  return {0., 0., 0.};
228  }
229  return elementary_charge *
230  charge_density.jmu_net().threevec().cross_product(dr) /
231  std::pow(dr.abs(), 3);
232  }
Here is the call graph for this function:
Here is the caller graph for this function:

◆ vdf_force() [1/2]

std::pair< ThreeVector, ThreeVector > smash::Potentials::vdf_force ( double  rhoB,
const double  drhoB_dt,
const ThreeVector  grad_rhoB,
const ThreeVector  gradrhoB_cross_vecjB,
const double  j0B,
const ThreeVector  grad_j0B,
const ThreeVector  vecjB,
const ThreeVector  dvecjB_dt,
const ThreeVector  curl_vecjB 
) const

Evaluates the electric and magnetic components of force in the VDF model given the derivatives of the baryon current \(j^{\mu}\).

Parameters
[in]rhoBrest frame baryon density in fm \(^{-3}\)
[in]drhoB_dttime derivative of the rest frame density
[in]grad_rhoBgradient of the rest frame density
[in]gradrhoB_cross_vecjBcross product of the gradient of the rest frame density and the 3-vector baryon current density
[in]j0Bcomputational frame baryon density in fm \(^{-3}\)
[in]grad_j0Bgradient of the computational frame baryon density
[in]vecjB3-vector baryon current
[in]dvecjB_dttime derivative of the computational frame 3-vector baryon current
[in]curl_vecjBcurl of the 3-vector baryon current
Returns
( \(E_{VDF}, B_{VDF}\)) [GeV/fm], where

\[\vec{E}_{VDF} = - F_1 \big((\vec{\nabla} \rho) j^0 + (\partial_t \rho) \vec{j}\big) - F_2 (\vec{\nabla} j^0 + \partial_t\vec{j})\]

is the electrical component of VDF force and

\[\vec{B}_{VDF} = F_1 (\vec{\nabla} \rho) \times \vec{j} + F_2 \vec{\nabla} \times \vec{j}\]

is the magnetic component of the VDF force, with

\[F_1 = \sum_i C_i (b_i - 2) \frac{\rho^{b_i - 3}}{\rho_0^{b_i - 1}} \,, \]

\[F_2 = \sum_i C_i \frac{\rho^{b_i - 2}}{\rho_0^{b_i - 1}} \,, \]

where \(\rho_0\) is the saturation density.

Definition at line 359 of file potentials.cc.

363  {
364  // this needs to be used to prevent trying to calculate something like
365  // (-rhoB)^{3.4}
366  const int sgn = rhoB > 0 ? 1 : -1;
367  ThreeVector E_component(0.0, 0.0, 0.0), B_component(0.0, 0.0, 0.0);
368  // to prevent NAN expressions
369  double abs_rhoB = std::abs(rhoB);
370  if (abs_rhoB < very_small_double) {
371  abs_rhoB = very_small_double;
372  }
373  if (use_vdf_) {
374  // F_1 and F_2 are multiplicative factors in front of the baryon current
375  // in the VDF potential
376  double F_1 = 0.0;
377  for (int i = 0; i < number_of_terms(); i++) {
378  F_1 += coeffs_[i] * (powers_[i] - 2.0) *
379  std::pow(abs_rhoB, powers_[i] - 3.0) /
380  std::pow(saturation_density_, powers_[i] - 1.0);
381  }
382  F_1 = F_1 * sgn;
383 
384  double F_2 = 0.0;
385  for (int i = 0; i < number_of_terms(); i++) {
386  F_2 += coeffs_[i] * std::pow(abs_rhoB, powers_[i] - 2.0) /
387  std::pow(saturation_density_, powers_[i] - 1.0);
388  }
389  F_2 = F_2 * sgn;
390 
391  E_component -= (F_1 * (grad_rhoB * j0B + drhoB_dt * vecjB) +
392  F_2 * (grad_j0B + dvecjB_dt));
393  B_component += F_1 * gradrhoB_cross_vecjB + F_2 * curl_vecjB;
394  }
395  return std::make_pair(E_component, B_component);
396 }
Here is the call graph for this function:

◆ vdf_force() [2/2]

std::pair< ThreeVector, ThreeVector > smash::Potentials::vdf_force ( const ThreeVector  grad_A_0,
const ThreeVector  dA_dt,
const ThreeVector  curl_vecA 
) const

Evaluates the electric and magnetic components of force in the VDF force given the derivatives of the VDF mean-field \(A^\mu\).

Parameters
[in]grad_A_0gradient of the zeroth component of the field A^mu
[in]dA_dttime derivative of the field A^mu
[in]curl_vecAcurl of the vector component of the field A^mu
Returns
( \(E_{VDF}, B_{VDF}\)) [GeV/fm], where

\[\vec{E}_{VDF} = - \vec{\nabla} A^0 - \partial_t\vec{A}\]

is the electrical component of VDF force and

\[\vec{B}_{VDF} = \vec{\nabla} \times \vec{A}\]

is the magnetic component of the VDF force.

Definition at line 399 of file potentials.cc.

401  {
402  ThreeVector E_component(0.0, 0.0, 0.0), B_component(0.0, 0.0, 0.0);
403  if (use_vdf_) {
404  E_component -= (grad_A_0 + dA_dt);
405  B_component += curl_A;
406  }
407  return std::make_pair(E_component, B_component);
408 }

◆ all_forces()

std::tuple< ThreeVector, ThreeVector, ThreeVector, ThreeVector > smash::Potentials::all_forces ( const ThreeVector r,
const ParticleList &  plist 
) const
virtual

Evaluates the electric and magnetic components of the forces at point r.

Point r is in the computational frame.

Parameters
[in]rArbitrary space point where potential gradient is calculated
[in]plistList of all particles to be used in \(j^{\mu}\) calculation. If the distance between particle and calculation point r, \( |r-r_i| > r_{cut} \) then particle input to density will be ignored.
Returns
( \(E_B, B_B, E_{I3}, B_{I3}\)) [GeV/fm], where \(E_B\): the electric component of the Skyrme or VDF force, \(B_B\): the magnetic component of the Skyrme or VDF force, \(E_{I3}\): the electric component of the symmetry force, \(B_{I3}\): the magnetic component of the symmetry force

Definition at line 434 of file potentials.cc.

434  {
435  const bool compute_gradient = true;
436  const bool smearing = true;
437  auto F_skyrme_or_VDF =
438  std::make_pair(ThreeVector(0., 0., 0.), ThreeVector(0., 0., 0.));
439  auto F_symmetry =
440  std::make_pair(ThreeVector(0., 0., 0.), ThreeVector(0., 0., 0.));
441 
442  const auto baryon_density_and_gradient = current_eckart(
443  r, plist, param_, DensityType::Baryon, compute_gradient, smearing);
444  double rhoB = std::get<0>(baryon_density_and_gradient);
445  const ThreeVector grad_j0B = std::get<2>(baryon_density_and_gradient);
446  const ThreeVector curl_vecjB = std::get<3>(baryon_density_and_gradient);
447  const FourVector djmuB_dt = std::get<4>(baryon_density_and_gradient);
448  if (use_skyrme_) {
449  F_skyrme_or_VDF =
450  skyrme_force(rhoB, grad_j0B, djmuB_dt.threevec(), curl_vecjB);
451  }
452 
453  if (use_symmetry_) {
454  const auto density_and_gradient =
456  compute_gradient, smearing);
457  const double rhoI3 = std::get<0>(density_and_gradient);
458  const ThreeVector grad_j0I3 = std::get<2>(density_and_gradient);
459  const ThreeVector curl_vecjI3 = std::get<3>(density_and_gradient);
460  const FourVector dvecjI3_dt = std::get<4>(density_and_gradient);
461  F_symmetry =
462  symmetry_force(rhoI3, grad_j0I3, dvecjI3_dt.threevec(), curl_vecjI3,
463  rhoB, grad_j0B, djmuB_dt.threevec(), curl_vecjB);
464  }
465 
466  if (use_vdf_) {
467  const FourVector jmuB = std::get<1>(baryon_density_and_gradient);
468  const FourVector djmuB_dx = std::get<5>(baryon_density_and_gradient);
469  const FourVector djmuB_dy = std::get<6>(baryon_density_and_gradient);
470  const FourVector djmuB_dz = std::get<7>(baryon_density_and_gradient);
471 
472  // safety check to not divide by zero
473  const int sgn = rhoB > 0 ? 1 : -1;
474  if (std::abs(rhoB) < very_small_double) {
475  rhoB = sgn * very_small_double;
476  }
477 
478  const double drhoB_dt =
479  (1 / rhoB) * (jmuB.x0() * djmuB_dt.x0() - jmuB.x1() * djmuB_dt.x1() -
480  jmuB.x2() * djmuB_dt.x2() - jmuB.x3() * djmuB_dt.x3());
481 
482  const double drhoB_dx =
483  (1 / rhoB) * (jmuB.x0() * djmuB_dx.x0() - jmuB.x1() * djmuB_dx.x1() -
484  jmuB.x2() * djmuB_dx.x2() - jmuB.x3() * djmuB_dx.x3());
485 
486  const double drhoB_dy =
487  (1 / rhoB) * (jmuB.x0() * djmuB_dy.x0() - jmuB.x1() * djmuB_dy.x1() -
488  jmuB.x2() * djmuB_dy.x2() - jmuB.x3() * djmuB_dy.x3());
489 
490  const double drhoB_dz =
491  (1 / rhoB) * (jmuB.x0() * djmuB_dz.x0() - jmuB.x1() * djmuB_dz.x1() -
492  jmuB.x2() * djmuB_dz.x2() - jmuB.x3() * djmuB_dz.x3());
493 
494  const FourVector drhoB_dxnu = {drhoB_dt, drhoB_dx, drhoB_dy, drhoB_dz};
495 
496  const ThreeVector grad_rhoB = drhoB_dxnu.threevec();
497  const ThreeVector vecjB = jmuB.threevec();
498  const ThreeVector Drho_cross_vecj = grad_rhoB.cross_product(vecjB);
499 
500  F_skyrme_or_VDF = vdf_force(
501  rhoB, drhoB_dt, drhoB_dxnu.threevec(), Drho_cross_vecj, jmuB.x0(),
502  grad_j0B, jmuB.threevec(), djmuB_dt.threevec(), curl_vecjB);
503  }
504 
505  return std::make_tuple(F_skyrme_or_VDF.first, F_skyrme_or_VDF.second,
506  F_symmetry.first, F_symmetry.second);
507 }
std::pair< ThreeVector, ThreeVector > symmetry_force(const double rhoI3, const ThreeVector grad_j0I3, const ThreeVector dvecjI3_dt, const ThreeVector curl_vecjI3, const double rhoB, const ThreeVector grad_j0B, const ThreeVector dvecjB_dt, const ThreeVector curl_vecjB) const
Evaluates the electric and magnetic components of the symmetry force.
Definition: potentials.cc:344
std::pair< ThreeVector, ThreeVector > vdf_force(double rhoB, const double drhoB_dt, const ThreeVector grad_rhoB, const ThreeVector gradrhoB_cross_vecjB, const double j0B, const ThreeVector grad_j0B, const ThreeVector vecjB, const ThreeVector dvecjB_dt, const ThreeVector curl_vecjB) const
Evaluates the electric and magnetic components of force in the VDF model given the derivatives of the...
Definition: potentials.cc:359
std::pair< ThreeVector, ThreeVector > skyrme_force(const double rhoB, const ThreeVector grad_j0B, const ThreeVector dvecjB_dt, const ThreeVector curl_vecjB) const
Evaluates the electric and magnetic components of the skyrme force.
Definition: potentials.cc:326
Here is the call graph for this function:
Here is the caller graph for this function:

◆ use_skyrme()

virtual bool smash::Potentials::use_skyrme ( ) const
inlinevirtual
Returns
Is Skyrme potential on?

Definition at line 305 of file potentials.h.

305 { return use_skyrme_; }
Here is the caller graph for this function:

◆ use_symmetry()

virtual bool smash::Potentials::use_symmetry ( ) const
inlinevirtual
Returns
Is symmetry potential on?

Definition at line 307 of file potentials.h.

307 { return use_symmetry_; }
Here is the caller graph for this function:

◆ use_coulomb()

virtual bool smash::Potentials::use_coulomb ( ) const
inlinevirtual
Returns
Is Coulomb potential on?

Definition at line 309 of file potentials.h.

309 { return use_coulomb_; }
bool use_coulomb_
Coulomb potential on/Off.
Definition: potentials.h:349
Here is the caller graph for this function:

◆ skyrme_a()

double smash::Potentials::skyrme_a ( ) const
inline
Returns
Skyrme parameter skyrme_a, in MeV

Definition at line 312 of file potentials.h.

312 { return skyrme_a_; }
Here is the caller graph for this function:

◆ skyrme_b()

double smash::Potentials::skyrme_b ( ) const
inline
Returns
Skyrme parameter skyrme_b, in MeV

Definition at line 314 of file potentials.h.

314 { return skyrme_b_; }
Here is the caller graph for this function:

◆ skyrme_tau()

double smash::Potentials::skyrme_tau ( ) const
inline
Returns
Skyrme parameter skyrme_tau

Definition at line 316 of file potentials.h.

316 { return skyrme_tau_; }
Here is the caller graph for this function:

◆ symmetry_S_pot()

double smash::Potentials::symmetry_S_pot ( ) const
inline
Returns
Skyrme parameter S_pot, in MeV

Definition at line 318 of file potentials.h.

318 { return symmetry_S_Pot_; }

◆ use_vdf()

virtual bool smash::Potentials::use_vdf ( ) const
inlinevirtual
Returns
Is VDF potential on?

Definition at line 321 of file potentials.h.

321 { return use_vdf_; }
Here is the caller graph for this function:

◆ saturation_density()

double smash::Potentials::saturation_density ( ) const
inline
Returns
Value of the saturation density used in the VDF potential

Definition at line 323 of file potentials.h.

323 { return saturation_density_; }
Here is the caller graph for this function:

◆ coeffs()

const std::vector<double>& smash::Potentials::coeffs ( ) const
inline
Returns
Vector of the VDF coefficients \(C_i\), coefficients_

Definition at line 325 of file potentials.h.

325 { return coeffs_; }
Here is the caller graph for this function:

◆ powers()

const std::vector<double>& smash::Potentials::powers ( ) const
inline
Returns
Vector of the VDF exponents \(b_i\), powers_

Definition at line 327 of file potentials.h.

327 { return powers_; }
Here is the caller graph for this function:

◆ number_of_terms()

int smash::Potentials::number_of_terms ( ) const
inline
Returns
Number of terms in the VDF potential

Definition at line 329 of file potentials.h.

329 { return powers_.size(); }
Here is the caller graph for this function:

◆ coulomb_r_cut()

double smash::Potentials::coulomb_r_cut ( ) const
inline
Returns
cutoff radius in ntegration for coulomb potential in fm

Definition at line 332 of file potentials.h.

332 { return coulomb_r_cut_; }
double coulomb_r_cut_
Cutoff in integration for coulomb potential.
Definition: potentials.h:389

◆ dVsym_drhoI3()

double smash::Potentials::dVsym_drhoI3 ( const double  rhoB,
const double  rhoI3 
) const
private

Calculate the derivative of the symmetry potential with respect to the isospin density in GeV * fm^3.

\[ \frac{\partial V_\mathrm{sym}}{\partial \rho_{I_3}} = 2\frac{S_\mathrm{Pot}}{\rho_0} + \frac{2\rho_{I_3}\left[12.3\left(\frac{\rho_B}{\rho_0}\right)^{2/3} + 20 \left(\frac{\rho_B}{\rho_0}\right)^\gamma\right]}{\rho_B^2} \]

Note
The isospin 3 density here is actually the density of I3 / I.
Parameters
[in]rhoBnet baryon density
[in]rhoI3isospin density
Returns
partial derivative of the symmetry potenital with respect to the isospin density.

Definition at line 410 of file potentials.cc.

410  {
411  double term1 = 2. * symmetry_S_Pot_ / nuclear_density;
413  double term2 = 2. * rhoI3 * symmetry_S(rhoB) / (rhoB * rhoB);
414  return mev_to_gev * (term1 + term2);
415  } else {
416  return mev_to_gev * term1;
417  }
418 }

◆ dVsym_drhoB()

double smash::Potentials::dVsym_drhoB ( const double  rhoB,
const double  rhoI3 
) const
private

Calculate the derivative of the symmetry potential with respect to the net baryon density in GeV * fm^3.

\[ \frac{\partial V_\mathrm{sym}}{\partial \rho_B} = \left(\frac{\rho_{I_3}}{\rho_B}\right)^2 \left[\frac{8.2}{\rho_0}\left(\frac{\rho_B}{\rho_0}\right)^{-1/3} + \frac{20\gamma}{\rho_B}\left(\frac{\rho_B}{\rho_0}\right)^\gamma\right] -2\frac{\rho_{I_3}^2}{\rho_B^3} \left[12.3\left(\frac{\rho_B}{\rho_0}\right)^{2/3} + 20\left(\frac{\rho_B}{\rho_0}\right)^\gamma\right]\]

Note
The isospin 3 density here is actually the density of I3 / I
Parameters
[in]rhoBnet baryon density
[in]rhoI3isospin density
Returns
partial derivative of the symmetry potenital with respect to the net baryon density.

Definition at line 420 of file potentials.cc.

420  {
422  double rhoB_over_rho0 = rhoB / nuclear_density;
423  double term1 = 8.2 * std::pow(rhoB_over_rho0, -1. / 3.) / nuclear_density +
424  20. * symmetry_gamma_ *
425  std::pow(rhoB_over_rho0, symmetry_gamma_) / rhoB;
426  double term2 = -2. * symmetry_S(rhoB) / rhoB;
427  return mev_to_gev * (term1 + term2) * rhoI3 * rhoI3 / (rhoB * rhoB);
428  } else {
429  return 0.;
430  }
431 }

Member Data Documentation

◆ param_

const DensityParameters smash::Potentials::param_
private

Struct that contains the gaussian smearing width \(\sigma\), the distance cutoff \(r_{\rm cut}\) and the testparticle number needed for the density calculation.

Definition at line 340 of file potentials.h.

◆ use_skyrme_

bool smash::Potentials::use_skyrme_
private

Skyrme potential on/off.

Definition at line 343 of file potentials.h.

◆ use_symmetry_

bool smash::Potentials::use_symmetry_
private

Symmetry potential on/off.

Definition at line 346 of file potentials.h.

◆ use_coulomb_

bool smash::Potentials::use_coulomb_
private

Coulomb potential on/Off.

Definition at line 349 of file potentials.h.

◆ use_vdf_

bool smash::Potentials::use_vdf_
private

VDF potential on/off.

Definition at line 352 of file potentials.h.

◆ skyrme_a_

double smash::Potentials::skyrme_a_
private

Parameter of skyrme potentials: the coefficient in front of \(\frac{\rho}{\rho_0}\) in GeV.

Definition at line 358 of file potentials.h.

◆ skyrme_b_

double smash::Potentials::skyrme_b_
private

Parameters of skyrme potentials: the coefficient in front of \((\frac{\rho}{\rho_0})^\tau\) in GeV.

Definition at line 364 of file potentials.h.

◆ skyrme_tau_

double smash::Potentials::skyrme_tau_
private

Parameters of skyrme potentials: the power index.

Definition at line 370 of file potentials.h.

◆ symmetry_S_Pot_

double smash::Potentials::symmetry_S_Pot_
private

Parameter S_Pot in the symmetry potential in MeV.

Definition at line 373 of file potentials.h.

◆ symmetry_is_rhoB_dependent_

bool smash::Potentials::symmetry_is_rhoB_dependent_ = false
private

Whether the baryon density dependence of the symmetry potential is included.

Definition at line 379 of file potentials.h.

◆ symmetry_gamma_

double smash::Potentials::symmetry_gamma_
private

Power \( \gamma \) in formula for \( S(\rho) \):

\[ S(\rho)=12.3\,\mathrm{MeV}\times \left(\frac{\rho}{\rho_0}\right)^{2/3}+20\,\mathrm{MeV}\times \left(\frac{\rho}{\rho_0}\right)^\gamma \]

Definition at line 386 of file potentials.h.

◆ coulomb_r_cut_

double smash::Potentials::coulomb_r_cut_
private

Cutoff in integration for coulomb potential.

Definition at line 389 of file potentials.h.

◆ saturation_density_

double smash::Potentials::saturation_density_
private

Saturation density of nuclear matter used in the VDF potential; it may vary between different parameterizations.

Definition at line 395 of file potentials.h.

◆ coeffs_

std::vector<double> smash::Potentials::coeffs_
private

Parameters of the VDF potential: coefficients \(C_i\), in GeV.

Definition at line 397 of file potentials.h.

◆ powers_

std::vector<double> smash::Potentials::powers_
private

Parameters of the VDF potential: exponents \(b_i\).

Definition at line 399 of file potentials.h.


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