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

Class implementing the quantum real gas HRG model. More...

#include <ThermalModelRealGas.h>

Inheritance diagram for thermalfist::ThermalModelRealGas:
thermalfist::ThermalModelBase

Public Member Functions

 ThermalModelRealGas (ThermalParticleSystem *TPS_, const ThermalModelParameters &params=ThermalModelParameters())
 Construct a new ThermalModelRealGas object.
 
virtual ~ThermalModelRealGas (void)
 Destroy the ThermalModelRealGas object.
 
void SetExcludedVolumeModel (ExcludedVolumeModelMultiBase *exvolmod)
 Set the excluded volume model for the real gas HRG model.
 
void SetMeanFieldModel (MeanFieldModelMultiBase *mfmod)
 Set the mean field model for the real gas HRG model.
 
ExcludedVolumeModelMultiBaseExcludedVolumeModel () const
 Get the excluded volume model used in the real gas HRG model.
 
MeanFieldModelMultiBaseMeanFieldModel () const
 Get the mean field model used in the real gas HRG model.
 
virtual void SetMultipleSolutionsMode (bool search)
 Whether to search for multiple solutions of the real gas model equations by considering different initial guesses in the Broyden's method.
 
bool UseMultipleSolutionsMode () const
 Whether to search for multiple solutions of the real gas model equations by considering different initial guesses in the Broyden's method.
 
double MuStar (int i) const
 Get the shifted chemical potential of a particle species.
 
std::vector< double > GetMuStar () const
 Get the vector of shifted chemical potentials for all particle species.
 
void SetMuStar (const std::vector< double > &MuStar)
 Set the vector of shifted chemical potentials for all particle species.
 
void FillChemicalPotentials ()
 Fill the chemical potentials of the particle species.
 
virtual void SetChemicalPotentials (const std::vector< double > &chem=std::vector< double >(0))
 Set the chemical potentials of the particle species.
 
virtual void CalculatePrimordialDensities ()
 Calculate the primordial densities of the particle species.
 
virtual std::vector< double > CalculateChargeFluctuations (const std::vector< double > &chgs, int order=4)
 Calculate the charge fluctuations of the particle species.
 
virtual std::vector< double > CalculateChargeFluctuationsOld (const std::vector< double > &chgs, int order=4)
 Calculate the charge fluctuations of the particle species (old method).
 
virtual std::vector< std::vector< double > > CalculateFluctuations (int order)
 Calculate the fluctuations of the particle species.
 
void CalculateTwoParticleCorrelations ()
 Calculate the two-particle correlations of the particle species.
 
void CalculateFluctuations ()
 Calculate the fluctuations.
 
virtual double CalculatePressure ()
 Calculate the pressure of the system.
 
virtual double CalculateEnergyDensity ()
 Calculate the energy density of the system.
 
virtual double CalculateEntropyDensity ()
 Calculate the entropy density of the system.
 
virtual double CalculateEnergyDensityDerivativeT ()
 Calculate the derivative of the energy density with respect to temperature.
 
virtual void CalculateTemperatureDerivatives ()
 Calculate the temperature derivatives of the system.
 
virtual double CalculateBaryonMatterEntropyDensity ()
 Calculate the baryon matter entropy density of the system.
 
virtual double CalculateMesonMatterEntropyDensity ()
 Calculate the meson matter entropy density of the system.
 
virtual double ParticleScalarDensity (int part)
 Calculate the scalar density of a particle species.
 
bool IsLastSolutionOK () const
 Check if the last solution was successful.
 
double DensityId (int index)
 Get the ideal gas density of a particle species.
 
const std::vector< std::vector< int > > & ComponentIndices () const
 Get the component indices of the particle species.
 
virtual double DeltaMu (int i) const
 Get the delta mu (chemical potential shift due to interactions) of a particle species.
 
double ChiArb (int charge)
 
