8#ifndef THERMALMODELBASE_H
9#define THERMALMODELBASE_H
79 virtual void FillVirial(
const std::vector<double> & ri = std::vector<double>(0));
112 void SetOMP(
bool openMP) { m_useOpenMP = openMP; }
495 double SoverB()
const {
return m_SBgoal; }
507 double QoverB()
const {
return m_QBgoal; }
517 double Volume()
const {
return m_Parameters.V; }
615 double muBinit = 0.,
double muQinit = 0.,
double muSinit = 0.,
double muCinit = 0.,
616 bool ConstrMuB =
true,
bool ConstrMuQ =
true,
bool ConstrMuS =
true,
bool ConstrMuC =
true);
645 double muBinit = 0.,
double muQinit = 0.,
double muSinit = 0.,
double muCinit = 0.,
646 bool ConstrMuB =
true,
bool ConstrMuQ =
true,
bool ConstrMuS =
true,
bool ConstrMuC =
true);
914 virtual std::vector<double>
CalculateChargeFluctuations(
const std::vector<double> &chgs,
int order = 4,
bool dimensionfull =
false);
933 virtual std::vector<double> CalculateGeneralizedSusceptibilities(
const std::vector<std::vector<double>> &chgs);
942 virtual void SetMultipleSolutionsMode(
bool search) { }
950 virtual bool UseMultipleSolutionsMode()
const {
return false; }
959 double Pressure() {
return CalculatePressure(); }
968 double HadronDensity() {
return CalculateHadronDensity(); }
971 double BaryonDensity() {
return CalculateBaryonDensity(); }
974 double ElectricChargeDensity() {
return CalculateChargeDensity(); }
977 double StrangenessDensity() {
return CalculateStrangenessDensity(); }
980 double CharmDensity() {
return CalculateCharmDensity(); }
983 double AbsoluteBaryonDensity() {
return CalculateAbsoluteBaryonDensity(); }
986 double AbsoluteElectricChargeDensity() {
return CalculateAbsoluteChargeDensity(); }
989 double AbsoluteStrangenessDensity() {
return CalculateAbsoluteStrangenessDensity(); }
992 double AbsoluteCharmDensity() {
return CalculateAbsoluteCharmDensity(); }
995 double HeatCapacityMu() {
return CalculateHeatCapacityMu(); }
998 [[deprecated(
"Use HeatCapacityMu() instead")]]
999 double SpecificHeatChem() {
return CalculateHeatCapacityMu(); }
1014 double cs2(
bool rhoBconst =
true,
bool rhoQconst =
true,
bool rhoSconst =
true,
bool rhoCconst =
true) {
1015 return CalculateAdiabaticSpeedOfSoundSquared(rhoBconst, rhoQconst, rhoSconst, rhoCconst);
1031 double cT2(
bool rhoBconst =
true,
bool rhoQconst =
true,
bool rhoSconst =
true,
bool rhoCconst =
true) {
1032 return CalculateIsothermalSpeedOfSoundSquared(rhoBconst, rhoQconst, rhoSconst, rhoCconst);
1042 double HeatCapacity(
bool rhoBconst =
true,
bool rhoQconst =
true,
bool rhoSconst =
true,
bool rhoCconst =
true) {
1043 return CalculateHeatCapacity(rhoBconst, rhoQconst, rhoSconst, rhoCconst);
1067 std::vector<std::vector<double>> PressureHessian();
1077 bool IsThermodynamicallyStable();
1081 virtual double CalculatePressure() = 0;
1082 virtual double CalculateEnergyDensity() = 0;
1083 virtual double CalculateEntropyDensity() = 0;
1084 virtual double CalculateHadronDensity();
1085 virtual double CalculateBaryonDensity();
1086 virtual double CalculateChargeDensity();
1087 virtual double CalculateStrangenessDensity();
1088 virtual double CalculateCharmDensity();
1089 virtual double CalculateAbsoluteBaryonDensity();
1090 virtual double CalculateAbsoluteChargeDensity();
1091 virtual double CalculateAbsoluteStrangenessDensity();
1092 virtual double CalculateAbsoluteCharmDensity();
1093 virtual double CalculateAbsoluteStrangenessDensityModulo();
1094 virtual double CalculateAbsoluteCharmDensityModulo();
1095 virtual double CalculateEnergyDensityDerivativeT() = 0;
1096 virtual double CalculateEntropyDensityDerivativeT() = 0;
1097 virtual double CalculateHeatCapacityMu();
1098 virtual double CalculateAdiabaticSpeedOfSoundSquared(
bool rhoBconst =
true,
bool rhoQconst =
true,
bool rhoSconst =
true,
bool rhoCconst =
true);
1099 virtual double CalculateIsothermalSpeedOfSoundSquared(
bool rhoBconst =
true,
bool rhoQconst =
true,
bool rhoSconst =
true,
bool rhoCconst =
true);
1100 virtual double CalculateHeatCapacity(
bool rhoBconst =
true,
bool rhoQconst =
true,
bool rhoSconst =
true,
bool rhoCconst =
true);
1104 virtual double CalculateArbitraryChargeDensity();
1107 virtual double CalculateBaryonMatterEntropyDensity() {
return 0.; }
1110 virtual double CalculateMesonMatterEntropyDensity() {
return 0.; }
1113 virtual double ParticleScaledVariance(
int) {
return 1.; }
1116 virtual double ParticleSkewness(
int) {
return 1.; }
1119 virtual double ParticleKurtosis(
int) {
return 1.; }
1123 virtual double CalculateEigenvolumeFraction() {
return 0.; }
1126 virtual double ParticleScalarDensity(
int i) = 0;
1128 virtual double GetMaxDiff()
const {
return m_MaxDiff; }
1133 virtual bool IsLastSolutionOK()
const {
return m_LastCalculationSuccessFlag; }
1139 double GetdndT(
int i)
const;
1146 double GetDensity(
long long PDGID,
int feeddown) {
return GetDensity(PDGID,
static_cast<Feeddown::Type>(feeddown)); }
1166 double GetYield(
long long PDGID,
Feeddown::Type feeddown) {
return GetDensity(PDGID, feeddown) *
Volume(); }
1168 std::vector<double> GetIdealGasDensities()
const;
1172 ThermalParticleSystem* TPS() {
return m_TPS; }
1178 const std::vector<double>& Densities()
const {
return m_densities; }
1179 std::vector<double>& Densities() {
return m_densities; }
1187 const std::vector<double>& TotalDensities()
const {
return m_densitiestotal; }
1192 const std::vector< std::vector<double> >& AllDensities()
const {
return m_densitiesbyfeeddown; }
1199 void SetTAG(
const std::string & tag) { m_TAG = tag; }
1202 const std::string& TAG()
const {
return m_TAG; }
1205 int PdgToId(
long long pdgid)
const {
return m_TPS->PdgToId(pdgid); }
1208 virtual void ResetCalculatedFlags();
1212 bool IsCalculated()
const {
return m_Calculated; }
1216 bool IsTwoParticleCorrelationsCalculated()
const {
return m_TwoParticleCorrelationsCalculated; }
1220 bool IsFluctuationsCalculated()
const {
return m_FluctuationsCalculated; }
1224 bool IsSusceptibilitiesCalculated()
const {
return m_SusceptibilitiesCalculated; }
1228 bool IsTemperatureDerivativesCalculated()
const {
return m_TemperatureDerivativesCalculated; }
1232 bool IsGCECalculated()
const {
return m_GCECalculated; }
1241 double ScaledVariancePrimordial(
int id)
const {
return (
id >= 0 &&
id <
static_cast<int>(m_wprim.size())) ? m_wprim[id] : 1.; }
1253 double ScaledVarianceTotal(
int id)
const {
return (
id >= 0 &&
id <
static_cast<int>(m_wtot.size())) ? m_wtot[id] : 1.; }
1262 double SkewnessPrimordial(
int id)
const {
return (
id >= 0 &&
id <
static_cast<int>(m_skewprim.size())) ? m_skewprim[id] : 1.; }
1274 double SkewnessTotal(
int id)
const {
return (
id >= 0 &&
id <
static_cast<int>(m_skewtot.size())) ? m_skewtot[id] : 1.; }
1283 double KurtosisPrimordial(
int id)
const {
return (
id >= 0 &&
id <
static_cast<int>(m_kurtprim.size())) ? m_kurtprim[id] : 1.; }
1295 double KurtosisTotal(
int id)
const {
return (
id >= 0 &&
id <
static_cast<int>(m_kurttot.size())) ? m_kurttot[id] : 1.; }
1376 double ChargedMultiplicity(
int type = 0);
1388 double ChargedScaledVariance(
int type = 0);
1397 double ChargedMultiplicityFinal(
int type = 0);
1406 double ChargedScaledVarianceFinal(
int type = 0);
1432 void SetDensityModelForParticleSpecies(
int i, GeneralizedDensity* density_model = NULL);
1437 void SetDensityModelForParticleSpeciesByPdg(
long long PDGID, GeneralizedDensity* density_model = NULL);
1442 void ClearDensityModels();
1444 const IdealGasFunctions::IdealGasFunctionsExtraConfig& GetIdealGasFunctionsExtraConfig()
const {
return m_IGFExtraConfig; }
1454 void SetMagneticField(
double B = 0.0,
int lmax = -1);
1460 void RecomputeThresholdsDueToMagneticField();
1468 std::vector<double> PartialPressures();
1471 void ClearMagneticField();
1474 ThermalModelParameters m_Parameters;
1475 ThermalParticleSystem* m_TPS;
1478 IdealGasFunctions::IdealGasFunctionsExtraConfig m_IGFExtraConfig;
1480 bool m_LastCalculationSuccessFlag;
1484 bool m_FeeddownCalculated;
1485 bool m_TwoParticleCorrelationsCalculated;
1486 bool m_FluctuationsCalculated;
1487 bool m_SusceptibilitiesCalculated;
1488 bool m_GCECalculated;
1491 bool m_QuantumStats;
1496 bool m_ConstrainMuB;
1497 bool m_ConstrainMuQ;
1498 bool m_ConstrainMuS;
1499 bool m_ConstrainMuC;
1505 std::vector<double> m_densities;
1506 std::vector<double> m_densitiestotal;
1508 std::vector< std::vector<double> > m_densitiesbyfeeddown;
1509 std::vector<double> m_Chem;
1512 std::vector<double> m_wprim;
1513 std::vector<double> m_wtot;
1516 std::vector<double> m_skewprim;
1517 std::vector<double> m_skewtot;
1520 std::vector<double> m_kurtprim;
1521 std::vector<double> m_kurttot;
1524 std::vector< std::vector<double> > m_PrimCorrel;
1525 std::vector< std::vector<double> > m_TotalCorrel;
1528 std::vector< std::vector<double> > m_dmusdmu;
1532 std::vector< std::vector<double> > m_PrimChargesCorrel;
1533 std::vector< std::vector<double> > m_FinalChargesCorrel;
1536 std::vector< std::vector<double> > m_Susc;
1539 std::vector< std::vector<double> > m_dSuscdT;
1542 std::vector< std::vector<double> > m_ProxySusc;
1545 bool m_TemperatureDerivativesCalculated;
1547 std::vector<double> m_dndT;
1549 std::vector<double> m_dmusdT;
1551 std::vector< std::vector<double> > m_PrimChi2sdT;
1557 std::string m_ValidityLog;
1567 virtual double MuShift(
int )
const {
return 0.; }
1571 void ResetChemicalPotentials();
1573 double GetDensity(
long long PDGID,
const std::vector<double> *dens);
1575 class BroydenEquationsChem :
public BroydenEquations
1578 BroydenEquationsChem(ThermalModelBase *model) : BroydenEquations(), m_THM(model) { m_N = 2; }
1579 std::vector<double> Equations(
const std::vector<double> &x);
1581 ThermalModelBase *m_THM;
1584 class BroydenJacobianChem :
public BroydenJacobian
1587 BroydenJacobianChem(ThermalModelBase *model) : BroydenJacobian(), m_THM(model) { }
1588 std::vector<double> Jacobian(
const std::vector<double> &x);
1590 ThermalModelBase *m_THM;
1593 class BroydenChem :
public Broyden
1596 BroydenChem(ThermalModelBase *THM, BroydenEquations *eqs = NULL, BroydenJacobian *jaco = NULL) : Broyden(eqs, jaco) { m_THM = THM; }
1597 ~BroydenChem(
void) { }
1598 std::vector<double> Solve(
const std::vector<double> &x0, BroydenSolutionCriterium *solcrit = NULL,
int max_iterations = MAX_ITERS);
1600 ThermalModelBase *m_THM;
1604 class BroydenEquationsChemTotals :
public BroydenEquations
1607 BroydenEquationsChemTotals(
const std::vector<int> & vConstr,
const std::vector<int> & vType,
const std::vector<double> & vTotals, ThermalModelBase *model);
1608 std::vector<double> Equations(
const std::vector<double> &x);
1610 std::vector<int> m_Constr;
1611 std::vector<int> m_Type;
1612 std::vector<double> m_Totals;
1613 ThermalModelBase *m_THM;
1616 class BroydenJacobianChemTotals :
public BroydenJacobian
1619 BroydenJacobianChemTotals(
const std::vector<int> & vConstr,
const std::vector<int> & vType,
const std::vector<double> & vTotals, ThermalModelBase *model) : BroydenJacobian(), m_Constr(vConstr), m_Type(vType), m_Totals(vTotals), m_THM(model) { }
1620 std::vector<double> Jacobian(
const std::vector<double> &x);
1622 std::vector<int> m_Constr;
1623 std::vector<int> m_Type;
1624 std::vector<double> m_Totals;
1625 ThermalModelBase *m_THM;
Implementation of the generic Broyden's method routines.
map< string, double > params
virtual void CalculateSusceptibilityMatrix()
Calculates the conserved charges susceptibility matrix.
virtual double FullIdealChemicalPotential(int i) const
Chemical potential entering the ideal gas expressions of particle species i.
virtual void FixParametersNoReset()
Method which actually implements ConstrainChemicalPotentialsNoReset() (for backward compatibility).
void SetSoverB(double SB)
The entropy per baryon ratio to be used to constrain the baryon chemical potential.
virtual void CalculateDensitiesGCE()
Calculates the particle densities in a grand-canonical ensemble.
virtual void ChangeTPS(ThermalParticleSystem *TPS)
Change the particle list.
virtual void SetGammaC(double gammaC)
Set the charm quark fugacity factor.
virtual void CalculateParticleChargeCorrelationMatrix()
Calculates the matrix of correlators between primordial (and also final) particle numbers and conserv...
virtual void SetTemperature(double T)
Set the temperature.
virtual void CalculateFluctuations()
Computes the fluctuation observables.
bool UsePartialChemicalEquilibrium()
Whether partial chemical equilibrium with additional chemical potentials is used.
void SetQoverB(double QB)
The electric-to-baryon charge ratio to be used to constrain the electric chemical potential.
virtual void SetRepulsion(int i, int j, double b)
Same as SetVirial() but with a more clear name on what is actually does.
void DisableBaryonAntiBaryonRepulsion()
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 which reproduce the specified total baryon,...
bool ConstrainMuS() const
virtual double TwoParticleSusceptibilityPrimordialByPdg(long long id1, long long id2)
Returns the computed primordial particle number (cross-)susceptibility for particles with pdg codes ...
virtual void DisableMesonMesonVirial()
void SetCanonicalVolumeRadius(double Rc)
Set the canonical correlation system radius.
bool ConstrainMuQ() const
double RepulsionCoefficient(int i, int j) const
void ConstrainMuB(bool constrain)
virtual std::vector< double > CalculateChargeFluctuations(const std::vector< double > &chgs, int order=4, bool dimensionfull=false)
Calculates fluctuations (diagonal susceptibilities) of an arbitrary "conserved" charge.
virtual void SetVirial(int, int, double)
Set the excluded volume coefficient .
virtual ~ThermalModelBase(void)
virtual void FillChemicalPotentials()
Sets the chemical potentials of all particles.
virtual void CalculateTemperatureDerivatives()
Computes the temperature derivatives of densities, shifted chemical potentials, and primordial hadron...
virtual void DisableMesonBaryonAttraction()
virtual void WriteInteractionParameters(const std::string &)
Write the QvdW interaction parameters to a file.
void ConstrainMuC(bool constrain)
void ConstrainChemicalPotentials(bool resetInitialValues=true)
Constrains the chemical potentials by the conservation laws imposed.
virtual void FixParameters()
Method which actually implements ConstrainChemicalPotentials() (for backward compatibility).
bool ConstrainMuC() const
virtual void DisableMesonBaryonVirial()
virtual double TwoParticleSusceptibilityFinalByPdg(long long id1, long long id2)
Returns the computed final particle number (cross-)susceptibility for particles with pdg codes id1 a...
virtual double TwoParticleSusceptibilityFinal(int i, int j) const
Returns the computed final particle number (cross-)susceptibility for particles with ids i and j....
virtual void SetGammaS(double gammaS)
Set the strange quark fugacity factor.
virtual void SetCharm(int C)
Set the total charm (for canonical ensemble only)
ThermalModelInteraction
Type of interactions included in the HRG model.
@ QvdW
Quantum van der Waals model.
@ MeanField
Mean field model. Not yet fully implemented.
@ DiagonalEV
Diagonal excluded volume model.
@ CrosstermsEV
Crossterms excluded volume model.
virtual double FinalNetParticleChargeSusceptibilityByPdg(long long id1, ConservedCharge::Name chg)
Returns the computed final (after decays) net particle vs conserved charge (cross-)susceptibility fo...
void ConstrainMuQ(bool constrain)
void SetResonanceWidthIntegrationType(ThermalParticle::ResonanceWidthIntegration type)
Set the ThermalParticle::ResonanceWidthIntegration scheme for all particles. Calls the corresponding ...
virtual void CalculateFeeddown()
Calculates the total densities which include feeddown contributions.
virtual double TwoParticleSusceptibilityTemperatureDerivativePrimordial(int i, int j) const
Returns the computed temperature derivative of the primordial particle number (cross-)susceptibility ...
virtual void ReadInteractionParameters(const std::string &)
Reads the QvdW interaction parameters from a file.
void SetUseWidth(bool useWidth)
Sets whether finite resonance widths are used. Deprecated.
virtual double AttractionCoefficient(int, int) const
QvdW mean field attraction coefficient .
virtual double FinalParticleChargeSusceptibilityByPdg(long long id1, ConservedCharge::Name chg)
Returns the computed final (after decays) particle vs conserved charge (cross-)susceptibility for pa...
virtual void DisableBaryonAntiBaryonAttraction()
void UsePartialChemicalEquilibrium(bool usePCE)
Sets whether partial chemical equilibrium with additional chemical potentials is used.
virtual double FinalParticleChargeSusceptibility(int i, ConservedCharge::Name chg) const
Returns the computed final (after decays) particle vs conserved charge (cross-)susceptibility for pa...
virtual void SetChemicalPotentials(const std::vector< double > &chem=std::vector< double >(0))
Sets the chemical potentials of all particles.
void SetOMP(bool openMP)
OpenMP support. Currently not used.
double ChemicalPotential(int i) const
Chemical potential of particle species i.
virtual void SetElectricCharge(int Q)
Set the total electric charge (for canonical ensemble only)
virtual void SetBaryonChemicalPotential(double muB)
Set the baryon chemical potential.
virtual void CalculateTwoParticleCorrelations()
Computes the fluctuations and correlations of the primordial particle numbers.
virtual double SusceptibilityDimensionfull(ConservedCharge::Name i, ConservedCharge::Name j) const
A 2nd order susceptibility of conserved charges.
void SetNormBratio(bool normBratio)
Whether branching ratios are renormalized to 100%.
void UpdateParameters()
Calls SetParameters() with current m_Parameters.
virtual void SetParameters(const ThermalModelParameters ¶ms)
The thermal parameters.
virtual void DisableMesonMesonAttraction()
virtual void SetRadius(double)
Set the same excluded volume radius parameter for all species.
virtual void SetElectricChemicalPotential(double muQ)
Set the electric chemical potential.
virtual double PrimordialNetParticleChargeSusceptibilityByPdg(long long id1, ConservedCharge::Name chg)
Returns the computed primordial net particle vs conserved charge (cross-)susceptibility for particle...
virtual double TwoParticleSusceptibilityTemperatureDerivativePrimordialByPdg(long long id1, long long id2)
Returns the computed temperature derivative of the primordial particle number (cross-)susceptibility ...
virtual double NetParticleSusceptibilityFinalByPdg(long long id1, long long id2)
Returns the computed (cross-)susceptibility between final net-particle numbers for pdg codes id1 and...
virtual void SetChemicalPotential(int i, double chem)
Sets the chemical potential of particle species i.
std::string ValidityCheckLog() const
All messaged which occured during the validation procedure in the ValidateCalculation() method.
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 which reproduce the specified baryon,...
void SetVolume(double Volume)
Sets the system volume.
virtual void ValidateCalculation()
Checks whether issues have occured during the calculation of particle densities in the CalculateDensi...
void SetResonanceWidthShape(ThermalParticle::ResonanceWidthShape shape)
Set the ThermalParticle::ResonanceWidthShape for all particles. Calls the corresponding method in TPS...
void DisableMesonMesonRepulsion()
int ComponentsNumber() const
Number of different particle species in the list.
virtual void SetBaryonCharge(int B)
Set the total baryon number (for canonical ensemble only)
virtual void DisableBaryonBaryonVirial()
virtual double VirialCoefficient(int, int) const
Excluded volume coefficient .
virtual void CalculateProxySusceptibilityMatrix()
Calculates the susceptibility matrix of conserved charges proxies.
virtual void DisableBaryonAntiBaryonVirial()
ThermalModelBase(ThermalParticleSystem *TPS, const ThermalModelParameters ¶ms=ThermalModelParameters())
Construct a new ThermalModelBase object.
virtual void DisableBaryonBaryonAttraction()
virtual double PrimordialParticleChargeSusceptibility(int i, ConservedCharge::Name chg) const
Returns the computed primordial particle vs conserved charge (cross-)susceptibility for particle wit...
virtual void CalculatePrimordialDensities()=0
Calculates the primordial densities of all species.
void SetCanonicalVolume(double Volume)
Set the canonical correlation volume V .
virtual void SetStrangenessChemicalPotential(double muS)
Set the strangeness chemical potential.
bool UseWidth() const
Whether finite resonance widths are considered.
void DisableMesonBaryonRepulsion()
void ConstrainMuS(bool constrain)
virtual double TwoParticleSusceptibilityPrimordial(int i, int j) const
Returns the computed primordial particle number (cross-)susceptibility for particles with ids i and ...
virtual void SetClusterExpansionOrder(int order)
Set the number of terms in the cluster expansion method. Calls the corresponding method in TPS().
virtual void CalculateDensities()
Calculates the primordial and total (after decays) densities of all species.
virtual void SetCalculationType(IdealGasFunctions::QStatsCalculationType type)
Sets the CalculationType() method to evaluate quantum statistics. Calls the corresponding method in T...
virtual void FillVirial(const std::vector< double > &ri=std::vector< double >(0))
Fills the excluded volume coefficients based on the provided radii parameters for all species.
const std::vector< double > & ChemicalPotentials() const
A vector of chemical potentials of all particles.
virtual void SetGammaq(double gammaq)
Set the light quark fugacity factor.
virtual void SetStrangeness(int S)
Set the total strangeness (for canonical ensemble only)
virtual void SetAttraction(int, int, double)
Set the vdW mean field attraction coefficient .
virtual double PrimordialParticleChargeSusceptibilityByPdg(long long id1, ConservedCharge::Name chg)
Returns the computed primordial particle vs conserved charge (cross-)susceptibility for particle wit...
const ThermalModelParameters & Parameters() const
virtual double NetParticleSusceptibilityPrimordialByPdg(long long id1, long long id2)
Returns the computed (cross-)susceptibility between primordial net-particle numbers for pdg codes id...
virtual void SetCharmChemicalPotential(double muC)
Set the charm chemical potential.
virtual void SetStatistics(bool stats)
double Volume() const
System volume (fm )
virtual void SetRadius(int, double)
Set the radius parameter for particle species i.
void DisableBaryonBaryonRepulsion()
ThermalModelEnsemble
The list of statistical ensembles.
@ SCE
Strangeness-canonical ensemble.
@ GCE
Grand canonical ensemble.
@ CCE
Charm-canonical ensemble.
virtual void CalculateTwoParticleFluctuationsDecays()
Computes particle number correlations and fluctuations for all final-state particles which are marked...
double CanonicalVolume() const
The canonical correlation volume V (fm )
void SetVolumeRadius(double R)
Sets the system radius.
bool QuantumStatistics() const
bool ConstrainMuB() const
ResonanceWidthIntegration
Treatment of finite resonance widths.
ResonanceWidthShape
Relativistic vs non-relativistic Breit-Wigner shape.
Class containing the particle list.
int PdgToId(long long pdg)
QStatsCalculationType
Identifies whether quantum statistics are to be computed using the cluster expansion or numerical int...
constexpr double Pi()
Pi constant.
The main namespace where all classes and functions of the Thermal-FIST library reside.
Name
Set of all conserved charges considered.
Structure containing all thermal parameters of the model.
Contains some extra mathematical functions used in the code.