Thermal-FIST 1.5
Package for hadron resonance gas model applications
Loading...
Searching...
No Matches
thermalfist::ThermalModelVDW Class Reference

Class implementing the quantum van der Waals HRG model. More...

#include <ThermalModelVDW.h>

Inheritance diagram for thermalfist::ThermalModelVDW:
thermalfist::ThermalModelBase thermalfist::ThermalModelEVCrossterms

Public Member Functions

 ThermalModelVDW (ThermalParticleSystem *TPS_, const ThermalModelParameters &params=ThermalModelParameters())
 Construct a new ThermalModelVDW object.
 
virtual ~ThermalModelVDW (void)
 Destroy the ThermalModelVDW object.
 
void FillVirialEV (const std::vector< std::vector< double > > &bij=std::vector< std::vector< double > >(0))
 Same as FillVirial() but uses the matrix of excluded-volume coefficients \( v_i \equiv b_{ii} \) as input instead of radii.
 
void SetVirialdT (int i, int j, double dbdT)
 Set the temperature derivative of the eigenvolume parameter \( \tilde{b}_{ij} \).
 
void SetAttractiondT (int i, int j, double dadT)
 Set the temperature derivative of the QvdW attraction parameter \( a_{ij} \).
 
double VirialCoefficientdT (int i, int j) const
 The temperature derivative of the eigenvolume parameter \( \tilde{b}_{ij} \).
 
double AttractionCoefficientdT (int i, int j) const
 The temperature derivative of the QvdW attraction parameter \( a_{ij} \).
 
void SetTemperatureDependentAB (bool Tdep)
 Sets whether temperature depedence of QvdW parameters should be considered.
 
bool TemperatureDependentAB () const
 Whether temperature depedence of QvdW parameters is considered.
 
virtual void SetMultipleSolutionsMode (bool search)
 Whether to search for multiple solutions of the QvdW equations by considering different initial guesses in the Broyden's method.
 
bool UseMultipleSolutionsMode () const
 Whether to search for multiple solutions of the QvdW equations by considering different initial guesses in the Broyden's method.
 
double MuStar (int i) const
 The shifted chemical potential of particle species i.
 
std::vector< double > GetMuStar () const
 
void SetMuStar (const std::vector< double > &MuStar)
 Set the vector of shifted chemical potentials.
 
void FillChemicalPotentials ()
 Sets the chemical potentials of all particles.
 
virtual void SetChemicalPotentials (const std::vector< double > &chem=std::vector< double >(0))
 Sets the chemical potentials of all particles.
 
virtual void FillVirial (const std::vector< double > &ri=std::vector< double >(0))
 Fills the excluded volume coefficients \( \tilde{b}_{ij} \) based on the provided radii parameters for all species.
 
virtual void FillAttraction (const std::vector< std::vector< double > > &aij=std::vector< std::vector< double > >(0))
 
virtual void ReadInteractionParameters (const std::string &filename)
 Reads the QvdW interaction parameters from a file.
 
virtual void WriteInteractionParameters (const std::string &filename)
 Write the QvdW interaction parameters to a file.
 
virtual void SetVirial (int i, int j, double b)
 Set the excluded volume coefficient \( \tilde{b}_{ij} \).
 
virtual void SetAttraction (int i, int j, double a)
 Set the vdW mean field attraction coefficient \( a_{ij} \).
 
double VirialCoefficient (int i, int j) const
 Excluded volume coefficient \( \tilde{b}_{ij} = 0 \).
 
double AttractionCoefficient (int i, int j) const
 QvdW mean field attraction coefficient \( a_{ij} \).
 
virtual void ChangeTPS (ThermalParticleSystem *TPS)
 Change the particle list.
 
virtual void CalculatePrimordialDensities ()
 Calculates the primordial densities of all species.
 
virtual std::vector< double > CalculateChargeFluctuations (const std::vector< double > &chgs, int order=4)
 Calculates fluctuations (diagonal susceptibilities) of an arbitrary "conserved" charge.
 
virtual std::vector< std::vector< double > > CalculateFluctuations (int order)
 
void CalculateTwoParticleCorrelations ()
 Computes the fluctuations and correlations of the primordial particle numbers.
 
void CalculateFluctuations ()
 Computes the fluctuation observables.
 
virtual double CalculatePressure ()
 
virtual double CalculateEnergyDensity ()
 
virtual double CalculateEntropyDensity ()
 
virtual double CalculateEnergyDensityDerivativeT ()
 
virtual void CalculateTemperatureDerivatives ()
 Computes the temperature derivatives of densities, shifted chemical potentials, and primordial hadron number susceptibilities.
 
virtual double CalculateBaryonMatterEntropyDensity ()
 
virtual double CalculateMesonMatterEntropyDensity ()
 
virtual double ParticleScalarDensity (int part)
 
bool IsLastSolutionOK () const
 
double DensityId (int index)
 
const std::vector< std::vector< int > > & VDWComponentIndices () const
 
virtual double DeltaMu (int i) const
 
const std::vector< std::vector< double > > & VirialMatrix () const
 
const std::vector< std::vector< double > > & AttractionMatrix () const
 
- Public Member Functions inherited from thermalfist::ThermalModelBase
 ThermalModelBase (ThermalParticleSystem *TPS, const ThermalModelParameters &params=ThermalModelParameters())
 Construct a new ThermalModelBase object.
 
virtual ~ThermalModelBase (void)
 
int ComponentsNumber () const
 Number of different particle species in the list.
 
