CoolProp 6.8.1dev
An open-source fluid property and humid air property database
Public Member Functions | List of all members
CoolProp::superancillary::SuperAncillary< ArrayType > Class Template Reference

Detailed Description

template<typename ArrayType = Eigen::ArrayXd>
class CoolProp::superancillary::SuperAncillary< ArrayType >

A superancillary object is formed of a number of one dimensional Chebyshev approximations, one for each phase, property pair.

Loaded from the file are density and pressure as functions of temperature, and a thermodynamic model can be used to build the

Definition at line 735 of file superancillary.h.

#include <superancillary.h>

Public Member Functions

 SuperAncillary (const nlohmann::json &j)
 
 SuperAncillary (const std::string &s)
 
const auto & get_approx1d (char k, short Q) const
 
const auto & get_invlnp ()
 Get a const reference to the inverse approximation for T(ln(p)) More...
 
const double get_pmin () const
 Get the minimum pressure in Pa. More...
 
const double get_pmax () const
 Get the maximum pressure in Pa. More...
 
const double get_Tmin () const
 Get the minimum temperature in K. More...
 
const double get_Tcrit_num () const
 Get the numerical critical temperature in K. More...
 
const double get_rhocrit_num () const
 Get the numerical critical density in mol/m^3. More...
 
void add_variable (char var, const std::function< double(double, double)> &caller)
 
double eval_sat (double T, char k, short Q) const
 
template<typename Container >
void eval_sat_many (const Container &T, char k, short Q, Container &y) const
 
template<typename Container >
void eval_sat_manyC (const Container T[], std::size_t N, char k, short Q, Container y[]) const
 
auto solve_for_T (double propval, char k, bool Q, unsigned int bits=64, unsigned int max_iter=100, double boundsftol=1e-13) const
 
auto get_vaporquality (double T, double propval, char k) const
 
auto get_T_from_p (double p)
 Use the inverted pressure superancillary to calculate temperature given pressure. More...
 
auto get_yval (double T, double q, char k) const
 Return the evaluated value of the thermodynamic variable, given the temperature and vapor quality. More...
 
template<typename Container >
void get_yval_many (const Container &T, char k, const Container &q, Container &y) const
 A vectorized version of get_yval for profiling in Python. More...
 
auto get_all_intersections (const char k, const double val, unsigned int bits, std::size_t max_iter, double boundsftol) const
 
std::optional< SuperAncillaryTwoPhaseSolutioniterate_for_Tq_XY (double Tmin, double Tmax, char ch1, double val1, char ch2, double val2, unsigned int bits, std::size_t max_iter, double boundsftol) const
 Iterate to find a value of temperature and vapor quality corresponding to the two given thermodynamic variables, if such a solution exists. This is the lower-level function used by the solve_XX methods. More...
 
std::optional< SuperAncillaryTwoPhaseSolutionsolve_for_Tq_DX (const double rho, const double propval, const char k, unsigned int bits, std::size_t max_iter, double boundsftol) const
 
template<typename Container >
void solve_for_Tq_DX_many (const Container &rho, const Container &propval, const char k, unsigned int bits, std::size_t max_iter, double boundsftol, Container &T, Container &q, Container &count)
 A vectorize version of solve_for_Tq_DX for use in the Python interface for profiling. More...
 

Constructor & Destructor Documentation

◆ SuperAncillary() [1/2]

template<typename ArrayType = Eigen::ArrayXd>
CoolProp::superancillary::SuperAncillary< ArrayType >::SuperAncillary ( const nlohmann::json &  j)
inline

Reading in a data structure in the JSON format of https://pubs.aip.org/aip/jpr/article/53/1/013102/3270194 which includes sets of Chebyshev expansions for rhoL, rhoV, and p

Definition at line 817 of file superancillary.h.

◆ SuperAncillary() [2/2]

template<typename ArrayType = Eigen::ArrayXd>
CoolProp::superancillary::SuperAncillary< ArrayType >::SuperAncillary ( const std::string &  s)
inline

Load the superancillary with the data passed in as a string blob. This constructor delegates directly to the the one that consumes JSON

Parameters
sThe string-encoded JSON data for the superancillaries

Definition at line 831 of file superancillary.h.

Member Function Documentation

◆ add_variable()

template<typename ArrayType = Eigen::ArrayXd>
void CoolProp::superancillary::SuperAncillary< ArrayType >::add_variable ( char  var,
const std::function< double(double, double)> &  caller 
)
inline

Using the provided function that gives y(T, rho), build the ancillaries for this variable based on the ancillaries for rhoL and rhoV

Parameters
varThe key for the property (H,S,U)
callerA function that takes temperature and molar density and returns the property of interest, molar enthalpy in the case of H, etc.