- 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.
 
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.
 
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 SetVirial (int, int, double)
 Set the excluded volume coefficient \( \tilde{b}_{ij} \).
 
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 SetAttraction (int, int, double)
 Set the vdW mean field attraction coefficient \( a_{ij} \).
 
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 ()
 
virtual void ReadInteractionParameters (const std::string &)
 Reads the QvdW interaction parameters from a file.
 
virtual void WriteInteractionParameters (const std::string &)
 Write the QvdW interaction parameters to a file.
 
virtual void ChangeTPS (ThermalParticleSystem *TPS)
 Change the particle list.
 
virtual double VirialCoefficient (int, int) const
 Excluded volume coefficient \( \tilde{b}_{ij} = 0 \).
 
double RepulsionCoefficient (int i, int j) const
 
virtual double AttractionCoefficient (int, int) const
 QvdW mean field attraction coefficient \( a_{ij} \).
 
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.
 

Public Attributes

std::vector< std::vector< double > > m_chi
 Vector of computed susceptibilities values.
 
std::vector< double > m_chiarb
 Vector of computed susceptibilities for a specified arbitraty charge.
 

Protected Member Functions

void CalculateComponentsMap ()
 Partitions particles species into sets that have identical pair interactions.
 
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.
 
virtual std::vector< double > SearchSingleSolutionUsingComponents (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, using the component mapping.
 
virtual std::vector< double > SearchSingleSolutionDirect (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 for all particles 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.
 
bool m_SearchMultipleSolutions
 Whether multiple solutions are considered.
 
bool m_LastBroydenSuccessFlag
 Whether Broyden's method was successfull.
 
bool m_ComponentMapCalculated
 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
 Mapping from particle species to dMuStar indices.
 
std::vector< int > m_MapFromdMuStar
 Mapping from dMuStar indices to particle species.
 
std::vector< std::vector< int > > m_dMuStarIndices
 Vector of component indices for each particle species.
 
ExcludedVolumeModelMultiBasem_exvolmod
 Excluded volume model used in the real gas HRG model.
 
MeanFieldModelMultiBasem_mfmod
 Mean field model used in the real gas HRG model.
 
ExcludedVolumeModelMultiBasem_exvolmodideal
 Excluded volume model object in the ideal gas limit.
 
MeanFieldModelMultiBasem_mfmodideal
 Mean field model object in the ideal gas limit.
 

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 real gas HRG model.

The model formulation can be found in

V. Vovchenko, Phys. Rev. C 96, 015206 (2017), https://arxiv.org/pdf/1701.06524

The system of transcendental equations for the "shifted" chemical potentials of hadrons is solved using the Broyden's method.

Examples
NeutronStars-CSHRG.cpp, SusceptibilitiesBQS.cpp, and ThermodynamicsBQS.cpp.

Definition at line 24 of file ThermalModelRealGas.h.

Constructor & Destructor Documentation

◆ ThermalModelRealGas()

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

Construct a new ThermalModelRealGas object.

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

Definition at line 25 of file ThermalModelRealGas.cpp.

◆ ~ThermalModelRealGas()

thermalfist::ThermalModelRealGas::~ThermalModelRealGas ( void )
virtual

Destroy the ThermalModelRealGas object.

Definition at line 45 of file ThermalModelRealGas.cpp.

Member Function Documentation

◆ CalculateBaryonMatterEntropyDensity()

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

Calculate the baryon matter entropy density of the system.

Returns
The baryon matter entropy density of the system (always returns 0).

Definition at line 210 of file ThermalModelRealGas.h.

◆ CalculateChargeFluctuations()

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

Calculate the charge fluctuations of the particle species.

Parameters
chgsVector of charges, one element per particle species.
orderOrder of the fluctuations (default is 4).
Returns
Vector of charge fluctuations, one element per particle species.

Reimplemented from thermalfist::ThermalModelBase.

Definition at line 284 of file ThermalModelRealGas.cpp.

◆ CalculateChargeFluctuationsOld()

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

Calculate the charge fluctuations of the particle species (old method).

Parameters
chgsVector of charges, one element per particle species.
orderOrder of the fluctuations (default is 4).
Returns
Vector of charge fluctuations, one element per particle species.

Definition at line 920 of file ThermalModelRealGas.cpp.

◆ CalculateComponentsMap()

void thermalfist::ThermalModelRealGas::CalculateComponentsMap ( )
protected

Partitions particles species into sets that have identical pair interactions.

This function is used to optimize the calculation of the real gas model equations.

Definition at line 82 of file ThermalModelRealGas.cpp.

◆ CalculateEnergyDensity()

double thermalfist::ThermalModelRealGas::CalculateEnergyDensity ( )
virtual

Calculate the energy density of the system.

Returns
The energy density of the system.

Definition at line 1652 of file ThermalModelRealGas.cpp.

◆ CalculateEnergyDensityDerivativeT()

double thermalfist::ThermalModelRealGas::CalculateEnergyDensityDerivativeT ( )
virtual

Calculate the derivative of the energy density with respect to temperature.

Returns
The derivative of the energy density with respect to temperature.

Definition at line 1619 of file ThermalModelRealGas.cpp.

◆ CalculateEntropyDensity()

double thermalfist::ThermalModelRealGas::CalculateEntropyDensity ( )
virtual

Calculate the entropy density of the system.

Returns
The entropy density of the system.

Definition at line 1670 of file ThermalModelRealGas.cpp.

◆ CalculateFluctuations() [1/2]

void thermalfist::ThermalModelRealGas::CalculateFluctuations ( )
virtual

Calculate the fluctuations.

Reimplemented from thermalfist::ThermalModelBase.

Definition at line 1596 of file ThermalModelRealGas.cpp.

◆ CalculateFluctuations() [2/2]

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

Calculate the fluctuations of the particle species.

Parameters
orderOrder of the fluctuations (default is not specified).
Returns
Vector of fluctuations, one element per particle species.
Examples
NeutronStars-CSHRG.cpp.

Definition at line 1197 of file ThermalModelRealGas.cpp.

◆ CalculateMesonMatterEntropyDensity()

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

Calculate the meson matter entropy density of the system.

Returns
The meson matter entropy density of the system (always returns 0).

Definition at line 217 of file ThermalModelRealGas.h.

◆ CalculatePressure()

double thermalfist::ThermalModelRealGas::CalculatePressure ( )
virtual

Calculate the pressure of the system.

Returns
The pressure of the system.

Definition at line 1687 of file ThermalModelRealGas.cpp.

◆ CalculatePrimordialDensities()

void thermalfist::ThermalModelRealGas::CalculatePrimordialDensities ( )
virtual

Calculate the primordial densities of the particle species.

Implements thermalfist::ThermalModelBase.

Examples
NeutronStars-CSHRG.cpp.

Definition at line 262 of file ThermalModelRealGas.cpp.

◆ CalculateTemperatureDerivatives()

void thermalfist::ThermalModelRealGas::CalculateTemperatureDerivatives ( )
virtual

Calculate the temperature derivatives of the system.

Reimplemented from thermalfist::ThermalModelBase.

Definition at line 1322 of file ThermalModelRealGas.cpp.

◆ CalculateTwoParticleCorrelations()

void thermalfist::ThermalModelRealGas::CalculateTwoParticleCorrelations ( )
virtual

Calculate the two-particle correlations of the particle species.

Reimplemented from thermalfist::ThermalModelBase.

Definition at line 1230 of file ThermalModelRealGas.cpp.

◆ ChiArb()

double thermalfist::ThermalModelRealGas::ChiArb ( int charge)

◆ ComponentIndices()

const std::vector< std::vector< int > > & thermalfist::ThermalModelRealGas::ComponentIndices ( ) const
inline

Get the component indices of the particle species.

Returns
Vector of component indices, one element per particle species.

Definition at line 249 of file ThermalModelRealGas.h.

◆ DeltaMu()

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

Get the delta mu (chemical potential shift due to interactions) of a particle species.

Parameters
i0-based particle species index.
Returns
The delta mu of the particle species.

Definition at line 257 of file ThermalModelRealGas.h.

◆ DensityId()

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

Get the ideal gas density of a particle species.

Parameters
index0-based particle species index.
Returns
The ideal gas density of the particle species.

Definition at line 240 of file ThermalModelRealGas.h.

◆ ExcludedVolumeModel()

ExcludedVolumeModelMultiBase * thermalfist::ThermalModelRealGas::ExcludedVolumeModel ( ) const
inline

Get the excluded volume model used in the real gas HRG model.

Returns
Pointer to the excluded volume model used.

Definition at line 60 of file ThermalModelRealGas.h.

◆ FillChemicalPotentials()

void thermalfist::ThermalModelRealGas::FillChemicalPotentials ( )
virtual

Fill the chemical potentials of the particle species.

Reimplemented from thermalfist::ThermalModelBase.

Definition at line 69 of file ThermalModelRealGas.cpp.

◆ GetMuStar()

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

Get the vector of shifted chemical potentials for all particle species.

Returns
Vector of shifted chemical potentials, one element per particle species.

Definition at line 106 of file ThermalModelRealGas.h.

◆ IsLastSolutionOK()

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

Check if the last solution was successful.

Returns
True if the last solution was successful, false otherwise.

Definition at line 232 of file ThermalModelRealGas.h.

◆ MeanFieldModel()

MeanFieldModelMultiBase * thermalfist::ThermalModelRealGas::MeanFieldModel ( ) const
inline

Get the mean field model used in the real gas HRG model.

Returns
Pointer to the mean field model used.

Definition at line 67 of file ThermalModelRealGas.h.

◆ MuShift()

double thermalfist::ThermalModelRealGas::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 1712 of file ThermalModelRealGas.cpp.

◆ MuStar()

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

Get the shifted chemical potential of a particle species.

Parameters
i0-based particle species index.
Returns
The shifted chemical potential of the particle species.

Definition at line 99 of file ThermalModelRealGas.h.

◆ ParticleScalarDensity()

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

Calculate the scalar density of a particle species.

Parameters
part0-based particle species index.
Returns
The scalar density of the particle species.

Definition at line 1706 of file ThermalModelRealGas.cpp.

◆ SearchMultipleSolutions()

vector< double > thermalfist::ThermalModelRealGas::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 188 of file ThermalModelRealGas.cpp.

◆ SearchSingleSolution()

std::vector< double > thermalfist::ThermalModelRealGas::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 121 of file ThermalModelRealGas.cpp.

◆ SearchSingleSolutionDirect()

vector< double > thermalfist::ThermalModelRealGas::SearchSingleSolutionDirect ( 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 for all particles 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 160 of file ThermalModelRealGas.cpp.

◆ SearchSingleSolutionUsingComponents()

vector< double > thermalfist::ThermalModelRealGas::SearchSingleSolutionUsingComponents ( 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, using the component mapping.

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

Definition at line 127 of file ThermalModelRealGas.cpp.

◆ SetChemicalPotentials()

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

Set the chemical potentials of the particle species.

Parameters
chemVector of chemical potentials, one element per particle species. If empty, the chemical potentials are set to zero.

Reimplemented from thermalfist::ThermalModelBase.

Definition at line 75 of file ThermalModelRealGas.cpp.

◆ SetExcludedVolumeModel()

void thermalfist::ThermalModelRealGas::SetExcludedVolumeModel ( ExcludedVolumeModelMultiBase * exvolmod)
inline

Set the excluded volume model for the real gas HRG model.

Parameters
exvolmodPointer to the excluded volume model to be used.
Examples
NeutronStars-CSHRG.cpp.

Definition at line 46 of file ThermalModelRealGas.h.

◆ SetMeanFieldModel()

void thermalfist::ThermalModelRealGas::SetMeanFieldModel ( MeanFieldModelMultiBase * mfmod)
inline

Set the mean field model for the real gas HRG model.

Parameters
mfmodPointer to the mean field model to be used.

Definition at line 53 of file ThermalModelRealGas.h.

◆ SetMultipleSolutionsMode()

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

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

Multiple solutions in the real gas 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.

Definition at line 82 of file ThermalModelRealGas.h.

◆ SetMuStar()

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

Set the vector of shifted chemical potentials for all particle species.

Parameters
MuStarVector of shifted chemical potentials, one element per particle species.

Definition at line 113 of file ThermalModelRealGas.h.

◆ SolveEquations()

void thermalfist::ThermalModelRealGas::SolveEquations ( )
protected

Solve the transcedental equations for the shifted chemical potentials

Definition at line 246 of file ThermalModelRealGas.cpp.

◆ UseMultipleSolutionsMode()

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

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

Returns
true Multiple solutions considered
false Multiple solutions not considered

Definition at line 91 of file ThermalModelRealGas.h.

Member Data Documentation

◆ m_chi

std::vector< std::vector<double> > thermalfist::ThermalModelRealGas::m_chi

Vector of computed susceptibilities values.

Definition at line 262 of file ThermalModelRealGas.h.

◆ m_chiarb

std::vector<double> thermalfist::ThermalModelRealGas::m_chiarb

Vector of computed susceptibilities for a specified arbitraty charge.

Definition at line 267 of file ThermalModelRealGas.h.

◆ m_ComponentMapCalculated

bool thermalfist::ThermalModelRealGas::m_ComponentMapCalculated
protected

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

Definition at line 352 of file ThermalModelRealGas.h.

◆ m_DensitiesId

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

Vector of ideal gas densities with shifted chemical potentials.

Definition at line 340 of file ThermalModelRealGas.h.

◆ m_dMuStarIndices

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

Vector of component indices for each particle species.

Definition at line 364 of file ThermalModelRealGas.h.

◆ m_exvolmod

ExcludedVolumeModelMultiBase* thermalfist::ThermalModelRealGas::m_exvolmod
protected

Excluded volume model used in the real gas HRG model.

Definition at line 368 of file ThermalModelRealGas.h.

◆ m_exvolmodideal

ExcludedVolumeModelMultiBase* thermalfist::ThermalModelRealGas::m_exvolmodideal
protected

Excluded volume model object in the ideal gas limit.

Definition at line 374 of file ThermalModelRealGas.h.

◆ m_LastBroydenSuccessFlag

bool thermalfist::ThermalModelRealGas::m_LastBroydenSuccessFlag
protected

Whether Broyden's method was successfull.

Definition at line 349 of file ThermalModelRealGas.h.

◆ m_MapFromdMuStar

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

Mapping from dMuStar indices to particle species.

Definition at line 361 of file ThermalModelRealGas.h.

◆ m_MapTodMuStar

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

Mapping from particle species to dMuStar indices.

Definition at line 358 of file ThermalModelRealGas.h.

◆ m_mfmod

MeanFieldModelMultiBase* thermalfist::ThermalModelRealGas::m_mfmod
protected

Mean field model used in the real gas HRG model.

Definition at line 371 of file ThermalModelRealGas.h.

◆ m_mfmodideal

MeanFieldModelMultiBase* thermalfist::ThermalModelRealGas::m_mfmodideal
protected

Mean field model object in the ideal gas limit.

Definition at line 377 of file ThermalModelRealGas.h.

◆ m_MuStar

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

Vector of the shifted chemical potentials.

Definition at line 355 of file ThermalModelRealGas.h.

◆ m_scaldens

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

Vector of scalar densities. Not used.

Definition at line 343 of file ThermalModelRealGas.h.

◆ m_SearchMultipleSolutions

bool thermalfist::ThermalModelRealGas::m_SearchMultipleSolutions
protected

Whether multiple solutions are considered.

Definition at line 346 of file ThermalModelRealGas.h.


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