bool UseWidth () const
 Whether finite resonance widths are considered.
 
void SetUseWidth (bool useWidth)
 Sets whether finite resonance widths are used. Deprecated.
 
void SetUseWidth (ThermalParticle::ResonanceWidthIntegration type)
 Sets the finite resonance widths scheme to use.
 
void SetNormBratio (bool normBratio)
 Whether branching ratios are renormalized to 100%.
 
bool NormBratio () const
 
void SetOMP (bool openMP)
 OpenMP support. Currently not used.
 
virtual void SetParameters (const ThermalModelParameters &params)
 The thermal parameters.
 
const ThermalModelParametersParameters () const
 
void UpdateParameters ()
 Calls SetParameters() with current m_Parameters.
 
virtual void SetTemperature (double T)
 Set the temperature.
 
virtual void SetBaryonChemicalPotential (double muB)
 Set the baryon chemical potential.
 
virtual void SetElectricChemicalPotential (double muQ)
 Set the electric chemical potential.
 
virtual void SetStrangenessChemicalPotential (double muS)
 Set the strangeness chemical potential.
 
virtual void SetCharmChemicalPotential (double muC)
 Set the charm chemical potential.
 
virtual void SetGammaq (double gammaq)
 Set the light quark fugacity factor.
 
virtual void SetGammaS (double gammaS)
 Set the strange quark fugacity factor.
 
virtual void SetGammaC (double gammaC)
 Set the charm quark fugacity factor.
 
virtual void SetBaryonCharge (int B)
 Set the total baryon number (for canonical ensemble only)
 
virtual void SetElectricCharge (int Q)
 Set the total electric charge (for canonical ensemble only)
 
virtual void SetStrangeness (int S)
 Set the total strangeness (for canonical ensemble only)
 
virtual void SetCharm (int C)
 Set the total charm (for canonical ensemble only)
 
virtual void SetRadius (double)
 Set the same excluded volume radius parameter for all species.
 
virtual void SetRadius (int, double)
 Set the radius parameter for particle species i.
 
virtual void SetRepulsion (int i, int j, double b)
 Same as SetVirial() but with a more clear name on what is actually does.
 
virtual void DisableMesonMesonVirial ()
 
void DisableMesonMesonRepulsion ()
 
virtual void DisableMesonMesonAttraction ()
 
virtual void DisableMesonBaryonVirial ()
 
void DisableMesonBaryonRepulsion ()
 
virtual void DisableMesonBaryonAttraction ()
 
virtual void DisableBaryonBaryonVirial ()
 
void DisableBaryonBaryonRepulsion ()
 
virtual void DisableBaryonBaryonAttraction ()
 
virtual void DisableBaryonAntiBaryonVirial ()
 
void DisableBaryonAntiBaryonRepulsion ()
 
virtual void DisableBaryonAntiBaryonAttraction ()
 
double RepulsionCoefficient (int i, int j) const
 
bool QuantumStatistics () const
 
virtual void SetStatistics (bool stats)
 
virtual void SetCalculationType (IdealGasFunctions::QStatsCalculationType type)
 Sets the CalculationType() method to evaluate quantum statistics. Calls the corresponding method in TPS().
 
virtual void SetClusterExpansionOrder (int order)
 Set the number of terms in the cluster expansion method. Calls the corresponding method in TPS().
 
void SetResonanceWidthShape (ThermalParticle::ResonanceWidthShape shape)
 Set the ThermalParticle::ResonanceWidthShape for all particles. Calls the corresponding method in TPS().
 
void SetResonanceWidthIntegrationType (ThermalParticle::ResonanceWidthIntegration type)
 Set the ThermalParticle::ResonanceWidthIntegration scheme for all particles. Calls the corresponding method in TPS().
 
const std::vector< double > & ChemicalPotentials () const
 A vector of chemical potentials of all particles.
 
double ChemicalPotential (int i) const
 Chemical potential of particle species i.
 
virtual void SetChemicalPotential (int i, double chem)
 Sets the chemical potential of particle species i.
 
virtual double FullIdealChemicalPotential (int i) const
 Chemical potential entering the ideal gas expressions of particle species i.
 
bool ConstrainMuB () const
 
void ConstrainMuB (bool constrain)
 
bool ConstrainMuQ () const
 
void ConstrainMuQ (bool constrain)
 
bool ConstrainMuS () const
 
void ConstrainMuS (bool constrain)
 
bool ConstrainMuC () const
 
void ConstrainMuC (bool constrain)
 
void UsePartialChemicalEquilibrium (bool usePCE)
 Sets whether partial chemical equilibrium with additional chemical potentials is used.
 
bool UsePartialChemicalEquilibrium ()
 Whether partial chemical equilibrium with additional chemical potentials is used.
 
void SetSoverB (double SB)
 The entropy per baryon ratio to be used to constrain the baryon chemical potential.
 
double SoverB () const
 
void SetQoverB (double QB)
 The electric-to-baryon charge ratio to be used to constrain the electric chemical potential.
 
double QoverB () const
 
void SetVolume (double Volume)
 Sets the system volume.
 
double Volume () const
 System volume (fm \(^3\))
 
void SetVolumeRadius (double R)
 Sets the system radius.
 
double CanonicalVolume () const
 The canonical correlation volume V \(_c\) (fm \(^3\))
 
void SetCanonicalVolume (double Volume)
 Set the canonical correlation volume V \(_c\).
 
void SetCanonicalVolumeRadius (double Rc)
 Set the canonical correlation system radius.
 
