Loading [MathJax]/extensions/tex2jax.js
Thermal-FIST 1.5
Package for hadron resonance gas model applications
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Modules Pages
thermalfist::ThermalModelBase Class Referenceabstract

Abstract base class for an HRG model implementation. More...

#include <ThermalModelBase.h>

Inheritance diagram for thermalfist::ThermalModelBase:
thermalfist::ThermalModelCanonical thermalfist::ThermalModelCanonicalCharm thermalfist::ThermalModelCanonicalStrangeness thermalfist::ThermalModelEVCrosstermsLegacy thermalfist::ThermalModelEVDiagonal thermalfist::ThermalModelIdeal thermalfist::ThermalModelRealGas thermalfist::ThermalModelVDW

Public Types

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...
 

Public Member Functions

 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().
 
virtual 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.
 
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 CalculatePrimordialDensities ()=0
 Calculates the primordial densities of all species.
 
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 CalculateFluctuations ()
 Computes the fluctuation observables.
 
virtual void CalculateTemperatureDerivatives ()
 Computes the temperature derivatives of densities, shifted chemical potentials, and primordial hadron number susceptibilities.
 
virtual void CalculateTwoParticleCorrelations ()
 Computes the fluctuations and correlations of the primordial particle numbers.
 
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.
 
virtual std::vector< double > CalculateChargeFluctuations (const std::vector< double > &chgs, int order=4)
 Calculates fluctuations (diagonal susceptibilities) of an arbitrary "conserved" charge.
 

Detailed Description

Abstract base class for an HRG model implementation.

Contains the base implementation of an arbitrary HRG model. The actual calculations of thermodynamic functions are implemented in derived classes.

Examples
CalculationTmu.cpp, SusceptibilitiesBQS.cpp, ThermodynamicsBQS.cpp, cpc1-HRG-TDep.cpp, cpc2-chi2-vs-T.cpp, cpc3-chi2NEQ.cpp, and cpc4-mcHRG.cpp.

Definition at line 28 of file ThermalModelBase.h.

Member Enumeration Documentation

◆ ThermalModelEnsemble

The list of statistical ensembles.

Enumerator
GCE 

Grand canonical ensemble.

CE 

Canonical ensemble.

SCE 

Strangeness-canonical ensemble.

CCE 

Charm-canonical ensemble.

Definition at line 35 of file ThermalModelBase.h.

◆ ThermalModelInteraction

Type of interactions included in the HRG model.

Enumerator
Ideal 

Ideal HRG model.

DiagonalEV 

Diagonal excluded volume model.

CrosstermsEV 

Crossterms excluded volume model.

QvdW 

Quantum van der Waals model.

RealGas 

Real gas model.

MeanField 

Mean field model. Not yet fully implemented.

Definition at line 46 of file ThermalModelBase.h.

Constructor & Destructor Documentation

◆ ThermalModelBase()

thermalfist::ThermalModelBase::ThermalModelBase ( ThermalParticleSystem * TPS,
const ThermalModelParameters & params = ThermalModelParameters() )

Construct a new ThermalModelBase object.

Parameters
TPSA pointer to the ThermalParticleSystem object containing the particle list
paramsThermalModelParameters object with current thermal parameters (default is ThermalModelParameters())

Definition at line 24 of file ThermalModelBase.cpp.

◆ ~ThermalModelBase()

virtual thermalfist::ThermalModelBase::~ThermalModelBase ( void )
inlinevirtual

Definition at line 63 of file ThermalModelBase.h.

Member Function Documentation

◆ AttractionCoefficient()

virtual double thermalfist::ThermalModelBase::AttractionCoefficient ( int ,
int  ) const
inlinevirtual

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 in thermalfist::ThermalModelVDW, and thermalfist::ThermalModelVDWCanonicalStrangeness.

Definition at line 344 of file ThermalModelBase.h.

◆ CalculateChargeFluctuations()

std::vector< double > thermalfist::ThermalModelBase::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 in thermalfist::ThermalModelEVCrosstermsLegacy, thermalfist::ThermalModelEVDiagonal, thermalfist::ThermalModelIdeal, thermalfist::ThermalModelRealGas, and thermalfist::ThermalModelVDW.

Examples
CalculationTmu.cpp, and cpc1-HRG-TDep.cpp.

Definition at line 652 of file ThermalModelBase.cpp.

◆ CalculateDensities()

void thermalfist::ThermalModelBase::CalculateDensities ( )
virtual

Calculates the primordial and total (after decays) densities of all species.

Examples
CalculationTmu.cpp, and cpc1-HRG-TDep.cpp.

Definition at line 626 of file ThermalModelBase.cpp.

◆ CalculateDensitiesGCE()

virtual void thermalfist::ThermalModelBase::CalculateDensitiesGCE ( )
inlinevirtual

Calculates the particle densities in a grand-canonical ensemble.

A non-GCE based derived class will override this method.

Reimplemented in thermalfist::ThermalModelCanonicalCharm, thermalfist::ThermalModelCanonicalStrangeness, thermalfist::ThermalModelEVCanonicalStrangeness, and thermalfist::ThermalModelVDWCanonicalStrangeness.

Definition at line 684 of file ThermalModelBase.h.

◆ CalculateFeeddown()

void thermalfist::ThermalModelBase::CalculateFeeddown ( )
virtual

Calculates the total densities which include feeddown contributions.

Calculation is based on the primordial densities computed by CalculateDensities().

Definition at line 407 of file ThermalModelBase.cpp.

◆ CalculateFluctuations()

