Thermal-FIST 1.6
Package for hadron resonance gas model applications
Loading...
Searching...
No Matches
ThermalModelBase.h
Go to the documentation of this file.
1/*
2 * Thermal-FIST package
3 *
4 * Copyright (c) 2014-2019 Volodymyr Vovchenko
5 *
6 * GNU General Public License (GPLv3 or later)
7 */
8#ifndef THERMALMODELBASE_H
9#define THERMALMODELBASE_H
10
11#include <string>
12
14#include "HRGBase/xMath.h"
15#include "HRGBase/Broyden.h"
16
17
18namespace thermalfist {
19
29 {
30 public:
36 GCE = 0,
37 CE = 1,
38 SCE = 2,
39 CCE = 3
40 };
41
54
62
63 virtual ~ThermalModelBase(void) { }
64
66 int ComponentsNumber() const { return static_cast<int>(m_densities.size()); }
67
79 virtual void FillVirial(const std::vector<double> & ri = std::vector<double>(0));
80
82 bool UseWidth() const { return m_UseWidth; }
83
92 void SetUseWidth(bool useWidth);
93
100
102
107 void SetNormBratio(bool normBratio);
108 bool NormBratio() const { return m_NormBratio; }
110
112 void SetOMP(bool openMP) { m_useOpenMP = openMP; }
113
115
120 virtual void SetParameters(const ThermalModelParameters& params);
121 const ThermalModelParameters& Parameters() const { return m_Parameters; }
123
128 void UpdateParameters() { SetParameters(m_Parameters); }
129
135 virtual void SetTemperature(double T);
136
142 virtual void SetBaryonChemicalPotential(double muB);
143
149 virtual void SetElectricChemicalPotential(double muQ);
150
156 virtual void SetStrangenessChemicalPotential(double muS);
157
163 virtual void SetCharmChemicalPotential(double muC);
164
170 virtual void SetGammaq(double gammaq);
171
177 virtual void SetGammaS(double gammaS);
178
184 virtual void SetGammaC(double gammaC);
185
191 virtual void SetBaryonCharge(int B);
192
198 virtual void SetElectricCharge(int Q);
199
205 virtual void SetStrangeness(int S);
206
212 virtual void SetCharm(int C);
213
219 virtual void SetRadius(double /*rad*/) { }
220
227 virtual void SetRadius(int /*i*/, double /*rad*/) { }
228
239 virtual void SetVirial(int /*i*/, int /*j*/, double /*b*/) { }
240
242 virtual void SetRepulsion(int i, int j, double b) { SetVirial(i,j,b); }
243
251 virtual void SetAttraction(int /*i*/, int /*j*/, double /*a*/) { }
252
254
256 virtual void DisableMesonMesonVirial();
259
262 virtual void DisableMesonMesonAttraction();
263
265
267 virtual void DisableMesonBaryonVirial();
270
273 virtual void DisableMesonBaryonAttraction();
274
276
278 virtual void DisableBaryonBaryonVirial();
281
284 virtual void DisableBaryonBaryonAttraction();
285
287
289 virtual void DisableBaryonAntiBaryonVirial();
292
296
304 virtual void ReadInteractionParameters(const std::string & /*filename*/) { }
305
313 virtual void WriteInteractionParameters(const std::string & /*filename*/) { }
314
320 virtual void ChangeTPS(ThermalParticleSystem *TPS);
321
323
333 virtual double VirialCoefficient(int /*i*/, int /*j*/) const { return 0.; }
334 double RepulsionCoefficient(int i, int j) const { return VirialCoefficient(i,j); }
336
344 virtual double AttractionCoefficient(int /*i*/, int /*j*/) const { return 0.; }
345
348 bool QuantumStatistics() const { return m_QuantumStats; }
349
352 virtual void SetStatistics(bool stats);
353
360 virtual void SetCalculationType(IdealGasFunctions::QStatsCalculationType type) { m_TPS->SetCalculationType(type); }
361
372 virtual void SetClusterExpansionOrder(int order) { m_TPS->SetClusterExpansionOrder(order); }
373
383 void SetResonanceWidthShape(ThermalParticle::ResonanceWidthShape shape) { m_TPS->SetResonanceWidthShape(shape); }
384
392 void SetResonanceWidthIntegrationType(ThermalParticle::ResonanceWidthIntegration type);// { m_TPS->SetResonanceWidthIntegrationType(type); }
393
401 virtual void FillChemicalPotentials();
402
411 virtual void SetChemicalPotentials(const std::vector<double> & chem = std::vector<double>(0));
412
418 const std::vector<double>& ChemicalPotentials() const { return m_Chem; }
419
426 double ChemicalPotential(int i) const;
427
435 virtual void SetChemicalPotential(int i, double chem);
436
445 virtual double FullIdealChemicalPotential(int i) const;
446
447
450 bool ConstrainMuB() const { return m_ConstrainMuB; }
451
454 void ConstrainMuB(bool constrain) { m_ConstrainMuB = constrain; }
455
458 bool ConstrainMuQ() const { return m_ConstrainMuQ; }
459
462 void ConstrainMuQ(bool constrain) { m_ConstrainMuQ = constrain; }
463
466 bool ConstrainMuS() const { return m_ConstrainMuS; }
467
470 void ConstrainMuS(bool constrain) { m_ConstrainMuS = constrain; }
471
474 bool ConstrainMuC() const { return m_ConstrainMuC; }
475
478 void ConstrainMuC(bool constrain) { m_ConstrainMuC = constrain; }
479
481 void UsePartialChemicalEquilibrium(bool usePCE) { m_PCE = usePCE; }
482
484 bool UsePartialChemicalEquilibrium() { return m_PCE; }
485
487
494 void SetSoverB(double SB) { m_SBgoal = SB; }
495 double SoverB() const { return m_SBgoal; }
497
499
506 void SetQoverB(double QB) { m_QBgoal = QB; }
507 double QoverB() const { return m_QBgoal; }
509
515 void SetVolume(double Volume) { m_Volume = Volume; m_Parameters.V = Volume; }
517 double Volume() const { return m_Parameters.V; }
518
527 void SetVolumeRadius(double R) { m_Volume = 4. / 3.*xMath::Pi() * R * R * R; m_Parameters.V = m_Volume; }
528
530 double CanonicalVolume() const { return m_Parameters.SVc; }
536
537 void SetCanonicalVolume(double Volume) { m_Parameters.SVc = Volume; }
538
547 void SetCanonicalVolumeRadius(double Rc) { m_Parameters.SVc = 4. / 3. * xMath::Pi() * Rc * Rc * Rc; }
548
549
550 //double StrangenessCanonicalVolume() const { return CanonicalVolume(); }
551 //void SetStrangenessCanonicalVolume(double Volume) { SetCanonicalVolume(Volume); }
552 //void SetStrangenessCanonicalVolumeRadius(double Radius) { SetCanonicalVolumeRadius(Radius); }
553
554 // Same as FixParameters but with a more clear name on what is actually does
570 void ConstrainChemicalPotentials(bool resetInitialValues = true);
571
572
578 virtual void FixParameters();
579
585 virtual void FixParametersNoReset();
586
614 virtual bool SolveChemicalPotentials(double totB = 0., double totQ = 0., double totS = 0., double totC = 0.,
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);
617
644 virtual bool FixChemicalPotentialsThroughDensities(double rhoB = 0., double rhoQ = 0., double rhoS = 0., double rhoC = 0.,
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);
647
654
660 virtual void CalculateDensities();
661
668 virtual void ValidateCalculation();
669
676 std::string ValidityCheckLog() const { return m_ValidityLog; }
677
684 virtual void CalculateDensitiesGCE() { CalculateDensities(); m_GCECalculated = true; }
685
693 virtual void CalculateFeeddown();
694
703 virtual void CalculateFluctuations();
704
709 virtual void CalculateTemperatureDerivatives();
710
721
740
746 virtual double TwoParticleSusceptibilityPrimordial(int i, int j) const;
747
753 virtual double TwoParticleSusceptibilityPrimordialByPdg(long long id1, long long id2);
754
760 virtual double TwoParticleSusceptibilityTemperatureDerivativePrimordial(int i, int j) const;
761
767 virtual double TwoParticleSusceptibilityTemperatureDerivativePrimordialByPdg(long long id1, long long id2);
768
774 virtual double NetParticleSusceptibilityPrimordialByPdg(long long id1, long long id2);
775
781 virtual double TwoParticleSusceptibilityFinal(int i, int j) const;
782
788 virtual double TwoParticleSusceptibilityFinalByPdg(long long id1, long long id2);
789
795 virtual double NetParticleSusceptibilityFinalByPdg(long long id1, long long id2);
796
802 virtual double PrimordialParticleChargeSusceptibility(int i, ConservedCharge::Name chg) const;
803
810
817
818
824 virtual double FinalParticleChargeSusceptibility(int i, ConservedCharge::Name chg) const;
825
831 virtual double FinalParticleChargeSusceptibilityByPdg(long long id1, ConservedCharge::Name chg);
832
838 virtual double FinalNetParticleChargeSusceptibilityByPdg(long long id1, ConservedCharge::Name chg);
839
840
851
867 virtual void CalculateSusceptibilityMatrix();
868
884
894
914 virtual std::vector<double> CalculateChargeFluctuations(const std::vector<double> &chgs, int order = 4, bool dimensionfull = false);
915
933 virtual std::vector<double> CalculateGeneralizedSusceptibilities(const std::vector<std::vector<double>> &chgs);
934
942 virtual void SetMultipleSolutionsMode(bool search) { }
943
950 virtual bool UseMultipleSolutionsMode() const { return false; }
951
952
953 //virtual double GetParticlePrimordialDensity(unsigned int);
954 //virtual double GetParticleTotalDensity(unsigned int);
955
956 // Equation of state etc.
957
959 double Pressure() { return CalculatePressure(); }
960
962 double EnergyDensity() { return CalculateEnergyDensity(); }
963
965 double EntropyDensity() { return CalculateEntropyDensity(); }
966
968 double HadronDensity() { return CalculateHadronDensity(); }
969
971 double BaryonDensity() { return CalculateBaryonDensity(); }
972
974 double ElectricChargeDensity() { return CalculateChargeDensity(); }
975
977 double StrangenessDensity() { return CalculateStrangenessDensity(); }
978
980 double CharmDensity() { return CalculateCharmDensity(); }
981
983 double AbsoluteBaryonDensity() { return CalculateAbsoluteBaryonDensity(); }
984
986 double AbsoluteElectricChargeDensity() { return CalculateAbsoluteChargeDensity(); }
987
989 double AbsoluteStrangenessDensity() { return CalculateAbsoluteStrangenessDensity(); }
990
992 double AbsoluteCharmDensity() { return CalculateAbsoluteCharmDensity(); }
993
995 double HeatCapacityMu() { return CalculateHeatCapacityMu(); }
996
998 [[deprecated("Use HeatCapacityMu() instead")]]
999 double SpecificHeatChem() { return CalculateHeatCapacityMu(); }
1000
1014 double cs2(bool rhoBconst = true, bool rhoQconst = true, bool rhoSconst = true, bool rhoCconst = true) {
1015 return CalculateAdiabaticSpeedOfSoundSquared(rhoBconst, rhoQconst, rhoSconst, rhoCconst);
1016 }
1017
1031 double cT2(bool rhoBconst = true, bool rhoQconst = true, bool rhoSconst = true, bool rhoCconst = true) {
1032 return CalculateIsothermalSpeedOfSoundSquared(rhoBconst, rhoQconst, rhoSconst, rhoCconst);
1033 }
1034
1042 double HeatCapacity(bool rhoBconst = true, bool rhoQconst = true, bool rhoSconst = true, bool rhoCconst = true) {
1043 return CalculateHeatCapacity(rhoBconst, rhoQconst, rhoSconst, rhoCconst);
1044 }
1045
1067 std::vector<std::vector<double>> PressureHessian();
1068
1077 bool IsThermodynamicallyStable();
1078
1080
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);
1102
1104 virtual double CalculateArbitraryChargeDensity();
1105
1107 virtual double CalculateBaryonMatterEntropyDensity() { return 0.; }
1108
1110 virtual double CalculateMesonMatterEntropyDensity() { return 0.; }
1111
1113 virtual double ParticleScaledVariance(int) { return 1.; }
1114
1116 virtual double ParticleSkewness(int) { return 1.; }
1117
1119 virtual double ParticleKurtosis(int) { return 1.; }
1120
1123 virtual double CalculateEigenvolumeFraction() { return 0.; }
1124
1126 virtual double ParticleScalarDensity(int i) = 0;
1127
1128 virtual double GetMaxDiff() const { return m_MaxDiff; }
1129
1133 virtual bool IsLastSolutionOK() const { return m_LastCalculationSuccessFlag; }
1134
1139 double GetdndT(int i) const;
1140
1146 double GetDensity(long long PDGID, int feeddown) { return GetDensity(PDGID, static_cast<Feeddown::Type>(feeddown)); }
1147
1156 double GetDensity(long long PDGID, Feeddown::Type feeddown);
1157
1166 double GetYield(long long PDGID, Feeddown::Type feeddown) { return GetDensity(PDGID, feeddown) * Volume(); }
1167
1168 std::vector<double> GetIdealGasDensities() const;
1169
1172 ThermalParticleSystem* TPS() { return m_TPS; }
1173
1178 const std::vector<double>& Densities() const { return m_densities; }
1179 std::vector<double>& Densities() { return m_densities; }
1180
1187 const std::vector<double>& TotalDensities() const { return m_densitiestotal; }
1188
1192 const std::vector< std::vector<double> >& AllDensities() const { return m_densitiesbyfeeddown; }
1193
1199 void SetTAG(const std::string & tag) { m_TAG = tag; }
1200
1202 const std::string& TAG() const { return m_TAG; }
1203
1205 int PdgToId(long long pdgid) const { return m_TPS->PdgToId(pdgid); }
1206
1208 virtual void ResetCalculatedFlags();
1209
1212 bool IsCalculated() const { return m_Calculated; }
1213
1216 bool IsTwoParticleCorrelationsCalculated() const { return m_TwoParticleCorrelationsCalculated; }
1217
1220 bool IsFluctuationsCalculated() const { return m_FluctuationsCalculated; }
1221
1224 bool IsSusceptibilitiesCalculated() const { return m_SusceptibilitiesCalculated; }
1225
1228 bool IsTemperatureDerivativesCalculated() const { return m_TemperatureDerivativesCalculated; }
1229
1232 bool IsGCECalculated() const { return m_GCECalculated; }
1233
1241 double ScaledVariancePrimordial(int id) const { return (id >= 0 && id < static_cast<int>(m_wprim.size())) ? m_wprim[id] : 1.; }
1242
1253 double ScaledVarianceTotal(int id) const { return (id >= 0 && id < static_cast<int>(m_wtot.size())) ? m_wtot[id] : 1.; }
1254
1262 double SkewnessPrimordial(int id) const { return (id >= 0 && id < static_cast<int>(m_skewprim.size())) ? m_skewprim[id] : 1.; }
1263
1274 double SkewnessTotal(int id) const { return (id >= 0 && id < static_cast<int>(m_skewtot.size())) ? m_skewtot[id] : 1.; }
1275
1283 double KurtosisPrimordial(int id) const { return (id >= 0 && id < static_cast<int>(m_kurtprim.size())) ? m_kurtprim[id] : 1.; }
1284
1295 double KurtosisTotal(int id) const { return (id >= 0 && id < static_cast<int>(m_kurttot.size())) ? m_kurttot[id] : 1.; }
1296
1305 double ConservedChargeDensity(ConservedCharge::Name chg);
1306
1315 double AbsoluteConservedChargeDensity(ConservedCharge::Name chg);
1316
1325 double ConservedChargeDensitydT(ConservedCharge::Name chg);
1326
1336 double Susc(ConservedCharge::Name i, ConservedCharge::Name j) const;// { return m_Susc[i][j]; }
1337
1347 double dSuscdT(ConservedCharge::Name i, ConservedCharge::Name j) const;
1348
1364 double ProxySusc(ConservedCharge::Name i, ConservedCharge::Name j) const;// { return m_ProxySusc[i][j]; }
1365
1376 double ChargedMultiplicity(int type = 0);
1377
1388 double ChargedScaledVariance(int type = 0);
1389
1397 double ChargedMultiplicityFinal(int type = 0);
1398
1406 double ChargedScaledVarianceFinal(int type = 0);
1407
1413 ThermalModelEnsemble Ensemble() { return m_Ensemble; }
1414
1415
1419 virtual bool IsConservedChargeCanonical(ConservedCharge::Name charge) const { return 0; }
1420
1426 ThermalModelInteraction InteractionModel() { return m_InteractionModel; }
1427
1428
1432 void SetDensityModelForParticleSpecies(int i, GeneralizedDensity* density_model = NULL);
1433
1437 void SetDensityModelForParticleSpeciesByPdg(long long PDGID, GeneralizedDensity* density_model = NULL);
1438
1442 void ClearDensityModels();
1443
1444 const IdealGasFunctions::IdealGasFunctionsExtraConfig& GetIdealGasFunctionsExtraConfig() const { return m_IGFExtraConfig; }
1445
1446
1454 void SetMagneticField(double B = 0.0, int lmax = -1);
1460 void RecomputeThresholdsDueToMagneticField();
1461
1468 std::vector<double> PartialPressures();
1469
1471 void ClearMagneticField();
1472
1473 protected:
1474 ThermalModelParameters m_Parameters;
1475 ThermalParticleSystem* m_TPS;
1476
1478 IdealGasFunctions::IdealGasFunctionsExtraConfig m_IGFExtraConfig;
1479
1480 bool m_LastCalculationSuccessFlag;
1481 double m_MaxDiff;
1482
1483 bool m_Calculated;
1484 bool m_FeeddownCalculated;
1485 bool m_TwoParticleCorrelationsCalculated;
1486 bool m_FluctuationsCalculated;
1487 bool m_SusceptibilitiesCalculated;
1488 bool m_GCECalculated;
1489 bool m_UseWidth;
1490 bool m_NormBratio;
1491 bool m_QuantumStats;
1492 double m_QBgoal;
1493 double m_SBgoal;
1494 double m_Volume;
1495
1496 bool m_ConstrainMuB;
1497 bool m_ConstrainMuQ;
1498 bool m_ConstrainMuS;
1499 bool m_ConstrainMuC;
1500
1501 bool m_PCE;
1502
1503 bool m_useOpenMP;
1504
1505 std::vector<double> m_densities;
1506 std::vector<double> m_densitiestotal;
1507 //std::vector<double> m_densitiestotalweak;
1508 std::vector< std::vector<double> > m_densitiesbyfeeddown;
1509 std::vector<double> m_Chem;
1510
1511 // Scaled variance
1512 std::vector<double> m_wprim;
1513 std::vector<double> m_wtot;
1514
1515 // Skewness
1516 std::vector<double> m_skewprim;
1517 std::vector<double> m_skewtot;
1518
1519 // Kurtosis
1520 std::vector<double> m_kurtprim;
1521 std::vector<double> m_kurttot;
1522
1523 // 2nd order correlations of primordial and total numbers
1524 std::vector< std::vector<double> > m_PrimCorrel;
1525 std::vector< std::vector<double> > m_TotalCorrel;
1526
1527 // Shifted chemical potential derivatives
1528 std::vector< std::vector<double> > m_dmusdmu;
1529
1530
1531 // Particle number-conserved charge correlators
1532 std::vector< std::vector<double> > m_PrimChargesCorrel;
1533 std::vector< std::vector<double> > m_FinalChargesCorrel;
1534
1535 // Conserved charges susceptibility matrix
1536 std::vector< std::vector<double> > m_Susc;
1537
1538 // Temperature derivative of the conserved charges susceptibility matrix (GeV^-1)
1539 std::vector< std::vector<double> > m_dSuscdT;
1540
1541 // Susceptibility matrix of net-p, net-Q, and net-K
1542 std::vector< std::vector<double> > m_ProxySusc;
1543
1544 // Temperature derivatives
1545 bool m_TemperatureDerivativesCalculated;
1546 // Densities
1547 std::vector<double> m_dndT;
1548 // (Shifted) chemical potentials
1549 std::vector<double> m_dmusdT;
1550 // Primordial mixed susceptibilities
1551 std::vector< std::vector<double> > m_PrimChi2sdT; // GeV^-1
1552
1553 // Cumulants of arbitrary charge calculation
1554 //std::vector< std::vector<double> > m_chi;
1555
1556 // Contains log of possible errors when checking the calculation
1557 std::string m_ValidityLog;
1558
1559 double m_wnSum;
1560
1561 std::string m_TAG;
1562
1563 ThermalModelEnsemble m_Ensemble;
1564 ThermalModelInteraction m_InteractionModel;
1565
1567 virtual double MuShift(int /*id*/) const { return 0.; }
1568
1569
1570 private:
1571 void ResetChemicalPotentials();
1572
1573 double GetDensity(long long PDGID, const std::vector<double> *dens);
1574
1575 class BroydenEquationsChem : public BroydenEquations
1576 {
1577 public:
1578 BroydenEquationsChem(ThermalModelBase *model) : BroydenEquations(), m_THM(model) { m_N = 2; }
1579 std::vector<double> Equations(const std::vector<double> &x);
1580 private:
1581 ThermalModelBase *m_THM;
1582 };
1583
1584 class BroydenJacobianChem : public BroydenJacobian
1585 {
1586 public:
1587 BroydenJacobianChem(ThermalModelBase *model) : BroydenJacobian(), m_THM(model) { }
1588 std::vector<double> Jacobian(const std::vector<double> &x);
1589 private:
1590 ThermalModelBase *m_THM;
1591 };
1592
1593 class BroydenChem : public Broyden
1594 {
1595 public:
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);
1599 private:
1600 ThermalModelBase *m_THM;
1601 };
1602
1603
1604 class BroydenEquationsChemTotals : public BroydenEquations
1605 {
1606 public:
1607 BroydenEquationsChemTotals(const std::vector<int> & vConstr, const std::vector<int> & vType, const std::vector<double> & vTotals, ThermalModelBase *model);// : BroydenEquations(), m_THM(model) { m_N = 3; }
1608 std::vector<double> Equations(const std::vector<double> &x);
1609 private:
1610 std::vector<int> m_Constr;
1611 std::vector<int> m_Type;
1612 std::vector<double> m_Totals;
1613 ThermalModelBase *m_THM;
1614 };
1615
1616 class BroydenJacobianChemTotals : public BroydenJacobian
1617 {
1618 public:
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);
1621 private:
1622 std::vector<int> m_Constr;
1623 std::vector<int> m_Type;
1624 std::vector<double> m_Totals;
1625 ThermalModelBase *m_THM;
1626 };
1627 };
1628
1629} // namespace thermalfist
1630
1631#endif // THERMALMODELBASE_H
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.
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 double TwoParticleSusceptibilityPrimordialByPdg(long long id1, long long id2)
Returns the computed primordial particle number (cross-)susceptibility for particles with pdg codes ...
void SetCanonicalVolumeRadius(double Rc)
Set the canonical correlation system radius.
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 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 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).
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 &params)
The thermal parameters.
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...
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 double VirialCoefficient(int, int) const
Excluded volume coefficient .
virtual void CalculateProxySusceptibilityMatrix()
Calculates the susceptibility matrix of conserved charges proxies.
ThermalModelBase(ThermalParticleSystem *TPS, const ThermalModelParameters &params=ThermalModelParameters())
Construct a new ThermalModelBase object.
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 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.
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.
ResonanceWidthIntegration
Treatment of finite resonance widths.
ResonanceWidthShape
Relativistic vs non-relativistic Breit-Wigner shape.
Class containing the particle list.
QStatsCalculationType
Identifies whether quantum statistics are to be computed using the cluster expansion or numerical int...
constexpr double Pi()
Pi constant.
Definition xMath.h:23
The main namespace where all classes and functions of the Thermal-FIST library reside.
Definition CosmicEoS.h:9
Name
Set of all conserved charges considered.
Structure containing all thermal parameters of the model.
Contains some extra mathematical functions used in the code.