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

Class implementing the diagonal excluded-volume model in the strangeness-canonical ensemble. More...

#include <ThermalModelEVCanonicalStrangeness.h>

Inheritance diagram for thermalfist::ThermalModelEVCanonicalStrangeness:
thermalfist::ThermalModelCanonicalStrangeness thermalfist::ThermalModelBase

Public Member Functions

 ThermalModelEVCanonicalStrangeness (ThermalParticleSystem *TPS, const ThermalModelParameters &params=ThermalModelParameters())
 Construct a new Thermal ModelEVCanonicalStrangeness object.
 
virtual ~ThermalModelEVCanonicalStrangeness (void)
 Destroy the ThermalModelEVCanonicalStrangeness 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
 
virtual void CalculateEnergyDensitiesGCE ()
 Calculates the grand-canonical energy densities.
 
virtual void CalculatePressuresGCE ()
 Calculates the grand-canonical pressures.
 
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 ()
 
virtual void CalculateDensitiesGCE ()
 Calculates the particle densities in a grand-canonical ensemble.
 
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
 
- Public Member Functions inherited from thermalfist::ThermalModelCanonicalStrangeness
 ThermalModelCanonicalStrangeness (ThermalParticleSystem *TPS, const ThermalModelParameters &params=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 &params)
 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 ()
 
- Public Member Functions inherited from thermalfist::ThermalModelBase
 ThermalModelBase (ThermalParticleSystem *TPS, const ThermalModelParameters &params=ThermalModelParameters())
 Construct a new ThermalModelBase object.
 
virtual ~ThermalModelBase (void)
 
int ComponentsNumber () const
 Number of different particle species in the list.
 
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 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 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 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 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 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 PrepareModelEV ()
 
void ClearModelEV ()
 Clears m_modelEV.
 
virtual double MuShift (int id) const
 The shift in the chemical potential of particle species i due to the excluded volume interactions.
 

Protected Attributes

ThermalModelEVDiagonalm_modelEV
 
std::vector< double > m_v
 
double m_PNS
 
double m_Suppression
 
double m_EVNS
 
double m_EVS
 
- Protected Attributes inherited from thermalfist::ThermalModelCanonicalStrangeness
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

- Public Types inherited from thermalfist::ThermalModelBase
enum  ThermalModelEnsemble { GCE = 0 , CE = 1 , SCE = 2 , CCE = 3 }
 The list of statistical ensembles. More...
 
enum  ThermalModelInteraction {
  Ideal = 0 , DiagonalEV = 1 , CrosstermsEV = 2 , QvdW = 3 ,
  RealGas = 4 , MeanField = 5
}
 Type of interactions included in the HRG model. More...
 

Detailed Description

Class implementing the diagonal excluded-volume 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 ThermalModelEVCanonicalStrangeness.h.

Constructor & Destructor Documentation

◆ ThermalModelEVCanonicalStrangeness()

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

Construct a new Thermal ModelEVCanonicalStrangeness object.

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

Definition at line 23 of file ThermalModelEVCanonicalStrangeness.cpp.

◆ ~ThermalModelEVCanonicalStrangeness()

thermalfist::ThermalModelEVCanonicalStrangeness::~ThermalModelEVCanonicalStrangeness ( void )
virtual

Destroy the ThermalModelEVCanonicalStrangeness object.

Definition at line 34 of file ThermalModelEVCanonicalStrangeness.cpp.

Member Function Documentation

◆ CalculateDensitiesGCE()

void thermalfist::ThermalModelEVCanonicalStrangeness::CalculateDensitiesGCE ( )
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 39 of file ThermalModelEVCanonicalStrangeness.cpp.

◆ CalculateEigenvolumeFraction()

double thermalfist::ThermalModelEVCanonicalStrangeness::CalculateEigenvolumeFraction ( )
virtual

Definition at line 269 of file ThermalModelEVCanonicalStrangeness.cpp.

◆ CalculateEnergyDensitiesGCE()

void thermalfist::ThermalModelEVCanonicalStrangeness::CalculateEnergyDensitiesGCE ( )
virtual

Calculates the grand-canonical energy densities.

Reimplemented from thermalfist::ThermalModelCanonicalStrangeness.

Definition at line 50 of file ThermalModelEVCanonicalStrangeness.cpp.

◆ CalculateEnergyDensity()

double thermalfist::ThermalModelEVCanonicalStrangeness::CalculateEnergyDensity ( )
virtual

◆ CalculateEntropyDensity()

double thermalfist::ThermalModelEVCanonicalStrangeness::CalculateEntropyDensity ( )
virtual

◆ CalculatePressure()

double thermalfist::ThermalModelEVCanonicalStrangeness::CalculatePressure ( )
virtual

◆ CalculatePressuresGCE()

void thermalfist::ThermalModelEVCanonicalStrangeness::CalculatePressuresGCE ( )
virtual

Calculates the grand-canonical pressures.

Reimplemented from thermalfist::ThermalModelCanonicalStrangeness.

Definition at line 60 of file ThermalModelEVCanonicalStrangeness.cpp.

◆ CalculatePrimordialDensities()