void thermalfist::ThermalModelBase::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 in thermalfist::ThermalModelCanonical, thermalfist::ThermalModelCanonicalCharm, thermalfist::ThermalModelCanonicalStrangeness, thermalfist::ThermalModelEVCrosstermsLegacy, thermalfist::ThermalModelEVDiagonal, thermalfist::ThermalModelIdeal, thermalfist::ThermalModelRealGas, and thermalfist::ThermalModelVDW.

Examples
CalculationTmu.cpp, SusceptibilitiesBQS.cpp, and cpc4-mcHRG.cpp.

Definition at line 1418 of file ThermalModelBase.cpp.

◆ CalculateParticleChargeCorrelationMatrix()

void thermalfist::ThermalModelBase::CalculateParticleChargeCorrelationMatrix ( )
virtual

Calculates the matrix of correlators between primordial (and also final) particle numbers and conserved charges.

Correlators are normalized by VT^3.

The calculation results are accessible through PrimordialParticleChargeCorrelation() and FinalParticleChargeCorrelation().

Definition at line 1366 of file ThermalModelBase.cpp.

◆ CalculatePrimordialDensities()

◆ CalculateProxySusceptibilityMatrix()

void thermalfist::ThermalModelBase::CalculateProxySusceptibilityMatrix ( )
virtual

Calculates the susceptibility matrix of conserved charges proxies.

The following proxies are used: final net proton number for baryon number, net charge as is, and net kaon number for net strangeness. Charm not yet considered.

Decay feeddown contributions are in accordance with the stability flags.

The calculation results are accessible through ProxySusc()

Definition at line 1326 of file ThermalModelBase.cpp.

◆ CalculateSusceptibilityMatrix()

void thermalfist::ThermalModelBase::CalculateSusceptibilityMatrix ( )
virtual

Calculates the conserved charges susceptibility matrix.

Computes \( \chi_{ijkl}^{BSQC}~=~\frac{\partial^{l+m+n+k}p/T^4}{\partial(\mu_B/T)^i \,\partial(\mu_S/T)^j \,\partial(\mu_Q/T)^k\,\partial(\mu_C/T)^l} \) where i+j+k+l = 2.

If temperature derivatives are available, calculates also the temperature derivative of the conserved charges susceptibility matrix.

\( \frac{\partial\chi_{ijkl}^{BSQC}}{\partial T} \)

The calculation results are accessible through Susc() and dSuscdT()

Definition at line 1292 of file ThermalModelBase.cpp.

◆ CalculateTemperatureDerivatives()

void thermalfist::ThermalModelBase::CalculateTemperatureDerivatives ( )
virtual

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

Reimplemented in thermalfist::ThermalModelEVDiagonal, thermalfist::ThermalModelIdeal, thermalfist::ThermalModelRealGas, and thermalfist::ThermalModelVDW.

Examples
SusceptibilitiesBQS.cpp.

Definition at line 1422 of file ThermalModelBase.cpp.

◆ CalculateTwoParticleCorrelations()

void thermalfist::ThermalModelBase::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 in thermalfist::ThermalModelCanonical, thermalfist::ThermalModelEVCrosstermsLegacy, thermalfist::ThermalModelEVDiagonal, thermalfist::ThermalModelIdeal, thermalfist::ThermalModelRealGas, and thermalfist::ThermalModelVDW.

Definition at line 964 of file ThermalModelBase.cpp.

◆ CalculateTwoParticleFluctuationsDecays()

void thermalfist::ThermalModelBase::CalculateTwoParticleFluctuationsDecays ( )
virtual

Computes particle number correlations and fluctuations for all final-state particles which are marked stable.

More specifically, computes the susceptibility matrix \( \frac{1}{VT^3} \, \langle \Delta N_i \Delta N_j \rangle \), where \( N_i \) is the final yield of species i, and indices i and j corresponds to stable particles only.

To incorporate probabilistic decays uses the formalism from paper V.V. Begun, M.I. Gorenstein, M. Hauer, V. Konchakovski, O.S. Zozulya, Phys.Rev. C 74, 044903 (2006), https://arxiv.org/pdf/nucl-th/0606036.pdf

Definition at line 969 of file ThermalModelBase.cpp.

◆ CanonicalVolume()

double thermalfist::ThermalModelBase::CanonicalVolume ( ) const
inline

The canonical correlation volume V \(_c\) (fm \(^3\))

Definition at line 530 of file ThermalModelBase.h.

◆ ChangeTPS()

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

◆ ChemicalPotential()

double thermalfist::ThermalModelBase::ChemicalPotential ( int i) const

Chemical potential of particle species i.

Parameters
i0-based index of particle species
Returns
double Chemical potential of particle species i

Definition at line 372 of file ThermalModelBase.cpp.

◆ ChemicalPotentials()

const std::vector< double > & thermalfist::ThermalModelBase::ChemicalPotentials ( ) const
inline

A vector of chemical potentials of all particles.

Returns
const std::vector<double>& A vector of chemical potentials of all particles

Definition at line 418 of file ThermalModelBase.h.

◆ ComponentsNumber()

int thermalfist::ThermalModelBase::ComponentsNumber ( ) const
inline

Number of different particle species in the list.

Definition at line 66 of file ThermalModelBase.h.

◆ ConstrainChemicalPotentials()

void thermalfist::ThermalModelBase::ConstrainChemicalPotentials ( bool resetInitialValues = true)

Constrains the chemical potentials \( \mu_B,\,\mu_Q,\,\mu_S,\,\mu_C \) by the conservation laws imposed.

This procedure uses the Broyden's method to solve the system of equations corresponding to the conservation laws. The actual implementation of the procedure is in FixParameters() and FixParametersNoReset() methods.

