![]() |
Thermal-FIST 1.5
Package for hadron resonance gas model applications
|
Class implementing the diagonal excluded-volume model. More...
#include <ThermalModelEVDiagonal.h>
Public Member Functions | |
ThermalModelEVDiagonal (ThermalParticleSystem *TPS, const ThermalModelParameters ¶ms=ThermalModelParameters()) | |
Construct a new ThermalModelEVDiagonal object. | |
virtual | ~ThermalModelEVDiagonal (void) |
Destroy the ThermalModelEVDiagonal object. | |
void | FillVirialEV (const std::vector< double > &vi=std::vector< double >(0)) |
Same as FillVirial() but uses the diagonal excluded-volume coefficients \( v_i \equiv b_{ii} \) as input instead of radii. | |
double | ExcludedVolume (int i) const |
The excluded volume parameter of particle species i. | |
double | CommonSuppressionFactor () |
The density suppression factor \( (1 - \sum_i v_i n_i) \), common for all species. | |
virtual void | DisableMesonMesonVirial () |
virtual void | DisableMesonMesonAttraction () |
virtual void | DisableMesonBaryonVirial () |
virtual void | DisableMesonBaryonAttraction () |
virtual void | DisableBaryonBaryonVirial () |
virtual void | DisableBaryonBaryonAttraction () |
virtual void | DisableBaryonAntiBaryonVirial () |
virtual void | DisableBaryonAntiBaryonAttraction () |
void | SetRadius (double rad) |
Set the same excluded volume radius parameter for all species. | |
void | SetRadius (int i, double rad) |
Set the radius parameter for particle species i. | |
void | FillVirial (const std::vector< double > &ri=std::vector< double >(0)) |
Fills the vector of particle eigenvolume parameters. | |
virtual void | ReadInteractionParameters (const std::string &filename) |
Read the set of eigenvolume parameters for all particles from an external file. | |
virtual void | WriteInteractionParameters (const std::string &filename) |
Write the set of eigenvolume parameters for all particles to an external file. | |
virtual double | CalculateEigenvolumeFraction () |
double | VirialCoefficient (int i, int j) const |
Return the eigenvolume parameter of particle species i. | |
void | SetVirial (int i, int j, double b) |
Sets the eigenvolume parameter of particle species i. | |
virtual void | ChangeTPS (ThermalParticleSystem *TPS) |
Change the particle list. | |
virtual void | CalculatePrimordialDensities () |
Calculates the primordial densities of all species. | |
virtual void | CalculateFluctuations () |
Computes the fluctuation observables. | |
virtual void | CalculateTwoParticleCorrelations () |
Computes the fluctuations and correlations of the primordial particle numbers. | |
virtual std::vector< double > | CalculateChargeFluctuations (const std::vector< double > &chgs, int order=4) |
Calculates fluctuations (diagonal susceptibilities) of an arbitrary "conserved" charge. | |
virtual double | CalculatePressure () |
virtual double | CalculateEnergyDensity () |
virtual double | CalculateEntropyDensity () |
virtual double | CalculateEnergyDensityDerivativeT () |
virtual void | CalculateTemperatureDerivatives () |
Computes the temperature derivatives of densities, shifted chemical potentials, and primordial hadron number susceptibilities. | |
virtual double | CalculateBaryonMatterEntropyDensity () |
virtual double | CalculateMesonMatterEntropyDensity () |
virtual double | ParticleScaledVariance (int) |
virtual double | ParticleSkewness (int) |
virtual double | ParticleKurtosis (int) |
virtual double | ParticleScalarDensity (int part) |
![]() | |
ThermalModelBase (ThermalParticleSystem *TPS, const ThermalModelParameters ¶ms=ThermalModelParameters()) | |
Construct a new ThermalModelBase object. | |
virtual | ~ThermalModelBase (void) |
int | ComponentsNumber () const |
Number of different particle species in the list. | |
bool | UseWidth () const |
Whether finite resonance widths are considered. | |
void | SetUseWidth (bool useWidth) |
Sets whether finite resonance widths are used. Deprecated. | |
void | SetUseWidth (ThermalParticle::ResonanceWidthIntegration type) |
Sets the finite resonance widths scheme to use. | |
void | SetNormBratio (bool normBratio) |
Whether branching ratios are renormalized to 100%. | |
bool | NormBratio () const |
void | SetOMP (bool openMP) |
OpenMP support. Currently not used. | |
virtual void | SetParameters (const ThermalModelParameters ¶ms) |
The thermal parameters. | |
const ThermalModelParameters & | Parameters () 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 | 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} \). | |
void | DisableMesonMesonRepulsion () |
void | DisableMesonBaryonRepulsion () |
void | DisableBaryonBaryonRepulsion () |
void | DisableBaryonAntiBaryonRepulsion () |
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 | CalculateDensities () |
Calculates the primordial and total (after decays) densities of all species. | |
virtual void | ValidateCalculation () |
Checks whether issues have occured during the calculation of particle densities in the CalculateDensities() method. | |
std::string | ValidityCheckLog () const |
All messaged which occured during the validation procedure in the ValidateCalculation() method. | |
virtual void | CalculateDensitiesGCE () |
Calculates the particle densities in a grand-canonical ensemble. | |
virtual void | CalculateFeeddown () |
Calculates the total densities which include feeddown contributions. | |
virtual void | CalculateTwoParticleFluctuationsDecays () |
Computes particle number correlations and fluctuations for all final-state particles which are marked stable. | |
virtual double | TwoParticleSusceptibilityPrimordial (int i, int j) const |
Returns the computed primordial particle number (cross-)susceptibility \( \frac{1}{VT^3} \, \langle \Delta N_i \Delta N_j \rangle \) for particles with ids i and j. CalculateFluctuations() must be called beforehand. | |
virtual double | TwoParticleSusceptibilityPrimordialByPdg (long long id1, long long id2) |
Returns the computed primordial particle number (cross-)susceptibility \( \frac{1}{VT^3} \, \langle \Delta N_i \Delta N_j \rangle \) for particles with pdg codes id1 and id2. CalculateFluctuations() must be called beforehand. | |
virtual double | TwoParticleSusceptibilityTemperatureDerivativePrimordial (int i, int j) const |
Returns the computed temperature derivative of the primordial particle number (cross-)susceptibility \( \frac{\partial}{\partial T} \frac{1}{VT^3} \, \langle \Delta N_i \Delta N_j \rangle \) for particles with ids i and j. CalculateFluctuations() must be called beforehand. | |
virtual double | TwoParticleSusceptibilityTemperatureDerivativePrimordialByPdg (long long id1, long long id2) |
Returns the computed temperature derivative of the primordial particle number (cross-)susceptibility \( \frac{\partial}{\partial T} \frac{1}{VT^3} \, \langle \Delta N_i \Delta N_j \rangle \) for particles with pdg codes id1 and id2. CalculateFluctuations() must be called beforehand. | |
virtual double | NetParticleSusceptibilityPrimordialByPdg (long long id1, long long id2) |
Returns the computed (cross-)susceptibility \( \frac{1}{VT^3} \, \langle \Delta N_i \Delta N_j \rangle \) between primordial net-particle numbers for pdg codes id1 and id2. CalculateFluctuations() must be called beforehand. | |
virtual double | TwoParticleSusceptibilityFinal (int i, int j) const |
Returns the computed final particle number (cross-)susceptibility \( \frac{1}{VT^3} \, \langle \Delta N_i \Delta N_j \rangle \) for particles with ids i and j. CalculateFluctuations() must be called beforehand. Both particle species must be those marked stable. | |
virtual double | TwoParticleSusceptibilityFinalByPdg (long long id1, long long id2) |
Returns the computed final particle number (cross-)susceptibility \( \frac{1}{VT^3} \, \langle \Delta N_i \Delta N_j \rangle \) for particles with pdg codes id1 and id2. CalculateFluctuations() must be called beforehand. Both particle species must be those marked stable. | |
virtual double | NetParticleSusceptibilityFinalByPdg (long long id1, long long id2) |
Returns the computed (cross-)susceptibility \( \frac{1}{VT^3} \, \langle \Delta N_i \Delta N_j \rangle \) between final net-particle numbers for pdg codes id1 and id2. CalculateFluctuations() must be called beforehand. | |
virtual double | PrimordialParticleChargeSusceptibility (int i, ConservedCharge::Name chg) const |
Returns the computed primordial particle vs conserved charge (cross-)susceptibility \( \frac{1}{VT^3} \, \langle \Delta N_i \Delta N_chg \rangle \) for particle with id i and conserved charge chg. CalculateFluctuations() must be called beforehand. | |
virtual double | PrimordialParticleChargeSusceptibilityByPdg (long long id1, ConservedCharge::Name chg) |
Returns the computed primordial particle vs conserved charge (cross-)susceptibility \( \frac{1}{VT^3} \, \langle \Delta N_i \Delta N_chg \rangle \) for particle with pdg code id1 and conserved charge chg. CalculateFluctuations() must be called beforehand. | |
virtual double | PrimordialNetParticleChargeSusceptibilityByPdg (long long id1, ConservedCharge::Name chg) |
Returns the computed primordial net particle vs conserved charge (cross-)susceptibility \( \frac{1}{VT^3} \, \langle \Delta N_i \Delta N_chg \rangle \) for particle with pdg code id1 and conserved charge chg. CalculateFluctuations() must be called beforehand. | |
virtual double | FinalParticleChargeSusceptibility (int i, ConservedCharge::Name chg) const |
Returns the computed final (after decays) particle vs conserved charge (cross-)susceptibility \( \frac{1}{VT^3} \, \langle \Delta N_i \Delta N_chg \rangle \) for particle with id i and conserved charge chg. CalculateFluctuations() must be called beforehand. | |
virtual double | FinalParticleChargeSusceptibilityByPdg (long long id1, ConservedCharge::Name chg) |
Returns the computed final (after decays) particle vs conserved charge (cross-)susceptibility \( \frac{1}{VT^3} \, \langle \Delta N_i \Delta N_chg \rangle \) for particle with pdg code id1 and conserved charge chg. CalculateFluctuations() must be called beforehand. | |
virtual double | FinalNetParticleChargeSusceptibilityByPdg (long long id1, ConservedCharge::Name chg) |
Returns the computed final (after decays) net particle vs conserved charge (cross-)susceptibility \( \frac{1}{VT^3} \, \langle \Delta N_i \Delta N_chg \rangle \) for particle with pdg code id1 and conserved charge chg. CalculateFluctuations() must be called beforehand. | |
virtual double | SusceptibilityDimensionfull (ConservedCharge::Name i, ConservedCharge::Name j) const |
A 2nd order susceptibility of conserved charges. | |
virtual void | CalculateSusceptibilityMatrix () |
Calculates the conserved charges susceptibility matrix. | |
virtual void | CalculateProxySusceptibilityMatrix () |
Calculates the susceptibility matrix of conserved charges proxies. | |
virtual void | CalculateParticleChargeCorrelationMatrix () |
Calculates the matrix of correlators between primordial (and also final) particle numbers and conserved charges. | |
Protected Member Functions | |
virtual void | SolvePressure () |
Solves the transcdental equation for the pressure. | |
double | Pressure (double P) |
Computes the l.h.s. of the transcdental equation for the pressure. | |
double | DensityId (int i, double Pressure) |
Calculate the ideal gas density of particle species i for the given pressure value. | |
double | PressureId (int i, double Pressure) |
Calculate the ideal gas pressure of particle species i for the given pressure value. | |
double | ScaledVarianceId (int i, double Pressure) |
Calculate the ideal gas scaled variance of particle species i number fluctuations for the given pressure value. | |
virtual double | MuShift (int i) const |
The shift in the chemical potential of particle species i due to the excluded volume interactions. | |
Protected Attributes | |
std::vector< double > | m_densitiesid |
std::vector< double > | m_densitiesidnoshift |
std::vector< double > | m_v |
double | m_Suppression |
double | m_Pressure |
double | m_Densityid |
double | m_TotalDensity |
EVSolution | m_sol |
Additional Inherited Members | |
![]() | |
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... | |
Class implementing the diagonal excluded-volume model.
The model formulation can be found in G.D. Yen, M.I. Gorenstein, W. Greiner, S.-N. Yang, Phys. Rev. C 56 2210, (1997), http://arxiv.org/pdf/nucl-th/9711062.pdf
The transcendental equation for the pressure is solved using the Broyden's method.
Definition at line 45 of file ThermalModelEVDiagonal.h.
thermalfist::ThermalModelEVDiagonal::ThermalModelEVDiagonal | ( | ThermalParticleSystem * | TPS, |
const ThermalModelParameters & | params = ThermalModelParameters() ) |
Construct a new ThermalModelEVDiagonal object.
TPS | A pointer to the ThermalParticleSystem object containing the particle list |
params | ThermalModelParameters object with current thermal parameters |
Definition at line 32 of file ThermalModelEVDiagonal.cpp.
|
virtual |
Destroy the ThermalModelEVDiagonal object.
Definition at line 45 of file ThermalModelEVDiagonal.cpp.
|
inlinevirtual |
Definition at line 181 of file ThermalModelEVDiagonal.h.
|
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.
chgs | A 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() |
order | Up to which order the susceptibilities are computed |
Reimplemented from thermalfist::ThermalModelBase.
Definition at line 409 of file ThermalModelEVDiagonal.cpp.
|
virtual |
Definition at line 652 of file ThermalModelEVDiagonal.cpp.
|
virtual |
Definition at line 602 of file ThermalModelEVDiagonal.cpp.
|
virtual |
Definition at line 284 of file ThermalModelEVDiagonal.cpp.
|
virtual |
Definition at line 613 of file ThermalModelEVDiagonal.cpp.
|
virtual |
Computes the fluctuation observables.
Includes the matrix of 2nd order susceptibilities of conserved charges, as well as particle number correlations and fluctuations, both for primordial yields and after decays.
Reimplemented from thermalfist::ThermalModelBase.
Definition at line 383 of file ThermalModelEVDiagonal.cpp.
|
inlinevirtual |
Definition at line 183 of file ThermalModelEVDiagonal.h.
|
virtual |
Definition at line 624 of file ThermalModelEVDiagonal.cpp.
|
virtual |
Calculates the primordial densities of all species.
Implements thermalfist::ThermalModelBase.
Definition at line 204 of file ThermalModelEVDiagonal.cpp.
|
virtual |
Computes the temperature derivatives of densities, shifted chemical potentials, and primordial hadron number susceptibilities.
Reimplemented from thermalfist::ThermalModelBase.
Definition at line 306 of file ThermalModelEVDiagonal.cpp.
|
virtual |
Computes the fluctuations and correlations of the primordial particle numbers.
More specifically, computes the susceptibility matrix \( \frac{1}{VT^3} \, \langle \Delta N_i^* \Delta N_j^* \rangle \), where \( N_i^* \) is the primordial yield of species i.
Reimplemented from thermalfist::ThermalModelBase.
Definition at line 244 of file ThermalModelEVDiagonal.cpp.
|
virtual |
Change the particle list.
TPS | A pointer to new particle list. |
Reimplemented from thermalfist::ThermalModelBase.
Definition at line 130 of file ThermalModelEVDiagonal.cpp.
double thermalfist::ThermalModelEVDiagonal::CommonSuppressionFactor | ( | ) |
The density suppression factor \( (1 - \sum_i v_i n_i) \), common for all species.
Definition at line 637 of file ThermalModelEVDiagonal.cpp.
|
protected |
Calculate the ideal gas density of particle species i for the given pressure value.
i | 0-based particle specie index |
Pressure | Input pressure (GeV fm \(^{-3}\)) |
Definition at line 136 of file ThermalModelEVDiagonal.cpp.
|
inlinevirtual |
Switches off QvdW attraction terms for all baryon-antibaryon pairs i.e. \( a_{ij} = 0 \) for baryon-antibaryon pairs
Reimplemented from thermalfist::ThermalModelBase.
Definition at line 102 of file ThermalModelEVDiagonal.h.
|
inlinevirtual |
Switches off eigenvolume terms for all baryon-antibaryon pairs i.e. \( \tilde{b}_{ij} = 0 \) for baryon-antibaryon pairs
Reimplemented from thermalfist::ThermalModelBase.
Definition at line 101 of file ThermalModelEVDiagonal.h.
|
inlinevirtual |
Switches off QvdW attraction terms for all baryon-baryon pairs i.e. \( a_{ij} = 0 \) for baryon-baryon pairs
Reimplemented from thermalfist::ThermalModelBase.
Definition at line 100 of file ThermalModelEVDiagonal.h.
|
inlinevirtual |
Switches off excluded volume terms for all baryon-baryon pairs i.e. \( \tilde{b}_{ij} = 0 \) for baryon-baryon pairs
Reimplemented from thermalfist::ThermalModelBase.
Definition at line 99 of file ThermalModelEVDiagonal.h.
|
inlinevirtual |
Switches off QvdW attraction terms for all meson-baryon pairs i.e. \( a_{ij} = 0 \) for meson-baryon pairs
Reimplemented from thermalfist::ThermalModelBase.
Definition at line 98 of file ThermalModelEVDiagonal.h.
|
inlinevirtual |
Switches off excluded volume terms for all meson-baryon pairs i.e. \( \tilde{b}_{ij} = 0 \) for meson-baryon pairs
Reimplemented from thermalfist::ThermalModelBase.
Definition at line 97 of file ThermalModelEVDiagonal.h.
|
inlinevirtual |
Switches off QvdW attraction terms for all meson-meson pairs i.e. \( a_{ij} = 0 \) for meson-meson pairs
Reimplemented from thermalfist::ThermalModelBase.
Definition at line 96 of file ThermalModelEVDiagonal.h.
|
inlinevirtual |
Differential treatment of pair interactions is not applicable for the diagonal excluded volume model
Reimplemented from thermalfist::ThermalModelBase.
Definition at line 95 of file ThermalModelEVDiagonal.h.
double thermalfist::ThermalModelEVDiagonal::ExcludedVolume | ( | int | i | ) | const |
The excluded volume parameter of particle species i.
i | 0-based particle species index. |
Definition at line 122 of file ThermalModelEVDiagonal.cpp.
|
virtual |
Fills the vector of particle eigenvolume parameters.
ri | A vector of radii parameters for all particles (in fm) |
Reimplemented from thermalfist::ThermalModelBase.
Definition at line 58 of file ThermalModelEVDiagonal.cpp.
void thermalfist::ThermalModelEVDiagonal::FillVirialEV | ( | const std::vector< double > & | vi = std::vector<double>(0) | ) |
Same as FillVirial() but uses the diagonal excluded-volume coefficients \( v_i \equiv b_{ii} \) as input instead of radii.
vi | A vector with diagonal excluded-volume coefficients for all species. 0-based indices of the vector must corresponds to the 0-based indices of the particle list TPS() |
Definition at line 69 of file ThermalModelEVDiagonal.cpp.
|
protectedvirtual |
The shift in the chemical potential of particle species i due to the excluded volume interactions.
Equal to \(-v_i p\)
i | 0-based particle specie index |
Definition at line 644 of file ThermalModelEVDiagonal.cpp.
|
inlinevirtual |
Definition at line 192 of file ThermalModelEVDiagonal.h.
|
virtual |
Definition at line 629 of file ThermalModelEVDiagonal.cpp.
|
inlinevirtual |
Definition at line 186 of file ThermalModelEVDiagonal.h.
|
inlinevirtual |
Definition at line 189 of file ThermalModelEVDiagonal.h.
|
protected |
Computes the l.h.s. of the transcdental equation for the pressure.
P | Input pressure (GeV fm \(^{-3}\)) |
Definition at line 154 of file ThermalModelEVDiagonal.cpp.
|
protected |
Calculate the ideal gas pressure of particle species i for the given pressure value.
i | 0-based particle specie index |
Pressure | Input pressure (GeV fm \(^{-3}\)) |
Definition at line 142 of file ThermalModelEVDiagonal.cpp.
|
virtual |
Read the set of eigenvolume parameters for all particles from an external file.
filename | Input file name |
Reimplemented from thermalfist::ThermalModelBase.
Definition at line 78 of file ThermalModelEVDiagonal.cpp.
|
protected |
Calculate the ideal gas scaled variance of particle species i number fluctuations for the given pressure value.
i | 0-based particle specie index |
Pressure | Input pressure (GeV fm \(^{-3}\)) |
Definition at line 148 of file ThermalModelEVDiagonal.cpp.
|
virtual |
Set the same excluded volume radius parameter for all species.
rad | Radius parameter (fm) |
Reimplemented from thermalfist::ThermalModelBase.
Definition at line 50 of file ThermalModelEVDiagonal.cpp.
|
virtual |
Set the radius parameter for particle species i.
i | 0-based index of particle species |
rad | Radius parameter (fm) |
Reimplemented from thermalfist::ThermalModelBase.
Definition at line 660 of file ThermalModelEVDiagonal.cpp.
|
virtual |
Sets the eigenvolume parameter of particle species i.
i | 0-based particle specie index |
j | Here must be equal to i |
b | The eigenvolume parameter (fm \(^{-3}\)) |
Reimplemented from thermalfist::ThermalModelBase.
Definition at line 671 of file ThermalModelEVDiagonal.cpp.
|
protectedvirtual |
Solves the transcdental equation for the pressure.
Solves \( p(T,\mu) = \sum_i p_i^{\rm id} (T, \mu_i - v_i p). \)
Definition at line 165 of file ThermalModelEVDiagonal.cpp.
|
virtual |
Return the eigenvolume parameter of particle species i.
i | 0-based particle specie index |
j | Irrelevant |
Reimplemented from thermalfist::ThermalModelBase.
Definition at line 666 of file ThermalModelEVDiagonal.cpp.
|
virtual |
Write the set of eigenvolume parameters for all particles to an external file.
One particle specie per line.
filename | Output file name |
Reimplemented from thermalfist::ThermalModelBase.
Definition at line 102 of file ThermalModelEVDiagonal.cpp.
|
protected |
Vector of ideal gas densities with shifted chemical potentials
Definition at line 258 of file ThermalModelEVDiagonal.h.
|
protected |
Vector of ideal gas densities without shifted chemical potentials
Definition at line 259 of file ThermalModelEVDiagonal.h.
|
protected |
Definition at line 263 of file ThermalModelEVDiagonal.h.
|
protected |
The solved pressure
Definition at line 262 of file ThermalModelEVDiagonal.h.
|
protected |
Axuiliary
Definition at line 265 of file ThermalModelEVDiagonal.h.
|
protected |
The common density suppression factor
Definition at line 261 of file ThermalModelEVDiagonal.h.
|
protected |
Definition at line 264 of file ThermalModelEVDiagonal.h.
|
protected |
Vector of eigenvolumes of all hadrons
Definition at line 260 of file ThermalModelEVDiagonal.h.