Thermal-FIST 1.6
Package for hadron resonance gas model applications
Loading...
Searching...
No Matches
ThermalModelEVCrosstermsLegacy.h
Go to the documentation of this file.
1/*
2 * Thermal-FIST package
3 *
4 * Copyright (c) 2015-2018 Volodymyr Vovchenko
5 *
6 * GNU General Public License (GPLv3 or later)
7 */
8#ifndef THERMALMODELEVCROSSTERMSLEGACY_H
9#define THERMALMODELEVCROSSTERMSLEGACY_H
10
12
13namespace thermalfist {
14
39 {
40 public:
48
54
55 // Override functions begin
56
57 virtual void FillVirial(const std::vector<double> & ri = std::vector<double>(0));
58
59 virtual void ReadInteractionParameters(const std::string &filename);
60
61 virtual void WriteInteractionParameters(const std::string &filename);
62
63 void SetRadius(double rad);
64
65 double VirialCoefficient(int i, int j) const;
66
67 void SetVirial(int i, int j, double b);
68
69 virtual void ChangeTPS(ThermalParticleSystem *TPS);
70
71 virtual void CalculatePrimordialDensities();
72
74
76
77 virtual std::vector<double> CalculateChargeFluctuations(const std::vector<double> &chgs, int order = 4, bool dimensionfull = false);
78
79 virtual double CalculatePressure();
80
81 virtual double CalculateEnergyDensity();
82
83 virtual double CalculateEntropyDensity();
84
85 // Dummy
86 virtual double CalculateBaryonMatterEntropyDensity() { return 0.; }
87 virtual double CalculateMesonMatterEntropyDensity() { return 0.; }
88
89 virtual double CalculateEnergyDensityDerivativeT() { exit(1); return 0.; } // Not implemented
90
91 virtual double CalculateEntropyDensityDerivativeT() { exit(1); return 0.; } // Not implemented
92
93 // TODO properly with excluded volume
94 virtual double ParticleScaledVariance(int /*part*/) { return 1.; }
95
96 // TODO properly with excluded volume
97 virtual double ParticleSkewness(int part) { return m_skewprim[part]; }
98
99 // TODO properly with excluded volume
100 virtual double ParticleKurtosis(int part) { return m_kurtprim[part]; }
101
102 // TODO properly with excluded volume
103 virtual double ParticleScalarDensity(int /*part*/) { return 0.; }
104
105 // Override functions end
106
107
108 const std::vector< std::vector<int> >& EVComponentIndices() const { return m_EVComponentIndices; }
109 virtual double DeltaMu(int i) const { return MuShift(i); }
110 const std::vector< std::vector<double> >& VirialMatrix() const { return m_Virial; }
111
112 protected:
120 virtual void SolvePressure(bool resetPartials = true); // Using Broyden's method
121
122 void SolvePressureIter(); // Using iteration method
123
125
127
138 virtual void SolveDiagonal();
139
148 double PartialPressureDiagonal(int i, double P);
149
157 double PressureDiagonalTotal(double P);
158
168 virtual double DensityId(int i, const std::vector<double>& pstars = std::vector<double>());
169
179 virtual double Pressure(int i, const std::vector<double>& pstars = std::vector<double>());
180
191 double ScaledVarianceId(int ind);
192
203 virtual double MuShift(int i) const;
204
205 std::vector<double> m_densitiesid;
206 std::vector<double> m_Ps;
207 std::vector< std::vector<double> > m_Virial;
208 double m_Pressure;
210
211
212 std::vector<int> m_MapToEVComponent;
213 std::vector<int> m_MapFromEVComponent;
214 std::vector< std::vector<int> > m_EVComponentIndices;
215
216 private:
217 class BroydenEquationsCRS : public BroydenEquations
218 {
219 public:
220 BroydenEquationsCRS(ThermalModelEVCrosstermsLegacy*model) : BroydenEquations(), m_THM(model) { m_N = model->Densities().size(); }
221 std::vector<double> Equations(const std::vector<double> &x);
222 private:
224 };
225
226 class BroydenEquationsCRSDEV : public BroydenEquations
227 {
228 public:
229 BroydenEquationsCRSDEV(ThermalModelEVCrosstermsLegacy*model) : BroydenEquations(), m_THM(model) { m_N = 1; }
230 std::vector<double> Equations(const std::vector<double> &x);
231 private:
233 };
234
235 class BroydenJacobianCRS : public BroydenJacobian
236 {
237 public:
238 BroydenJacobianCRS(ThermalModelEVCrosstermsLegacy*model) : BroydenJacobian(), m_THM(model) { }
239 std::vector<double> Jacobian(const std::vector<double> &x);
240 private:
242 };
243
244 class BroydenSolutionCriteriumCRS : public Broyden::BroydenSolutionCriterium
245 {
246 public:
247 BroydenSolutionCriteriumCRS(ThermalModelEVCrosstermsLegacy*model, double relative_error = Broyden::TOL) : Broyden::BroydenSolutionCriterium(relative_error), m_THM(model) { }
248 virtual bool IsSolved(const std::vector<double>& x, const std::vector<double>& f, const std::vector<double>& xdelta = std::vector<double>()) const;
249 protected:
251 };
252 };
253
254} // namespace thermalfist
255
256#endif
map< string, double > params
BroydenSolutionCriterium(double maximum_error=TOL)
Definition Broyden.h:153
Abstract class which defines the system of non-linear equations to be solved by the Broyden's method.
Definition Broyden.h:32
int m_N
The number of equations.
Definition Broyden.h:66
BroydenEquations()=default
Default constructor. Does nothing.
static const double TOL
Default desired solution accuracy.
Definition Broyden.h:183
BroydenJacobian(BroydenEquations *eqs=NULL)
Construct a new BroydenJacobian object.
Definition Broyden.h:89
ThermalModelBase(ThermalParticleSystem *TPS, const ThermalModelParameters &params=ThermalModelParameters())
Construct a new ThermalModelBase object.
Class implementing the crossterms excluded-volume model.
const std::vector< std::vector< double > > & VirialMatrix() const
virtual double Pressure(int i, const std::vector< double > &pstars=std::vector< double >())
Calculate the ideal gas pressure of particle species i for the given values of partial pressures.
virtual void WriteInteractionParameters(const std::string &filename)
Write the QvdW interaction parameters to a file.
double VirialCoefficient(int i, int j) const
Excluded volume coefficient .
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.
void SetVirial(int i, int j, double b)
Set the excluded volume coefficient .
void CalculateFluctuations()
Computes the fluctuation observables.
double PressureDiagonalTotal(double P)
The total pressure for the given input pressure in the diagonal model.
virtual double MuShift(int i) const
The shift in the chemical potential of particle species i due to the excluded volume interactions.
virtual ~ThermalModelEVCrosstermsLegacy(void)
Destroy the ThermalModelEVCrossterms object.
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.
virtual double DensityId(int i, const std::vector< double > &pstars=std::vector< double >())
Calculate the ideal gas density of particle species i for the given values of partial pressures.
void CalculateTwoParticleCorrelations()
Computes the fluctuations and correlations of the primordial particle numbers.
virtual void SolvePressure(bool resetPartials=true)
Solves the system of transcdental equations for the pressure using the Broyden's method.
virtual void ReadInteractionParameters(const std::string &filename)
Reads the QvdW interaction parameters from a file.
virtual void SolveDiagonal()
Solves the transcendental equation of the corresponding diagonal EV model.
void SetRadius(double rad)
Set the same excluded volume radius parameter for all species.
virtual void ChangeTPS(ThermalParticleSystem *TPS)
Change the particle list.
const std::vector< std::vector< int > > & EVComponentIndices() const
ThermalModelEVCrosstermsLegacy(ThermalParticleSystem *TPS, const ThermalModelParameters &params=ThermalModelParameters())
Construct a new ThermalModelEVCrossterms object.
double PartialPressureDiagonal(int i, double P)
The "partial pressure" of hadron species i for the given total pressure in the diagonal model.
virtual void CalculatePrimordialDensities()
Calculates the primordial densities of all species.
double ScaledVarianceId(int ind)
Calculate the ideal gas scaled variance of particle species i number fluctuations for the given value...
Class containing the particle list.
The main namespace where all classes and functions of the Thermal-FIST library reside.
Definition CosmicEoS.h:9
Structure containing all thermal parameters of the model.