Parameters
resetInitialValuesWhether initial guess values for \( \mu_B,\,\mu_Q,\,\mu_S,\,\mu_C \) are reset or current values will be used
Examples
NeutronStars-CSHRG.cpp, PCE-Saha-LHC.cpp, and cpc4-mcHRG.cpp.

Definition at line 447 of file ThermalModelBase.cpp.

◆ ConstrainMuB() [1/2]

bool thermalfist::ThermalModelBase::ConstrainMuB ( ) const
inline

Whether the baryon chemical potential is to be constrained by a fixed entropy per baryon ratio SoverB().

Definition at line 450 of file ThermalModelBase.h.

◆ ConstrainMuB() [2/2]

void thermalfist::ThermalModelBase::ConstrainMuB ( bool constrain)
inline

Sets whether the baryon chemical potential is to be constrained by a fixed entropy per baryon ratio SoverB().

Definition at line 454 of file ThermalModelBase.h.

◆ ConstrainMuC() [1/2]

bool thermalfist::ThermalModelBase::ConstrainMuC ( ) const
inline

Whether the charm chemical potential is to be constrained by the condition of charm neutrality

Examples
cpc4-mcHRG.cpp.

Definition at line 474 of file ThermalModelBase.h.

◆ ConstrainMuC() [2/2]

void thermalfist::ThermalModelBase::ConstrainMuC ( bool constrain)
inline

Sets whether the charm chemical potential is to be constrained by the condition of charm neutrality

Definition at line 478 of file ThermalModelBase.h.

◆ ConstrainMuQ() [1/2]

bool thermalfist::ThermalModelBase::ConstrainMuQ ( ) const
inline

Whether the electric chemical potential is to be constrained by a fixed electric-to-baryon ratio QoverB().

Examples
CalculationTmu.cpp, and cpc4-mcHRG.cpp.

Definition at line 458 of file ThermalModelBase.h.

◆ ConstrainMuQ() [2/2]

void thermalfist::ThermalModelBase::ConstrainMuQ ( bool constrain)
inline

Sets whether the electric chemical potential is to be constrained by a fixed electric-to-baryon ratio QoverB().

Definition at line 462 of file ThermalModelBase.h.

◆ ConstrainMuS() [1/2]

bool thermalfist::ThermalModelBase::ConstrainMuS ( ) const
inline

Whether the strangeness chemical potential is to be constrained by the condition of strangeness neutrality

Examples
CalculationTmu.cpp, NeutronStars-CSHRG.cpp, and cpc4-mcHRG.cpp.

Definition at line 466 of file ThermalModelBase.h.

◆ ConstrainMuS() [2/2]

void thermalfist::ThermalModelBase::ConstrainMuS ( bool constrain)
inline

Sets whether the strangeness chemical potential is to be constrained by the condition of strangeness neutrality

Definition at line 470 of file ThermalModelBase.h.

◆ DisableBaryonAntiBaryonAttraction()

void thermalfist::ThermalModelBase::DisableBaryonAntiBaryonAttraction ( )
virtual

Switches off QvdW attraction terms for all baryon-antibaryon pairs i.e. \( a_{ij} = 0 \) for baryon-antibaryon pairs

Reimplemented in thermalfist::ThermalModelEVDiagonal.

Definition at line 302 of file ThermalModelBase.cpp.

◆ DisableBaryonAntiBaryonRepulsion()

void thermalfist::ThermalModelBase::DisableBaryonAntiBaryonRepulsion ( )
inline

Definition at line 290 of file ThermalModelBase.h.

◆ DisableBaryonAntiBaryonVirial()

void thermalfist::ThermalModelBase::DisableBaryonAntiBaryonVirial ( )
virtual

Switches off eigenvolume terms for all baryon-antibaryon pairs i.e. \( \tilde{b}_{ij} = 0 \) for baryon-antibaryon pairs

Reimplemented in thermalfist::ThermalModelEVDiagonal.

Definition at line 289 of file ThermalModelBase.cpp.

◆ DisableBaryonBaryonAttraction()

void thermalfist::ThermalModelBase::DisableBaryonBaryonAttraction ( )
virtual

Switches off QvdW attraction terms for all baryon-baryon pairs i.e. \( a_{ij} = 0 \) for baryon-baryon pairs

Reimplemented in thermalfist::ThermalModelEVDiagonal.

Definition at line 276 of file ThermalModelBase.cpp.

◆ DisableBaryonBaryonRepulsion()

void thermalfist::ThermalModelBase::DisableBaryonBaryonRepulsion ( )
inline

Definition at line 279 of file ThermalModelBase.h.

◆ DisableBaryonBaryonVirial()

void thermalfist::ThermalModelBase::DisableBaryonBaryonVirial ( )
virtual

Switches off excluded volume terms for all baryon-baryon pairs i.e. \( \tilde{b}_{ij} = 0 \) for baryon-baryon pairs

Reimplemented in thermalfist::ThermalModelEVDiagonal.

Definition at line 263 of file ThermalModelBase.cpp.

◆ DisableMesonBaryonAttraction()

void thermalfist::ThermalModelBase::DisableMesonBaryonAttraction ( )
virtual

Switches off QvdW attraction terms for all meson-baryon pairs i.e. \( a_{ij} = 0 \) for meson-baryon pairs

Reimplemented in thermalfist::ThermalModelEVDiagonal.

Definition at line 250 of file ThermalModelBase.cpp.

◆ DisableMesonBaryonRepulsion()

void thermalfist::ThermalModelBase::DisableMesonBaryonRepulsion ( )
inline

Definition at line 268 of file ThermalModelBase.h.

◆ DisableMesonBaryonVirial()

void thermalfist::ThermalModelBase::DisableMesonBaryonVirial ( )
virtual

