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

Class implementing the ideal HRG model in the canonical ensemble. More...

#include <ThermalModelCanonical.h>

Inheritance diagram for thermalfist::ThermalModelCanonical:
thermalfist::ThermalModelBase

Public Member Functions

 ThermalModelCanonical (ThermalParticleSystem *TPS, const ThermalModelParameters &params=ThermalModelParameters())
 Construct a new ThermalModelCanonical object.
 
virtual ~ThermalModelCanonical (void)
 Destroy the ThermalModelCanonical object.
 
virtual void CalculateQuantumNumbersRange (bool computeFluctuations=false)
 Calculates the range of quantum numbers values for which it is necessary to compute the canonical partition functions.
 
virtual void CalculatePartitionFunctions (double Vc=-1.)
 Calculates all necessary canonical partition functions.
 
virtual bool IsParticleCanonical (const ThermalParticle &part)
 Determines whether the specified ThermalParticle is treat canonically or grand-canonically in the present setup.
 
virtual void ConserveBaryonCharge (bool conserve=true)
 Specifies whether the baryon number is treated canonically.
 
virtual void ConserveElectricCharge (bool conserve=true)
 Specifies whether the electric charge is treated canonically.
 
virtual void ConserveStrangeness (bool conserve=true)
 Specifies whether the strangeness charge is treated canonically.
 
virtual void ConserveCharm (bool conserve=true)
 Specifies whether the charm charge is treated canonically.
 
virtual bool IsConservedChargeCanonical (ConservedCharge::Name charge) const
 
virtual double GetGCEDensity (int i) const
 Density of particle species i in the grand-canonical ensemble.
 
int IntegrationIterationsMultiplier () const
 The multiplier of the number of iterations in the numerical integration.
 
void SetIntegrationIterationsMultiplier (int multiplier)
 Assigns the multiplier of the number of iterations in the numerical integration.
 
void ChangeTPS (ThermalParticleSystem *TPS)
 Change the particle list.
 
virtual void SetStatistics (bool stats)
 
virtual void FixParameters ()
 Method which actually implements ConstrainChemicalPotentials() (for backward compatibility).
 
virtual void FixParametersNoReset ()
 Method which actually implements ConstrainChemicalPotentialsNoReset() (for backward compatibility).
 
virtual void CalculatePrimordialDensities ()
 Calculates the primordial densities of all species.
 
virtual void ValidateCalculation ()
 Checks whether issues have occured during the calculation of particle densities in the CalculateDensities() method.
 
virtual double ParticleScaledVariance (int part)
 
virtual void CalculateTwoParticleCorrelations ()
 Computes the fluctuations and correlations of the primordial particle numbers.
 
virtual void CalculateFluctuations ()
 Computes the fluctuation observables.
 
virtual double CalculateEnergyDensity ()
 
virtual double CalculatePressure ()
 
virtual double CalculateEntropyDensity ()
 
virtual double CalculateEnergyDensityDerivativeT ()
 
virtual double ParticleScalarDensity (int)
 Returns the scalar density of a particle species.
 
- 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.
 
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.
 
bool UseWidth () const
 Whether finite resonance widths are considered.
 
void SetUseWidth (bool useWidth)
 Sets whether finite resonance widths are used. Deprecated.
 
void SetUseWidth (ThermalParticle::ResonanceWidthIntegration type)
 Sets the finite resonance widths scheme to use.
 
void SetNormBratio (bool normBratio)
 Whether branching ratios are renormalized to 100%.
 
bool NormBratio () const
 
void SetOMP (bool openMP)
 OpenMP support. Currently not used.
 
virtual void SetParameters (const ThermalModelParameters &params)
 The thermal parameters.
 
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 SetStrangenessChemicalPotential (double muS)
 Set the strangeness chemical potential.
 
virtual void SetCharmChemicalPotential (double muC)
 Set the charm chemical potential.
 
virtual void SetGammaq (double gammaq)
 Set the light quark fugacity factor.
 
virtual void SetGammaS (double gammaS)
 Set the strange quark fugacity factor.
 
virtual void SetGammaC (double gammaC)
 Set the charm quark fugacity factor.
 