Definition at line 882 of file superancillary.h.

◆ eval_sat()

template<typename ArrayType = Eigen::ArrayXd>
double CoolProp::superancillary::SuperAncillary< ArrayType >::eval_sat ( double  T,
char  k,
short  Q 
) const
inline

Given the value of Q in {0,1}, evaluate one of the the ChebyshevApproximation1D

Parameters
TTemperature, in K
kProperty key, in {D,P,H,S,U}
QVapor quality, in {0,1}

Definition at line 928 of file superancillary.h.

◆ eval_sat_many()

template<typename ArrayType = Eigen::ArrayXd>
template<typename Container >
void CoolProp::superancillary::SuperAncillary< ArrayType >::eval_sat_many ( const Container &  T,
char  k,
short  Q,
Container &  y 
) const
inline

A vectorized version of eval_sat for wrapping in Python interface and profiling

Definition at line 941 of file superancillary.h.

◆ eval_sat_manyC()

template<typename ArrayType = Eigen::ArrayXd>
template<typename Container >
void CoolProp::superancillary::SuperAncillary< ArrayType >::eval_sat_manyC ( const Container  T[],
std::size_t  N,
char  k,
short  Q,
Container  y[] 
) const
inline

A vectorized version of eval_sat for wrapping in Python interface and profiling

Definition at line 953 of file superancillary.h.

◆ get_all_intersections()

template<typename ArrayType = Eigen::ArrayXd>
auto CoolProp::superancillary::SuperAncillary< ArrayType >::get_all_intersections ( const char  k,
const double  val,
unsigned int  bits,
std::size_t  max_iter,
double  boundsftol 
) const
inline

Determine all the values of temperature that correspond to intersections with the superancillary function, for both the vapor and liquid phases

Parameters
kProperty key, in {D,P,H,S,U}
valValue of the thermodynamic variable
bitspassed to toms748 algorithm
max_iterMaximum allowed number of function calls
boundsftolA functional value stopping condition to test on the endpoints

Definition at line 1056 of file superancillary.h.

◆ get_approx1d()

template<typename ArrayType = Eigen::ArrayXd>
const auto & CoolProp::superancillary::SuperAncillary< ArrayType >::get_approx1d ( char  k,
short  Q 
) const
inline

Get a const reference to a ChebyshevApproximation1D

Parameters
kThe key for the property (D,S,H,P,U)
QThe vapor quality, either 0 or 1

Definition at line 838 of file superancillary.h.

◆ get_invlnp()

template<typename ArrayType = Eigen::ArrayXd>
const auto & CoolProp::superancillary::SuperAncillary< ArrayType >::get_invlnp ( )
inline

Get a const reference to the inverse approximation for T(ln(p))

Definition at line 857 of file superancillary.h.

◆ get_pmax()

template<typename ArrayType = Eigen::ArrayXd>
const double CoolProp::superancillary::SuperAncillary< ArrayType >::get_pmax ( ) const
inline

Get the maximum pressure in Pa.

Definition at line 869 of file superancillary.h.

◆ get_pmin()

template<typename ArrayType = Eigen::ArrayXd>
const double CoolProp::superancillary::SuperAncillary< ArrayType >::get_pmin ( ) const
inline

Get the minimum pressure in Pa.

Definition at line 867 of file superancillary.h.

◆ get_rhocrit_num()

template<typename ArrayType = Eigen::ArrayXd>
const double CoolProp::superancillary::SuperAncillary< ArrayType >::get_rhocrit_num ( ) const
inline

Get the numerical critical density in mol/m^3.

Definition at line 875 of file superancillary.h.

◆ get_T_from_p()

template<typename ArrayType = Eigen::ArrayXd>
auto CoolProp::superancillary::SuperAncillary< ArrayType >::get_T_from_p ( double  p)
inline

Use the inverted pressure superancillary to calculate temperature given pressure.

Parameters
pThe pressure (not its logarithm!), in Pa
Returns
T The temperature, in K

Definition at line 999 of file superancillary.h.

◆ get_Tcrit_num()

template<typename ArrayType = Eigen::ArrayXd>
const double CoolProp::superancillary::SuperAncillary< ArrayType >::get_Tcrit_num ( ) const
inline

Get the numerical critical temperature in K.

Definition at line 873 of file superancillary.h.

◆ get_Tmin()

template<typename ArrayType = Eigen::ArrayXd>
const double CoolProp::superancillary::SuperAncillary< ArrayType >::get_Tmin ( ) const
inline

Get the minimum temperature in K.

Definition at line 871 of file superancillary.h.

◆ get_vaporquality()