Switches off excluded volume terms for all meson-baryon pairs i.e. \( \tilde{b}_{ij} = 0 \) for meson-baryon pairs

Reimplemented in thermalfist::ThermalModelEVDiagonal.

Definition at line 237 of file ThermalModelBase.cpp.

◆ DisableMesonMesonAttraction()

void thermalfist::ThermalModelBase::DisableMesonMesonAttraction ( )
virtual

Switches off QvdW attraction terms for all meson-meson pairs i.e. \( a_{ij} = 0 \) for meson-meson pairs

Reimplemented in thermalfist::ThermalModelEVDiagonal.

Definition at line 225 of file ThermalModelBase.cpp.

◆ DisableMesonMesonRepulsion()

void thermalfist::ThermalModelBase::DisableMesonMesonRepulsion ( )
inline

Definition at line 257 of file ThermalModelBase.h.

◆ DisableMesonMesonVirial()

void thermalfist::ThermalModelBase::DisableMesonMesonVirial ( )
virtual

Switches off excluded volume terms for all meson-meson pairs i.e. \( \tilde{b}_{ij} = 0 \) for meson-meson pairs

Reimplemented in thermalfist::ThermalModelEVDiagonal.

Definition at line 213 of file ThermalModelBase.cpp.

◆ FillChemicalPotentials()

void thermalfist::ThermalModelBase::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 in thermalfist::ThermalModelRealGas, and thermalfist::ThermalModelVDW.

Examples
PCE-Saha-LHC.cpp, cpc1-HRG-TDep.cpp, and cpc2-chi2-vs-T.cpp.

Definition at line 357 of file ThermalModelBase.cpp.

◆ FillVirial()

void thermalfist::ThermalModelBase::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 in thermalfist::ThermalModelEVCanonicalStrangeness, thermalfist::ThermalModelEVCrosstermsLegacy, thermalfist::ThermalModelEVDiagonal, thermalfist::ThermalModelVDW, and thermalfist::ThermalModelVDWCanonicalStrangeness.

Definition at line 94 of file ThermalModelBase.cpp.

◆ FinalNetParticleChargeSusceptibilityByPdg()

double thermalfist::ThermalModelBase::FinalNetParticleChargeSusceptibilityByPdg ( long long id1,
ConservedCharge::Name chg )
virtual

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.

Definition at line 1277 of file ThermalModelBase.cpp.

◆ FinalParticleChargeSusceptibility()

double thermalfist::ThermalModelBase::FinalParticleChargeSusceptibility ( int i,
ConservedCharge::Name chg ) const
virtual

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.

Definition at line 1257 of file ThermalModelBase.cpp.

◆ FinalParticleChargeSusceptibilityByPdg()

double thermalfist::ThermalModelBase::FinalParticleChargeSusceptibilityByPdg ( long long id1,
ConservedCharge::Name chg )
virtual

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.

Definition at line 1266 of file ThermalModelBase.cpp.

◆ FixChemicalPotentialsThroughDensities()

bool thermalfist::ThermalModelBase::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 )
virtual

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.

Uses the Broyden's method to determine the chemical potentials. The resulting chemical potentials are stored in Parameters().

Parameters
rhoBThe total desired net baryon charge of the system.
rhoQThe total desired electric charge of the system.
rhoSThe total desired net strangeness of the system.
rhoCThe total desired net charm of the system.
muBinitThe initial guess for the baryon chemical potential.
muQinitThe initial guess for the electric chemical potential.
muSinitThe initial guess for the strangeness chemical potential.
muCinitThe initial guess for the charm chemical potential.
ConstrMuBWhether the baryon chemical potential should be constrained.
ConstrMuQWhether the electric chemical potential should be constrained.
ConstrMuSWhether the strangeness chemical potential should be constrained.
ConstrMuCWhether the charm chemical potential should be constrained.
Returns
true is chemical potentials were contrained successfully, false otherwise

Definition at line 529 of file ThermalModelBase.cpp.

◆ FixParameters()

void thermalfist::ThermalModelBase::FixParameters ( )
virtual

Method which actually implements ConstrainChemicalPotentials() (for backward compatibility).

Reimplemented in thermalfist::ThermalModelCanonical, thermalfist::ThermalModelCanonicalCharm, and thermalfist::ThermalModelCanonicalStrangeness.

Examples
CalculationTmu.cpp.

Definition at line 455 of file ThermalModelBase.cpp.

◆ FixParametersNoReset()

void thermalfist::ThermalModelBase::FixParametersNoReset ( )
virtual

Method which actually implements ConstrainChemicalPotentialsNoReset() (for backward compatibility).

Reimplemented in thermalfist::ThermalModelCanonical.

Definition at line 488 of file ThermalModelBase.cpp.

◆ FullIdealChemicalPotential()

double thermalfist::ThermalModelBase::FullIdealChemicalPotential ( int i) const
virtual

Chemical potential entering the ideal gas expressions of particle species i.

Includes chemical non-equilibrium and EV/vdW effects

Parameters
i0-based index of particle species
Returns
double Chemical potential of particle species i

Definition at line 388 of file ThermalModelBase.cpp.

◆ NetParticleSusceptibilityFinalByPdg()

double thermalfist::ThermalModelBase::NetParticleSusceptibilityFinalByPdg ( long long id1,
long long id2 )
virtual

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.

Definition at line 1186 of file ThermalModelBase.cpp.

◆ NetParticleSusceptibilityPrimordialByPdg()

double thermalfist::ThermalModelBase::NetParticleSusceptibilityPrimordialByPdg ( long long id1,
long long id2 )
virtual

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.

Definition at line 1117 of file ThermalModelBase.cpp.

◆ NormBratio()

bool thermalfist::ThermalModelBase::NormBratio ( ) const
inline