void thermalfist::ThermalModelEVCanonicalStrangeness::CalculatePrimordialDensities ( )
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 72 of file ThermalModelEVCanonicalStrangeness.cpp.

◆ ClearModelEV()

void thermalfist::ThermalModelEVCanonicalStrangeness::ClearModelEV ( )
protected

Clears m_modelEV.

Definition at line 189 of file ThermalModelEVCanonicalStrangeness.cpp.

◆ ExcludedVolume()

double thermalfist::ThermalModelEVCanonicalStrangeness::ExcludedVolume ( int i) const

Definition at line 262 of file ThermalModelEVCanonicalStrangeness.cpp.

◆ FillVirial()

void thermalfist::ThermalModelEVCanonicalStrangeness::FillVirial ( const std::vector< double > & ri = std::vector<double>(0))
virtual

Fills the vector of particle eigenvolume parameters.

Parameters
riA vector of radii parameters for all particles (in fm)

Reimplemented from thermalfist::ThermalModelBase.

Definition at line 209 of file ThermalModelEVCanonicalStrangeness.cpp.

◆ FillVirialEV()

void thermalfist::ThermalModelEVCanonicalStrangeness::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.

Parameters
viA 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 219 of file ThermalModelEVCanonicalStrangeness.cpp.

◆ IsConservedChargeCanonical()

bool thermalfist::ThermalModelEVCanonicalStrangeness::IsConservedChargeCanonical ( ConservedCharge::Name charge) const
virtual

◆ MuShift()

double thermalfist::ThermalModelEVCanonicalStrangeness::MuShift ( int id) const
protectedvirtual

The shift in the chemical potential of particle species i due to the excluded volume interactions.

Equal to \(-v_i p\)

Parameters
i0-based particle specie index
Returns
The shift in the chemical potential

Definition at line 150 of file ThermalModelEVCanonicalStrangeness.cpp.

◆ PrepareModelEV()

void thermalfist::ThermalModelEVCanonicalStrangeness::PrepareModelEV ( )
protected

Calculates the necessary auxiliary quantities in the Diagonal EV model in the GCE consisting of non-strange particles

Definition at line 159 of file ThermalModelEVCanonicalStrangeness.cpp.

◆ ReadInteractionParameters()

void thermalfist::ThermalModelEVCanonicalStrangeness::ReadInteractionParameters ( const std::string & filename)
virtual

Read the set of eigenvolume parameters for all particles from an external file.

Parameters
filenameInput file name

Reimplemented from thermalfist::ThermalModelBase.

Definition at line 227 of file ThermalModelEVCanonicalStrangeness.cpp.

◆ SetRadius() [1/2]

void thermalfist::ThermalModelEVCanonicalStrangeness::SetRadius ( double )
virtual

Set the same excluded volume radius parameter for all species.

Parameters
radRadius parameter (fm)

Reimplemented from thermalfist::ThermalModelBase.

Definition at line 200 of file ThermalModelEVCanonicalStrangeness.cpp.

◆ SetRadius() [2/2]

void thermalfist::ThermalModelEVCanonicalStrangeness::SetRadius ( int ,
double  )
virtual

Set the radius parameter for particle species i.

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

Reimplemented from thermalfist::ThermalModelBase.

Definition at line 278 of file ThermalModelEVCanonicalStrangeness.cpp.

◆ WriteInteractionParameters()

void thermalfist::ThermalModelEVCanonicalStrangeness::WriteInteractionParameters ( const std::string & filename)
virtual

Write the set of eigenvolume parameters for all particles to an external file.

One particle specie per line.

Parameters
filenameOutput file name

Reimplemented from thermalfist::ThermalModelBase.

Definition at line 251 of file ThermalModelEVCanonicalStrangeness.cpp.

Member Data Documentation

◆ m_EVNS

double thermalfist::ThermalModelEVCanonicalStrangeness::m_EVNS
protected

Total eigenvolume of all non-strange hadrons

Definition at line 110 of file ThermalModelEVCanonicalStrangeness.h.

◆ m_EVS

double thermalfist::ThermalModelEVCanonicalStrangeness::m_EVS
protected

Total eigenvolume of all strange hadrons

Definition at line 111 of file ThermalModelEVCanonicalStrangeness.h.

◆ m_modelEV

ThermalModelEVDiagonal* thermalfist::ThermalModelEVCanonicalStrangeness::m_modelEV
protected

Pointer to the diagonal EV model in the GCE with non-strange particles only

Definition at line 106 of file ThermalModelEVCanonicalStrangeness.h.

◆ m_PNS

double thermalfist::ThermalModelEVCanonicalStrangeness::m_PNS
protected

Pressure of all non-strange hadrons

Definition at line 108 of file ThermalModelEVCanonicalStrangeness.h.

◆ m_Suppression

double thermalfist::ThermalModelEVCanonicalStrangeness::m_Suppression
protected

Common suppression factor, from non-strange hadrons

Definition at line 109 of file ThermalModelEVCanonicalStrangeness.h.

◆ m_v

std::vector<double> thermalfist::ThermalModelEVCanonicalStrangeness::m_v
protected

Vector of eigenvolumes of all hadrons

Definition at line 107 of file ThermalModelEVCanonicalStrangeness.h.


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