virtual void SetBaryonCharge (int B)
 Set the total baryon number (for canonical ensemble only)
 
virtual void SetElectricCharge (int Q)
 Set the total electric charge (for canonical ensemble only)
 
virtual void SetStrangeness (int S)
 Set the total strangeness (for canonical ensemble only)
 
virtual void SetCharm (int C)
 Set the total charm (for canonical ensemble only)
 
virtual void 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 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 void ReadInteractionParameters (const std::string &)
 Reads the QvdW interaction parameters from a file.
 
virtual void WriteInteractionParameters (const std::string &)
 Write the QvdW interaction parameters to a file.
 
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 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.
 
std::string ValidityCheckLog () const
 All messaged which occured during the validation procedure in the ValidateCalculation() method.
 
virtual void CalculateDensitiesGCE ()
 Calculates the particle densities in a grand-canonical ensemble.
 
virtual void CalculateFeeddown ()
 Calculates the total densities which include feeddown contributions.
 
virtual void CalculateTemperatureDerivatives ()
 Computes the temperature derivatives of densities, shifted chemical potentials, and primordial hadron number susceptibilities.
 
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 Attributes

std::vector< QuantumNumbersm_QNvec
 A set of QuantumNumbers combinations for which it is necessary to compute the canonical partition function.
 
std::map< QuantumNumbers, int > m_QNMap
 Maps QuantumNumbers combinations to a 1-dimensional index.
 
std::vector< double > m_Corr
 A vector of chemical factors.
 
std::vector< double > m_PartialZ
 The computed canonical partition function.
 
int m_IntegrationIterationsMultiplier
 A multiplier to increase the number of iterations during the numerical integration used to calculate the partition functions.
 
int m_BMAX
 Maximum baryon number.
 
int m_QMAX
 Maximum electric charge.
 
int m_SMAX
 Maximum strangeness.
 
int m_CMAX
 Maximum charm.
 
int m_BMAX_list
 
int m_QMAX_list
 
int m_SMAX_list
 
int m_CMAX_list
 
double m_MultExp
 Exponential multiplier for canonical partition function calculations.
 
double m_MultExpBanalyt
 Exponential multiplier for analytical baryon fugacity calculations.
 
ThermalModelIdealm_modelgce
 Pointer to a ThermalModelIdeal object used for GCE calculations.
 
int m_BCE
 Flag indicating if baryon charge is conserved canonically.
 
int m_QCE
 Flag indicating if electric charge is conserved canonically.
 
int m_SCE
 Flag indicating if strangeness is conserved canonically.
 
int m_CCE
 Flag indicating if charm is conserved canonically.
 
bool m_Banalyt
 Flag indicating whether the analytical calculation of baryon fugacity is used.
 

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 ideal HRG model in the canonical ensemble.

Calculation of particle densities proceeds as described in F. Becattini, U. Heinz, Z. Phys. C76, 269 (1997), https://arxiv.org/pdf/hep-ph/9702274.pdf

Particle number fluctuations are determined through the derivatives with respect to the fictitous fugacities, as is usually done in HRG, see e.g. https://arxiv.org/pdf/nucl-th/0404056.pdf

Definition at line 64 of file ThermalModelCanonical.h.

Constructor & Destructor Documentation

◆ ThermalModelCanonical()

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

Construct a new ThermalModelCanonical object.

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

Definition at line 23 of file ThermalModelCanonical.cpp.

◆ ~ThermalModelCanonical()

thermalfist::ThermalModelCanonical::~ThermalModelCanonical ( void )
virtual

Destroy the ThermalModelCanonical object.

Definition at line 44 of file ThermalModelCanonical.cpp.

Member Function Documentation

◆ CalculateEnergyDensity()

double thermalfist::ThermalModelCanonical::CalculateEnergyDensity ( )
virtual

Definition at line 755 of file ThermalModelCanonical.cpp.

◆ CalculateEnergyDensityDerivativeT()

virtual double thermalfist::ThermalModelCanonical::CalculateEnergyDensityDerivativeT ( )
inlinevirtual

Definition at line 218 of file ThermalModelCanonical.h.

