CoolProp
6.6.0
An open-source fluid property and humid air property database
|
A class to do newton raphson solver mixture bubble point and dew point calculations.
A class is used rather than a function so that it is easier to store iteration histories, additional output values, etc. This class is used in PhaseEnvelopeRoutines for the construction of the phase envelope
This class only handles bubble and dew lines. The independent variables are the first N-1 mole fractions in the incipient phase along with one of T, p, or \(\rho''\).
These methods are based on the work of Gernert, FPE, 2014, the thesis of Gernert, as well as much help from Andreas Jaeger of Uni. Bochum.
There are N residuals from
\(F_i = \ln f_i(T, p, \mathbf{x}) - \ln f_i(T, p, \mathbf{y})\) for \(i = 1, ... N\)
if either T or p are imposed. In this case a solver is used to find \(\rho\) given \(T\) and \(p\).
In the case that \(\rho''\) is imposed, the case is much nicer and the first \(N\) residuals are
\(F_i = \ln f_i(T, \rho', \mathbf{x}) - \ln f_i(T, \rho'', \mathbf{y})\) for \(i = 1, ... N\)
which requires no iteration. A final residual is set up to ensure the same pressures in the phases:
\(p(T,\rho',\mathbf{x})-p(T,\rho'',\mathbf{y}) = 0\)
Documentation of the derivatives needed can be found in the work of Gernert, FPE, 2014
Definition at line 476 of file VLERoutines.h.
#include <VLERoutines.h>
Public Member Functions | |
newton_raphson_saturation () | |
void | resize (std::size_t N) |
void | pre_call () |
void | call (HelmholtzEOSMixtureBackend &HEOS, const std::vector< CoolPropDbl > &z, std::vector< CoolPropDbl > &z_incipient, newton_raphson_saturation_options &IO) |
Call the Newton-Raphson VLE Solver. More... | |
void | build_arrays () |
Build the arrays for the Newton-Raphson solve. More... | |
void | check_Jacobian () |
Check the derivatives in the Jacobian using numerical derivatives. More... | |
Public Attributes | |
newton_raphson_saturation_options::imposed_variable_options | imposed_variable |
CoolPropDbl | error_rms |
CoolPropDbl | rhomolar_liq |
CoolPropDbl | rhomolar_vap |
CoolPropDbl | T |
CoolPropDbl | p |
CoolPropDbl | min_rel_change |
std::size_t | N |
bool | logging |
bool | bubble_point |
int | Nsteps |
Eigen::MatrixXd | J |
HelmholtzEOSMixtureBackend * | HEOS |
CoolPropDbl | dTsat_dPsat |
CoolPropDbl | dPsat_dTsat |
std::vector< CoolPropDbl > | K |
std::vector< CoolPropDbl > | x |
std::vector< CoolPropDbl > | y |
Eigen::VectorXd | r |
Eigen::VectorXd | err_rel |
std::vector< SuccessiveSubstitutionStep > | step_logger |
|
inline |
Definition at line 492 of file VLERoutines.h.
void CoolProp::SaturationSolvers::newton_raphson_saturation::build_arrays | ( | ) |
Build the arrays for the Newton-Raphson solve.
This method builds the Jacobian matrix, the sensitivity matrix, etc.
Definition at line 1384 of file VLERoutines.cpp.
void CoolProp::SaturationSolvers::newton_raphson_saturation::call | ( | HelmholtzEOSMixtureBackend & | HEOS, |
const std::vector< CoolPropDbl > & | z, | ||
std::vector< CoolPropDbl > & | z_incipient, | ||
newton_raphson_saturation_options & | IO | ||
) |
Call the Newton-Raphson VLE Solver.
This solver must be passed reasonable guess values for the mole fractions, densities, etc. You may want to take a few steps of successive substitution before you start with Newton Raphson.
HEOS | HelmholtzEOSMixtureBackend instance |
z | Bulk mole fractions [-] |
z_incipient | Initial guesses for the mole fractions of the incipient phase [-] |
IO | The input/output data structure |
Definition at line 1284 of file VLERoutines.cpp.
void CoolProp::SaturationSolvers::newton_raphson_saturation::check_Jacobian | ( | ) |
Check the derivatives in the Jacobian using numerical derivatives.
Definition at line 1204 of file VLERoutines.cpp.
|
inline |
Definition at line 497 of file VLERoutines.h.
void CoolProp::SaturationSolvers::newton_raphson_saturation::resize | ( | std::size_t | N | ) |
Definition at line 1187 of file VLERoutines.cpp.
bool CoolProp::SaturationSolvers::newton_raphson_saturation::bubble_point |
Definition at line 483 of file VLERoutines.h.
CoolPropDbl CoolProp::SaturationSolvers::newton_raphson_saturation::dPsat_dTsat |
Definition at line 487 of file VLERoutines.h.
CoolPropDbl CoolProp::SaturationSolvers::newton_raphson_saturation::dTsat_dPsat |
Definition at line 487 of file VLERoutines.h.
Eigen::VectorXd CoolProp::SaturationSolvers::newton_raphson_saturation::err_rel |
Definition at line 489 of file VLERoutines.h.
CoolPropDbl CoolProp::SaturationSolvers::newton_raphson_saturation::error_rms |
Definition at line 480 of file VLERoutines.h.
HelmholtzEOSMixtureBackend* CoolProp::SaturationSolvers::newton_raphson_saturation::HEOS |
Definition at line 486 of file VLERoutines.h.
newton_raphson_saturation_options::imposed_variable_options CoolProp::SaturationSolvers::newton_raphson_saturation::imposed_variable |
Definition at line 479 of file VLERoutines.h.
Eigen::MatrixXd CoolProp::SaturationSolvers::newton_raphson_saturation::J |
Definition at line 485 of file VLERoutines.h.
std::vector<CoolPropDbl> CoolProp::SaturationSolvers::newton_raphson_saturation::K |
Definition at line 488 of file VLERoutines.h.
bool CoolProp::SaturationSolvers::newton_raphson_saturation::logging |
Definition at line 482 of file VLERoutines.h.
CoolPropDbl CoolProp::SaturationSolvers::newton_raphson_saturation::min_rel_change |
Definition at line 480 of file VLERoutines.h.
std::size_t CoolProp::SaturationSolvers::newton_raphson_saturation::N |
Definition at line 481 of file VLERoutines.h.
int CoolProp::SaturationSolvers::newton_raphson_saturation::Nsteps |
Definition at line 484 of file VLERoutines.h.
CoolPropDbl CoolProp::SaturationSolvers::newton_raphson_saturation::p |
Definition at line 480 of file VLERoutines.h.
Eigen::VectorXd CoolProp::SaturationSolvers::newton_raphson_saturation::r |
Definition at line 489 of file VLERoutines.h.
CoolPropDbl CoolProp::SaturationSolvers::newton_raphson_saturation::rhomolar_liq |
Definition at line 480 of file VLERoutines.h.
CoolPropDbl CoolProp::SaturationSolvers::newton_raphson_saturation::rhomolar_vap |
Definition at line 480 of file VLERoutines.h.
std::vector<SuccessiveSubstitutionStep> CoolProp::SaturationSolvers::newton_raphson_saturation::step_logger |
Definition at line 490 of file VLERoutines.h.
CoolPropDbl CoolProp::SaturationSolvers::newton_raphson_saturation::T |
Definition at line 480 of file VLERoutines.h.
std::vector<CoolPropDbl> CoolProp::SaturationSolvers::newton_raphson_saturation::x |
Definition at line 488 of file VLERoutines.h.
std::vector<CoolPropDbl> CoolProp::SaturationSolvers::newton_raphson_saturation::y |
Definition at line 488 of file VLERoutines.h.