Loading [MathJax]/extensions/tex2jax.js
Thermal-FIST 1.5
Package for hadron resonance gas model applications
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Modules Pages
ExcludedVolumeModelsMulti.cpp
Go to the documentation of this file.
2#include <cmath>
3#include <map>
4
5#include <Eigen/Dense>
6
7using namespace Eigen;
8
9namespace thermalfist {
10
12 const std::vector<double>& n = m_densities;
13 double ret = 1.;
14 for (int ii = 0; ii < n.size(); ++ii)
15 ret -= m_b[ii] * n[ii];
16 return ret;
17 }
18
19 double ExcludedVolumeModelDiagonalVDW::df(int i, int j) const {
20 const std::vector<double>& n = m_densities;
21 return -m_b[j];
22 }
23
25 {
26 const std::vector<double>& n = m_densities;
27 double ret = 0.;
28 for (int ii = 0; ii < n.size(); ++ii)
29 ret -= m_dbdT[ii] * n[ii];
30 return ret;
31 }
32
33 std::vector<double> ExcludedVolumeModelDiagonalVDW::nsol(const std::vector<double>& nid) {
34 std::vector<double> ret = nid;
35 double denom = 1.;
36 for (int ii = 0; ii < nid.size(); ++ii)
37 denom += m_b[ii] * nid[ii];
38 for (int ii = 0; ii < nid.size(); ++ii)
39 ret[ii] = nid[ii] / denom;
40 return ret;
41 }
42
44 {
45 m_components = std::vector<int>(m_N, 0);
46 m_componentsFrom = std::vector<int>();
47 std::map<std::pair<double,double>, int> MapEV;
49 for (int i = 0; i < m_N; ++i) {
50 std::pair<double, double> param = std::make_pair(m_b[i], m_dbdT[i]);
51 if (MapEV.count(param) == 0) {
52 MapEV[param] = m_componentsNumber;
53 m_componentsFrom.push_back(i);
55 }
56 m_components[i] = MapEV[param];
57 }
58 }
59
61 const std::vector<double>& n = m_densities;
62 double ret = 1.;
63 for (int j = 0; j < n.size(); ++j)
64 ret -= m_b[j][i] * n[j];
65 return ret;
66 }
67
68 double ExcludedVolumeModelCrosstermsVDW::df(int i, int j) const {
69 const std::vector<double>& n = m_densities;
70 return -m_b[j][i];
71 }
72
74 const std::vector<double>& n = m_densities;
75 double ret = 0.;
76 for (int j = 0; j < n.size(); ++j)
77 ret -= m_dbdT[j][i] * n[j];
78 return ret;
79 }
80
81 std::vector<double> ExcludedVolumeModelCrosstermsVDW::nsol(const std::vector<double>& nid) {
82 std::vector<double> ret = nid;
83 int NN = nid.size();
84
85 MatrixXd densMatrix(NN, NN);
86 VectorXd solVector(NN), xVector(NN);
87
88 for (int i = 0; i < NN; ++i)
89 for (int j = 0; j < NN; ++j) {
90 densMatrix(i, j) = m_b[j][i] * nid[i];
91 if (i == j) densMatrix(i, j) += 1.;
92 }
93
94 PartialPivLU<MatrixXd> decomp(densMatrix);
95
96 for (int i = 0; i < NN; ++i) xVector[i] = nid[i];
97 solVector = decomp.solve(xVector);
98
99 for (int i = 0; i < NN; ++i) ret[i] = solVector[i];
100
101 return ret;
102 }
103
105 {
106 m_components = std::vector<int>(m_N, 0);
107 m_componentsFrom = std::vector<int>();
108 std::map<std::vector<double>, int> MapEV;
110 for (int i = 0; i < m_N; ++i) {
111 std::vector<double> param(0);
112 for (int j = 0; j < m_N; ++j)
113 param.push_back(m_b[i][j]);
114 for (int j = 0; j < m_N; ++j)
115 param.push_back(m_b[j][i]);
116 for (int j = 0; j < m_N; ++j)
117 param.push_back(m_dbdT[i][j]);
118 for (int j = 0; j < m_N; ++j)
119 param.push_back(m_dbdT[j][i]);
120
121 if (MapEV.count(param) == 0) {
122 MapEV[param] = m_componentsNumber;
123 m_componentsFrom.push_back(i);
125 }
126 m_components[i] = MapEV[param];
127 }
128 }
129
137
139 {
140 return m_evmodelsingle->f(m_eta);
141 }
142
144 {
145 return m_evmodelsingle->df(1, m_eta) * (m_b[j] / 4.);
146 }
147
148 double ExcludedVolumeModelDiagonalGeneralized::d2f(int i, int j, int k) const
149 {
150 return m_evmodelsingle->df(2, m_eta) * (m_b[j] / 4.) * (m_b[k] / 4.);
151 }
152
153 double ExcludedVolumeModelDiagonalGeneralized::d3f(int i, int j, int k, int l) const
154 {
155 return m_evmodelsingle->df(3, m_eta) * (m_b[j] / 4.) * (m_b[k] / 4.) * (m_b[l] / 4.);
156 }
157
158 double ExcludedVolumeModelDiagonalGeneralized::d4f(int i, int j, int k, int l, int m) const
159 {
160 return m_evmodelsingle->df(4, m_eta) * (m_b[j] / 4.) * (m_b[k] / 4.) * (m_b[l] / 4.) * (m_b[m] / 4.);
161 }
162
164 {
165 double ret = 0.;
166 for (int i = 0; i < m_N; ++i)
167 ret += m_dbdT[i] * m_densities[i] / 4.;
168 ret *= m_evmodelsingle->df(1, m_eta);
169 return ret;
170 }
171
172 std::vector<double> ExcludedVolumeModelDiagonalGeneralized::nsol(const std::vector<double>& nid)
173 {
174 double etaid = GetEta(nid);
175 double eta = m_evmodelsingle->etasol(etaid);
176 double fval = m_evmodelsingle->f(eta);
177 std::vector<double> ret(m_N, 0.);
178 for (int i = 0; i < m_N; ++i)
179 ret[i] = fval * nid[i];
180 return ret;
181 }
182
188
190 {
191 m_components = std::vector<int>(m_N, 0);
192 m_componentsFrom = std::vector<int>();
193 std::map<std::pair<double, double>, int> MapEV;
195 for (int i = 0; i < m_N; ++i) {
196 std::pair<double, double> param = std::make_pair(m_b[i], m_dbdT[i]);
197 if (MapEV.count(param) == 0) {
198 MapEV[param] = m_componentsNumber;
199 m_componentsFrom.push_back(i);
201 }
202 m_components[i] = MapEV[param];
203 }
204 }
205
206 double ExcludedVolumeModelDiagonalGeneralized::GetEta(const std::vector<double>& n) const
207 {
208 double ret = 0.0;
209 for (int i = 0; i < n.size(); ++i)
210 ret += m_b[i] * n[i];
211 ret /= 4.;
212 return ret;
213 }
214
215 std::vector<double> ExcludedVolumeModelMultiBase::BroydenEquationsEVMulti::Equations(const std::vector<double>& x)
216 {
217 std::vector<double> ret(m_N, 0.);
218
219 std::vector<double> dens(m_EVM->m_N, 0.);
220 for (int i = 0; i < m_EVM->m_N; ++i) {
221 int ti = i;
222 if (m_componentsMode)
223 ti = m_EVM->ComponentIndices()[i];
224 dens[i] = x[ti] * m_ntil->operator[](i);
225 }
226 m_EVM->SetDensities(dens);
227
228 for (int i = 0; i < m_N; ++i) {
229 int fi = i;
230 if (m_componentsMode)
231 fi = m_EVM->ComponentIndicesFrom()[i];
232 ret[i] = x[i] - m_EVM->f(fi);
233 }
234
235 return ret;
236 }
237
238 std::vector<double> ExcludedVolumeModelMultiBase::BroydenJacobianEVMulti::Jacobian(const std::vector<double>& x)
239 {
240 int N = x.size();
241 std::vector<double> ret(N * N);
242
243 std::vector<double> dens(m_EVM->m_N, 0.);
244 for (int i = 0; i < m_EVM->m_N; ++i) {
245 int ti = i;
246 if (m_componentsMode)
247 ti = m_EVM->ComponentIndices()[i];
248 dens[i] = x[ti] * m_ntil->operator[](i);
249 }
250 m_EVM->SetDensities(dens);
251
252 std::vector<double> ntilcomp(N, 0.);
253 for (int i = 0; i < m_EVM->m_N; ++i) {
254 int ti = i;
255 if (m_componentsMode)
256 ti = m_EVM->ComponentIndices()[i];
257 ntilcomp[ti] = m_ntil->operator[](i);
258 }
259
260 for (int i = 0; i < N; ++i) {
261 int fi = i;
262 if (m_componentsMode)
263 fi = m_EVM->ComponentIndicesFrom()[i];
264 for (int j = 0; j < N; ++j) {
265 int fj = j;
266 if (m_componentsMode)
267 fj = m_EVM->ComponentIndicesFrom()[j];
268 ret[i * N + j] = 0.;
269 if (i == j)
270 ret[i * N + j] += 1.;
271 ret[i * N + j] += -m_EVM->df(fi, fj) * ntilcomp[j];
272 }
273 }
274
275 return ret;;
276 }
277
278 std::vector<double> ExcludedVolumeModelMultiBase::nsolBroyden(const std::vector<double>& ntil)
279 {
280 BroydenEquationsEVMulti eqs(this, &ntil);
281 BroydenJacobianEVMulti jac(this, &ntil);
282 Broyden broydn(&eqs, &jac);
284
285 std::vector<double> ret(m_N, 0.);
286
287 ret = broydn.Solve(ret, &crit);
288 for (int i = 0; i < m_N; ++i)
289 ret[i] *= ntil[i];
290
291 return ret;
292 }
293
294 std::vector<double> ExcludedVolumeModelMultiBase::nsolBroydenComponents(const std::vector<double>& ntil)
295 {
296 BroydenEquationsEVMulti eqs(this, &ntil, true);
297 BroydenJacobianEVMulti jac(this, &ntil, true);
298 Broyden broydn(&eqs, &jac);
300
301 std::vector<double> ret(m_N, 0.);
302 std::vector<double> sol(m_componentsNumber, 0.);
303
304 sol = broydn.Solve(sol, &crit);
305 for (int i = 0; i < m_N; ++i)
306 ret[i] = ntil[i] * sol[m_components[i]];
307
308 return ret;
309 }
310
312 {
313 m_components = std::vector<int>(m_N, 0);
315 m_componentsFrom = std::vector<int>(1, 0);
316 }
317
318
320 {
322 for (int i = 0; i < m_componentsNumber; ++i) {
324 }
325 //for (int i = 0; i < m_N; ++i) {
326 // m_etas[i] = GetEta(i, m_densities);
327 //}
328 }
329
331 {
332 m_components = std::vector<int>(m_N, 0);
333 m_componentsFrom = std::vector<int>();
334 std::map<std::vector<double>, int> MapEV;
336 for (int i = 0; i < m_N; ++i) {
337 std::vector<double> param(0);
338 for (int j = 0; j < m_N; ++j)
339 param.push_back(m_b[i][j]);
340 for (int j = 0; j < m_N; ++j)
341 param.push_back(m_b[j][i]);
342 for (int j = 0; j < m_N; ++j)
343 param.push_back(m_dbdT[i][j]);
344 for (int j = 0; j < m_N; ++j)
345 param.push_back(m_dbdT[j][i]);
346
347 if (MapEV.count(param) == 0) {
348 MapEV[param] = m_componentsNumber;
349 m_componentsFrom.push_back(i);
351 }
352 m_components[i] = MapEV[param];
353 }
354
356 for (int i = 0; i < m_componentsNumber; ++i) {
357 for (int j = 0; j < m_componentsNumber; ++j) {
358 if (i != j) {
361 }
362 }
363 }
364
365 printf("EV Components: %d Are Disconnected: %d\n", m_componentsNumber, (int)m_componentsDisconnected);
366 }
367
368 double ExcludedVolumeModelCrosstermsGeneralized::GetEta(int i, const std::vector<double>& n) const
369 {
370 double ret = 0.0;
371 for (int j = 0; j < n.size(); ++j)
372 ret += m_b[j][i] * n[j];
373 ret /= 4.;
374 return ret;
375 }
376
384
386 {
387 //return m_evmodelsingle->f(m_etas[i]);
388 return m_evmodelsingle->f(m_etas[ComponentIndices()[i]]);
389 }
390
392 {
393 //return m_evmodelsingle->df(1, m_etas[i]) * (m_b[j][i] / 4.);
394 return m_evmodelsingle->df(1, m_etas[ComponentIndices()[i]]) * (m_b[j][i] / 4.);
395 }
396
397 double ExcludedVolumeModelCrosstermsGeneralized::d2f(int i, int j, int k) const
398 {
399 //return m_evmodelsingle->df(2, m_etas[i]) * (m_b[j][i] / 4.) * (m_b[k][i] / 4.);
400 return m_evmodelsingle->df(2, m_etas[ComponentIndices()[i]]) * (m_b[j][i] / 4.) * (m_b[k][i] / 4.);
401 }
402
403 double ExcludedVolumeModelCrosstermsGeneralized::d3f(int i, int j, int k, int l) const
404 {
405 //return m_evmodelsingle->df(3, m_etas[i]) * (m_b[j][i] / 4.) * (m_b[k][i] / 4.) * (m_b[l][i] / 4.);
406 return m_evmodelsingle->df(3, m_etas[ComponentIndices()[i]]) * (m_b[j][i] / 4.) * (m_b[k][i] / 4.) * (m_b[l][i] / 4.);
407 }
408
409 double ExcludedVolumeModelCrosstermsGeneralized::d4f(int i, int j, int k, int l, int m) const
410 {
411 //return m_evmodelsingle->df(4, m_etas[i]) * (m_b[j][i] / 4.) * (m_b[k][i] / 4.) * (m_b[l][i] / 4.) * (m_b[m][i] / 4.);
412 return m_evmodelsingle->df(4, m_etas[ComponentIndices()[i]]) * (m_b[j][i] / 4.) * (m_b[k][i] / 4.) * (m_b[l][i] / 4.) * (m_b[m][i] / 4.);
413 }
414
416 {
417 double ret = 0.;
418 for (int j = 0; j < m_N; ++j)
419 ret += m_dbdT[j][i] * m_densities[j] / 4.;
420 //ret *= m_evmodelsingle->df(1, m_etas[i]);
421 ret *= m_evmodelsingle->df(1, m_etas[ComponentIndices()[i]]);
422 return ret;
423 }
424
425 std::vector<double> ExcludedVolumeModelCrosstermsGeneralized::nsol(const std::vector<double>& nid)
426 {
428 return nsolBroydenComponents(nid);
429 return nsolBroyden(nid);
430 }
431
432 // Disconnected components, use binary search separately for each component
433 std::vector<double> nids_comp(m_componentsNumber, 0.);
434 for (int i = 0; i < nid.size(); ++i)
435 nids_comp[m_components[i]] += nid[i];
436
437 std::vector<double> fvals(m_componentsNumber, 0.);
438 for (int ic = 0; ic < m_componentsNumber; ++ic) {
439 int ti = m_componentsFrom[ic];
440 double etaid = m_b[ti][ti] * nids_comp[ic] / 4.;
441 double eta = m_evmodelsingle->etasol(etaid);
442 fvals[ic] = m_evmodelsingle->f(eta);
443 }
444 std::vector<double> ret(m_N, 0.);
445 for (int i = 0; i < nid.size(); ++i) {
446 ret[i] = fvals[m_components[i]] * nid[i];
447 }
448 return ret;
449
450
451
452 //std::vector<double> etasid(m_N);
453 //for (int i = 0; i < m_N; ++i)
454 // etasid[i] = GetEta(i, nid);
455 //for (int i = 0; i < m_N; ++i) {
456 // double etaid = etasid[i];
457 // double eta = m_evmodelsingle->etasol(etaid);
458 // double fval = m_evmodelsingle->f(eta);
459 // ret[i] = fval * nid[i];
460 //}
461 //return ret;
462 }
463
464 void ExcludedVolumeModelComponents::SetDensities(const std::vector<double>& n)
465 {
467
468 for (int i = 0; i < m_densities_components.size(); ++i)
470
471 for (int i = 0; i < n.size(); ++i)
473 }
474
476 {
478 for (int i = m_N - 1; i >= 0; --i) {
480 }
481 }
482
484 {
485 for (int i = 0; i < m_componentsNumber; ++i) {
486 if (m_evmodels[i] != NULL) {
487 delete m_evmodels[i];
488 m_evmodels[i] = NULL;
489 }
490 }
491 }
492
494 {
495 int tind = m_components[i];
496 double eta = m_b[tind] * m_densities_components[tind] / 4.;
497 return m_evmodels[tind]->f(eta);
498 }
499
500 double ExcludedVolumeModelComponents::df(int i, int j) const
501 {
502 bool fl = 1;
503 fl &= (m_components[i] == m_components[j]);
504 if (!fl) return 0.;
505
506 int tind = m_components[i];
507 double eta = m_b[tind] * m_densities_components[tind] / 4.;
508 return m_evmodels[tind]->df(1, eta) * (m_b[tind] / 4.);
509 }
510
511 double ExcludedVolumeModelComponents::d2f(int i, int j, int k) const
512 {
513 bool fl = 1;
514 fl &= (m_components[i] == m_components[j]);
515 fl &= (m_components[i] == m_components[k]);
516 if (!fl) return 0.;
517
518 int tind = m_components[i];
519 double eta = m_b[tind] * m_densities_components[tind] / 4.;
520 return m_evmodels[tind]->df(2, eta) * (m_b[tind] / 4.) * (m_b[tind] / 4.);
521 }
522
523 double ExcludedVolumeModelComponents::d3f(int i, int j, int k, int l) const
524 {
525 bool fl = 1;
526 fl &= (m_components[i] == m_components[j]);
527 fl &= (m_components[i] == m_components[k]);
528 fl &= (m_components[i] == m_components[l]);
529 if (!fl) return 0.;
530
531 int tind = m_components[i];
532 double eta = m_b[tind] * m_densities_components[tind] / 4.;
533 return m_evmodels[tind]->df(3, eta) * (m_b[tind] / 4.) * (m_b[tind] / 4.) * (m_b[tind] / 4.);
534 }
535
536 double ExcludedVolumeModelComponents::d4f(int i, int j, int k, int l, int m) const
537 {
538 bool fl = 1;
539 fl &= (m_components[i] == m_components[j]);
540 fl &= (m_components[i] == m_components[k]);
541 fl &= (m_components[i] == m_components[l]);
542 fl &= (m_components[i] == m_components[m]);
543 if (!fl) return 0.;
544
545 int tind = m_components[i];
546 double eta = m_b[tind] * m_densities_components[tind] / 4.;
547 return m_evmodels[tind]->df(4, eta) * (m_b[tind] / 4.) * (m_b[tind] / 4.) * (m_b[tind] / 4.) * (m_b[tind] / 4.);
548 }
549
551 {
552 int tind = m_components[i];
553 double eta = m_b[tind] * m_densities_components[tind] / 4.;
554 return m_evmodels[tind]->df(1, eta) * (m_dbdT[tind] * m_densities_components[tind] / 4.);
555 }
556
557 std::vector<double> ExcludedVolumeModelComponents::nsol(const std::vector<double>& nid)
558 {
559 std::vector<double> nids_comp(m_componentsNumber, 0.);
560 for (int i = 0; i < nid.size(); ++i)
561 nids_comp[m_components[i]] += nid[i];
562
563 std::vector<double> fvals(m_componentsNumber, 0.);
564 for (int ic = 0; ic < m_componentsNumber; ++ic) {
565 double etaid = m_b[ic] * nids_comp[ic] / 4.;
566 double eta = m_evmodels[ic]->etasol(etaid);
567 fvals[ic] = m_evmodels[ic]->f(eta);
568 }
569 std::vector<double> ret(m_N, 0.);
570 for (int i = 0; i < nid.size(); ++i) {
571 ret[i] = fvals[m_components[i]] * nid[i];
572 }
573 return ret;
574 }
575
576} // namespace thermalfist
Sub-class where it is determined whether the required accuracy is achieved in the Broyden's method.
Definition Broyden.h:144
int m_N
The number of equations.
Definition Broyden.h:66
Class implementing the Broyden method to solve a system of non-linear equations.
Definition Broyden.h:132
virtual std::vector< double > Solve(const std::vector< double > &x0, BroydenSolutionCriterium *solcrit=NULL, int max_iterations=MAX_ITERS)
Definition Broyden.cpp:83
virtual ~ExcludedVolumeModelComponents()
Destructor for the ExcludedVolumeModelComponents class.
virtual double d2f(int i, int j, int k) const
Calculates the second derivative of the suppression factor.
virtual double d4f(int i, int j, int k, int l, int m) const
Calculates the fourth derivative of the suppression factor.
virtual double d3f(int i, int j, int k, int l) const
Calculates the third derivative of the suppression factor.
virtual void SetDensities(const std::vector< double > &n)
Sets the densities of particle species.
virtual void ComputeComponents()
Computes the components based on the excluded volume parameters.
virtual double dfdT(int i) const
Calculates the temperature derivative of the suppression factor.
virtual std::vector< double > nsol(const std::vector< double > &nid)
Solves for the actual densities given the ideal gas densities.
virtual double f(int i) const
Calculates the suppression factor for species i.
std::vector< ExcludedVolumeModelBase * > m_evmodels
virtual double df(int i, int j) const
Calculates the first derivative of the suppression factor.
virtual double d3f(int i, int j, int k, int l) const
Calculates the third derivative of the suppression factor.
virtual double d2f(int i, int j, int k) const
Calculates the second derivative of the suppression factor.
virtual void SetDensities(const std::vector< double > &n)
Sets the densities of particle species.
virtual std::vector< double > nsol(const std::vector< double > &nid)
Solves for the actual densities given the ideal gas densities.
virtual double dfdT(int i) const
Calculates the temperature derivative of the suppression factor.
double GetEta(int i, const std::vector< double > &n) const
Calculates the eta parameter for the given densities and species index.
virtual double df(int i, int j) const
Calculates the first derivative of the suppression factor.
virtual double d4f(int i, int j, int k, int l, int m) const
Calculates the fourth derivative of the suppression factor.
virtual void ComputeComponents()
Computes the components based on the excluded volume parameters.
virtual double f(int i) const
Calculates the suppression factor for species i.
virtual ~ExcludedVolumeModelCrosstermsGeneralized()
Destructor for the ExcludedVolumeModelCrosstermsGeneralized class.
virtual double df(int i, int j) const
Calculates the first derivative of the suppression factor.
virtual void ComputeComponents()
Computes the components based on the excluded volume parameters.
virtual std::vector< double > nsol(const std::vector< double > &nid)
Solves for the actual densities given the ideal gas densities.
virtual double dfdT(int i) const
Calculates the temperature derivative of the suppression factor.
virtual double f(int i) const
Calculates the suppression factor for species i.
virtual ~ExcludedVolumeModelDiagonalGeneralized()
Destructor for the ExcludedVolumeModelDiagonalGeneralized class.
virtual double f(int i) const
Calculates the suppression factor for species i.
virtual double df(int i, int j) const
Calculates the first derivative of the suppression factor.
virtual double d2f(int i, int j, int k) const
Calculates the second derivative of the suppression factor.
double GetEta(const std::vector< double > &n) const
Calculates the eta parameter for the given densities.
virtual double dfdT(int i) const
Calculates the temperature derivative of the suppression factor.
virtual void ComputeComponents()
Computes the components based on the excluded volume parameters.
virtual std::vector< double > nsol(const std::vector< double > &nid)
Solves for the actual densities given the ideal gas densities.
virtual double d3f(int i, int j, int k, int l) const
Calculates the third derivative of the suppression factor.
virtual double d4f(int i, int j, int k, int l, int m) const
Calculates the fourth derivative of the suppression factor.
virtual void SetDensities(const std::vector< double > &n)
Sets the densities of particle species.
virtual void ComputeComponents()
Computes the components based on the excluded volume parameters.
virtual double f(int i) const
Calculates the suppression factor for species i.
virtual double df(int i, int j) const
Calculates the first derivative of the suppression factor.
virtual double dfdT(int i) const
Calculates the temperature derivative of the suppression factor.
virtual std::vector< double > nsol(const std::vector< double > &nid)
Solves for the actual densities given the ideal gas densities.
std::vector< double > Equations(const std::vector< double > &x)
Calculates the equations for Broyden's method.
std::vector< double > Jacobian(const std::vector< double > &x)
Calculates the Jacobian matrix for Broyden's method.
virtual std::vector< double > nsolBroydenComponents(const std::vector< double > &ntil)
Solves for the actual densities using Broyden's method, considering components.
virtual const std::vector< int > & ComponentIndices() const
Gets the component indices.
virtual void ComputeComponents()
Computes the components based on the excluded volume parameters.
virtual void SetDensities(const std::vector< double > &n)
Sets the densities of particle species.
virtual std::vector< double > nsolBroyden(const std::vector< double > &ntil)
Solves for the actual densities using Broyden's method.
The main namespace where all classes and functions of the Thermal-FIST library reside.
Definition CosmicEoS.h:9