Definition at line 108 of file ThermalModelBase.h.

◆ Parameters()

const ThermalModelParameters & thermalfist::ThermalModelBase::Parameters ( ) const
inline
Examples
NeutronStars-CSHRG.cpp, PCE-Saha-LHC.cpp, and cpc4-mcHRG.cpp.

Definition at line 121 of file ThermalModelBase.h.

◆ PrimordialNetParticleChargeSusceptibilityByPdg()

double thermalfist::ThermalModelBase::PrimordialNetParticleChargeSusceptibilityByPdg ( long long id1,
ConservedCharge::Name chg )
virtual

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.

Definition at line 1242 of file ThermalModelBase.cpp.

◆ PrimordialParticleChargeSusceptibility()

double thermalfist::ThermalModelBase::PrimordialParticleChargeSusceptibility ( int i,
ConservedCharge::Name chg ) const
virtual

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.

Definition at line 1222 of file ThermalModelBase.cpp.

◆ PrimordialParticleChargeSusceptibilityByPdg()

double thermalfist::ThermalModelBase::PrimordialParticleChargeSusceptibilityByPdg ( long long id1,
ConservedCharge::Name chg )
virtual

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.

Definition at line 1231 of file ThermalModelBase.cpp.

◆ QoverB()

double thermalfist::ThermalModelBase::QoverB ( ) const
inline

Definition at line 507 of file ThermalModelBase.h.

◆ QuantumStatistics()

bool thermalfist::ThermalModelBase::QuantumStatistics ( ) const
inline

Whether quantum statistics is used, 0 - Boltzmann, 1 - Quantum

Definition at line 348 of file ThermalModelBase.h.

◆ ReadInteractionParameters()

virtual void thermalfist::ThermalModelBase::ReadInteractionParameters ( const std::string & )
inlinevirtual

Reads the QvdW interaction parameters from a file.

Actual implementation is in a derived class.

Parameters
filenameFile with interaction parameters.

Reimplemented in thermalfist::ThermalModelEVCanonicalStrangeness, thermalfist::ThermalModelEVCrossterms, thermalfist::ThermalModelEVCrosstermsLegacy, thermalfist::ThermalModelEVDiagonal, thermalfist::ThermalModelVDW, and thermalfist::ThermalModelVDWCanonicalStrangeness.

Definition at line 304 of file ThermalModelBase.h.

◆ RepulsionCoefficient()

double thermalfist::ThermalModelBase::RepulsionCoefficient ( int i,
int j ) const
inline

Definition at line 334 of file ThermalModelBase.h.

◆ SetAttraction()

virtual void thermalfist::ThermalModelBase::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 in thermalfist::ThermalModelEVCrossterms, thermalfist::ThermalModelVDW, and thermalfist::ThermalModelVDWCanonicalStrangeness.

Examples
CalculationTmu.cpp, cpc1-HRG-TDep.cpp, and cpc2-chi2-vs-T.cpp.

Definition at line 251 of file ThermalModelBase.h.

◆ SetBaryonCharge()

void thermalfist::ThermalModelBase::SetBaryonCharge ( int B)
virtual

Set the total baryon number (for canonical ensemble only)

Parameters
BTotal baryon number

Definition at line 189 of file ThermalModelBase.cpp.

◆ SetBaryonChemicalPotential()

void thermalfist::ThermalModelBase::SetBaryonChemicalPotential ( double muB)
virtual

Set the baryon chemical potential.

Parameters
muBBaryon chemical potential (GeV)
Examples
BagModelFit.cpp, CalculationTmu.cpp, NeutronStars-CSHRG.cpp, SusceptibilitiesBQS.cpp, ThermodynamicsBQS.cpp, cpc1-HRG-TDep.cpp, cpc2-chi2-vs-T.cpp, and cpc4-mcHRG.cpp.

Definition at line 149 of file ThermalModelBase.cpp.

◆ SetCalculationType()

virtual void thermalfist::ThermalModelBase::SetCalculationType ( IdealGasFunctions::QStatsCalculationType type)
inlinevirtual

Sets the CalculationType() method to evaluate quantum statistics. Calls the corresponding method in TPS().

Parameters
typeMethod to evaluate quantum statistics.
Examples
NeutronStars-CSHRG.cpp, SusceptibilitiesBQS.cpp, ThermodynamicsBQS.cpp, and cpc3-chi2NEQ.cpp.

Definition at line 360 of file ThermalModelBase.h.

◆ SetCanonicalVolume()

void thermalfist::ThermalModelBase::SetCanonicalVolume ( double Volume)
inline

Set the canonical correlation volume V \(_c\).

Parameters
VolumeCanonical correlation volume (fm \(^3\))

Definition at line 537 of file ThermalModelBase.h.

◆ SetCanonicalVolumeRadius()

void thermalfist::ThermalModelBase::SetCanonicalVolumeRadius ( double Rc)
inline

Set the canonical correlation system radius.

The canonical correlation volume is computed as \(V_c = \frac{4\pi}{3} \, R_c^3\)

Parameters
Radius

Definition at line 547 of file ThermalModelBase.h.

◆ SetCharm()

void thermalfist::ThermalModelBase::SetCharm ( int C)
virtual

Set the total charm (for canonical ensemble only)

Parameters
BTotal charm

Definition at line 207 of file ThermalModelBase.cpp.

◆ SetCharmChemicalPotential()

void thermalfist::ThermalModelBase::SetCharmChemicalPotential ( double muC)
virtual

Set the charm chemical potential.

Parameters
muCCharm chemical potential (GeV)

Reimplemented in thermalfist::ThermalModelCanonicalCharm.