void ConstrainChemicalPotentials (bool resetInitialValues=true)
 Constrains the chemical potentials \( \mu_B,\,\mu_Q,\,\mu_S,\,\mu_C \) by the conservation laws imposed.
 
virtual void FixParameters ()
 Method which actually implements ConstrainChemicalPotentials() (for backward compatibility).
 
virtual void FixParametersNoReset ()
 Method which actually implements ConstrainChemicalPotentialsNoReset() (for backward compatibility).
 
virtual bool SolveChemicalPotentials (double totB=0., double totQ=0., double totS=0., double totC=0., double muBinit=0., double muQinit=0., double muSinit=0., double muCinit=0., bool ConstrMuB=true, bool ConstrMuQ=true, bool ConstrMuS=true, bool ConstrMuC=true)
 The procedure which calculates the chemical potentials \( \mu_B,\,\mu_Q,\,\mu_S,\,\mu_C \) which reproduce the specified total baryon, electric, strangeness, and charm charges of the system.
 
virtual bool FixChemicalPotentialsThroughDensities (double rhoB=0., double rhoQ=0., double rhoS=0., double rhoC=0., double muBinit=0., double muQinit=0., double muSinit=0., double muCinit=0., bool ConstrMuB=true, bool ConstrMuQ=true, bool ConstrMuS=true, bool ConstrMuC=true)
 The procedure which calculates the chemical potentials \( \mu_B,\,\mu_Q,\,\mu_S,\,\mu_C \) which reproduce the specified baryon, electric, strangeness, and charm densities.
 
virtual void CalculateDensities ()
 Calculates the primordial and total (after decays) densities of all species.
 
virtual void ValidateCalculation ()
 Checks whether issues have occured during the calculation of particle densities in the CalculateDensities() method.
 
std::string ValidityCheckLog () const
 All messaged which occured during the validation procedure in the ValidateCalculation() method.
 
virtual void CalculateDensitiesGCE ()
 Calculates the particle densities in a grand-canonical ensemble.
 
virtual void CalculateFeeddown ()
 Calculates the total densities which include feeddown contributions.
 
virtual void CalculateTwoParticleFluctuationsDecays ()
 Computes particle number correlations and fluctuations for all final-state particles which are marked stable.
 
virtual double TwoParticleSusceptibilityPrimordial (int i, int j) const
 Returns the computed primordial particle number (cross-)susceptibility \( \frac{1}{VT^3} \, \langle \Delta N_i \Delta N_j \rangle \) for particles with ids i and j. CalculateFluctuations() must be called beforehand.
 
virtual double TwoParticleSusceptibilityPrimordialByPdg (long long id1, long long id2)
 Returns the computed primordial particle number (cross-)susceptibility \( \frac{1}{VT^3} \, \langle \Delta N_i \Delta N_j \rangle \) for particles with pdg codes id1 and id2. CalculateFluctuations() must be called beforehand.
 
virtual double TwoParticleSusceptibilityTemperatureDerivativePrimordial (int i, int j) const
 Returns the computed temperature derivative of the primordial particle number (cross-)susceptibility \( \frac{\partial}{\partial T} \frac{1}{VT^3} \, \langle \Delta N_i \Delta N_j \rangle \) for particles with ids i and j. CalculateFluctuations() must be called beforehand.
 
virtual double TwoParticleSusceptibilityTemperatureDerivativePrimordialByPdg (long long id1, long long id2)
 Returns the computed temperature derivative of the primordial particle number (cross-)susceptibility \( \frac{\partial}{\partial T} \frac{1}{VT^3} \, \langle \Delta N_i \Delta N_j \rangle \) for particles with pdg codes id1 and id2. CalculateFluctuations() must be called beforehand.
 
virtual double NetParticleSusceptibilityPrimordialByPdg (long long id1, long long id2)
 Returns the computed (cross-)susceptibility \( \frac{1}{VT^3} \, \langle \Delta N_i \Delta N_j \rangle \) between primordial net-particle numbers for pdg codes id1 and id2. CalculateFluctuations() must be called beforehand.
 
virtual double TwoParticleSusceptibilityFinal (int i, int j) const
 Returns the computed final particle number (cross-)susceptibility \( \frac{1}{VT^3} \, \langle \Delta N_i \Delta N_j \rangle \) for particles with ids i and j. CalculateFluctuations() must be called beforehand. Both particle species must be those marked stable.
 
virtual double TwoParticleSusceptibilityFinalByPdg (long long id1, long long id2)
 Returns the computed final particle number (cross-)susceptibility \( \frac{1}{VT^3} \, \langle \Delta N_i \Delta N_j \rangle \) for particles with pdg codes id1 and id2. CalculateFluctuations() must be called beforehand. Both particle species must be those marked stable.
 
virtual double NetParticleSusceptibilityFinalByPdg (long long id1, long long id2)
 Returns the computed (cross-)susceptibility \( \frac{1}{VT^3} \, \langle \Delta N_i \Delta N_j \rangle \) between final net-particle numbers for pdg codes id1 and id2. CalculateFluctuations() must be called beforehand.
 
virtual double PrimordialParticleChargeSusceptibility (int i, ConservedCharge::Name chg) const
 Returns the computed primordial particle vs conserved charge (cross-)susceptibility \( \frac{1}{VT^3} \, \langle \Delta N_i \Delta N_chg \rangle \) for particle with id i and conserved charge chg. CalculateFluctuations() must be called beforehand.
 
