Thermal-FIST 1.5
Package for hadron resonance gas model applications
Loading...
Searching...
No Matches
ExcludedVolumeModelsMulti.h
Go to the documentation of this file.
1/*
2 * Thermal-FIST package
3 *
4 * Copyright (c) 2017-2022 Volodymyr Vovchenko
5 *
6 * GNU General Public License (GPLv3 or later)
7 */
8#ifndef EXCLUDEVOLUMEMODELSMULTI_H
9#define EXCLUDEVOLUMEMODELSMULTI_H
10#include <cmath>
11#include <vector>
12
14#include "HRGBase/Broyden.h"
15
16namespace thermalfist {
17
26 public:
33 m_N(N),
34 m_densities(std::vector<double>(m_N,0.)),
35 m_components(std::vector<int>(m_N, 0))
36 {
38 }
40
47 virtual double f(int i) const { return 1.; }
48
56 virtual double df(int i, int j) const { return 0.; }
57
66 virtual double d2f(int i, int j, int k) const { return 0.; }
67
77 virtual double d3f(int i, int j, int k, int l) const { return 0.; }
78
89 virtual double d4f(int i, int j, int k, int l, int m) const { return 0.; }
90
97 virtual double dfdT(int i) const { return 0.; }
98
105 virtual std::vector<double> nsol(const std::vector<double>& ntil) { return ntil; }
106
113 virtual std::vector<double> nsolBroyden(const std::vector<double>& ntil);
114
121 virtual std::vector<double> nsolBroydenComponents(const std::vector<double>& ntil);
122
128 virtual void SetDensities(const std::vector<double>& n) { m_densities = n; }
129
135 virtual const std::vector<int>& ComponentIndices() const { return m_components; }
136
142 virtual const std::vector<int>& ComponentIndicesFrom() const { return m_componentsFrom; }
143
149 virtual const int ComponentsNumber() const { return m_componentsNumber; }
150 protected:
154 virtual void ComputeComponents();
155 int m_N;
156 std::vector<double> m_densities;
157 std::vector<int> m_components;
158 std::vector<int> m_componentsFrom;
160
162 {
163 public:
171 BroydenEquationsEVMulti(ExcludedVolumeModelMultiBase* evmodel, const std::vector<double>* ntil, bool componentsMode = false) :
173 m_EVM(evmodel),
174 m_ntil(ntil),
175 m_componentsMode(componentsMode)
176 {
177 m_N = evmodel->m_N;
178 if (m_componentsMode)
179 m_N = evmodel->m_componentsNumber;
180 }
181
188 std::vector<double> Equations(const std::vector<double>& x);
189 private:
191 const std::vector<double>* m_ntil;
192 bool m_componentsMode;
193 };
194
196 {
197 public:
206 const std::vector<double>* ntil,
207 bool componentsMode = false) :
209 m_EVM(evmodel),
210 m_ntil(ntil),
211 m_componentsMode(componentsMode)
212 {
213 }
214
221 std::vector<double> Jacobian(const std::vector<double>& x);
222 private:
224 const std::vector<double>* m_ntil;
225 bool m_componentsMode;
226 };
227
228
229 };
230
237 public:
244 ExcludedVolumeModelDiagonalVDW(const std::vector<double>& b, const std::vector<double>& dbdT = std::vector<double>())
245 : ExcludedVolumeModelMultiBase(b.size()), m_b(b), m_dbdT(dbdT) {
246 if (m_dbdT.size() == 0)
247 m_dbdT = std::vector<double>(m_b.size(), 0.);
249 }
250
256 virtual double f(int i) const;
264 virtual double df(int i, int j) const;
273 virtual double d2f(int i, int j, int k) const { return 0.; }
283 virtual double d3f(int i, int j, int k, int l) const { return 0.; }
294 virtual double d4f(int i, int j, int k, int l, int m) const { return 0.; }
301 virtual double dfdT(int i) const;
308 virtual std::vector<double> nsol(const std::vector<double>& nid);
309 protected:
313 virtual void ComputeComponents();
314 std::vector<double> m_b;
315 std::vector<double> m_dbdT;
316 };
317
325 public:
334 ExcludedVolumeModelBase *evmodelsingle,
335 const std::vector<double>& b,
336 const std::vector<double>& dbdT = std::vector<double>())
337 : ExcludedVolumeModelMultiBase(b.size()), m_evmodelsingle(evmodelsingle), m_b(b), m_dbdT(dbdT), m_eta(0.) {
338 if (m_dbdT.size() == 0)
339 m_dbdT = std::vector<double>(m_b.size(), 0.);
341 }
342
352 virtual double f(int i) const;
360 virtual double df(int i, int j) const;
369 virtual double d2f(int i, int j, int k) const;
379 virtual double d3f(int i, int j, int k, int l) const;
390 virtual double d4f(int i, int j, int k, int l, int m) const;
397 virtual double dfdT(int i) const;
404 virtual std::vector<double> nsol(const std::vector<double>& nid);
410 virtual void SetDensities(const std::vector<double>& n);
411 protected:
415 virtual void ComputeComponents();
422 double GetEta(const std::vector<double>& n) const;
424 std::vector<double> m_b;
425 std::vector<double> m_dbdT;
426 double m_eta;
427 };
428
429
438 public:
445 ExcludedVolumeModelCrosstermsVDW(const std::vector< std::vector<double> >& b, const std::vector< std::vector<double> >& dbdT = std::vector< std::vector<double> >())
446 : ExcludedVolumeModelMultiBase(b.size()), m_b(b), m_dbdT(dbdT) {
447 if (m_dbdT.size() == 0)
448 m_dbdT = std::vector< std::vector<double> >(m_b.size(), std::vector<double>(m_b.size(), 0.));
450 }
451
457 virtual double f(int i) const;
465 virtual double df(int i, int j) const;
474 virtual double d2f(int i, int j, int k) const { return 0.; }
484 virtual double d3f(int i, int j, int k, int l) const { return 0.; }
495 virtual double d4f(int i, int j, int k, int l, int m) const { return 0.; }
502 virtual double dfdT(int i) const;
509 virtual std::vector<double> nsol(const std::vector<double>& nid);
510 protected:
514 virtual void ComputeComponents();
515 std::vector< std::vector<double> > m_b;
516 std::vector< std::vector<double> > m_dbdT;
517 };
518
529 public:
538 ExcludedVolumeModelBase* evmodelsingle,
539 const std::vector< std::vector<double> >& b,
540 const std::vector< std::vector<double> >& dbdT = std::vector< std::vector<double> >()
541 )
542 : ExcludedVolumeModelMultiBase(b.size()), m_evmodelsingle(evmodelsingle), m_b(b), m_dbdT(dbdT), m_componentsDisconnected(false) {
543 if (m_dbdT.size() == 0)
544 m_dbdT = std::vector< std::vector<double> >(m_b.size(), std::vector<double>(m_b.size(), 0.));
546 m_etas = std::vector<double>(m_componentsNumber, 0.);
547 }
548
558 virtual double f(int i) const;
566 virtual double df(int i, int j) const;
575 virtual double d2f(int i, int j, int k) const;
585 virtual double d3f(int i, int j, int k, int l) const;
596 virtual double d4f(int i, int j, int k, int l, int m) const;
603 virtual double dfdT(int i) const;
610 virtual std::vector<double> nsol(const std::vector<double>& nid);
616 virtual void SetDensities(const std::vector<double>& n);
617 protected:
621 virtual void ComputeComponents();
629 double GetEta(int i, const std::vector<double>& n) const;
631 std::vector< std::vector<double> > m_b;
632 std::vector< std::vector<double> > m_dbdT;
633 std::vector<double> m_etas;
635 };
636
645 public:
656 int components,
657 const std::vector<ExcludedVolumeModelBase*>& evmods,
658 const std::vector<int>& ind,
659 const std::vector<double>& b,
660 const std::vector<double>& dbdT = std::vector<double>()
661 )
662 : ExcludedVolumeModelMultiBase(ind.size()),
663 //m_N_components(components),
664 m_evmodels(evmods),
665 //m_indices(ind),
666 m_b(b),
667 m_dbdT(dbdT),
668 m_densities_components(std::vector<double>(b.size(), 0.))
669 {
670 if (m_dbdT.size() == 0)
671 m_dbdT.resize(m_b.size(), 0.);
672 m_components = ind;
673 m_componentsNumber = components;
675 }
676
686 virtual double f(int i) const;
694 virtual double df(int i, int j) const;
703 virtual double d2f(int i, int j, int k) const;
713 virtual double d3f(int i, int j, int k, int l) const;
724 virtual double d4f(int i, int j, int k, int l, int m) const;
731 virtual double dfdT(int i) const;
738 virtual std::vector<double> nsol(const std::vector<double>& nid);
744 virtual void SetDensities(const std::vector<double>& n);
745 protected:
749 virtual void ComputeComponents();
750 //int m_N_components;
751 std::vector<ExcludedVolumeModelBase*> m_evmodels;
752 //std::vector<int> m_indices;
753 std::vector<double> m_b;
754 std::vector<double> m_dbdT;
755 std::vector<double> m_densities_components;
756 };
757
758} // namespace thermalfist
759
760#endif // REALGASMODELS_H
Implementation of the generic Broyden's method routines.
int m_N
The number of equations.
Definition Broyden.h:66
BroydenEquations()=default
Default constructor. Does nothing.
BroydenJacobian(BroydenEquations *eqs=NULL)
Construct a new BroydenJacobian object.
Definition Broyden.h:89
Base class implementing the ideal gas.
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.
ExcludedVolumeModelComponents(int components, const std::vector< ExcludedVolumeModelBase * > &evmods, const std::vector< int > &ind, const std::vector< double > &b, const std::vector< double > &dbdT=std::vector< double >())
Constructor for the ExcludedVolumeModelComponents class.
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.
ExcludedVolumeModelCrosstermsGeneralized(ExcludedVolumeModelBase *evmodelsingle, const std::vector< std::vector< double > > &b, const std::vector< std::vector< double > > &dbdT=std::vector< std::vector< double > >())
Constructor for the ExcludedVolumeModelCrosstermsGeneralized class.
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 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 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.
ExcludedVolumeModelCrosstermsVDW(const std::vector< std::vector< double > > &b, const std::vector< std::vector< double > > &dbdT=std::vector< std::vector< double > >())
Constructor for the ExcludedVolumeModelCrosstermsVDW class.
virtual double d3f(int i, int j, int k, int l) const
Calculates the third 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 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.
ExcludedVolumeModelDiagonalGeneralized(ExcludedVolumeModelBase *evmodelsingle, const std::vector< double > &b, const std::vector< double > &dbdT=std::vector< double >())
Constructor for the ExcludedVolumeModelDiagonalGeneralized class.
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 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.
ExcludedVolumeModelDiagonalVDW(const std::vector< double > &b, const std::vector< double > &dbdT=std::vector< double >())
Constructor for the ExcludedVolumeModelDiagonalVDW class.
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 d2f(int i, int j, int k) const
Calculates the second derivative of the suppression factor.
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.
BroydenEquationsEVMulti(ExcludedVolumeModelMultiBase *evmodel, const std::vector< double > *ntil, bool componentsMode=false)
Constructor for the BroydenEquationsEVMulti class.
BroydenJacobianEVMulti(ExcludedVolumeModelMultiBase *evmodel, const std::vector< double > *ntil, bool componentsMode=false)
Constructor for the BroydenJacobianEVMulti class.
std::vector< double > Jacobian(const std::vector< double > &x)
Calculates the Jacobian matrix for Broyden's method.
virtual double f(int i) const
Calculates the suppression factor for species i.
virtual std::vector< double > nsolBroydenComponents(const std::vector< double > &ntil)
Solves for the actual densities using Broyden's method, considering components.
virtual const int ComponentsNumber() const
Gets the number of 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 double dfdT(int i) const
Calculates the temperature 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 double df(int i, int j) const
Calculates the first derivative of the suppression factor.
virtual const std::vector< int > & ComponentIndicesFrom() const
Gets the component indices from.
ExcludedVolumeModelMultiBase(int N)
Constructor for the ExcludedVolumeModelMultiBase class.
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.
virtual double d2f(int i, int j, int k) const
Calculates the second derivative of the suppression factor.
virtual std::vector< double > nsol(const std::vector< double > &ntil)
Solves for the actual densities given the ideal gas densities.
The main namespace where all classes and functions of the Thermal-FIST library reside.
Definition CosmicEoS.h:9