Examples
BagModelFit.cpp, cpc1-HRG-TDep.cpp, cpc2-chi2-vs-T.cpp, and cpc4-mcHRG.cpp.

Definition at line 170 of file ThermalModelBase.cpp.

◆ SetChemicalPotential()

void thermalfist::ThermalModelBase::SetChemicalPotential ( int i,
double chem )
virtual

Sets the chemical potential of particle species i.

Parameters
i0-based index of particle species
chemvalue of the chemical potential

Definition at line 380 of file ThermalModelBase.cpp.

◆ SetChemicalPotentials()

void thermalfist::ThermalModelBase::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 in thermalfist::ThermalModelRealGas, and thermalfist::ThermalModelVDW.

Definition at line 363 of file ThermalModelBase.cpp.

◆ SetClusterExpansionOrder()

virtual void thermalfist::ThermalModelBase::SetClusterExpansionOrder ( int order)
inlinevirtual

Set the number of terms in the cluster expansion method. Calls the corresponding method in TPS().

Sets the same value for all particles. To set individually for each particle use ThermalParticle::SetClusterExpansionOrder().

Parameters
orderNumber of terms.

Definition at line 372 of file ThermalModelBase.h.

◆ SetElectricCharge()

void thermalfist::ThermalModelBase::SetElectricCharge ( int Q)
virtual

Set the total electric charge (for canonical ensemble only)

Parameters
BTotal electric charge

Definition at line 195 of file ThermalModelBase.cpp.

◆ SetElectricChemicalPotential()

void thermalfist::ThermalModelBase::SetElectricChemicalPotential ( double muQ)
virtual

Set the electric chemical potential.

Parameters
muQElectric chemical potential (GeV)
Examples
BagModelFit.cpp, NeutronStars-CSHRG.cpp, SusceptibilitiesBQS.cpp, ThermodynamicsBQS.cpp, cpc1-HRG-TDep.cpp, cpc2-chi2-vs-T.cpp, and cpc4-mcHRG.cpp.

Definition at line 156 of file ThermalModelBase.cpp.

◆ SetGammaC()

void thermalfist::ThermalModelBase::SetGammaC ( double gammaC)
virtual

Set the charm quark fugacity factor.

Parameters
gammaCCharm quark fugacity factor

Definition at line 183 of file ThermalModelBase.cpp.

◆ SetGammaq()

void thermalfist::ThermalModelBase::SetGammaq ( double gammaq)
virtual

Set the light quark fugacity factor.

Parameters
gammaqLight quark fugacity factor
Examples
CalculationTmu.cpp.

Definition at line 315 of file ThermalModelBase.cpp.

◆ SetGammaS()

void thermalfist::ThermalModelBase::SetGammaS ( double gammaS)
virtual

Set the strange quark fugacity factor.

Parameters
gammaSStrange quark fugacity factor
Examples
CalculationTmu.cpp.

Definition at line 177 of file ThermalModelBase.cpp.

◆ SetNormBratio()

void thermalfist::ThermalModelBase::SetNormBratio ( bool normBratio)

Whether branching ratios are renormalized to 100%.

Parameters
normBratioWhether branching ratios shoul be renormalized to 100%.
Examples
cpc4-mcHRG.cpp.

Definition at line 117 of file ThermalModelBase.cpp.

◆ SetOMP()

void thermalfist::ThermalModelBase::SetOMP ( bool openMP)
inline

OpenMP support. Currently not used.

Definition at line 112 of file ThermalModelBase.h.

◆ SetParameters()

void thermalfist::ThermalModelBase::SetParameters ( const ThermalModelParameters & params)
virtual

The thermal parameters.

Parameters
paramsThe thermal parameters to be used

Reimplemented in thermalfist::ThermalModelCanonicalCharm, and thermalfist::ThermalModelCanonicalStrangeness.

Examples
PCE-Saha-LHC.cpp, and cpc4-mcHRG.cpp.

Definition at line 137 of file ThermalModelBase.cpp.

◆ SetQoverB()

void thermalfist::ThermalModelBase::SetQoverB ( double QB)
inline

The electric-to-baryon charge ratio to be used to constrain the electric chemical potential.

Parameters
QBThe electric-to-baryon charge ratio
Examples
CalculationTmu.cpp, and NeutronStars-CSHRG.cpp.

Definition at line 506 of file ThermalModelBase.h.

◆ SetRadius() [1/2]

virtual void thermalfist::ThermalModelBase::SetRadius ( double )
inlinevirtual

Set the same excluded volume radius parameter for all species.

Parameters
radRadius parameter (fm)

Reimplemented in thermalfist::ThermalModelEVCanonicalStrangeness, thermalfist::ThermalModelEVCrossterms, thermalfist::ThermalModelEVCrosstermsLegacy, and thermalfist::ThermalModelEVDiagonal.

Examples
cpc1-HRG-TDep.cpp, and cpc2-chi2-vs-T.cpp.

Definition at line 219 of file ThermalModelBase.h.

◆ SetRadius() [2/2]

virtual void thermalfist::ThermalModelBase::SetRadius ( int ,
double  )
inlinevirtual

Set the radius parameter for particle species i.

Parameters
i0-based index of particle species
radRadius parameter (fm)

Reimplemented in thermalfist::ThermalModelEVCanonicalStrangeness, and thermalfist::ThermalModelEVDiagonal.

Definition at line 227 of file ThermalModelBase.h.

◆ SetRepulsion()

virtual void thermalfist::ThermalModelBase::SetRepulsion ( int i,
int j,
double b )
inlinevirtual

Same as SetVirial() but with a more clear name on what is actually does.

Definition at line 242 of file ThermalModelBase.h.

◆ SetResonanceWidthIntegrationType()