virtual double PrimordialParticleChargeSusceptibilityByPdg (long long id1, ConservedCharge::Name chg)
 Returns the computed primordial particle vs conserved charge (cross-)susceptibility \( \frac{1}{VT^3} \, \langle \Delta N_i \Delta N_chg \rangle \) for particle with pdg code id1 and conserved charge chg. CalculateFluctuations() must be called beforehand.
 
virtual double PrimordialNetParticleChargeSusceptibilityByPdg (long long id1, ConservedCharge::Name chg)
 Returns the computed primordial net particle vs conserved charge (cross-)susceptibility \( \frac{1}{VT^3} \, \langle \Delta N_i \Delta N_chg \rangle \) for particle with pdg code id1 and conserved charge chg. CalculateFluctuations() must be called beforehand.
 
virtual double FinalParticleChargeSusceptibility (int i, ConservedCharge::Name chg) const
 Returns the computed final (after decays) particle vs conserved charge (cross-)susceptibility \( \frac{1}{VT^3} \, \langle \Delta N_i \Delta N_chg \rangle \) for particle with id i and conserved charge chg. CalculateFluctuations() must be called beforehand.
 
virtual double FinalParticleChargeSusceptibilityByPdg (long long id1, ConservedCharge::Name chg)
 Returns the computed final (after decays) particle vs conserved charge (cross-)susceptibility \( \frac{1}{VT^3} \, \langle \Delta N_i \Delta N_chg \rangle \) for particle with pdg code id1 and conserved charge chg. CalculateFluctuations() must be called beforehand.
 
virtual double FinalNetParticleChargeSusceptibilityByPdg (long long id1, ConservedCharge::Name chg)
 Returns the computed final (after decays) net particle vs conserved charge (cross-)susceptibility \( \frac{1}{VT^3} \, \langle \Delta N_i \Delta N_chg \rangle \) for particle with pdg code id1 and conserved charge chg. CalculateFluctuations() must be called beforehand.
 
virtual double SusceptibilityDimensionfull (ConservedCharge::Name i, ConservedCharge::Name j) const
 A 2nd order susceptibility of conserved charges.
 
virtual void CalculateSusceptibilityMatrix ()
 Calculates the conserved charges susceptibility matrix.
 
virtual void CalculateProxySusceptibilityMatrix ()
 Calculates the susceptibility matrix of conserved charges proxies.
 
virtual void CalculateParticleChargeCorrelationMatrix ()
 Calculates the matrix of correlators between primordial (and also final) particle numbers and conserved charges.
 

Protected Member Functions

std::vector< double > ComputeNp (const std::vector< double > &dmustar)
 
std::vector< double > ComputeNp (const std::vector< double > &dmustar, const std::vector< double > &ns)
 
void CalculateVDWComponentsMap ()
 Partitions particles species into sets that have identical VDW parameters.
 
virtual std::vector< double > SearchSingleSolution (const std::vector< double > &muStarInit)
 Uses the Broyden method with a provided initial guess to determine the shifted chemical potentials by solving the transcendental equations with the Broyden's method.
 
std::vector< double > SearchMultipleSolutions (int iters=300)
 Uses the Broyden method with different initial guesses to look for different possible solutions of the transcendental equations for shifted chemical potentials.
 
void SolveEquations ()
 
virtual double MuShift (int id) const
 The shift in the chemical potential of particle species i due to the QvdW interactions.
 

Protected Attributes

std::vector< double > m_DensitiesId
 Vector of ideal gas densities with shifted chemical potentials.
 
std::vector< double > m_scaldens
 Vector of scalar densities. Not used.
 
std::vector< std::vector< double > > m_Virial
 
std::vector< std::vector< double > > m_Attr
 Matrix of the attractive QvdW coefficients \( a_{ij} \).
 
std::vector< std::vector< double > > m_VirialdT
 
std::vector< std::vector< double > > m_AttrdT
 
bool m_TemperatureDependentAB
 
bool m_SearchMultipleSolutions
 Whether multiple solutions are considered.
 
bool m_LastBroydenSuccessFlag
 Whether Broyden's method was successfull.
 
bool m_VDWComponentMapCalculated
 Whether the mapping to components with the same VDW parameters has been calculated.
 
std::vector< double > m_MuStar
 Vector of the shifted chemical potentials.
 
std::vector< int > m_MapTodMuStar
 
std::vector< int > m_MapFromdMuStar
 
std::vector< std::vector< int > > m_dMuStarIndices
 
std::vector< std::vector< double > > m_chi
 
std::vector< double > m_chiarb
 

Additional Inherited Members

- Public Types inherited from thermalfist::ThermalModelBase
enum  ThermalModelEnsemble { GCE = 0 , CE = 1 , SCE = 2 , CCE = 3 }
 The list of statistical ensembles. More...
 
enum  ThermalModelInteraction {
  Ideal = 0 , DiagonalEV = 1 , CrosstermsEV = 2 , QvdW = 3 ,
  RealGas = 4 , MeanField = 5
}
 Type of interactions included in the HRG model. More...
 

Detailed Description

Class implementing the quantum van der Waals HRG model.

The model formulation can be found in

V. Vovchenko, M.I. Gorenstein, H. Stoecker, Phys. Rev. Lett. 118, 182301 (2017), http://arxiv.org/pdf/1609.03975.pdf

and in, in more detail, in

V. Vovchenko, A. Motornenko, P. Alba, M.I. Gorenstein, L.M. Satarov, H. Stoecker, Phys. Rev. C 96, 045202 (2017), http://arxiv.org/pdf/1707.09215.pdf

