32 m_TAG =
"ThermalModelCanonical";
35 m_InteractionModel =
Ideal;
63 for (
int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
87 if (computeFluctuations) {
100 std::cout <<
"BMAX = " <<
m_BMAX <<
"\tQMAX = " <<
m_QMAX <<
"\tSMAX = " <<
m_SMAX <<
"\tCMAX = " <<
m_CMAX << std::endl;
128 m_QuantumStats = stats;
129 for (
int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
130 m_TPS->Particle(i).UseStatistics(stats);
142 m_ConstrainMuB &= !
m_BCE;
143 m_ConstrainMuQ &= !
m_QCE;
144 m_ConstrainMuS &= !
m_SCE;
145 m_ConstrainMuC &= !
m_CCE;
151 m_ConstrainMuB &= !
m_BCE;
152 m_ConstrainMuQ &= !
m_QCE;
153 m_ConstrainMuS &= !
m_SCE;
154 m_ConstrainMuC &= !
m_CCE;
160 assert(m_IGFExtraConfig.MagneticField.B == 0.);
162 m_FluctuationsCalculated =
false;
169 m_Parameters.muB = 0.0;
170 m_Parameters.muQ = 0.0;
171 m_Parameters.muS = 0.0;
172 m_Parameters.muC = 0.0;
177 m_Parameters.muB = 0.0;
179 m_Parameters.muQ = 0.0;
181 m_Parameters.muS = 0.0;
183 m_Parameters.muC = 0.0;
190 for (
size_t i = 0; i < m_densities.size(); ++i) {
202 if (ind <
static_cast<int>(
m_Corr.size()))
205 cout <<
"ThermalModelCanonical::CalculatePrimordialDensities: Warning! No canonical partition function for this particle!" << endl;
212 if (ind <
static_cast<int>(
m_Corr.size()))
232 double totB = CalculateBaryonDensity() * m_Parameters.SVc;
233 if (fabs(m_Parameters.B - totB) > TOL) {
235 sprintf(cc,
"**WARNING** ThermalModelCanonical: Inaccurate calculation of total baryon number.\
239\n", m_Parameters.B, totB);
243 m_ValidityLog.append(cc);
245 m_LastCalculationSuccessFlag =
false;
251 double totQ = CalculateChargeDensity() * m_Parameters.SVc;
252 if (fabs(m_Parameters.Q - totQ) > TOL) {
253 sprintf(cc,
"**WARNING** ThermalModelCanonical: Inaccurate calculation of total electric charge.\
257\n", m_Parameters.Q, totQ);
261 m_ValidityLog.append(cc);
263 m_LastCalculationSuccessFlag =
false;
269 double totS = CalculateStrangenessDensity() * m_Parameters.SVc;
270 if (fabs(m_Parameters.S - totS) > TOL) {
271 sprintf(cc,
"**WARNING** ThermalModelCanonical: Inaccurate calculation of total strangeness.\
275\n", m_Parameters.S, totS);
279 m_ValidityLog.append(cc);
281 m_LastCalculationSuccessFlag =
false;
287 double totC = CalculateCharmDensity() * m_Parameters.SVc;
288 if (fabs(m_Parameters.C - totC) > TOL) {
289 sprintf(cc,
"**WARNING** ThermalModelCanonical: Inaccurate calculation of total charm.\
293\n", m_Parameters.C, totC);
297 m_ValidityLog.append(cc);
299 m_LastCalculationSuccessFlag =
false;
307 Vc = m_Parameters.SVc;
313 for (
int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
314 int i2 = m_TPS->PdgToId(-m_TPS->Particle(i).PdgId());
316 if (std::abs(m_Chem[i] - m_Chem[i2]) > 1.e-8) {
317 throw std::runtime_error(
"ThermalModelCanonical::CalculatePartitionFunctions: Partial chemical equilibrium canonical ensemble only supported if particle-antiparticle fugacities are symmetric!");
323 bool AllMuZero =
true;
324 for (
int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
336 for (
int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
342 throw std::invalid_argument(
"ThermalModelCanonical: neutral particle cannot have non-zero ce charges");
344 if (ind <
static_cast<int>(Nsx.size()))
351 if (ind <
static_cast<int>(Nsx.size())) {
354 Nsx[ind] += tdens * cosh(m_Chem[i] / m_Parameters.T);
355 Nsy[ind] += tdens * sinh(m_Chem[i] / m_Parameters.T);
367 if (ind <
static_cast<int>(Nsx.size())) {
370 Nsx[ind] += tdens * cosh(n * m_Chem[i] / m_Parameters.T);
371 Nsy[ind] += tdens * sinh(n * m_Chem[i] / m_Parameters.T);
383 for (
int i = 0; i < static_cast<int>(Nsx.size()); ++i) {
389 int nmax = max(3, (
int)sqrt(m_Parameters.B*m_Parameters.B + m_Parameters.Q*m_Parameters.Q + m_Parameters.S*m_Parameters.S + m_Parameters.C*m_Parameters.C));
390 if (m_Parameters.B == 0 && m_Parameters.Q == 0 && m_Parameters.S == 0 && m_Parameters.C == 0)
395 int nmaxB = max(4, m_Parameters.B);
396 int nmaxQ = max(4, m_Parameters.Q);
397 int nmaxS = max(4, m_Parameters.S);
398 int nmaxC = max(4, m_Parameters.C);
402 nmaxB = nmaxQ = nmaxS = nmaxC = nmax;
410 for (
size_t i = 0; i <
m_PartialZ.size(); ++i) {
418 int maxB = 2 * nmaxB;
422 for (
int iB = 0; iB < maxB; ++iB) {
424 vector<double> xlegB, wlegB;
427 double aB = iB * dphiB;
428 if (iB >= nmaxB) aB =
xMath::Pi() - (iB + 1) * dphiB;
429 double bB = aB + dphiB;
441 int maxS = 2 * nmaxS;
445 for (
int iS = 0; iS < maxS; ++iS) {
446 vector<double> xlegS, wlegS;
449 double aS = iS * dphiS;
450 if (iS >= nmaxS) aS =
xMath::Pi() - (iS + 1) * dphiS;
451 double bS = aS + dphiS;
466 for (
int iQ = 0; iQ < maxQ; ++iQ) {
467 vector<double> xlegQ, wlegQ;
470 double aQ = iQ * dphiQ;
471 double bQ = aQ + dphiQ;
482 int maxC = 2 * nmaxC;
486 for (
int iC = 0; iC < maxC; ++iC) {
487 vector<double> xlegC, wlegC;
489 double aC = iC * dphiC;
490 if (iC >= nmaxC) aC =
xMath::Pi() - (iC + 1) * dphiC;
491 double bC = aC + dphiC;
501 for (
size_t iBt = 0; iBt < xlegB.size(); ++iBt) {
502 for (
size_t iSt = 0; iSt < xlegS.size(); ++iSt) {
503 for (
size_t iQt = 0; iQt < xlegQ.size(); ++iQt) {
504 for (
size_t iCt = 0; iCt < xlegC.size(); ++iCt) {
506 double wx = 0., wy = 0., mx = 0., my = 0.;
507 for (
size_t i = 0; i <
m_PartialZ.size(); ++i) {
514 cosph[i] = cos(tS*xlegS[iSt] + tQ * xlegQ[iQt] + tC * xlegC[iCt]);
515 sinph[i] = sin(tS*xlegS[iSt] + tQ * xlegQ[iQt] + tC * xlegC[iCt]);
518 wx += Nsx[i] * cosph[i];
519 wy += Nsx[i] * sinph[i];
522 mx += Nsx[i] * (cosph[i] - 1.);
526 cosph[i] = cos(tB*xlegB[iBt] + tS * xlegS[iSt] + tQ * xlegQ[iQt] + tC * xlegC[iCt]);
527 mx += Nsx[i] * (cosph[i] - 1.);
530 sinph[i] = sin(tB*xlegB[iBt] + tS * xlegS[iSt] + tQ * xlegQ[iQt] + tC * xlegC[iCt]);
531 my += Nsy[i] * sinph[i];
540 wmod = sqrt(wx*wx + wy * wy);
541 warg = atan2(wy, wx);
544 for (
size_t iN = 0; iN <
m_PartialZ.size(); ++iN) {
545 int tBg = m_Parameters.B -
m_QNvec[iN].B;
546 int tQg = m_Parameters.Q -
m_QNvec[iN].Q;
547 int tSg = m_Parameters.S -
m_QNvec[iN].S;
548 int tCg = m_Parameters.C -
m_QNvec[iN].C;
552 m_PartialZ[iN] += wlegB[iBt] * wlegS[iSt] * wlegQ[iQt] * wlegC[iCt] *
553 cos(tBg * xlegB[iBt] + tSg * xlegS[iSt] + tQg * xlegQ[iQt] + tCg * xlegC[iCt] - tBg * warg) *
559 m_PartialZ[iN] += wlegB[iBt] * wlegS[iSt] * wlegQ[iQt] * wlegC[iCt] * cos(tBg*xlegB[iBt] + tSg * xlegS[iSt] + tQg * xlegQ[iQt] + tCg * xlegC[iCt]) * exp(mx);
561 m_PartialZ[iN] += wlegB[iBt] * wlegS[iSt] * wlegQ[iQt] * wlegC[iCt] * exp(mx)
562 * (cos(tBg*xlegB[iBt] + tSg * xlegS[iSt] + tQg * xlegQ[iQt] + tCg * xlegC[iCt]) * cos(my)
563 + sin(tBg*xlegB[iBt] + tSg * xlegS[iSt] + tQg * xlegQ[iQt] + tCg * xlegC[iCt]) * sin(my));
578 for (
size_t iN = 0; iN <
m_PartialZ.size(); ++iN) {
593 for (
size_t iN = 0; iN <
m_PartialZ.size(); ++iN) {
603 double ret1 = 0., ret2 = 0., ret3 = 0.;
606 return tpart.
ScaledVariance(m_Parameters, m_UseWidth, m_Chem[part]);
615 if (ind <
static_cast<int>(
m_Corr.size()) && ind2 <
static_cast<int>(
m_Corr.size()))
618 if (ind <
static_cast<int>(
m_Corr.size()))
622 double ret1num = 0., ret1zn = 0.;
628 if (ind <
static_cast<int>(
m_Corr.size())) {
629 ret1num +=
m_Corr[ind] * n * densityClusterN;
630 ret1zn +=
m_Corr[ind] * densityClusterN;
636 if (ind <
static_cast<int>(
m_Corr.size()) && ind2 <
static_cast<int>(
m_Corr.size()))
645 ret1 = ret1num / ret1zn;
646 ret2 = ret2 / ret1zn;
647 ret3 = -ret1zn * m_Parameters.SVc;
649 return ret1 + ret2 + ret3;
653 int NN = m_densities.size();
655 vector<double> yld(NN, 0);
656 vector<double> ret1num(NN, 0);
657 vector< vector<double> > ret2num(NN, vector<double>(NN, 0.));
659 for (
int i = 0; i < NN; ++i)
660 yld[i] = m_densities[i] * m_Parameters.SVc;
662 for (
int i = 0; i < NN; ++i) {
665 ret1num[i] = tpart.
ScaledVariance(m_Parameters, m_UseWidth, m_Chem[i]) * yld[i];
678 if (ind <
static_cast<int>(
m_Corr.size()))
679 ret1num[i] +=
m_Corr[ind] * n * densityClusterN * m_Parameters.SVc;
684 for (
int i = 0; i < NN; ++i) {
685 for (
int j = 0; j < NN; ++j) {
698 ret2num[i][j] = yld[i] * yld[j];
701 for (
int n1 = 1; n1 <= n1max; ++n1) {
702 for (
int n2 = 1; n2 <= n2max; ++n2) {
712 if (ind <
static_cast<int>(
m_Corr.size()))
713 ret2num[i][j] +=
m_Corr[ind] * densityClusterN1 * densityClusterN2 * m_Parameters.SVc * m_Parameters.SVc;
721 m_PrimCorrel.resize(NN);
722 for (
int i = 0; i < NN; ++i)
723 m_PrimCorrel[i].resize(NN);
724 m_TotalCorrel = m_PrimCorrel;
726 for (
int i = 0; i < NN; ++i) {
727 for (
int j = 0; j < NN; ++j) {
728 m_PrimCorrel[i][j] = 0.;
730 m_PrimCorrel[i][j] += ret1num[i] / yld[i];
731 m_PrimCorrel[i][j] += ret2num[i][j] / yld[i];
732 m_PrimCorrel[i][j] += -yld[j];
733 m_PrimCorrel[i][j] *= yld[i] / m_Parameters.SVc / m_Parameters.T;
735 m_PrimCorrel[i][j] = 0.;
739 for (
int i = 0; i < NN; ++i) {
740 m_wprim[i] = m_PrimCorrel[i][i];
741 if (m_densities[i] > 0.) m_wprim[i] *= m_Parameters.T / m_densities[i];
742 else m_wprim[i] = 1.;
743 if (m_wprim[i] != m_wprim[i]) m_wprim[i] = 1.;
746 m_TwoParticleCorrelationsCalculated =
true;
761 m_FluctuationsCalculated =
true;
763 for (
size_t i = 0; i < m_wprim.size(); ++i) {
775 for (
int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
785 if (ind <
static_cast<int>(
m_Corr.size()))
791 if (ind <
static_cast<int>(
m_Corr.size()))
805 for (
int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
815 if (ind <
static_cast<int>(
m_Corr.size()))
822 if (ind <
static_cast<int>(
m_Corr.size()))
837 ret += -m_Parameters.muB / m_Parameters.T * m_Parameters.B / m_Parameters.SVc;
839 ret += -m_Parameters.muB * CalculateBaryonDensity() / m_Parameters.T;
842 ret += -m_Parameters.muQ / m_Parameters.T * m_Parameters.Q / m_Parameters.SVc;
844 ret += -m_Parameters.muQ * CalculateChargeDensity() / m_Parameters.T;
847 ret += -m_Parameters.muS / m_Parameters.T * m_Parameters.S / m_Parameters.SVc;
849 ret += -m_Parameters.muS * CalculateStrangenessDensity() / m_Parameters.T;
852 ret += -m_Parameters.muC / m_Parameters.T * m_Parameters.C / m_Parameters.SVc;
854 ret += -m_Parameters.muC * CalculateCharmDensity() / m_Parameters.T;
887 void ThermalModelCanonical::PrepareModelGCE()
896 m_Parameters.muB = 0.0;
899 m_Parameters.muS = m_Parameters.muB / 3.;
902 m_Parameters.muQ = -m_Parameters.muB / 30.;
904 m_Parameters.muC = 0.;
907 m_Parameters.muB, m_Parameters.muQ, m_Parameters.muS, m_Parameters.muC,
908 static_cast<bool>(
m_BCE),
909 static_cast<bool>(
m_QCE),
910 static_cast<bool>(
m_SCE),
911 static_cast<bool>(
m_CCE));
931 void ThermalModelCanonical::CleanModelGCE()
map< string, double > params
virtual void CalculateSusceptibilityMatrix()
Calculates the conserved charges susceptibility matrix.
virtual void FixParametersNoReset()
Method which actually implements ConstrainChemicalPotentialsNoReset() (for backward compatibility).
virtual void ChangeTPS(ThermalParticleSystem *TPS)
Change the particle list.
virtual void CalculateParticleChargeCorrelationMatrix()
Calculates the matrix of correlators between primordial (and also final) particle numbers and conserv...
bool UsePartialChemicalEquilibrium()
Whether partial chemical equilibrium with additional chemical potentials is used.
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,...
virtual void FillChemicalPotentials()
Sets the chemical potentials of all particles.
virtual void FixParameters()
Method which actually implements ConstrainChemicalPotentials() (for backward compatibility).
void SetUseWidth(bool useWidth)
Sets whether finite resonance widths are used. Deprecated.
virtual void SetChemicalPotentials(const std::vector< double > &chem=std::vector< double >(0))
Sets the chemical potentials of all particles.
virtual void ValidateCalculation()
Checks whether issues have occured during the calculation of particle densities in the CalculateDensi...
virtual void CalculateProxySusceptibilityMatrix()
Calculates the susceptibility matrix of conserved charges proxies.
ThermalModelBase(ThermalParticleSystem *TPS, const ThermalModelParameters ¶ms=ThermalModelParameters())
Construct a new ThermalModelBase object.
virtual void CalculateDensities()
Calculates the primordial and total (after decays) densities of all species.
const ThermalModelParameters & Parameters() const
virtual void CalculateTwoParticleFluctuationsDecays()
Computes particle number correlations and fluctuations for all final-state particles which are marked...
std::map< QuantumNumbers, int > m_QNMap
Maps QuantumNumbers combinations to a 1-dimensional index.
double m_MultExpBanalyt
Exponential multiplier for analytical baryon fugacity calculations.
int m_IntegrationIterationsMultiplier
A multiplier to increase the number of iterations during the numerical integration used to calculate ...
int m_BMAX
Maximum baryon number.
virtual void SetStatistics(bool stats)
std::vector< double > m_PartialZ
The computed canonical partition function.
int m_QCE
Flag indicating if electric charge is conserved canonically.
void ChangeTPS(ThermalParticleSystem *TPS)
Change the particle list.
virtual void CalculateTwoParticleCorrelations()
Computes the fluctuations and correlations of the primordial particle numbers.
int m_CCE
Flag indicating if charm is conserved canonically.
virtual void CalculateFluctuations()
Computes the fluctuation observables.
bool m_PartialZCalculatedFlucts
double m_MultExp
Exponential multiplier for canonical partition function calculations.
virtual double ParticleScaledVariance(int part)
std::vector< QuantumNumbers > m_QNvec
A set of QuantumNumbers combinations for which it is necessary to compute the canonical partition fun...
std::vector< double > m_Corr
A vector of chemical factors.
virtual double GetGCEDensity(int i) const
Density of particle species i in the grand-canonical ensemble.
virtual bool IsParticleCanonical(const ThermalParticle &part)
Determines whether the specified ThermalParticle is treat canonically or grand-canonically in the pre...
virtual double CalculateEntropyDensity()
virtual void CalculatePartitionFunctions(double Vc=-1.)
Calculates all necessary canonical partition functions.
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 par...
virtual void FixParametersNoReset()
Method which actually implements ConstrainChemicalPotentialsNoReset() (for backward compatibility).
virtual double CalculatePressure()
int m_QMAX
Maximum electric charge.
int m_SMAX
Maximum strangeness.
int m_BCE
Flag indicating if baryon charge is conserved canonically.
ThermalModelIdeal * m_modelgce
Pointer to a ThermalModelIdeal object used for GCE calculations.
virtual void ValidateCalculation()
Checks whether issues have occured during the calculation of particle densities in the CalculateDensi...
virtual bool IsConservedChargeCanonical(ConservedCharge::Name charge) const
bool m_PartialZCalculated
bool m_Banalyt
Flag indicating whether the analytical calculation of baryon fugacity is used.
int m_SCE
Flag indicating if strangeness is conserved canonically.
virtual double CalculateEnergyDensity()
virtual void FixParameters()
Method which actually implements ConstrainChemicalPotentials() (for backward compatibility).
virtual void CalculatePrimordialDensities()
Calculates the primordial densities of all species.
ThermalModelCanonical(ThermalParticleSystem *TPS, const ThermalModelParameters ¶ms=ThermalModelParameters())
Construct a new ThermalModelCanonical object.
Class implementing the Ideal HRG model.
Class containing all information about a particle specie.
int BaryonCharge() const
Particle's baryon number.
int Statistics() const
Particle's statistics.
int Strangeness() const
Particle's strangeness.
void SetCalculationType(IdealGasFunctions::QStatsCalculationType type)
Sets the CalculationType() method to evaluate quantum statistics.
double Density(const ThermalModelParameters ¶ms, IdealGasFunctions::Quantity type=IdealGasFunctions::ParticleDensity, bool useWidth=0, double mu=0.) const
Computes a specified ideal gas thermodynamic function.
IdealGasFunctions::QStatsCalculationType CalculationType() const
Method to evaluate quantum statistics.
int ElectricCharge() const
Particle's electric charge.
int Charm() const
Particle's charm.
int ClusterExpansionOrder() const
Number of terms in the cluster expansion method.
double DensityCluster(int n, const ThermalModelParameters ¶ms, IdealGasFunctions::Quantity type=IdealGasFunctions::ParticleDensity, bool useWidth=0, double mu=0.) const
double ScaledVariance(const ThermalModelParameters ¶ms, bool useWidth=0, double mu=0.) const
Computes the scaled variance of particle number fluctuations in the ideal gas. Computes the scaled va...
Class containing the particle list.
void GetCoefsIntegrateLegendre10(double a, double b, std::vector< double > *x, std::vector< double > *w)
constexpr double Pi()
Pi constant.
double BesselIexp(int n, double x)
Modified Bessel function I_n(x), divided by exponential factor.
The main namespace where all classes and functions of the Thermal-FIST library reside.
Name
Set of all conserved charges considered.
@ BaryonCharge
Baryon number.
@ StrangenessCharge
Strangeness.
@ ElectricCharge
Electric charge.
Struct containing a set of quantum numbers: Baryon number, electric charge, strangeness,...
Structure containing all thermal parameters of the model.
Contains some extra mathematical functions used in the code.