◆ CalculateEntropyDensity()

double thermalfist::ThermalModelCanonical::CalculateEntropyDensity ( )
virtual

Definition at line 816 of file ThermalModelCanonical.cpp.

◆ CalculateFluctuations()

void thermalfist::ThermalModelCanonical::CalculateFluctuations ( )
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. Restricted to 2nd moments. TODO: Higher moments

Reimplemented from thermalfist::ThermalModelBase.

Definition at line 738 of file ThermalModelCanonical.cpp.

◆ CalculatePartitionFunctions()

void thermalfist::ThermalModelCanonical::CalculatePartitionFunctions ( double Vc = -1.)
virtual

Calculates all necessary canonical partition functions.

This corresponds to Eq. (8) in https://arxiv.org/pdf/hep-ph/9702274.pdf

Integrals are performed numerically using quadratures.

If multi-baryon states (light nuclei) are not included in the list, and quantum statistics for baryons is neglected, the integral over the baryon fugacity is performed analytically as described in (https://arxiv.org/pdf/nucl-th/0112021.pdf)[https://arxiv.org/pdf/nucl-th/0112021.pdf].

Parameters
Vc

TODO: Extra cross-check the factor 2, can be important for entropy density, irrelevant for everything else

Definition at line 297 of file ThermalModelCanonical.cpp.

◆ CalculatePressure()

double thermalfist::ThermalModelCanonical::CalculatePressure ( )
virtual

Definition at line 785 of file ThermalModelCanonical.cpp.

◆ CalculatePrimordialDensities()

void thermalfist::ThermalModelCanonical::CalculatePrimordialDensities ( )
virtual

Calculates the primordial densities of all species.

Implements thermalfist::ThermalModelBase.

Definition at line 156 of file ThermalModelCanonical.cpp.

◆ CalculateQuantumNumbersRange()

void thermalfist::ThermalModelCanonical::CalculateQuantumNumbersRange ( bool computeFluctuations = false)
virtual

Calculates the range of quantum numbers values for which it is necessary to compute the canonical partition functions.

Parameters
computeFluctuationsWhether it will be necessary to compute fluctuations. If that is the case the range is doubled for every quantum number.

Definition at line 53 of file ThermalModelCanonical.cpp.

◆ CalculateTwoParticleCorrelations()

void thermalfist::ThermalModelCanonical::CalculateTwoParticleCorrelations ( )
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 643 of file ThermalModelCanonical.cpp.

◆ ChangeTPS()

void thermalfist::ThermalModelCanonical::ChangeTPS ( ThermalParticleSystem * TPS)
virtual

Change the particle list.

Parameters
TPSA pointer to new particle list.

Reimplemented from thermalfist::ThermalModelBase.

Definition at line 49 of file ThermalModelCanonical.cpp.

◆ ConserveBaryonCharge()

virtual void thermalfist::ThermalModelCanonical::ConserveBaryonCharge ( bool conserve = true)
inlinevirtual

Specifies whether the baryon number is treated canonically.

By default the baryon number is treated canonically.

Parameters
conservetrue – canonically, false – grand-canonically

Definition at line 130 of file ThermalModelCanonical.h.

◆ ConserveCharm()

virtual void thermalfist::ThermalModelCanonical::ConserveCharm ( bool conserve = true)
inlinevirtual

Specifies whether the charm charge is treated canonically.

By default the charm charge is treated canonically.

Parameters
conservetrue – canonically, false – grand-canonically

Definition at line 157 of file ThermalModelCanonical.h.

◆ ConserveElectricCharge()

virtual void thermalfist::ThermalModelCanonical::ConserveElectricCharge ( bool conserve = true)
inlinevirtual

Specifies whether the electric charge is treated canonically.

By default the electric charge is treated canonically.

Parameters
conservetrue – canonically, false – grand-canonically

Definition at line 139 of file ThermalModelCanonical.h.

◆ ConserveStrangeness()

virtual void thermalfist::ThermalModelCanonical::ConserveStrangeness ( bool conserve = true)
inlinevirtual

Specifies whether the strangeness charge is treated canonically.

By default the strangeness charge is treated canonically.

Parameters
conservetrue – canonically, false – grand-canonically

Definition at line 148 of file ThermalModelCanonical.h.

◆ FixParameters()

void thermalfist::ThermalModelCanonical::FixParameters ( )
virtual

Method which actually implements ConstrainChemicalPotentials() (for backward compatibility).

Reimplemented from thermalfist::ThermalModelBase.

Definition at line 137 of file ThermalModelCanonical.cpp.

◆ FixParametersNoReset()

void thermalfist::ThermalModelCanonical::FixParametersNoReset ( )
virtual

Method which actually implements ConstrainChemicalPotentialsNoReset() (for backward compatibility).

Reimplemented from thermalfist::ThermalModelBase.

Definition at line 146 of file ThermalModelCanonical.cpp.

◆ GetGCEDensity()

double thermalfist::ThermalModelCanonical::GetGCEDensity ( int i) const
virtual

Density of particle species i in the grand-canonical ensemble.

Parameters
i0-based index of particle species
Returns
double The grand-canonical density

Definition at line 843 of file ThermalModelCanonical.cpp.

◆ IntegrationIterationsMultiplier()

int thermalfist::ThermalModelCanonical::IntegrationIterationsMultiplier ( ) const
inline

The multiplier of the number of iterations in the numerical integration.

Returns
The multiplier

Definition at line 174 of file ThermalModelCanonical.h.

◆ IsConservedChargeCanonical()

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

Definition at line 858 of file ThermalModelCanonical.cpp.

◆ IsParticleCanonical()

bool thermalfist::ThermalModelCanonical::IsParticleCanonical ( const ThermalParticle & part)
virtual

Determines whether the specified ThermalParticle is treat canonically or grand-canonically in the present setup.

This depends on whether the particle carries any of the exactly conserved charges.

Parameters
partThe particle species.
Returns
true Canonically.
false Grand-canonically.

Definition at line 848 of file ThermalModelCanonical.cpp.

◆ ParticleScalarDensity()

virtual double thermalfist::ThermalModelCanonical::ParticleScalarDensity ( int )
inlinevirtual

Returns the scalar density of a particle species.

This method is not implemented for the canonical ensemble model, hence it returns a constant value of 0.0.

Parameters
partThe particle species index.
Returns
double The scalar density, which is 0.0 in this implementation.

Definition at line 229 of file ThermalModelCanonical.h.

◆ ParticleScaledVariance()

double thermalfist::ThermalModelCanonical::ParticleScaledVariance ( int part)
virtual

Definition at line 591 of file ThermalModelCanonical.cpp.

◆ SetIntegrationIterationsMultiplier()

void thermalfist::ThermalModelCanonical::SetIntegrationIterationsMultiplier ( int multiplier)
inline

Assigns the multiplier of the number of iterations in the numerical integration.

The minimum value of multiplier is 1. Increase to improve the numerical accuracy of the canonical ensemble calculations.

Parameters
Themultiplier

Definition at line 183 of file ThermalModelCanonical.h.

◆ SetStatistics()

void thermalfist::ThermalModelCanonical::SetStatistics ( bool stats)
virtual

Set whether quantum statistics is used, 0 - Boltzmann, 1 - Quantum

Reimplemented from thermalfist::ThermalModelBase.

Definition at line 124 of file ThermalModelCanonical.cpp.

◆ ValidateCalculation()

void thermalfist::ThermalModelCanonical::ValidateCalculation ( )
virtual

Checks whether issues have occured during the calculation of particle densities in the CalculateDensities() method.

Reimplemented from thermalfist::ThermalModelBase.

Definition at line 215 of file ThermalModelCanonical.cpp.

Member Data Documentation

◆ m_Banalyt

bool thermalfist::ThermalModelCanonical::m_Banalyt
protected

Flag indicating whether the analytical calculation of baryon fugacity is used.

Definition at line 312 of file ThermalModelCanonical.h.

◆ m_BCE

int thermalfist::ThermalModelCanonical::m_BCE
protected

Flag indicating if baryon charge is conserved canonically.

Definition at line 304 of file ThermalModelCanonical.h.

◆ m_BMAX

int thermalfist::ThermalModelCanonical::m_BMAX
protected

Maximum baryon number.

Definition at line 292 of file ThermalModelCanonical.h.

◆ m_BMAX_list

int thermalfist::ThermalModelCanonical::m_BMAX_list
protected

Definition at line 296 of file ThermalModelCanonical.h.

◆ m_CCE

int thermalfist::ThermalModelCanonical::m_CCE
protected

Flag indicating if charm is conserved canonically.

Definition at line 307 of file ThermalModelCanonical.h.

◆ m_CMAX

int thermalfist::ThermalModelCanonical::m_CMAX
protected

Maximum charm.

Definition at line 295 of file ThermalModelCanonical.h.

◆ m_CMAX_list

int thermalfist::ThermalModelCanonical::m_CMAX_list
protected

Definition at line 296 of file ThermalModelCanonical.h.

◆ m_Corr

std::vector<double> thermalfist::ThermalModelCanonical::m_Corr
protected

A vector of chemical factors.

Chemical factors define the canonical corrections to the grand canonical thermodynamic functions.

Definition at line 271 of file ThermalModelCanonical.h.

◆ m_IntegrationIterationsMultiplier

int thermalfist::ThermalModelCanonical::m_IntegrationIterationsMultiplier
protected

A multiplier to increase the number of iterations during the numerical integration used to calculate the partition functions.

Set with SetIntegrationIterationsMultiplier()

Definition at line 290 of file ThermalModelCanonical.h.

◆ m_modelgce

ThermalModelIdeal* thermalfist::ThermalModelCanonical::m_modelgce
protected

Pointer to a ThermalModelIdeal object used for GCE calculations.

Definition at line 302 of file ThermalModelCanonical.h.

◆ m_MultExp

double thermalfist::ThermalModelCanonical::m_MultExp
protected

Exponential multiplier for canonical partition function calculations.

Definition at line 298 of file ThermalModelCanonical.h.

◆ m_MultExpBanalyt

double thermalfist::ThermalModelCanonical::m_MultExpBanalyt
protected

Exponential multiplier for analytical baryon fugacity calculations.

Definition at line 299 of file ThermalModelCanonical.h.

◆ m_PartialZ

std::vector<double> thermalfist::ThermalModelCanonical::m_PartialZ
protected

The computed canonical partition function.

The partition functions are computed up to a factor which does not depend of the quantum numbers since only ratios of the partition functions – the chemical factors – are of relevance.

Definition at line 282 of file ThermalModelCanonical.h.

◆ m_QCE

int thermalfist::ThermalModelCanonical::m_QCE
protected

Flag indicating if electric charge is conserved canonically.

Definition at line 305 of file ThermalModelCanonical.h.

◆ m_QMAX

int thermalfist::ThermalModelCanonical::m_QMAX
protected

Maximum electric charge.

Definition at line 293 of file ThermalModelCanonical.h.

◆ m_QMAX_list

int thermalfist::ThermalModelCanonical::m_QMAX_list
protected

Definition at line 296 of file ThermalModelCanonical.h.

◆ m_QNMap

std::map<QuantumNumbers, int> thermalfist::ThermalModelCanonical::m_QNMap
protected

Maps QuantumNumbers combinations to a 1-dimensional index.

Definition at line 262 of file ThermalModelCanonical.h.

◆ m_QNvec

std::vector<QuantumNumbers> thermalfist::ThermalModelCanonical::m_QNvec
protected

A set of QuantumNumbers combinations for which it is necessary to compute the canonical partition function.

Definition at line 255 of file ThermalModelCanonical.h.

◆ m_SCE

int thermalfist::ThermalModelCanonical::m_SCE
protected

Flag indicating if strangeness is conserved canonically.

Definition at line 306 of file ThermalModelCanonical.h.

◆ m_SMAX

int thermalfist::ThermalModelCanonical::m_SMAX
protected

Maximum strangeness.

Definition at line 294 of file ThermalModelCanonical.h.

◆ m_SMAX_list

int thermalfist::ThermalModelCanonical::m_SMAX_list
protected

Definition at line 296 of file ThermalModelCanonical.h.


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