The system of transcendental equations for the "shifted" chemical potentials of hadrons, Eq. (15) in the latter reference, is solved using the Broyden's method.

Examples
SusceptibilitiesBQS.cpp, and ThermodynamicsBQS.cpp.

Definition at line 36 of file ThermalModelVDW.h.

Constructor & Destructor Documentation

◆ ThermalModelVDW()

thermalfist::ThermalModelVDW::ThermalModelVDW ( ThermalParticleSystem * TPS_,
const ThermalModelParameters & params = ThermalModelParameters() )

Construct a new ThermalModelVDW object.

Parameters
TPSA pointer to the ThermalParticleSystem object containing the particle list
paramsThermalModelParameters object with current thermal parameters

Definition at line 34 of file ThermalModelVDW.cpp.

◆ ~ThermalModelVDW()

thermalfist::ThermalModelVDW::~ThermalModelVDW ( void )
virtual

Destroy the ThermalModelVDW object.

Definition at line 54 of file ThermalModelVDW.cpp.

Member Function Documentation

◆ AttractionCoefficient()

double thermalfist::ThermalModelVDW::AttractionCoefficient ( int ,
int  ) const
virtual

QvdW mean field attraction coefficient \( a_{ij} \).

Parameters
i0-based index of the first particle species
j0-based index of the second particle species
Returns
double Coefficient \( a_{ij} = 0 \)

Reimplemented from thermalfist::ThermalModelBase.

Definition at line 1155 of file ThermalModelVDW.cpp.

◆ AttractionCoefficientdT()

double thermalfist::ThermalModelVDW::AttractionCoefficientdT ( int i,
int j ) const

The temperature derivative of the QvdW attraction parameter \( a_{ij} \).

Parameters
i0-based index of the first particle species
j0-based index of the second particle species
Returns
\( d a_{ij} / dT \) in the units of fm \(^3\)

Definition at line 1169 of file ThermalModelVDW.cpp.

◆ AttractionMatrix()

const std::vector< std::vector< double > > & thermalfist::ThermalModelVDW::AttractionMatrix ( ) const
inline

Definition at line 210 of file ThermalModelVDW.h.

◆ CalculateBaryonMatterEntropyDensity()

virtual double thermalfist::ThermalModelVDW::CalculateBaryonMatterEntropyDensity ( )
inlinevirtual

Definition at line 195 of file ThermalModelVDW.h.

◆ CalculateChargeFluctuations()

vector< double > thermalfist::ThermalModelVDW::CalculateChargeFluctuations ( const std::vector< double > & chgs,
int order = 4 )
virtual

Calculates fluctuations (diagonal susceptibilities) of an arbitrary "conserved" charge.

Each particle specie is assumed to carry a conserved charge with a value provided by an input vector.

Restricted to the grand canonical ensemble.

Parameters
chgsA vector with conserved charge values for all species. 0-based indices of the vector must correspond to the 0-based indices of the particle list TPS()
orderUp to which order the susceptibilities are computed
Returns
std::vector<double> A vector with computed values of diagonal susceptibilities

Reimplemented from thermalfist::ThermalModelBase.

Definition at line 562 of file ThermalModelVDW.cpp.

◆ CalculateEnergyDensity()

double thermalfist::ThermalModelVDW::CalculateEnergyDensity ( )
virtual

Definition at line 1079 of file ThermalModelVDW.cpp.

◆ CalculateEnergyDensityDerivativeT()

double thermalfist::ThermalModelVDW::CalculateEnergyDensityDerivativeT ( )
virtual

Definition at line 882 of file ThermalModelVDW.cpp.

◆ CalculateEntropyDensity()

double thermalfist::ThermalModelVDW::CalculateEntropyDensity ( )
virtual

Definition at line 1102 of file ThermalModelVDW.cpp.

◆ CalculateFluctuations() [1/2]

void thermalfist::ThermalModelVDW::CalculateFluctuations ( )
virtual

Computes the fluctuation observables.

Includes the matrix of 2nd order susceptibilities of conserved charges, as well as particle number correlations and fluctuations, both for primordial yields and after decays.

Reimplemented from thermalfist::ThermalModelBase.

Definition at line 859 of file ThermalModelVDW.cpp.

◆ CalculateFluctuations() [2/2]

vector< vector< double > > thermalfist::ThermalModelVDW::CalculateFluctuations ( int order)
virtual

Definition at line 747 of file ThermalModelVDW.cpp.

◆ CalculateMesonMatterEntropyDensity()

virtual double thermalfist::ThermalModelVDW::CalculateMesonMatterEntropyDensity ( )
inlinevirtual

Definition at line 197 of file ThermalModelVDW.h.

◆ CalculatePressure()

double thermalfist::ThermalModelVDW::CalculatePressure ( )
virtual

Definition at line 1122 of file ThermalModelVDW.cpp.

◆ CalculatePrimordialDensities()

void thermalfist::ThermalModelVDW::CalculatePrimordialDensities ( )
virtual

Calculates the primordial densities of all species.

Implements thermalfist::ThermalModelBase.

Definition at line 401 of file ThermalModelVDW.cpp.

◆ CalculateTemperatureDerivatives()

void thermalfist::ThermalModelVDW::CalculateTemperatureDerivatives ( )
virtual

Computes the temperature derivatives of densities, shifted chemical potentials, and primordial hadron number susceptibilities.

Reimplemented from thermalfist::ThermalModelBase.

Definition at line 913 of file ThermalModelVDW.cpp.

◆ CalculateTwoParticleCorrelations()