template<typename ArrayType = Eigen::ArrayXd>
auto CoolProp::superancillary::SuperAncillary< ArrayType >::get_vaporquality ( double  T,
double  propval,
char  k 
) const
inline

Get the non-iterative vapor quality q given the temperature T and the value of the thermodynamic variable

Parameters
TTemperature, in K
propvalThe value of the given thermodynamic variable
kProperty key, in {D,P,H,S,U}
Returns
The non-iterative vapor quality based on the values from the superancillary functions

Definition at line 980 of file superancillary.h.

◆ get_yval()

template<typename ArrayType = Eigen::ArrayXd>
auto CoolProp::superancillary::SuperAncillary< ArrayType >::get_yval ( double  T,
double  q,
char  k 
) const
inline

Return the evaluated value of the thermodynamic variable, given the temperature and vapor quality.

Parameters
TTemperature, in K
qVapor quality, in [0,1]
kProperty key, in {D,P,H,S,U}

Definition at line 1010 of file superancillary.h.

◆ get_yval_many()

template<typename ArrayType = Eigen::ArrayXd>
template<typename Container >
void CoolProp::superancillary::SuperAncillary< ArrayType >::get_yval_many ( const Container &  T,
char  k,
const Container &  q,
Container &  y 
) const
inline

A vectorized version of get_yval for profiling in Python.

Definition at line 1028 of file superancillary.h.

◆ iterate_for_Tq_XY()

template<typename ArrayType = Eigen::ArrayXd>
std::optional< SuperAncillaryTwoPhaseSolution > CoolProp::superancillary::SuperAncillary< ArrayType >::iterate_for_Tq_XY ( double  Tmin,
double  Tmax,
char  ch1,
double  val1,
char  ch2,
double  val2,
unsigned int  bits,
std::size_t  max_iter,
double  boundsftol 
) const
inline

Iterate to find a value of temperature and vapor quality corresponding to the two given thermodynamic variables, if such a solution exists. This is the lower-level function used by the solve_XX methods.

Note
The temperature range must bound the solution, you might need to call get_all_intersections and parse its solutions to construct bounded intervals
Parameters
TminMinimum temperature, in K
TmaxMaximum temperature, in K
ch1The key for the first variable, in {T,D,P,H,S,U}
val1The value for the first variable
ch2The key for the second variable, in {T,D,P,H,S,U}
val2The value for the second variable
bitspassed to toms748 algorithm
max_iterMaximum allowed number of function calls
boundsftolA functional value stopping condition to test on the endpoints

Definition at line 1084 of file superancillary.h.

◆ solve_for_T()

template<typename ArrayType = Eigen::ArrayXd>
auto CoolProp::superancillary::SuperAncillary< ArrayType >::solve_for_T ( double  propval,
char  k,
bool  Q,
unsigned int  bits = 64,
unsigned int  max_iter = 100,
double  boundsftol = 1e-13 
) const
inline

A convenience function to pass off to the ChebyshevApproximation1D and do an inversion calculation for a value of the variable for a saturated state

Parameters
propvalThe value of the property
kProperty key, in {D,P,H,S,U}
QVapor quality, in {0,1}
bitspassed to toms748 algorithm
max_iterMaximum allowed number of function calls
boundsftolA functional value stopping condition to test on the endpoints

Definition at line 970 of file superancillary.h.

◆ solve_for_Tq_DX()

template<typename ArrayType = Eigen::ArrayXd>
std::optional< SuperAncillaryTwoPhaseSolution > CoolProp::superancillary::SuperAncillary< ArrayType >::solve_for_Tq_DX ( const double  rho,
const double  propval,
const char  k,
unsigned int  bits,
std::size_t  max_iter,
double  boundsftol 
) const
inline

Given a saturated density and another property other than T, solve for the temperature and vapor quality

Parameters
rhoThe molar density
propvalThe value of the other property
kProperty key, in {D,P,H,S,U}
bitspassed to toms748 algorithm
max_iterMaximum allowed number of function calls
boundsftolA functional value stopping condition to test on the endpoints

Definition at line 1142 of file superancillary.h.

◆ solve_for_Tq_DX_many()

template<typename ArrayType = Eigen::ArrayXd>
template<typename Container >
void CoolProp::superancillary::SuperAncillary< ArrayType >::solve_for_Tq_DX_many ( const Container &  rho,
const Container &  propval,
const char  k,
unsigned int  bits,
std::size_t  max_iter,
double  boundsftol,
Container &  T,
Container &  q,
Container &  count 
)
inline

A vectorize version of solve_for_Tq_DX for use in the Python interface for profiling.

Definition at line 1168 of file superancillary.h.


The documentation for this class was generated from the following file: