![]() |
Thermal-FIST 1.5
Package for hadron resonance gas model applications
|
Class implementing the quantum van der Waals model in the strangeness-canonical ensemble. More...
#include <ThermalModelVDWCanonicalStrangeness.h>
Public Member Functions | |
ThermalModelVDWCanonicalStrangeness (ThermalParticleSystem *TPS, const ThermalModelParameters ¶ms=ThermalModelParameters()) | |
Construct a new Thermal ModelEVCanonicalStrangeness object. | |
virtual | ~ThermalModelVDWCanonicalStrangeness (void) |
Destroy the ThermalModelEVCanonicalStrangeness object. | |
void | FillVirialEV (const std::vector< std::vector< double > > &bij=std::vector< std::vector< double > >(0)) |
Same as FillVirial() but uses the matrix of excluded-volume coefficients \( v_i \equiv b_{ii} \) as input instead of radii. | |
void | 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. | |
void | FillAttraction (const std::vector< std::vector< double > > &aij=std::vector< std::vector< double > >(0)) |
void | SetVirial (int i, int j, double b) |
Set the excluded volume coefficient \( \tilde{b}_{ij} \). | |
void | SetAttraction (int i, int j, double a) |
Set the vdW mean field attraction coefficient \( a_{ij} \). | |
double | VirialCoefficient (int i, int j) const |
Excluded volume coefficient \( \tilde{b}_{ij} = 0 \). | |
double | AttractionCoefficient (int i, int j) const |
QvdW mean field attraction coefficient \( a_{ij} \). | |
virtual void | ReadInteractionParameters (const std::string &filename) |
Reads the QvdW interaction parameters from a file. | |
virtual void | WriteInteractionParameters (const std::string &filename) |
Write the QvdW interaction parameters to a file. | |
virtual void | CalculateDensitiesGCE () |
Calculates the particle densities in a grand-canonical ensemble. | |
virtual void | CalculateEnergyDensitiesGCE () |
Calculates the grand-canonical energy densities. | |
virtual void | CalculatePressuresGCE () |
Calculates the grand-canonical pressures. | |
virtual void | CalculatePrimordialDensities () |
Calculates the primordial densities of all species. | |
virtual double | CalculateEnergyDensity () |
virtual double | CalculateEntropyDensity () |
virtual double | CalculatePressure () |
virtual bool | IsConservedChargeCanonical (ConservedCharge::Name charge) const |
![]() | |
ThermalModelCanonicalStrangeness (ThermalParticleSystem *TPS, const ThermalModelParameters ¶ms=ThermalModelParameters()) | |
Construct a new ThermalModelCanonicalStrangeness object. | |
virtual | ~ThermalModelCanonicalStrangeness (void) |
Destroy the ThermalModelCanonicalStrangeness object. | |
virtual void | CalculateSums (double Vc) |
Calculates the strangeness-canonical partition functions. | |
const std::vector< double > & | DensitiesGCE () const |
A vector of the grand-canonical particle number densities. | |
virtual void | SetParameters (const ThermalModelParameters ¶ms) |
The thermal parameters. | |
virtual void | SetStrangenessChemicalPotential (double muS) |
Override the base class method to always set \( \mu_S \) to zero. | |
virtual void | ChangeTPS (ThermalParticleSystem *TPS) |
Change the particle list. | |
virtual void | FixParameters () |
Method which actually implements ConstrainChemicalPotentials() (for backward compatibility). | |
virtual void | SetStatistics (bool stats) |
void | CalculateFluctuations () |
Dummy function. Fluctuations not yet supported. | |
virtual double | CalculateBaryonMatterEntropyDensity () |
virtual double | CalculateMesonMatterEntropyDensity () |
virtual double | ParticleScaledVariance (int) |
virtual double | ParticleScalarDensity (int) |
virtual double | CalculateEnergyDensityDerivativeT () |
![]() | |
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. | |
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 | SetCharmChemicalPotential (double muC) |
Set the charm chemical potential. | |
virtual void | SetGammaq (double gammaq) |
Set the light quark fugacity factor. | |
virtual void | SetGammaS (double gammaS) |
Set the strange quark fugacity factor. | |
virtual void | SetGammaC (double gammaC) |
Set the charm quark fugacity factor. | |
virtual void | SetBaryonCharge (int B) |
Set the total baryon number (for canonical ensemble only) | |
virtual void | SetElectricCharge (int Q) |
Set the total electric charge (for canonical ensemble only) | |
virtual void | SetStrangeness (int S) |
Set the total strangeness (for canonical ensemble only) | |
virtual void | SetCharm (int C) |
Set the total charm (for canonical ensemble only) | |
virtual void | SetRadius (double) |
Set the same excluded volume radius parameter for all species. | |
virtual void | SetRadius (int, double) |
Set the radius parameter for particle species i. | |
virtual void | SetRepulsion (int i, int j, double b) |
Same as SetVirial() but with a more clear name on what is actually does. | |
virtual void | DisableMesonMesonVirial () |
void | DisableMesonMesonRepulsion () |
virtual void | DisableMesonMesonAttraction () |
virtual void | DisableMesonBaryonVirial () |
void | DisableMesonBaryonRepulsion () |
virtual void | DisableMesonBaryonAttraction () |
virtual void | DisableBaryonBaryonVirial () |
void | DisableBaryonBaryonRepulsion () |
virtual void | DisableBaryonBaryonAttraction () |
virtual void | DisableBaryonAntiBaryonVirial () |
void | DisableBaryonAntiBaryonRepulsion () |
virtual void | DisableBaryonAntiBaryonAttraction () |
double | RepulsionCoefficient (int i, int j) const |
bool | QuantumStatistics () const |
virtual void | 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 | 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 | CalculateFeeddown () |
Calculates the total densities which include feeddown contributions. | |
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. | |
Protected Member Functions | |
void | PrepareModelVDW () |
void | ClearModelVDW () |
Clears m_modelVDW. | |
virtual void | CalculateSums (const std::vector< double > &Vcs) |
Calculates the necessary strangeness-canonical partition functions. | |
virtual double | MuShift (int id) const |
The shift in the chemical potential of particle species i due to the QvdW interactions. | |
Protected Attributes | |
ThermalModelVDWFull * | m_modelVDW |
std::vector< std::vector< double > > | m_Virial |
std::vector< std::vector< double > > | m_Attr |
double | m_PNS |
std::vector< double > | m_MuStar |
std::vector< double > | m_Suppression |
![]() | |
std::vector< double > | m_densitiesGCE |
std::vector< double > | m_energydensitiesGCE |
std::vector< double > | m_pressuresGCE |
std::vector< int > | m_StrVals |
std::map< int, int > | m_StrMap |
std::vector< double > | m_Zsum |
std::vector< double > | m_partialS |
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 quantum van der Waals model in the strangeness-canonical ensemble.
NOTE
The calculations are approximate and assume that strange particles form a small part of the total system, i.e. their contribution to the total density/pressure is close to negligible. Calculations may not be accurate if this condition is not fulfilled.
Definition at line 36 of file ThermalModelVDWCanonicalStrangeness.h.
thermalfist::ThermalModelVDWCanonicalStrangeness::ThermalModelVDWCanonicalStrangeness | ( | ThermalParticleSystem * | TPS, |
const ThermalModelParameters & | params = ThermalModelParameters() ) |
Construct a new Thermal ModelEVCanonicalStrangeness object.
TPS | A pointer to the ThermalParticleSystem object containing the particle list |
params | ThermalModelParameters object with current thermal parameters |
Definition at line 24 of file ThermalModelVDWCanonicalStrangeness.cpp.
|
virtual |
Destroy the ThermalModelEVCanonicalStrangeness object.
Definition at line 37 of file ThermalModelVDWCanonicalStrangeness.cpp.
|
virtual |
QvdW mean field attraction coefficient \( a_{ij} \).
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 360 of file ThermalModelVDWCanonicalStrangeness.cpp.
|
virtual |
Calculates the particle densities in a grand-canonical ensemble.
A non-GCE based derived class will override this method.
Reimplemented from thermalfist::ThermalModelCanonicalStrangeness.
Definition at line 42 of file ThermalModelVDWCanonicalStrangeness.cpp.
|
virtual |
Calculates the grand-canonical energy densities.
Reimplemented from thermalfist::ThermalModelCanonicalStrangeness.
Definition at line 53 of file ThermalModelVDWCanonicalStrangeness.cpp.
|
virtual |
Reimplemented from thermalfist::ThermalModelCanonicalStrangeness.
Definition at line 153 of file ThermalModelVDWCanonicalStrangeness.cpp.
|
virtual |
Reimplemented from thermalfist::ThermalModelCanonicalStrangeness.
Definition at line 172 of file ThermalModelVDWCanonicalStrangeness.cpp.
|
virtual |
Reimplemented from thermalfist::ThermalModelCanonicalStrangeness.
Definition at line 191 of file ThermalModelVDWCanonicalStrangeness.cpp.
|
virtual |
Calculates the grand-canonical pressures.
Reimplemented from thermalfist::ThermalModelCanonicalStrangeness.
Definition at line 63 of file ThermalModelVDWCanonicalStrangeness.cpp.
|
virtual |
Calculates the primordial densities of all species.
Calculates the primordial densities by applying the canonical corrections factors to the grand-canonical densities.
Reimplemented from thermalfist::ThermalModelCanonicalStrangeness.
Definition at line 118 of file ThermalModelVDWCanonicalStrangeness.cpp.
|
protectedvirtual |
Calculates the necessary strangeness-canonical partition functions.
Calculations are performed in a reduced strangeness correlation volume due to the eigenvolume corrections from non-strange particles
Vcs | A vector of effective correlation volumes which is seen by each particle species |
Definition at line 74 of file ThermalModelVDWCanonicalStrangeness.cpp.
|
protected |
Clears m_modelVDW.
Definition at line 260 of file ThermalModelVDWCanonicalStrangeness.cpp.
void thermalfist::ThermalModelVDWCanonicalStrangeness::FillAttraction | ( | const std::vector< std::vector< double > > & | aij = std::vector< std::vector<double> >(0) | ) |
Definition at line 300 of file ThermalModelVDWCanonicalStrangeness.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 271 of file ThermalModelVDWCanonicalStrangeness.cpp.
void thermalfist::ThermalModelVDWCanonicalStrangeness::FillVirialEV | ( | const std::vector< std::vector< double > > & | bij = std::vector< std::vector<double> >(0) | ) |
Same as FillVirial() but uses the matrix of excluded-volume coefficients \( v_i \equiv b_{ii} \) as input instead of radii.
bij | A vector with crossterms excluded-volume coefficients for all pairs of particle species |
Definition at line 292 of file ThermalModelVDWCanonicalStrangeness.cpp.
|
virtual |
Reimplemented from thermalfist::ThermalModelCanonicalStrangeness.
Definition at line 367 of file ThermalModelVDWCanonicalStrangeness.cpp.
|
protectedvirtual |
The shift in the chemical potential of particle species i due to the QvdW interactions.
i | 0-based particle specie index |
Definition at line 205 of file ThermalModelVDWCanonicalStrangeness.cpp.
|
protected |
Calculates the necessary auxiliary quantities in the QvdW model in the GCE consisting of non-strange particles
Definition at line 214 of file ThermalModelVDWCanonicalStrangeness.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 308 of file ThermalModelVDWCanonicalStrangeness.cpp.
|
inlinevirtual |
Set the vdW mean field attraction coefficient \( a_{ij} \).
i | 0-based index of the first particle species |
j | 0-based index of the second particle species |
a | vdW mean field attraction parameter \( a_{ij} \) (GeV fm \(^3\)) |
Reimplemented from thermalfist::ThermalModelBase.
Definition at line 64 of file ThermalModelVDWCanonicalStrangeness.h.
|
inlinevirtual |
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 62 of file ThermalModelVDWCanonicalStrangeness.h.
|
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 353 of file ThermalModelVDWCanonicalStrangeness.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 338 of file ThermalModelVDWCanonicalStrangeness.cpp.
|
protected |
Matrix of the excluded volume coefficients \( \tilde{b}_{ij} \)
Definition at line 118 of file ThermalModelVDWCanonicalStrangeness.h.
|
protected |
Pointer to the QvdW model in the GCE with non-strange particles only
Definition at line 116 of file ThermalModelVDWCanonicalStrangeness.h.
|
protected |
Vector of the shifted chemical potentials
Definition at line 120 of file ThermalModelVDWCanonicalStrangeness.h.
|
protected |
Matrix of the attractive QvdW coefficients \( a_{ij} \) Pressure of all non-strange hadrons
Definition at line 119 of file ThermalModelVDWCanonicalStrangeness.h.
|
protected |
Common suppression factor, from non-strange hadrons
Definition at line 121 of file ThermalModelVDWCanonicalStrangeness.h.
|
protected |
Definition at line 117 of file ThermalModelVDWCanonicalStrangeness.h.