void thermalfist::ThermalModelVDW::CalculateTwoParticleCorrelations ( )
virtual

Computes the fluctuations and correlations of the primordial particle numbers.

More specifically, computes the susceptibility matrix \( \frac{1}{VT^3} \, \langle \Delta N_i^* \Delta N_j^* \rangle \), where \( N_i^* \) is the primordial yield of species i.

Reimplemented from thermalfist::ThermalModelBase.

Definition at line 780 of file ThermalModelVDW.cpp.

◆ CalculateVDWComponentsMap()

void thermalfist::ThermalModelVDW::CalculateVDWComponentsMap ( )
protected

Partitions particles species into sets that have identical VDW parameters.

Definition at line 253 of file ThermalModelVDW.cpp.

◆ ChangeTPS()

void thermalfist::ThermalModelVDW::ChangeTPS ( ThermalParticleSystem * TPS)
virtual

Change the particle list.

Parameters
TPSA pointer to new particle list.

Reimplemented from thermalfist::ThermalModelBase.

Definition at line 187 of file ThermalModelVDW.cpp.

◆ ComputeNp() [1/2]

std::vector< double > thermalfist::ThermalModelVDW::ComputeNp ( const std::vector< double > & dmustar)
protected

Returns vector of particle densities for given values of shifted chemical potentials

Parameters
dmustarA vector of shifted chemical potentials
Returns
std::vector<double> Vector of particle densities

Definition at line 196 of file ThermalModelVDW.cpp.

◆ ComputeNp() [2/2]

std::vector< double > thermalfist::ThermalModelVDW::ComputeNp ( const std::vector< double > & dmustar,
const std::vector< double > & ns )
protected

Same as ComputeNp(const std::vector<double>&) but using the vector of ideal gas densities as input instead of calculating it

Definition at line 205 of file ThermalModelVDW.cpp.

◆ DeltaMu()

virtual double thermalfist::ThermalModelVDW::DeltaMu ( int i) const
inlinevirtual

Definition at line 208 of file ThermalModelVDW.h.

◆ DensityId()

double thermalfist::ThermalModelVDW::DensityId ( int index)
inline

Definition at line 203 of file ThermalModelVDW.h.

◆ FillAttraction()

void thermalfist::ThermalModelVDW::FillAttraction ( const std::vector< std::vector< double > > & aij = std::vector< std::vector<double> >(0))
virtual

Definition at line 108 of file ThermalModelVDW.cpp.

◆ FillChemicalPotentials()

void thermalfist::ThermalModelVDW::FillChemicalPotentials ( )
virtual

Sets the chemical potentials of all particles.

Uses the current values of \( \mu_B,\,\mu_Q,\,\mu_S,\,\mu_Q \) to set \( \mu_i \).

Reimplemented from thermalfist::ThermalModelBase.

Definition at line 58 of file ThermalModelVDW.cpp.

◆ FillVirial()

void thermalfist::ThermalModelVDW::FillVirial ( const std::vector< double > & ri = std::vector<double>(0))
virtual

Fills the excluded volume coefficients \( \tilde{b}_{ij} \) based on the provided radii parameters for all species.

Fills the coefficients in accordance with Eqs. (5) and (7) here https://arxiv.org/pdf/1606.06218.pdf

Parameters
riA vector with radii parameters for all species. 0-based indices of the vector must corresponds to the 0-based indices of the particle list TPS()

Reimplemented from thermalfist::ThermalModelBase.

Definition at line 71 of file ThermalModelVDW.cpp.

◆ FillVirialEV()

void thermalfist::ThermalModelVDW::FillVirialEV ( const std::vector< std::vector< double > > & bij = std::vector< std::vector<double> >(0))

Same as FillVirial() but uses the matrix of excluded-volume coefficients \( v_i \equiv b_{ii} \) as input instead of radii.

Parameters
bijA vector with crossterms excluded-volume coefficients for all pairs of particle species

Definition at line 97 of file ThermalModelVDW.cpp.

◆ GetMuStar()

std::vector< double > thermalfist::ThermalModelVDW::GetMuStar ( ) const
inline

Returns vector of shifted chemical potentials, one element per each species

Definition at line 145 of file ThermalModelVDW.h.

◆ IsLastSolutionOK()

bool thermalfist::ThermalModelVDW::IsLastSolutionOK ( ) const
inline

Definition at line 201 of file ThermalModelVDW.h.

◆ MuShift()

double thermalfist::ThermalModelVDW::MuShift ( int id) const
protectedvirtual

The shift in the chemical potential of particle species i due to the QvdW interactions.

Parameters
i0-based particle specie index
Returns
The shift in the chemical potential

Definition at line 1140 of file ThermalModelVDW.cpp.

◆ MuStar()

double thermalfist::ThermalModelVDW::MuStar ( int i) const
inline

The shifted chemical potential of particle species i.

Definition at line 141 of file ThermalModelVDW.h.

◆ ParticleScalarDensity()

double thermalfist::ThermalModelVDW::ParticleScalarDensity ( int part)
virtual

Definition at line 1134 of file ThermalModelVDW.cpp.

◆ ReadInteractionParameters()

void thermalfist::ThermalModelVDW::ReadInteractionParameters ( const std::string & )
virtual

Reads the QvdW interaction parameters from a file.

Actual implementation is in a derived class.

Parameters
filenameFile with interaction parameters.

Reimplemented from thermalfist::ThermalModelBase.

Reimplemented in thermalfist::ThermalModelEVCrossterms.

