8#ifndef EXCLUDEVOLUMEMODELSMULTI_H
9#define EXCLUDEVOLUMEMODELSMULTI_H
47 virtual double f(
int i)
const {
return 1.; }
56 virtual double df(
int i,
int j)
const {
return 0.; }
66 virtual double d2f(
int i,
int j,
int k)
const {
return 0.; }
77 virtual double d3f(
int i,
int j,
int k,
int l)
const {
return 0.; }
89 virtual double d4f(
int i,
int j,
int k,
int l,
int m)
const {
return 0.; }
97 virtual double dfdT(
int i)
const {
return 0.; }
105 virtual std::vector<double>
nsol(
const std::vector<double>& ntil) {
return ntil; }
113 virtual std::vector<double>
nsolBroyden(
const std::vector<double>& ntil);
175 m_componentsMode(componentsMode)
178 if (m_componentsMode)
188 std::vector<double>
Equations(
const std::vector<double>& x);
191 const std::vector<double>* m_ntil;
192 bool m_componentsMode;
206 const std::vector<double>* ntil,
207 bool componentsMode =
false) :
211 m_componentsMode(componentsMode)
221 std::vector<double>
Jacobian(
const std::vector<double>& x);
224 const std::vector<double>* m_ntil;
225 bool m_componentsMode;
247 m_dbdT = std::vector<double>(
m_b.size(), 0.);
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);
335 const std::vector<double>& b,
336 const std::vector<double>& dbdT = std::vector<double>())
339 m_dbdT = std::vector<double>(
m_b.size(), 0.);
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);
422 double GetEta(
const std::vector<double>& n)
const;
448 m_dbdT = std::vector< std::vector<double> >(
m_b.size(), std::vector<double>(
m_b.size(), 0.));
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);
515 std::vector< std::vector<double> >
m_b;
516 std::vector< std::vector<double> >
m_dbdT;
539 const std::vector< std::vector<double> >& b,
540 const std::vector< std::vector<double> >& dbdT = std::vector< std::vector<double> >()
544 m_dbdT = std::vector< std::vector<double> >(
m_b.size(), std::vector<double>(
m_b.size(), 0.));
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);
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;
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>()
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);
Implementation of the generic Broyden's method routines.
int m_N
The number of equations.
BroydenEquations()=default
Default constructor. Does nothing.
BroydenJacobian(BroydenEquations *eqs=NULL)
Construct a new BroydenJacobian object.
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.
std::vector< double > m_b
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
std::vector< double > m_densities_components
std::vector< double > m_dbdT
virtual double df(int i, int j) const
Calculates the first derivative of the suppression factor.
ExcludedVolumeModelBase * m_evmodelsingle
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.
std::vector< std::vector< double > > m_dbdT
bool m_componentsDisconnected
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.
std::vector< std::vector< double > > m_b
virtual void ComputeComponents()
Computes the components based on the excluded volume parameters.
std::vector< double > m_etas
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.
std::vector< std::vector< double > > m_b
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.
std::vector< std::vector< double > > m_dbdT
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.
ExcludedVolumeModelBase * m_evmodelsingle
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.
std::vector< double > m_dbdT
virtual void ComputeComponents()
Computes the components based on the excluded volume parameters.
std::vector< double > m_b
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.
std::vector< double > m_dbdT
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.
std::vector< double > m_b
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.
std::vector< int > m_componentsFrom
virtual std::vector< double > nsolBroydenComponents(const std::vector< double > &ntil)
Solves for the actual densities using Broyden's method, considering components.
std::vector< double > m_densities
std::vector< int > m_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 ~ExcludedVolumeModelMultiBase()
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.