![]() |
Thermal-FIST
1.3
Package for hadron resonance gas model applications
|
Class implementing the crossterms excluded-volume model. More...
#include <ThermalModelEVCrossterms.h>
Public Member Functions | |
ThermalModelEVCrossterms (ThermalParticleSystem *TPS, const ThermalModelParameters ¶ms=ThermalModelParameters()) | |
Construct a new ThermalModelEVCrossterms object. More... | |
virtual | ~ThermalModelEVCrossterms (void) |
Destroy the ThermalModelEVCrossterms object. More... | |
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. More... | |
virtual void | ReadInteractionParameters (const std::string &filename) |
Reads the QvdW interaction parameters from a file. More... | |
virtual void | WriteInteractionParameters (const std::string &filename) |
Write the QvdW interaction parameters to a file. More... | |
void | SetRadius (double rad) |
Set the same excluded volume radius parameter for all species. More... | |
double | VirialCoefficient (int i, int j) const |
Excluded volume coefficient \( \tilde{b}_{ij} = 0 \). More... | |
void | SetVirial (int i, int j, double b) |
Set the excluded volume coefficient \( \tilde{b}_{ij} \). More... | |
virtual void | ChangeTPS (ThermalParticleSystem *TPS) |
Change the particle list. More... | |
virtual void | CalculatePrimordialDensities () |
Calculates the primordial densities of all species. More... | |
void | CalculateTwoParticleCorrelations () |
Computes the fluctuations and correlations of the primordial particle numbers. More... | |
void | CalculateFluctuations () |
Computes the fluctuation observables. More... | |
virtual std::vector< double > | CalculateChargeFluctuations (const std::vector< double > &chgs, int order=4) |
Calculates fluctuations (diagonal susceptibilities) of an arbitrary "conserved" charge. More... | |
virtual double | CalculatePressure () |
Implementation of the equation of state functions. More... | |
virtual double | CalculateEnergyDensity () |
virtual double | CalculateEntropyDensity () |
virtual double | CalculateBaryonMatterEntropyDensity () |
The fraction of entropy carried by baryons (Ideal GCE only) More... | |
virtual double | CalculateMesonMatterEntropyDensity () |
The fraction of entropy carried by mesons (Ideal GCE only) More... | |
virtual double | ParticleScaledVariance (int) |
Scaled variance of primordial particle number fluctuations for species i. More... | |
virtual double | ParticleSkewness (int part) |
Skewness of primordial particle number fluctuations for species i. More... | |
virtual double | ParticleKurtosis (int part) |
Kurtosis of primordial particle number fluctuations for species i. More... | |
virtual double | ParticleScalarDensity (int) |
The scalar density of the particle species i. More... | |
![]() | |
ThermalModelBase (ThermalParticleSystem *TPS, const ThermalModelParameters ¶ms=ThermalModelParameters()) | |
Construct a new ThermalModelBase object. More... | |
virtual | ~ThermalModelBase (void) |
int | ComponentsNumber () const |
Number of different particle species in the list. More... | |
bool | UseWidth () const |
Whether finite resonance widths are considered. More... | |
void | SetUseWidth (bool useWidth) |
Sets whether finite resonance widths are used. Deprecated. More... | |
void | SetUseWidth (ThermalParticle::ResonanceWidthIntegration type) |
Sets the finite resonance widths scheme to use. More... | |
void | SetOMP (bool openMP) |
OpenMP support. Currently not used. More... | |
void | UpdateParameters () |
Calls SetParameters() with current m_Parameters. More... | |
virtual void | SetTemperature (double T) |
Set the temperature. More... | |
virtual void | SetBaryonChemicalPotential (double muB) |
Set the baryon chemical potential. More... | |
virtual void | SetElectricChemicalPotential (double muQ) |
Set the electric chemical potential. More... | |
virtual void | SetStrangenessChemicalPotential (double muS) |
Set the strangeness chemical potential. More... | |
virtual void | SetCharmChemicalPotential (double muC) |
Set the charm chemical potential. More... | |
virtual void | SetGammaq (double gammaq) |
Set the light quark fugacity factor. More... | |
virtual void | SetGammaS (double gammaS) |
Set the strange quark fugacity factor. More... | |
virtual void | SetGammaC (double gammaC) |
Set the charm quark fugacity factor. More... | |
virtual void | SetBaryonCharge (int B) |
Set the total baryon number (for canonical ensemble only) More... | |
virtual void | SetElectricCharge (int Q) |
Set the total electric charge (for canonical ensemble only) More... | |
virtual void | SetStrangeness (int S) |
Set the total strangeness (for canonical ensemble only) More... | |
virtual void | SetCharm (int C) |
Set the total charm (for canonical ensemble only) More... | |
virtual void | SetRadius (int, double) |
Set the radius parameter for particle species i. More... | |
virtual void | SetAttraction (int, int, double) |
Set the vdW mean field attraction coefficient \( a_{ij} \). More... | |
virtual void | DisableMesonMesonAttraction () |
virtual void | DisableMesonBaryonAttraction () |
virtual void | DisableBaryonBaryonAttraction () |
virtual void | DisableBaryonAntiBaryonAttraction () |
virtual double | AttractionCoefficient (int, int) const |
QvdW mean field attraction coefficient \( a_{ij} \). More... | |
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(). More... | |
virtual void | SetClusterExpansionOrder (int order) |
Set the number of terms in the cluster expansion method. Calls the corresponding method in TPS(). More... | |
void | SetResonanceWidthShape (ThermalParticle::ResonanceWidthShape shape) |
Set the ThermalParticle::ResonanceWidthShape for all particles. Calls the corresponding method in TPS(). More... | |
void | SetResonanceWidthIntegrationType (ThermalParticle::ResonanceWidthIntegration type) |
Set the ThermalParticle::ResonanceWidthIntegration scheme for all particles. Calls the corresponding method in TPS(). More... | |
virtual void | FillChemicalPotentials () |
Sets the chemical potentials of all particles. More... | |
virtual void | SetChemicalPotentials (const std::vector< double > &chem=std::vector< double >(0)) |
Sets the chemical potentials of all particles. More... | |
const std::vector< double > & | ChemicalPotentials () const |
A vector of chemical potentials of all particles. More... | |
double | ChemicalPotential (int i) const |
Chemical potential of particle species i. More... | |
virtual double | FullIdealChemicalPotential (int i) const |
Chemical potential entering the ideal gas expressions of particle species i. More... | |
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. More... | |
bool | UsePartialChemicalEquilibrium () |
Whether partial chemical equilibrium with additional chemical potentials is used. More... | |
void | SetVolume (double Volume) |
Sets the system volume. More... | |
double | Volume () const |
System volume (fm \(^3\)) More... | |
void | SetVolumeRadius (double R) |
Sets the system radius. More... | |
double | CanonicalVolume () const |
The canonical correlation volume V \(_c\) (fm \(^3\)) More... | |
void | SetCanonicalVolume (double Volume) |
Set the canonical correlation volume V \(_c\). More... | |
void | SetCanonicalVolumeRadius (double Rc) |
Set the canonical correlation system radius. More... | |
void | ConstrainChemicalPotentials (bool resetInitialValues=true) |
Constrains the chemical potentials \( \mu_B,\,\mu_Q,\,\mu_S,\,\mu_C \) by the conservation laws imposed. More... | |
virtual void | FixParameters () |
Method which actually implements ConstrainChemicalPotentials() (for backward compatibility). More... | |
virtual void | FixParametersNoReset () |
Method which actually implements ConstrainChemicalPotentialsNoReset() (for backward compatibility). More... | |
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_Q \) which reproduce the specified total baryon, electric, strangeness, and charm charges of the system. More... | |
virtual void | CalculateDensities () |
Calculates the primordial and total (after decays) densities of all species. More... | |
virtual void | ValidateCalculation () |
Checks whether issues have occured during the calculation of particle densities in the CalculateDensities() method. More... | |
std::string | ValidityCheckLog () const |
All messaged which occured during the validation procedure in the ValidateCalculation() method. More... | |
virtual void | CalculateDensitiesGCE () |
Calculates the particle densities in a grand-canonical ensemble. More... | |
virtual void | CalculateFeeddown () |
Calculates the total densities which include feeddown contributions. More... | |
virtual void | CalculateTwoParticleFluctuationsDecays () |
Computes particle number correlations and fluctuations for all final-state particles which are marked stable. More... | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
virtual void | CalculateSusceptibilityMatrix () |
Calculates the conserved charges susceptibility matrix. More... | |
virtual void | CalculateProxySusceptibilityMatrix () |
Calculates the susceptibility matrix of conserved charges proxies. More... | |
virtual void | CalculateParticleChargeCorrelationMatrix () |
Calculates the matrix of correlators between primordial (and also final) particle numbers and conserved charges. More... | |
double | Pressure () |
System pressure (GeV fm \(^{-3}\)) More... | |
double | EnergyDensity () |
Energy density (GeV fm \(^{-3}\)) More... | |
double | EntropyDensity () |
Entropy density (fm \(^{-3}\)) More... | |
double | HadronDensity () |
Total number density of all particle (fm \(^{-3}\)) More... | |
double | BaryonDensity () |
Net baryon density (fm \(^{-3}\)) More... | |
double | ElectricChargeDensity () |
Electric charge density (fm \(^{-3}\)) More... | |
double | StrangenessDensity () |
Net strangeness density (fm \(^{-3}\)) More... | |
double | CharmDensity () |
Net charm density (fm \(^{-3}\)) More... | |
double | AbsoluteBaryonDensity () |
Absolute baryon number density (baryons + antibaryons) (fm \(^{-3}\)) More... | |
double | AbsoluteElectricChargeDensity () |
Absolute electric charge density (Q+ + Q-) (fm \(^{-3}\)) More... | |
double | AbsoluteStrangenessDensity () |
Absolute strange quark content density (fm \(^{-3}\)) More... | |
double | AbsoluteCharmDensity () |
Absolute charm quark content density (fm \(^{-3}\)) More... | |
virtual double | CalculateArbitraryChargeDensity () |
Computes the density of the auxiliary ArbitraryCharge() More... | |
virtual double | CalculateEigenvolumeFraction () |
virtual double | GetMaxDiff () const |
virtual bool | IsLastSolutionOK () const |
double | GetDensity (long long PDGID, int feeddown) |
Same as GetDensity(int,Feeddown::Type) More... | |
double | GetDensity (long long PDGID, Feeddown::Type feeddown) |
Particle number density of species with a specified PDG ID and feeddown. More... | |
double | GetYield (long long PDGID, Feeddown::Type feeddown) |
Particle number yield of species with a specified PDG ID and feeddown. More... | |
std::vector< double > | GetIdealGasDensities () const |
ThermalParticleSystem * | TPS () |
const std::vector< double > & | Densities () const |
const std::vector< double > & | TotalDensities () const |
const std::vector< std::vector< double > > & | AllDensities () const |
void | SetTAG (const std::string &tag) |
Set the tag for this ThermalModelBase object. More... | |
const std::string & | TAG () const |
The tag of this ThermalModelBase object. More... | |
void | ResetCalculatedFlags () |
Reset all flags which correspond to a calculation status. More... | |
bool | IsCalculated () const |
bool | IsFluctuationsCalculated () const |
bool | IsGCECalculated () const |
double | ScaledVariancePrimordial (int id) const |
Scaled variance of primordial particle number fluctuations for particle species id. More... | |
double | ScaledVarianceTotal (int id) const |
Scaled variance of final particle number fluctuations for particle species id. More... | |
double | SkewnessPrimordial (int id) const |
Normalized skewness of primordial particle number fluctuations for particle species id. More... | |
double | SkewnessTotal (int id) const |
Normalized skewness of final particle number fluctuations for particle species id. More... | |
double | KurtosisPrimordial (int id) const |
Normalized excess kurtosis of primordial particle number fluctuations for particle species id. More... | |
double | KurtosisTotal (int id) const |
Normalized excess kurtosis of final particle number fluctuations for particle species id. More... | |
double | ConservedChargeDensity (ConservedCharge::Name chg) |
A density of a conserved charge (in fm^-3) More... | |
double | Susc (ConservedCharge::Name i, ConservedCharge::Name j) const |
A 2nd order susceptibility of conserved charges. More... | |
double | ProxySusc (ConservedCharge::Name i, ConservedCharge::Name j) const |
A 2nd order susceptibility of conserved charges proxies. More... | |
double | ChargedMultiplicity (int type=0) |
Multiplicity of charged particles. More... | |
double | ChargedScaledVariance (int type=0) |
Scaled variance of charged particles. More... | |
double | ChargedMultiplicityFinal (int type=0) |
Multiplicity of charged particles including the feeddown contributions in accordance with the stability flags. More... | |
double | ChargedScaledVarianceFinal (int type=0) |
Scaled variance of charged particles including the feeddown contributions in accordance with the stability flags. More... | |
ThermalModelEnsemble | Ensemble () |
The statistical ensemble of the current HRG model. More... | |
virtual bool | IsConservedChargeCanonical (ConservedCharge::Name charge) const |
Whether the given conserved charge is treated canonically. More... | |
ThermalModelInteraction | InteractionModel () |
The interactions present in the current HRG model. More... | |
void | SetNormBratio (bool normBratio) |
Whether branching ratios are renormalized to 100%. More... | |
bool | NormBratio () const |
virtual void | SetParameters (const ThermalModelParameters ¶ms) |
The thermal parameters. More... | |
const ThermalModelParameters & | Parameters () const |
virtual void | DisableMesonMesonVirial () |
void | DisableMesonMesonRepulsion () |
virtual void | DisableMesonBaryonVirial () |
void | DisableMesonBaryonRepulsion () |
virtual void | DisableBaryonBaryonVirial () |
void | DisableBaryonBaryonRepulsion () |
virtual void | DisableBaryonAntiBaryonVirial () |
void | DisableBaryonAntiBaryonRepulsion () |
double | RepulsionCoefficient (int i, int j) const |
void | SetSoverB (double SB) |
The entropy per baryon ratio to be used to constrain the baryon chemical potential. More... | |
double | SoverB () const |
void | SetQoverB (double QB) |
The electric-to-baryon charge ratio to be used to constrain the electric chemical potential. More... | |
double | QoverB () const |
virtual double | CalculateHadronDensity () |
virtual double | CalculateBaryonDensity () |
virtual double | CalculateChargeDensity () |
virtual double | CalculateStrangenessDensity () |
virtual double | CalculateCharmDensity () |
virtual double | CalculateAbsoluteBaryonDensity () |
virtual double | CalculateAbsoluteChargeDensity () |
virtual double | CalculateAbsoluteStrangenessDensity () |
virtual double | CalculateAbsoluteCharmDensity () |
Protected Member Functions | |
virtual void | SolvePressure (bool resetPartials=true) |
Solves the system of transcdental equations for the pressure using the Broyden's method. More... | |
void | SolvePressureIter () |
virtual void | CalculatePrimordialDensitiesNoReset () |
virtual void | CalculatePrimordialDensitiesIter () |
virtual void | SolveDiagonal () |
Solves the transcendental equation of the corresponding diagonal EV model. More... | |
double | PartialPressureDiagonal (int i, double P) |
The "partial pressure" of hadron species i for the given total pressure in the diagonal model. More... | |
double | PressureDiagonalTotal (double P) |
The total pressure for the given input pressure in the diagonal model. More... | |
virtual double | DensityId (int i, const std::vector< double > &pstars=std::vector< double >()) |
Calculate the ideal gas density of particle species i for the given values of partial pressures. More... | |
virtual double | Pressure (int i, const std::vector< double > &pstars=std::vector< double >()) |
Calculate the ideal gas pressure of particle species i for the given values of partial pressures. More... | |
double | ScaledVarianceId (int ind) |
Calculate the ideal gas scaled variance of particle species i number fluctuations for the given values of partial pressures. More... | |
virtual double | MuShift (int i) const |
The shift in the chemical potential of particle species i due to the excluded volume interactions. More... | |
Protected Attributes | |
std::vector< double > | m_densitiesid |
std::vector< double > | m_Ps |
std::vector< std::vector< double > > | m_Virial |
double | m_Pressure |
double | m_TotalEntropyDensity |
![]() | |
ThermalModelParameters | m_Parameters |
ThermalParticleSystem * | m_TPS |
bool | m_LastCalculationSuccessFlag |
double | m_MaxDiff |
bool | m_Calculated |
bool | m_FeeddownCalculated |
bool | m_FluctuationsCalculated |
bool | m_GCECalculated |
bool | m_UseWidth |
bool | m_NormBratio |
bool | m_QuantumStats |
double | m_QBgoal |
double | m_SBgoal |
double | m_Volume |
bool | m_ConstrainMuB |
bool | m_ConstrainMuQ |
bool | m_ConstrainMuS |
bool | m_ConstrainMuC |
bool | m_PCE |
bool | m_useOpenMP |
std::vector< double > | m_densities |
std::vector< double > | m_densitiestotal |
std::vector< std::vector< double > > | m_densitiesbyfeeddown |
std::vector< double > | m_Chem |
std::vector< double > | m_wprim |
std::vector< double > | m_wtot |
std::vector< double > | m_skewprim |
std::vector< double > | m_skewtot |
std::vector< double > | m_kurtprim |
std::vector< double > | m_kurttot |
std::vector< std::vector< double > > | m_PrimCorrel |
std::vector< std::vector< double > > | m_TotalCorrel |
std::vector< std::vector< double > > | m_PrimChargesCorrel |
std::vector< std::vector< double > > | m_FinalChargesCorrel |
std::vector< std::vector< double > > | m_Susc |
std::vector< std::vector< double > > | m_ProxySusc |
std::string | m_ValidityLog |
double | m_wnSum |
std::string | m_TAG |
ThermalModelEnsemble | m_Ensemble |
ThermalModelInteraction | m_InteractionModel |
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 crossterms excluded-volume model.
The model formulation can be found in
M.I. Gorenstein, A.P. Kostyuk, Y.D. Krivenko, J.Phys. G 25, L75 (1999), http://arxiv.org/pdf/nucl-th/9906068.pdf
and in
V. Vovchenko, H. Stoecker, Phys. Rev. C 95, 044904 (2017), http://arxiv.org/pdf/1606.06218.pdf
The system of transcendental equations for "partial pressures" of hadrons is solved using the Broyden's method.
Definition at line 36 of file ThermalModelEVCrossterms.h.
thermalfist::ThermalModelEVCrossterms::ThermalModelEVCrossterms | ( | ThermalParticleSystem * | TPS, |
const ThermalModelParameters & | params = ThermalModelParameters() |
||
) |
Construct a new ThermalModelEVCrossterms object.
TPS | A pointer to the ThermalParticleSystem object containing the particle list |
params | ThermalModelParameters object with current thermal parameters |
Definition at line 33 of file ThermalModelEVCrossterms.cpp.
|
virtual |
Destroy the ThermalModelEVCrossterms object.
Definition at line 46 of file ThermalModelEVCrossterms.cpp.
|
inlinevirtual |
The fraction of entropy carried by baryons (Ideal GCE only)
Reimplemented from thermalfist::ThermalModelBase.
Definition at line 84 of file ThermalModelEVCrossterms.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 561 of file ThermalModelEVCrossterms.cpp.
|
virtual |
Implements thermalfist::ThermalModelBase.
Definition at line 748 of file ThermalModelEVCrossterms.cpp.
|
virtual |
Implements thermalfist::ThermalModelBase.
Definition at line 757 of file ThermalModelEVCrossterms.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 544 of file ThermalModelEVCrossterms.cpp.
|
inlinevirtual |
The fraction of entropy carried by mesons (Ideal GCE only)
Reimplemented from thermalfist::ThermalModelBase.
Definition at line 85 of file ThermalModelEVCrossterms.h.
|
virtual |
Implementation of the equation of state functions.
Implements thermalfist::ThermalModelBase.
Definition at line 762 of file ThermalModelEVCrossterms.cpp.
|
virtual |
Calculates the primordial densities of all species.
Implements thermalfist::ThermalModelBase.
Definition at line 240 of file ThermalModelEVCrossterms.cpp.
|
protectedvirtual |
Definition at line 396 of file ThermalModelEVCrossterms.cpp.
|
protectedvirtual |
Definition at line 302 of file ThermalModelEVCrossterms.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 459 of file ThermalModelEVCrossterms.cpp.
|
virtual |
Change the particle list.
TPS | A pointer to new particle list. |
Reimplemented from thermalfist::ThermalModelBase.
Definition at line 143 of file ThermalModelEVCrossterms.cpp.
|
protectedvirtual |
Calculate the ideal gas density of particle species i for the given values of partial pressures.
i | 0-based particle specie index |
Pressure | Input vector of partial pressures (GeV fm \(^{-3}\)) |
Definition at line 149 of file ThermalModelEVCrossterms.cpp.
|
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
ri | A vector with radii parameters for all species. 0-based indices of the vector must corresponds to the 0-based indices of the particle list TPS() |
Reimplemented from thermalfist::ThermalModelBase.
Definition at line 50 of file ThermalModelEVCrossterms.cpp.
|
protectedvirtual |
The shift in the chemical potential of particle species i due to the excluded volume interactions.
Equal to \(- \sum_j \tilde{b_{ij}} p_j\)
i | 0-based particle specie index |
Reimplemented from thermalfist::ThermalModelBase.
Definition at line 768 of file ThermalModelEVCrossterms.cpp.
|
protected |
The "partial pressure" of hadron species i for the given total pressure in the diagonal model.
i | 0-based particle species index. |
P | Input pressure (GeV fm \(^{-3}\)) |
Definition at line 182 of file ThermalModelEVCrossterms.cpp.
|
inlinevirtual |
Kurtosis of primordial particle number fluctuations for species i.
Reimplemented from thermalfist::ThermalModelBase.
Definition at line 94 of file ThermalModelEVCrossterms.h.
|
inlinevirtual |
The scalar density of the particle species i.
Implements thermalfist::ThermalModelBase.
Definition at line 97 of file ThermalModelEVCrossterms.h.
|
inlinevirtual |
Scaled variance of primordial particle number fluctuations for species i.
Reimplemented from thermalfist::ThermalModelBase.
Definition at line 88 of file ThermalModelEVCrossterms.h.
|
inlinevirtual |
Skewness of primordial particle number fluctuations for species i.
Reimplemented from thermalfist::ThermalModelBase.
Definition at line 91 of file ThermalModelEVCrossterms.h.
|
protectedvirtual |
Calculate the ideal gas pressure of particle species i for the given values of partial pressures.
i | 0-based particle specie index |
Pressure | Input vector of partial pressures (GeV fm \(^{-3}\)) |
Definition at line 162 of file ThermalModelEVCrossterms.cpp.
|
protected |
The total pressure for the given input pressure in the diagonal model.
Input | pressure (GeV fm \(^{-3}\)) |
Definition at line 190 of file ThermalModelEVCrossterms.cpp.
|
virtual |
Reads the QvdW interaction parameters from a file.
Actual implementation is in a derived class.
filename | File with interaction parameters. |
Reimplemented from thermalfist::ThermalModelBase.
Definition at line 71 of file ThermalModelEVCrossterms.cpp.
|
protected |
Calculate the ideal gas scaled variance of particle species i number fluctuations for the given values of partial pressures.
i | 0-based particle specie index |
Pressure | Input vector of partial pressures (GeV fm \(^{-3}\)) |
Definition at line 175 of file ThermalModelEVCrossterms.cpp.
|
virtual |
Set the same excluded volume radius parameter for all species.
rad | Radius parameter (fm) |
Reimplemented from thermalfist::ThermalModelBase.
Definition at line 127 of file ThermalModelEVCrossterms.cpp.
|
virtual |
Set the excluded volume coefficient \( \tilde{b}_{ij} \).
Excluded parameter for repulsive interaction between particle species i and j.
i | 0-based index of the first particle species |
j | 0-based index of the second particle species |
b | Excluded volume parameter \( \tilde{b}_{ij} \) (fm \(^3\)) |
Reimplemented from thermalfist::ThermalModelBase.
Definition at line 137 of file ThermalModelEVCrossterms.cpp.
|
protectedvirtual |
Solves the transcendental equation of the corresponding diagonal EV model.
The diagonal EV model here has \(v_i = \tilde{b}_{ii}\).
The partial pressures of the diagonal model will be recorded into m_Ps.
Definition at line 197 of file ThermalModelEVCrossterms.cpp.
|
protectedvirtual |
Solves the system of transcdental equations for the pressure using the Broyden's method.
Solves \( p_i(T,\mu) = p_i^{\rm id} (T, \mu_i - \sum_j \tilde{b_{ij}} p_j). \)
Definition at line 216 of file ThermalModelEVCrossterms.cpp.
|
protected |
Definition at line 371 of file ThermalModelEVCrossterms.cpp.
|
virtual |
Excluded volume coefficient \( \tilde{b}_{ij} = 0 \).
Excluded parameter for repulsive interaction between particle species i and j.
i | 0-based index of the first particle species |
j | 0-based index of the second particle species |
Reimplemented from thermalfist::ThermalModelBase.
Definition at line 131 of file ThermalModelEVCrossterms.cpp.
|
virtual |
Write the QvdW interaction parameters to a file.
Actual implementation is in a derived class.
filename | Output file. |
Reimplemented from thermalfist::ThermalModelBase.
Definition at line 96 of file ThermalModelEVCrossterms.cpp.
|
protected |
Vector of ideal gas densities with shifted chemical potentials
Definition at line 194 of file ThermalModelEVCrossterms.h.
|
protected |
The (solved) total pressure
Definition at line 197 of file ThermalModelEVCrossterms.h.
|
protected |
Vector of (solved) partial pressures
Definition at line 195 of file ThermalModelEVCrossterms.h.
|
protected |
The (solved) entropy pressure
Definition at line 198 of file ThermalModelEVCrossterms.h.
|
protected |
Matrix of virial (excluded-volume) coefficients \( \tilde{b}_{ij} \)
Definition at line 196 of file ThermalModelEVCrossterms.h.