void thermalfist::ThermalModelBase::SetResonanceWidthIntegrationType ( ThermalParticle::ResonanceWidthIntegration type)

Set the ThermalParticle::ResonanceWidthIntegration scheme for all particles. Calls the corresponding method in TPS().

Parameters
typeThermalParticle::ResonanceWidthIntegration

Definition at line 347 of file ThermalModelBase.cpp.

◆ SetResonanceWidthShape()

void thermalfist::ThermalModelBase::SetResonanceWidthShape ( ThermalParticle::ResonanceWidthShape shape)
inline

Set the ThermalParticle::ResonanceWidthShape for all particles. Calls the corresponding method in TPS().

This global setting can be overriden for individual particles by ThermalParticle::SetResonanceWidthShape().

Parameters
shapeThermalParticle::ResonanceWidthShape

Definition at line 383 of file ThermalModelBase.h.

◆ SetSoverB()

void thermalfist::ThermalModelBase::SetSoverB ( double SB)
inline

The entropy per baryon ratio to be used to constrain the baryon chemical potential.

Parameters
SBThe entropy per baryon

Definition at line 494 of file ThermalModelBase.h.

◆ SetStatistics()

void thermalfist::ThermalModelBase::SetStatistics ( bool stats)
virtual

◆ SetStrangeness()

void thermalfist::ThermalModelBase::SetStrangeness ( int S)
virtual

Set the total strangeness (for canonical ensemble only)

Parameters
BTotal strangeness

Definition at line 201 of file ThermalModelBase.cpp.

◆ SetStrangenessChemicalPotential()

void thermalfist::ThermalModelBase::SetStrangenessChemicalPotential ( double muS)
virtual

Set the strangeness chemical potential.

Parameters
muSStrangeness chemical potential (GeV)

Reimplemented in thermalfist::ThermalModelCanonicalStrangeness.

Examples
BagModelFit.cpp, NeutronStars-CSHRG.cpp, SusceptibilitiesBQS.cpp, ThermodynamicsBQS.cpp, cpc1-HRG-TDep.cpp, cpc2-chi2-vs-T.cpp, and cpc4-mcHRG.cpp.

Definition at line 163 of file ThermalModelBase.cpp.

◆ SetTemperature()

void thermalfist::ThermalModelBase::SetTemperature ( double T)
virtual

Set the temperature.

Parameters
TTemperature (GeV)
Examples
CalculationTmu.cpp, NeutronStars-CSHRG.cpp, SusceptibilitiesBQS.cpp, ThermodynamicsBQS.cpp, cpc1-HRG-TDep.cpp, and cpc4-mcHRG.cpp.

Definition at line 143 of file ThermalModelBase.cpp.

◆ SetUseWidth() [1/2]

void thermalfist::ThermalModelBase::SetUseWidth ( bool useWidth)

Sets whether finite resonance widths are used. Deprecated.

If the widths are to be used, the ThermalParticle::ResonanceWidthIntegration::BWTwoGamma scheme will be used.

Parameters
useWidthWhether to use the finite resonance widths.
Examples
BagModelFit.cpp, CalculationTmu.cpp, CosmicTrajectory.cpp, NeutronStars-CSHRG.cpp, SusceptibilitiesBQS.cpp, ThermodynamicsBQS.cpp, cpc1-HRG-TDep.cpp, cpc2-chi2-vs-T.cpp, cpc3-chi2NEQ.cpp, and cpc4-mcHRG.cpp.

Definition at line 98 of file ThermalModelBase.cpp.

◆ SetUseWidth() [2/2]

void thermalfist::ThermalModelBase::SetUseWidth ( ThermalParticle::ResonanceWidthIntegration type)

Sets the finite resonance widths scheme to use.

Parameters
typeThermalParticle::ResonanceWidthIntegration scheme

Definition at line 110 of file ThermalModelBase.cpp.

◆ SetVirial()

virtual void thermalfist::ThermalModelBase::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 in thermalfist::ThermalModelEVCrosstermsLegacy, thermalfist::ThermalModelEVDiagonal, thermalfist::ThermalModelVDW, and thermalfist::ThermalModelVDWCanonicalStrangeness.

Examples
CalculationTmu.cpp, cpc1-HRG-TDep.cpp, and cpc2-chi2-vs-T.cpp.

Definition at line 239 of file ThermalModelBase.h.

◆ SetVolume()

void thermalfist::ThermalModelBase::SetVolume ( double Volume)
inline

Sets the system volume.

Parameters
VolumeSystem volume (fm \(^3\))

Definition at line 515 of file ThermalModelBase.h.

◆ SetVolumeRadius()

void thermalfist::ThermalModelBase::SetVolumeRadius ( double R)
inline

Sets the system radius.

The system volume is computed as \(V = \frac{4\pi}{3} \, R^3\)

Parameters
VolumeSystem radius (fm)
Examples
CalculationTmu.cpp.

Definition at line 527 of file ThermalModelBase.h.

◆ SolveChemicalPotentials()

bool thermalfist::ThermalModelBase::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 )
virtual

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.

Uses the Broyden's method to determine the chemical potentials. The resulting chemical potentials are stored in Parameters(). Useful for canonical ensemble applications.

Parameters
totBThe total desired net baryon charge of the system.
totQThe total desired electric charge of the system.
totSThe total desired net strangeness of the system.
totCThe total desired net charm of the system.
muBinitThe initial guess for the baryon chemical potential.
muQinitThe initial guess for the electric chemical potential.
muSinitThe initial guess for the strangeness chemical potential.
muCinitThe initial guess for the charm chemical potential.
ConstrMuBWhether the baryon chemical potential should be constrained.
ConstrMuQWhether the electric chemical potential should be constrained.
ConstrMuSWhether the strangeness chemical potential should be constrained.
ConstrMuCWhether the charm chemical potential should be constrained.
Returns
true is chemical potentials were contrained successfully, false otherwise