Definition at line 119 of file ThermalModelVDW.cpp.

◆ SearchMultipleSolutions()

vector< double > thermalfist::ThermalModelVDW::SearchMultipleSolutions ( int iters = 300)
protected

Uses the Broyden method with different initial guesses to look for different possible solutions of the transcendental equations for shifted chemical potentials.

Looks for the solution with the largest pressure.

Parameters
itersNumber of different initial guesses to try
Returns
std::vector<double> The solution with the largest pressure among those which were found

Definition at line 319 of file ThermalModelVDW.cpp.

◆ SearchSingleSolution()

vector< double > thermalfist::ThermalModelVDW::SearchSingleSolution ( const std::vector< double > & muStarInit)
protectedvirtual

Uses the Broyden method with a provided initial guess to determine the shifted chemical potentials by solving the transcendental equations with the Broyden's method.

Parameters
muStarInitInitial guess for the shifted chemical potentials
Returns
std::vector<double> The solved shifted chemical potentials

Definition at line 290 of file ThermalModelVDW.cpp.

◆ SetAttraction()

virtual void thermalfist::ThermalModelVDW::SetAttraction ( int ,
int ,
double  )
inlinevirtual

Set the vdW mean field attraction coefficient \( a_{ij} \).

Parameters
i0-based index of the first particle species
j0-based index of the second particle species
avdW mean field attraction parameter \( a_{ij} \) (GeV fm \(^3\))

Reimplemented from thermalfist::ThermalModelBase.

Reimplemented in thermalfist::ThermalModelEVCrossterms.

Definition at line 166 of file ThermalModelVDW.h.

◆ SetAttractiondT()

void thermalfist::ThermalModelVDW::SetAttractiondT ( int i,
int j,
double dadT )
inline

Set the temperature derivative of the QvdW attraction parameter \( a_{ij} \).

Parameters
i0-based index of the first particle species
j0-based index of the second particle species
dadT\( d a_{ij} / dT \) in the units of fm \(^3\)

Definition at line 79 of file ThermalModelVDW.h.

◆ SetChemicalPotentials()

void thermalfist::ThermalModelVDW::SetChemicalPotentials ( const std::vector< double > & chem = std::vector<double>(0))
virtual

Sets the chemical potentials of all particles.

Parameters
chemA vector with chemical potentials of all species. 0-based indices of the vector must corresponds to the 0-based indices of the particle list TPS()

Reimplemented from thermalfist::ThermalModelBase.

Definition at line 64 of file ThermalModelVDW.cpp.

◆ SetMultipleSolutionsMode()

virtual void thermalfist::ThermalModelVDW::SetMultipleSolutionsMode ( bool search)
inlinevirtual

Whether to search for multiple solutions of the QvdW equations by considering different initial guesses in the Broyden's method.

Multiple solutions in the QvdW model appear e.g. below the critical temperature of the liquid-gas phase transition.

Parameters
searchWhether multiple solutions of the QvdW equations should be considered. False by default.

Reimplemented in thermalfist::ThermalModelEVCrossterms.

Definition at line 129 of file ThermalModelVDW.h.

◆ SetMuStar()

void thermalfist::ThermalModelVDW::SetMuStar ( const std::vector< double > & MuStar)
inline

Set the vector of shifted chemical potentials.

Definition at line 148 of file ThermalModelVDW.h.

◆ SetTemperatureDependentAB()

void thermalfist::ThermalModelVDW::SetTemperatureDependentAB ( bool Tdep)
inline

Sets whether temperature depedence of QvdW parameters should be considered.

Parameters
Tdeptrue – considered, false – not considered

Definition at line 107 of file ThermalModelVDW.h.

◆ SetVirial()

virtual void thermalfist::ThermalModelVDW::SetVirial ( int ,
int ,
double  )
inlinevirtual

Set the excluded volume coefficient \( \tilde{b}_{ij} \).

Excluded parameter for repulsive interaction between particle species i and j.

Parameters
i0-based index of the first particle species
j0-based index of the second particle species
bExcluded volume parameter \( \tilde{b}_{ij} \) (fm \(^3\))

Reimplemented from thermalfist::ThermalModelBase.

Definition at line 164 of file ThermalModelVDW.h.

◆ SetVirialdT()

void thermalfist::ThermalModelVDW::SetVirialdT ( int i,
int j,
double dbdT )
inline

Set the temperature derivative of the eigenvolume parameter \( \tilde{b}_{ij} \).

Parameters
i0-based index of the first particle species
j0-based index of the second particle species
dbdT\( d \tilde{b}_{ij} / dT \) in the units of fm \(^3\) GeV \(^{-1}\)

Definition at line 69 of file ThermalModelVDW.h.

◆ SolveEquations()

void thermalfist::ThermalModelVDW::SolveEquations ( )
protected

Solve the transcedental equations for the shifted chemical potentials

Definition at line 385 of file ThermalModelVDW.cpp.

◆ TemperatureDependentAB()

bool thermalfist::ThermalModelVDW::TemperatureDependentAB ( ) const
inline

Whether temperature depedence of QvdW parameters is considered.

Returns
true considered
false not considered

Definition at line 116 of file ThermalModelVDW.h.

◆ UseMultipleSolutionsMode()

bool thermalfist::ThermalModelVDW::UseMultipleSolutionsMode ( ) const
inline

Whether to search for multiple solutions of the QvdW equations by considering different initial guesses in the Broyden's method.

Returns
true Multiple solutions considered
false Multiple solutions not considered