Definition at line 539 of file ThermalModelBase.cpp.

◆ SoverB()

double thermalfist::ThermalModelBase::SoverB ( ) const
inline

Definition at line 495 of file ThermalModelBase.h.

◆ SusceptibilityDimensionfull()

double thermalfist::ThermalModelBase::SusceptibilityDimensionfull ( ConservedCharge::Name i,
ConservedCharge::Name j ) const
virtual

A 2nd order susceptibility of conserved charges.

\( \chi_{11}^{c_i c_j} \)

Parameters
iFirst conserved charge
jSecond conserved charge
Returns
Susceptibility \( d P^2 / d \mu_i \mu_j \) (GeV \(^{-2}\))
Examples
NeutronStars-CSHRG.cpp.

Definition at line 2120 of file ThermalModelBase.cpp.

◆ TwoParticleSusceptibilityFinal()

double thermalfist::ThermalModelBase::TwoParticleSusceptibilityFinal ( int i,
int j ) const
virtual

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.

Definition at line 1153 of file ThermalModelBase.cpp.

◆ TwoParticleSusceptibilityFinalByPdg()

double thermalfist::ThermalModelBase::TwoParticleSusceptibilityFinalByPdg ( long long id1,
long long id2 )
virtual

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.

Definition at line 1169 of file ThermalModelBase.cpp.

◆ TwoParticleSusceptibilityPrimordial()

double thermalfist::ThermalModelBase::TwoParticleSusceptibilityPrimordial ( int i,
int j ) const
virtual

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.

Definition at line 1066 of file ThermalModelBase.cpp.

◆ TwoParticleSusceptibilityPrimordialByPdg()

double thermalfist::ThermalModelBase::TwoParticleSusceptibilityPrimordialByPdg ( long long id1,
long long id2 )
virtual

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.

Definition at line 1075 of file ThermalModelBase.cpp.

◆ TwoParticleSusceptibilityTemperatureDerivativePrimordial()

double thermalfist::ThermalModelBase::TwoParticleSusceptibilityTemperatureDerivativePrimordial ( int i,
int j ) const
virtual

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.

Definition at line 1092 of file ThermalModelBase.cpp.

◆ TwoParticleSusceptibilityTemperatureDerivativePrimordialByPdg()

double thermalfist::ThermalModelBase::TwoParticleSusceptibilityTemperatureDerivativePrimordialByPdg ( long long id1,
long long id2 )
virtual

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.

Definition at line 1100 of file ThermalModelBase.cpp.

◆ UpdateParameters()

void thermalfist::ThermalModelBase::UpdateParameters ( )
inline

Calls SetParameters() with current m_Parameters.

Definition at line 128 of file ThermalModelBase.h.

◆ UsePartialChemicalEquilibrium() [1/2]

bool thermalfist::ThermalModelBase::UsePartialChemicalEquilibrium ( )
inline

Whether partial chemical equilibrium with additional chemical potentials is used.

Definition at line 484 of file ThermalModelBase.h.

◆ UsePartialChemicalEquilibrium() [2/2]

void thermalfist::ThermalModelBase::UsePartialChemicalEquilibrium ( bool usePCE)
inline

Sets whether partial chemical equilibrium with additional chemical potentials is used.

Definition at line 481 of file ThermalModelBase.h.

◆ UseWidth()

bool thermalfist::ThermalModelBase::UseWidth ( ) const
inline

Whether finite resonance widths are considered.

Definition at line 82 of file ThermalModelBase.h.

◆ ValidateCalculation()

void thermalfist::ThermalModelBase::ValidateCalculation ( )
virtual

Checks whether issues have occured during the calculation of particle densities in the CalculateDensities() method.

Reimplemented in thermalfist::ThermalModelCanonical.

Definition at line 633 of file ThermalModelBase.cpp.

◆ ValidityCheckLog()

std::string thermalfist::ThermalModelBase::ValidityCheckLog ( ) const
inline

All messaged which occured during the validation procedure in the ValidateCalculation() method.

Returns
std::string List of messages, one per line

Definition at line 676 of file ThermalModelBase.h.

◆ VirialCoefficient()

virtual double thermalfist::ThermalModelBase::VirialCoefficient ( int ,
int  ) const
inlinevirtual

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 in thermalfist::ThermalModelEVCrosstermsLegacy, thermalfist::ThermalModelEVDiagonal, thermalfist::ThermalModelVDW, and thermalfist::ThermalModelVDWCanonicalStrangeness.

Definition at line 333 of file ThermalModelBase.h.

◆ Volume()

double thermalfist::ThermalModelBase::Volume ( ) const
inline

System volume (fm \(^3\))

Examples
CalculationTmu.cpp, and PCE-Saha-LHC.cpp.

Definition at line 517 of file ThermalModelBase.h.

◆ WriteInteractionParameters()

virtual void thermalfist::ThermalModelBase::WriteInteractionParameters ( const std::string & )
inlinevirtual

Write the QvdW interaction parameters to a file.

Actual implementation is in a derived class.

Parameters
filenameOutput file.

Reimplemented in thermalfist::ThermalModelEVCanonicalStrangeness, thermalfist::ThermalModelEVCrossterms, thermalfist::ThermalModelEVCrosstermsLegacy, thermalfist::ThermalModelEVDiagonal, thermalfist::ThermalModelVDW, and thermalfist::ThermalModelVDWCanonicalStrangeness.

Definition at line 313 of file ThermalModelBase.h.


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