Definition at line 138 of file ThermalModelVDW.h.

◆ VDWComponentIndices()

const std::vector< std::vector< int > > & thermalfist::ThermalModelVDW::VDWComponentIndices ( ) const
inline

Definition at line 207 of file ThermalModelVDW.h.

◆ VirialCoefficient()

double thermalfist::ThermalModelVDW::VirialCoefficient ( int ,
int  ) const
virtual

Excluded volume coefficient \( \tilde{b}_{ij} = 0 \).

Excluded parameter for repulsive interaction between particle species i and j.

Parameters
i0-based index of the first particle species
j0-based index of the second particle species
Returns
double Coefficient \( \tilde{b}_{ij} = 0 \)

Reimplemented from thermalfist::ThermalModelBase.

Definition at line 1148 of file ThermalModelVDW.cpp.

◆ VirialCoefficientdT()

double thermalfist::ThermalModelVDW::VirialCoefficientdT ( int i,
int j ) const

The temperature derivative of the eigenvolume parameter \( \tilde{b}_{ij} \).

Parameters
i0-based index of the first particle species
j0-based index of the second particle species
Returns
dbdT \( d \tilde{b}_{ij} / dT \) in the units of fm \(^3\) GeV \(^{-1}\)

Definition at line 1162 of file ThermalModelVDW.cpp.

◆ VirialMatrix()

const std::vector< std::vector< double > > & thermalfist::ThermalModelVDW::VirialMatrix ( ) const
inline

Definition at line 209 of file ThermalModelVDW.h.

◆ WriteInteractionParameters()

void thermalfist::ThermalModelVDW::WriteInteractionParameters ( const std::string & )
virtual

Write the QvdW interaction parameters to a file.

Actual implementation is in a derived class.

Parameters
filenameOutput file.

Reimplemented from thermalfist::ThermalModelBase.

Reimplemented in thermalfist::ThermalModelEVCrossterms.

Definition at line 151 of file ThermalModelVDW.cpp.

Member Data Documentation

◆ m_Attr

std::vector< std::vector<double> > thermalfist::ThermalModelVDW::m_Attr
protected

Matrix of the attractive QvdW coefficients \( a_{ij} \).

Definition at line 275 of file ThermalModelVDW.h.

◆ m_AttrdT

std::vector< std::vector<double> > thermalfist::ThermalModelVDW::m_AttrdT
protected

Matrix of the temperature derivatives of the attractive QvdW coefficients \( d a_{ij} / dT \)

Definition at line 283 of file ThermalModelVDW.h.

◆ m_chi

std::vector< std::vector<double> > thermalfist::ThermalModelVDW::m_chi
protected

Definition at line 307 of file ThermalModelVDW.h.

◆ m_chiarb

std::vector<double> thermalfist::ThermalModelVDW::m_chiarb
protected

Definition at line 309 of file ThermalModelVDW.h.

◆ m_DensitiesId

std::vector<double> thermalfist::ThermalModelVDW::m_DensitiesId
protected

Vector of ideal gas densities with shifted chemical potentials.

Definition at line 266 of file ThermalModelVDW.h.

◆ m_dMuStarIndices

std::vector< std::vector<int> > thermalfist::ThermalModelVDW::m_dMuStarIndices
protected

Definition at line 305 of file ThermalModelVDW.h.

◆ m_LastBroydenSuccessFlag

bool thermalfist::ThermalModelVDW::m_LastBroydenSuccessFlag
protected

Whether Broyden's method was successfull.

Definition at line 293 of file ThermalModelVDW.h.

◆ m_MapFromdMuStar

std::vector<int> thermalfist::ThermalModelVDW::m_MapFromdMuStar
protected

Definition at line 303 of file ThermalModelVDW.h.

◆ m_MapTodMuStar

std::vector<int> thermalfist::ThermalModelVDW::m_MapTodMuStar
protected

Definition at line 301 of file ThermalModelVDW.h.

◆ m_MuStar

std::vector<double> thermalfist::ThermalModelVDW::m_MuStar
protected

Vector of the shifted chemical potentials.

Definition at line 299 of file ThermalModelVDW.h.

◆ m_scaldens

std::vector<double> thermalfist::ThermalModelVDW::m_scaldens
protected

Vector of scalar densities. Not used.

Definition at line 269 of file ThermalModelVDW.h.

◆ m_SearchMultipleSolutions

bool thermalfist::ThermalModelVDW::m_SearchMultipleSolutions
protected

Whether multiple solutions are considered.

Definition at line 290 of file ThermalModelVDW.h.

◆ m_TemperatureDependentAB

bool thermalfist::ThermalModelVDW::m_TemperatureDependentAB
protected

Whether temperature depedence of QvdW parameters is considered.

Definition at line 287 of file ThermalModelVDW.h.

◆ m_VDWComponentMapCalculated

bool thermalfist::ThermalModelVDW::m_VDWComponentMapCalculated
protected

Whether the mapping to components with the same VDW parameters has been calculated.

Definition at line 296 of file ThermalModelVDW.h.

◆ m_Virial

std::vector< std::vector<double> > thermalfist::ThermalModelVDW::m_Virial
protected

Definition at line 272 of file ThermalModelVDW.h.

◆ m_VirialdT

std::vector< std::vector<double> > thermalfist::ThermalModelVDW::m_VirialdT
protected

Matrix of the temperature derivatives of the virial (excluded-volume) coefficients \( d \tilde{b}_{ij} / dT \)

Definition at line 279 of file ThermalModelVDW.h.


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