CoolProp
6.6.1dev
An open-source fluid property and humid air property database
|
Definition at line 14 of file IncompressibleBackend.h.
#include <IncompressibleBackend.h>
Public Member Functions | |
IncompressibleBackend () | |
virtual | ~IncompressibleBackend () |
std::string | backend_name (void) |
IncompressibleBackend (IncompressibleFluid *fluid) | |
IncompressibleBackend (const std::string &fluid_name) | |
IncompressibleBackend (const std::vector< std::string > &component_names) | |
bool | using_mole_fractions (void) |
bool | using_mass_fractions (void) |
bool | using_volu_fractions (void) |
void | update (CoolProp::input_pairs input_pair, double value1, double value2) |
Updating function for incompressible fluid. More... | |
std::string | fluid_param_string (const std::string &ParamName) |
Return a string from the backend for the mixture/fluid - backend dependent - could be CAS #, name, etc. More... | |
bool | clear () |
Clear all the cached values. More... | |
void | set_reference_state (double T0=20+273.15, double p0=101325, double x0=0.0, double h0=0.0, double s0=0.0) |
Update the reference values and clear the state. More... | |
void | set_mole_fractions (const std::vector< CoolPropDbl > &mole_fractions) |
Set the mole fractions. More... | |
const std::vector< CoolPropDbl > & | get_mole_fractions (void) |
Get the mole fractions of the fluid. More... | |
void | set_mass_fractions (const std::vector< CoolPropDbl > &mass_fractions) |
Set the mass fractions. More... | |
void | set_volu_fractions (const std::vector< CoolPropDbl > &volu_fractions) |
Set the volume fractions. More... | |
void | check_status () |
Check if the mole fractions have been set, etc. More... | |
double | rhomass (void) |
Return the mass density in kg/m^3. More... | |
double | hmass (void) |
Return the mass enthalpy in J/kg. More... | |
double | smass (void) |
Return the molar entropy in J/mol/K. More... | |
double | umass (void) |
Return the molar internal energy in J/mol. More... | |
double | cmass (void) |
Return the mass constant pressure specific heat in J/kg/K. More... | |
double | drhodTatPx (void) |
double | dsdTatPx (void) |
double | dhdTatPx (void) |
double | dsdTatPxdT (void) |
double | dhdTatPxdT (void) |
double | dsdpatTx (void) |
double | dhdpatTx (void) |
double | T_ref (void) |
Return the temperature in K. More... | |
double | p_ref (void) |
Return the pressure in Pa. More... | |
double | x_ref (void) |
Return the composition. More... | |
double | h_ref (void) |
Return the mass enthalpy in J/kg. More... | |
double | s_ref (void) |
Return the molar entropy in J/mol/K. More... | |
double | hmass_ref (void) |
Return the mass enthalpy in J/kg. More... | |
double | smass_ref (void) |
Return the molar entropy in J/mol/K. More... | |
CoolPropDbl | DmassP_flash (CoolPropDbl rhomass, CoolPropDbl p) |
Calculate T given pressure and density. More... | |
CoolPropDbl | HmassP_flash (CoolPropDbl hmass, CoolPropDbl p) |
Calculate T given pressure and enthalpy. More... | |
CoolPropDbl | PSmass_flash (CoolPropDbl p, CoolPropDbl smass) |
Calculate T given pressure and entropy. More... | |
CoolPropDbl | calc_rhomass (void) |
We start with the functions that do not need a reference state. More... | |
CoolPropDbl | calc_cmass (void) |
CoolPropDbl | calc_cpmass (void) |
CoolPropDbl | calc_cvmass (void) |
CoolPropDbl | calc_viscosity (void) |
Using this backend, calculate the viscosity in Pa-s. More... | |
CoolPropDbl | calc_conductivity (void) |
Using this backend, calculate the thermal conductivity in W/m/K. More... | |
CoolPropDbl | calc_T_freeze (void) |
CoolPropDbl | calc_melting_line (int param, int given, CoolPropDbl value) |
We start with the functions that do not need a reference state. More... | |
CoolPropDbl | calc_umass (void) |
CoolPropDbl | calc_hmass (void) |
... and continue with the ones that depend on reference conditions. More... | |
CoolPropDbl | calc_smass (void) |
CoolPropDbl | raw_calc_hmass (double T, double p, double x) |
Functions that can be used with the solver, they miss the reference values! More... | |
CoolPropDbl | raw_calc_smass (double T, double p, double x) |
CoolPropDbl | calc_Tmax (void) |
Constants from the fluid object. More... | |
CoolPropDbl | calc_Tmin (void) |
Using this backend, calculate the minimum temperature in K. More... | |
CoolPropDbl | calc_fraction_min (void) |
Get the minimum fraction (mole, mass, volume) for incompressible fluid. More... | |
CoolPropDbl | calc_fraction_max (void) |
Get the maximum fraction (mole, mass, volume) for incompressible fluid. More... | |
std::string | calc_name (void) |
Using this backend, get the name of the fluid. More... | |
std::string | calc_description (void) |
Using this backend, get the description of the fluid. More... | |
Public Member Functions inherited from CoolProp::AbstractState | |
AbstractState () | |
virtual | ~AbstractState () |
void | set_T (CoolPropDbl T) |
Set the internal variable T without a flash call (expert use only!) More... | |
virtual void | set_reference_stateS (const std::string &reference_state) |
Set the reference state based on a string representation. More... | |
virtual void | set_reference_stateD (double T, double rhomolar, double hmolar0, double smolar0) |
std::vector< CoolPropDbl > | mole_fractions_liquid (void) |
Get the mole fractions of the equilibrium liquid phase. More... | |
std::vector< double > | mole_fractions_liquid_double (void) |
Get the mole fractions of the equilibrium liquid phase (but as a double for use in SWIG wrapper) More... | |
std::vector< CoolPropDbl > | mole_fractions_vapor (void) |
Get the mole fractions of the equilibrium vapor phase. More... | |
std::vector< double > | mole_fractions_vapor_double (void) |
Get the mole fractions of the equilibrium vapor phase (but as a double for use in SWIG wrapper) More... | |
virtual const std::vector< CoolPropDbl > | get_mass_fractions (void) |
Get the mass fractions of the fluid. More... | |
virtual void | update_with_guesses (CoolProp::input_pairs input_pair, double Value1, double Value2, const GuessesStructure &guesses) |
virtual bool | available_in_high_level (void) |
std::vector< std::string > | fluid_names (void) |
Return a vector of strings of the fluid names that are in use. More... | |
virtual const double | get_fluid_constant (std::size_t i, parameters param) const |
virtual void | set_binary_interaction_double (const std::string &CAS1, const std::string &CAS2, const std::string ¶meter, const double value) |
Set binary mixture floating point parameter (EXPERT USE ONLY!!!) More... | |
virtual void | set_binary_interaction_double (const std::size_t i, const std::size_t j, const std::string ¶meter, const double value) |
Set binary mixture floating point parameter (EXPERT USE ONLY!!!) More... | |
virtual void | set_binary_interaction_string (const std::string &CAS1, const std::string &CAS2, const std::string ¶meter, const std::string &value) |
Set binary mixture string parameter (EXPERT USE ONLY!!!) More... | |
virtual void | set_binary_interaction_string (const std::size_t i, const std::size_t j, const std::string ¶meter, const std::string &value) |
Set binary mixture string parameter (EXPERT USE ONLY!!!) More... | |
virtual double | get_binary_interaction_double (const std::string &CAS1, const std::string &CAS2, const std::string ¶meter) |
Get binary mixture double value (EXPERT USE ONLY!!!) More... | |
virtual double | get_binary_interaction_double (const std::size_t i, const std::size_t j, const std::string ¶meter) |
Get binary mixture double value (EXPERT USE ONLY!!!) More... | |
virtual std::string | get_binary_interaction_string (const std::string &CAS1, const std::string &CAS2, const std::string ¶meter) |
Get binary mixture string value (EXPERT USE ONLY!!!) More... | |
virtual void | apply_simple_mixing_rule (std::size_t i, std::size_t j, const std::string &model) |
Apply a simple mixing rule (EXPERT USE ONLY!!!) More... | |
virtual void | set_cubic_alpha_C (const size_t i, const std::string ¶meter, const double c1, const double c2, const double c3) |
Set the cubic alpha function's constants: More... | |
virtual void | set_fluid_parameter_double (const size_t i, const std::string ¶meter, const double value) |
Set fluid parameter (currently the volume translation parameter for cubic) More... | |
virtual double | get_fluid_parameter_double (const size_t i, const std::string ¶meter) |
Double fluid parameter (currently the volume translation parameter for cubic) More... | |
virtual bool | clear_comp_change () |
When the composition changes, clear all cached values that are only dependent on composition, but not the thermodynamic state. More... | |
virtual const CoolProp::SimpleState & | get_reducing_state () |
const CoolProp::SimpleState & | get_state (const std::string &state) |
Get a desired state point - backend dependent. More... | |
double | Tmin (void) |
Get the minimum temperature in K. More... | |
double | Tmax (void) |
Get the maximum temperature in K. More... | |
double | pmax (void) |
Get the maximum pressure in Pa. More... | |
double | Ttriple (void) |
Get the triple point temperature in K. More... | |
phases | phase (void) |
Get the phase of the state. More... | |
void | specify_phase (phases phase) |
Specify the phase for all further calculations with this state class. More... | |
void | unspecify_phase (void) |
Unspecify the phase and go back to calculating it based on the inputs. More... | |
double | T_critical (void) |
Return the critical temperature in K. More... | |
double | p_critical (void) |
Return the critical pressure in Pa. More... | |
double | rhomolar_critical (void) |
Return the critical molar density in mol/m^3. More... | |
double | rhomass_critical (void) |
Return the critical mass density in kg/m^3. More... | |
std::vector< CriticalState > | all_critical_points (void) |
Return the vector of critical points, including points that are unstable or correspond to negative pressure. More... | |
void | build_spinodal () |
Construct the spinodal curve for the mixture (or pure fluid) More... | |
SpinodalData | get_spinodal_data () |
Get the data from the spinodal curve constructed in the call to build_spinodal() More... | |
void | criticality_contour_values (double &L1star, double &M1star) |
Calculate the criticality contour values \(\mathcal{L}_1^*\) and \(\mathcal{M}_1^*\). More... | |
double | tangent_plane_distance (const double T, const double p, const std::vector< double > &w, const double rhomolar_guess=-1) |
double | T_reducing (void) |
Return the reducing point temperature in K. More... | |
double | rhomolar_reducing (void) |
Return the molar density at the reducing point in mol/m^3. More... | |
double | rhomass_reducing (void) |
Return the mass density at the reducing point in kg/m^3. More... | |
double | p_triple (void) |
Return the triple point pressure in Pa. More... | |
std::string | name () |
Return the name - backend dependent. More... | |
std::string | description () |
Return the description - backend dependent. More... | |
double | dipole_moment () |
Return the dipole moment in C-m (1 D = 3.33564e-30 C-m) More... | |
double | keyed_output (parameters key) |
Retrieve a value by key. More... | |
double | trivial_keyed_output (parameters key) |
A trivial keyed output like molar mass that does not depend on the state. More... | |
double | saturated_liquid_keyed_output (parameters key) |
Get an output from the saturated liquid state by key. More... | |
double | saturated_vapor_keyed_output (parameters key) |
Get an output from the saturated vapor state by key. More... | |
double | T (void) |
Return the temperature in K. More... | |
double | rhomolar (void) |
Return the molar density in mol/m^3. More... | |
double | rhomass (void) |
Return the mass density in kg/m^3. More... | |
double | p (void) |
Return the pressure in Pa. More... | |
double | Q (void) |
Return the vapor quality (mol/mol); Q = 0 for saturated liquid. More... | |
double | tau (void) |
Return the reciprocal of the reduced temperature ( \(\tau = T_c/T\)) More... | |
double | delta (void) |
Return the reduced density ( \(\delta = \rho/\rho_c\)) More... | |
double | molar_mass (void) |
Return the molar mass in kg/mol. More... | |
double | acentric_factor (void) |
Return the acentric factor. More... | |
double | gas_constant (void) |
Return the mole-fraction weighted gas constant in J/mol/K. More... | |
double | Bvirial (void) |
Return the B virial coefficient. More... | |
double | dBvirial_dT (void) |
Return the derivative of the B virial coefficient with respect to temperature. More... | |
double | Cvirial (void) |
Return the C virial coefficient. More... | |
double | dCvirial_dT (void) |
Return the derivative of the C virial coefficient with respect to temperature. More... | |
double | compressibility_factor (void) |
Return the compressibility factor \( Z = p/(rho R T) \). More... | |
double | hmolar (void) |
Return the molar enthalpy in J/mol. More... | |
double | hmolar_residual (void) |
Return the residual molar enthalpy in J/mol. More... | |
double | hmass (void) |
Return the mass enthalpy in J/kg. More... | |
double | hmolar_excess (void) |
Return the excess molar enthalpy in J/mol. More... | |
double | hmass_excess (void) |
Return the excess mass enthalpy in J/kg. More... | |
double | smolar (void) |
Return the molar entropy in J/mol/K. More... | |
double | smolar_residual (void) |
Return the residual molar entropy (as a function of temperature and density) in J/mol/K. More... | |
double | neff (void) |
Return the effective hardness of interaction. More... | |
double | smass (void) |
Return the molar entropy in J/kg/K. More... | |
double | smolar_excess (void) |
Return the molar entropy in J/mol/K. More... | |
double | smass_excess (void) |
Return the molar entropy in J/kg/K. More... | |
double | umolar (void) |
Return the molar internal energy in J/mol. More... | |
double | umass (void) |
Return the mass internal energy in J/kg. More... | |
double | umolar_excess (void) |
Return the excess internal energy in J/mol. More... | |
double | umass_excess (void) |
Return the excess internal energy in J/kg. More... | |
double | cpmolar (void) |
Return the molar constant pressure specific heat in J/mol/K. More... | |
double | cpmass (void) |
Return the mass constant pressure specific heat in J/kg/K. More... | |
double | cp0molar (void) |
Return the molar constant pressure specific heat for ideal gas part only in J/mol/K. More... | |
double | cp0mass (void) |
Return the mass constant pressure specific heat for ideal gas part only in J/kg/K. More... | |
double | cvmolar (void) |
Return the molar constant volume specific heat in J/mol/K. More... | |
double | cvmass (void) |
Return the mass constant volume specific heat in J/kg/K. More... | |
double | gibbsmolar (void) |
Return the Gibbs energy in J/mol. More... | |
double | gibbsmolar_residual (void) |
Return the residual Gibbs energy in J/mol. More... | |
double | gibbsmass (void) |
Return the Gibbs energy in J/kg. More... | |
double | gibbsmolar_excess (void) |
Return the excess Gibbs energy in J/mol. More... | |
double | gibbsmass_excess (void) |
Return the excess Gibbs energy in J/kg. More... | |
double | helmholtzmolar (void) |
Return the Helmholtz energy in J/mol. More... | |
double | helmholtzmass (void) |
Return the Helmholtz energy in J/kg. More... | |
double | helmholtzmolar_excess (void) |
Return the excess Helmholtz energy in J/mol. More... | |
double | helmholtzmass_excess (void) |
Return the excess Helmholtz energy in J/kg. More... | |
double | volumemolar_excess (void) |
Return the excess volume in m^3/mol. More... | |
double | volumemass_excess (void) |
Return the excess volume in m^3/kg. More... | |
double | speed_sound (void) |
Return the speed of sound in m/s. More... | |
double | isothermal_compressibility (void) |
Return the isothermal compressibility \( \kappa = -\frac{1}{v}\left.\frac{\partial v}{\partial p}\right|_T=\frac{1}{\rho}\left.\frac{\partial \rho}{\partial p}\right|_T\) in 1/Pa. More... | |
double | isobaric_expansion_coefficient (void) |
Return the isobaric expansion coefficient \( \beta = \frac{1}{v}\left.\frac{\partial v}{\partial T}\right|_p = -\frac{1}{\rho}\left.\frac{\partial \rho}{\partial T}\right|_p\) in 1/K. More... | |
double | isentropic_expansion_coefficient (void) |
Return the isentropic expansion coefficient \( \kappa_s = -\frac{c_p}{c_v}\frac{v}{p}\left.\frac{\partial p}{\partial v}\right|_T = \frac{\rho}{p}\left.\frac{\partial p}{\partial \rho}\right|_s\). More... | |
double | fugacity_coefficient (std::size_t i) |
Return the fugacity coefficient of the i-th component of the mixture. More... | |
std::vector< double > | fugacity_coefficients () |
Return a vector of the fugacity coefficients for all components in the mixture. More... | |
double | fugacity (std::size_t i) |
Return the fugacity of the i-th component of the mixture. More... | |
double | chemical_potential (std::size_t i) |
Return the chemical potential of the i-th component of the mixture. More... | |
double | fundamental_derivative_of_gas_dynamics (void) |
Return the fundamental derivative of gas dynamics \( \Gamma \). More... | |
double | PIP () |
Return the phase identification parameter (PIP) of G. Venkatarathnam and L.R. Oellrich, "Identification of the phase of a fluid using partial derivatives of pressure, volume, and temperature without reference to saturation properties: Applications in phase equilibria calculations". More... | |
void | true_critical_point (double &T, double &rho) |
Calculate the "true" critical point for pure fluids where dpdrho|T and d2p/drho2|T are equal to zero. More... | |
void | ideal_curve (const std::string &type, std::vector< double > &T, std::vector< double > &p) |
Calculate an ideal curve for a pure fluid. More... | |
CoolPropDbl | first_partial_deriv (parameters Of, parameters Wrt, parameters Constant) |
The first partial derivative in homogeneous phases. More... | |
CoolPropDbl | second_partial_deriv (parameters Of1, parameters Wrt1, parameters Constant1, parameters Wrt2, parameters Constant2) |
The second partial derivative in homogeneous phases. More... | |
CoolPropDbl | first_saturation_deriv (parameters Of1, parameters Wrt1) |
The first partial derivative along the saturation curve. More... | |
CoolPropDbl | second_saturation_deriv (parameters Of1, parameters Wrt1, parameters Wrt2) |
The second partial derivative along the saturation curve. More... | |
double | first_two_phase_deriv (parameters Of, parameters Wrt, parameters Constant) |
Calculate the first "two-phase" derivative as described by Thorade and Sadaat, EAS, 2013. More... | |
double | second_two_phase_deriv (parameters Of, parameters Wrt1, parameters Constant1, parameters Wrt2, parameters Constant2) |
Calculate the second "two-phase" derivative as described by Thorade and Sadaat, EAS, 2013. More... | |
double | first_two_phase_deriv_splined (parameters Of, parameters Wrt, parameters Constant, double x_end) |
Calculate the first "two-phase" derivative as described by Thorade and Sadaat, EAS, 2013. More... | |
void | build_phase_envelope (const std::string &type="") |
Construct the phase envelope for a mixture. More... | |
const CoolProp::PhaseEnvelopeData & | get_phase_envelope_data () |
After having calculated the phase envelope, return the phase envelope data. More... | |
virtual bool | has_melting_line (void) |
Return true if the fluid has a melting line - default is false, but can be re-implemented by derived class. More... | |
double | melting_line (int param, int given, double value) |
double | saturation_ancillary (parameters param, int Q, parameters given, double value) |
double | viscosity (void) |
Return the viscosity in Pa-s. More... | |
void | viscosity_contributions (CoolPropDbl &dilute, CoolPropDbl &initial_density, CoolPropDbl &residual, CoolPropDbl &critical) |
Return the viscosity contributions, each in Pa-s. More... | |
double | conductivity (void) |
Return the thermal conductivity in W/m/K. More... | |
void | conductivity_contributions (CoolPropDbl &dilute, CoolPropDbl &initial_density, CoolPropDbl &residual, CoolPropDbl &critical) |
Return the thermal conductivity contributions, each in W/m/K. More... | |
double | surface_tension (void) |
Return the surface tension in N/m. More... | |
double | Prandtl (void) |
Return the Prandtl number (dimensionless) More... | |
void | conformal_state (const std::string &reference_fluid, CoolPropDbl &T, CoolPropDbl &rhomolar) |
Find the conformal state needed for ECS. More... | |
void | change_EOS (const std::size_t i, const std::string &EOS_name) |
Change the equation of state for a given component to a specified EOS. More... | |
CoolPropDbl | alpha0 (void) |
Return the term \( \alpha^0 \). More... | |
CoolPropDbl | dalpha0_dDelta (void) |
Return the term \( \alpha^0_{\delta} \). More... | |
CoolPropDbl | dalpha0_dTau (void) |
Return the term \( \alpha^0_{\tau} \). More... | |
CoolPropDbl | d2alpha0_dDelta2 (void) |
Return the term \( \alpha^0_{\delta\delta} \). More... | |
CoolPropDbl | d2alpha0_dDelta_dTau (void) |
Return the term \( \alpha^0_{\delta\tau} \). More... | |
CoolPropDbl | d2alpha0_dTau2 (void) |
Return the term \( \alpha^0_{\tau\tau} \). More... | |
CoolPropDbl | d3alpha0_dTau3 (void) |
Return the term \( \alpha^0_{\tau\tau\tau} \). More... | |
CoolPropDbl | d3alpha0_dDelta_dTau2 (void) |
Return the term \( \alpha^0_{\delta\tau\tau} \). More... | |
CoolPropDbl | d3alpha0_dDelta2_dTau (void) |
Return the term \( \alpha^0_{\delta\delta\tau} \). More... | |
CoolPropDbl | d3alpha0_dDelta3 (void) |
Return the term \( \alpha^0_{\delta\delta\delta} \). More... | |
CoolPropDbl | alphar (void) |
Return the term \( \alpha^r \). More... | |
CoolPropDbl | dalphar_dDelta (void) |
Return the term \( \alpha^r_{\delta} \). More... | |
CoolPropDbl | dalphar_dTau (void) |
Return the term \( \alpha^r_{\tau} \). More... | |
CoolPropDbl | d2alphar_dDelta2 (void) |
Return the term \( \alpha^r_{\delta\delta} \). More... | |
CoolPropDbl | d2alphar_dDelta_dTau (void) |
Return the term \( \alpha^r_{\delta\tau} \). More... | |
CoolPropDbl | d2alphar_dTau2 (void) |
Return the term \( \alpha^r_{\tau\tau} \). More... | |
CoolPropDbl | d3alphar_dDelta3 (void) |
Return the term \( \alpha^r_{\delta\delta\delta} \). More... | |
CoolPropDbl | d3alphar_dDelta2_dTau (void) |
Return the term \( \alpha^r_{\delta\delta\tau} \). More... | |
CoolPropDbl | d3alphar_dDelta_dTau2 (void) |
Return the term \( \alpha^r_{\delta\tau\tau} \). More... | |
CoolPropDbl | d3alphar_dTau3 (void) |
Return the term \( \alpha^r_{\tau\tau\tau} \). More... | |
CoolPropDbl | d4alphar_dDelta4 (void) |
Return the term \( \alpha^r_{\delta\delta\delta\delta} \). More... | |
CoolPropDbl | d4alphar_dDelta3_dTau (void) |
Return the term \( \alpha^r_{\delta\delta\delta\tau} \). More... | |
CoolPropDbl | d4alphar_dDelta2_dTau2 (void) |
Return the term \( \alpha^r_{\delta\delta\tau\tau} \). More... | |
CoolPropDbl | d4alphar_dDelta_dTau3 (void) |
Return the term \( \alpha^r_{\delta\tau\tau\tau} \). More... | |
CoolPropDbl | d4alphar_dTau4 (void) |
Return the term \( \alpha^r_{\tau\tau\tau\tau} \). More... | |
Protected Member Functions | |
void | set_fractions (const std::vector< CoolPropDbl > &fractions) |
Set the fractions. More... | |
CoolPropDbl | calc_first_partial_deriv (parameters Of, parameters Wrt, parameters Constant) |
Calculate the first partial derivative for the desired derivative. More... | |
double | calc_drhodTatPx (double T, double p, double x) |
Partial derivative of density with respect to temperature at constant pressure and composition. More... | |
double | calc_dsdTatPx (double T, double p, double x) |
Partial derivative of entropy with respect to temperature at constant pressure and composition. More... | |
double | calc_dhdTatPx (double T, double p, double x) |
Partial derivative of enthalpy with respect to temperature at constant pressure and composition. More... | |
double | calc_dsdTatPxdT (double T, double p, double x) |
double | calc_dhdTatPxdT (double T, double p, double x) |
double | calc_dsdpatTx (double rho, double drhodTatPx) |
Partial derivative of entropy with respect to pressure at constant temperature and composition. More... | |
double | calc_dhdpatTx (double T, double rho, double drhodTatPx) |
Partial derivative of enthalpy with respect to pressure at constant temperature and composition. More... | |
Protected Member Functions inherited from CoolProp::AbstractState | |
bool | isSupercriticalPhase (void) |
bool | isHomogeneousPhase (void) |
bool | isTwoPhase (void) |
virtual CoolPropDbl | calc_hmolar (void) |
Using this backend, calculate the molar enthalpy in J/mol. More... | |
virtual CoolPropDbl | calc_hmolar_residual (void) |
Using this backend, calculate the residual molar enthalpy in J/mol. More... | |
virtual CoolPropDbl | calc_smolar (void) |
Using this backend, calculate the molar entropy in J/mol/K. More... | |
virtual CoolPropDbl | calc_smolar_residual (void) |
Using this backend, calculate the residual molar entropy in J/mol/K. More... | |
virtual CoolPropDbl | calc_neff (void) |
Using this backend, calculate effective hardness of interaction. More... | |
virtual CoolPropDbl | calc_umolar (void) |
Using this backend, calculate the molar internal energy in J/mol. More... | |
virtual CoolPropDbl | calc_cpmolar (void) |
Using this backend, calculate the molar constant-pressure specific heat in J/mol/K. More... | |
virtual CoolPropDbl | calc_cpmolar_idealgas (void) |
Using this backend, calculate the ideal gas molar constant-pressure specific heat in J/mol/K. More... | |
virtual CoolPropDbl | calc_cvmolar (void) |
Using this backend, calculate the molar constant-volume specific heat in J/mol/K. More... | |
virtual CoolPropDbl | calc_gibbsmolar (void) |
Using this backend, calculate the molar Gibbs function in J/mol. More... | |
virtual CoolPropDbl | calc_gibbsmolar_residual (void) |
Using this backend, calculate the residual molar Gibbs function in J/mol. More... | |
virtual CoolPropDbl | calc_helmholtzmolar (void) |
Using this backend, calculate the molar Helmholtz energy in J/mol. More... | |
virtual CoolPropDbl | calc_speed_sound (void) |
Using this backend, calculate the speed of sound in m/s. More... | |
virtual CoolPropDbl | calc_isothermal_compressibility (void) |
Using this backend, calculate the isothermal compressibility \( \kappa = -\frac{1}{v}\left.\frac{\partial v}{\partial p}\right|_T=\frac{1}{\rho}\left.\frac{\partial \rho}{\partial p}\right|_T\) in 1/Pa. More... | |
virtual CoolPropDbl | calc_isobaric_expansion_coefficient (void) |
Using this backend, calculate the isobaric expansion coefficient \( \beta = \frac{1}{v}\left.\frac{\partial v}{\partial T}\right|_p = -\frac{1}{\rho}\left.\frac{\partial \rho}{\partial T}\right|_p\) in 1/K. More... | |
virtual CoolPropDbl | calc_isentropic_expansion_coefficient (void) |
Using this backend, calculate the isentropic expansion coefficient \( \kappa_s = -\frac{c_p}{c_v}\frac{v}{p}\left.\frac{\partial p}{\partial v}\right|_T = \frac{\rho}{p}\left.\frac{\partial p}{\partial \rho}\right|_s\). More... | |
virtual CoolPropDbl | calc_surface_tension (void) |
Using this backend, calculate the surface tension in N/m. More... | |
virtual CoolPropDbl | calc_molar_mass (void) |
Using this backend, calculate the molar mass in kg/mol. More... | |
virtual CoolPropDbl | calc_acentric_factor (void) |
Using this backend, calculate the acentric factor. More... | |
virtual CoolPropDbl | calc_pressure (void) |
Using this backend, calculate the pressure in Pa. More... | |
virtual CoolPropDbl | calc_gas_constant (void) |
Using this backend, calculate the universal gas constant \(R_u\) in J/mol/K. More... | |
virtual CoolPropDbl | calc_fugacity_coefficient (std::size_t i) |
Using this backend, calculate the fugacity coefficient (dimensionless) More... | |
virtual std::vector< CoolPropDbl > | calc_fugacity_coefficients () |
Using this backend, calculate the fugacity in Pa. More... | |
virtual CoolPropDbl | calc_fugacity (std::size_t i) |
Using this backend, calculate the fugacity in Pa. More... | |
virtual CoolPropDbl | calc_chemical_potential (std::size_t i) |
Using this backend, calculate the chemical potential in J/mol. More... | |
virtual CoolPropDbl | calc_PIP (void) |
Using this backend, calculate the phase identification parameter (PIP) More... | |
virtual void | calc_excess_properties (void) |
Using this backend, calculate and cache the excess properties. More... | |
virtual CoolPropDbl | calc_alphar (void) |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r\) (dimensionless) More... | |
virtual CoolPropDbl | calc_dalphar_dDelta (void) |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta}\) (dimensionless) More... | |
virtual CoolPropDbl | calc_dalphar_dTau (void) |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\tau}\) (dimensionless) More... | |
virtual CoolPropDbl | calc_d2alphar_dDelta2 (void) |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\delta}\) (dimensionless) More... | |
virtual CoolPropDbl | calc_d2alphar_dDelta_dTau (void) |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\tau}\) (dimensionless) More... | |
virtual CoolPropDbl | calc_d2alphar_dTau2 (void) |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\tau\tau}\) (dimensionless) More... | |
virtual CoolPropDbl | calc_d3alphar_dDelta3 (void) |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\delta\delta}\) (dimensionless) More... | |
virtual CoolPropDbl | calc_d3alphar_dDelta2_dTau (void) |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\delta\tau}\) (dimensionless) More... | |
virtual CoolPropDbl | calc_d3alphar_dDelta_dTau2 (void) |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\tau\tau}\) (dimensionless) More... | |
virtual CoolPropDbl | calc_d3alphar_dTau3 (void) |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\tau\tau\tau}\) (dimensionless) More... | |
virtual CoolPropDbl | calc_d4alphar_dDelta4 (void) |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\delta\delta\delta}\) (dimensionless) More... | |
virtual CoolPropDbl | calc_d4alphar_dDelta3_dTau (void) |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\delta\delta\tau}\) (dimensionless) More... | |
virtual CoolPropDbl | calc_d4alphar_dDelta2_dTau2 (void) |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\delta\tau\tau}\) (dimensionless) More... | |
virtual CoolPropDbl | calc_d4alphar_dDelta_dTau3 (void) |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\delta\tau\tau\tau}\) (dimensionless) More... | |
virtual CoolPropDbl | calc_d4alphar_dTau4 (void) |
Using this backend, calculate the residual Helmholtz energy term \(\alpha^r_{\tau\tau\tau\tau}\) (dimensionless) More... | |
virtual CoolPropDbl | calc_alpha0 (void) |
Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0\) (dimensionless) More... | |
virtual CoolPropDbl | calc_dalpha0_dDelta (void) |
Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\delta}\) (dimensionless) More... | |
virtual CoolPropDbl | calc_dalpha0_dTau (void) |
Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\tau}\) (dimensionless) More... | |
virtual CoolPropDbl | calc_d2alpha0_dDelta_dTau (void) |
Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\delta\delta}\) (dimensionless) More... | |
virtual CoolPropDbl | calc_d2alpha0_dDelta2 (void) |
Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\delta\tau}\) (dimensionless) More... | |
virtual CoolPropDbl | calc_d2alpha0_dTau2 (void) |
Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\tau\tau}\) (dimensionless) More... | |
virtual CoolPropDbl | calc_d3alpha0_dDelta3 (void) |
Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\delta\delta\delta}\) (dimensionless) More... | |
virtual CoolPropDbl | calc_d3alpha0_dDelta2_dTau (void) |
Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\delta\delta\tau}\) (dimensionless) More... | |
virtual CoolPropDbl | calc_d3alpha0_dDelta_dTau2 (void) |
Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\delta\tau\tau}\) (dimensionless) More... | |
virtual CoolPropDbl | calc_d3alpha0_dTau3 (void) |
Using this backend, calculate the ideal-gas Helmholtz energy term \(\alpha^0_{\tau\tau\tau}\) (dimensionless) More... | |
virtual void | calc_reducing_state (void) |
virtual CoolPropDbl | calc_pmax (void) |
Using this backend, calculate the maximum pressure in Pa. More... | |
virtual CoolPropDbl | calc_GWP20 (void) |
Using this backend, calculate the 20-year global warming potential (GWP) More... | |
virtual CoolPropDbl | calc_GWP100 (void) |
Using this backend, calculate the 100-year global warming potential (GWP) More... | |
virtual CoolPropDbl | calc_GWP500 (void) |
Using this backend, calculate the 500-year global warming potential (GWP) More... | |
virtual CoolPropDbl | calc_ODP (void) |
Using this backend, calculate the ozone depletion potential (ODP) More... | |
virtual CoolPropDbl | calc_flame_hazard (void) |
Using this backend, calculate the flame hazard. More... | |
virtual CoolPropDbl | calc_health_hazard (void) |
Using this backend, calculate the health hazard. More... | |
virtual CoolPropDbl | calc_physical_hazard (void) |
Using this backend, calculate the physical hazard. More... | |
virtual CoolPropDbl | calc_dipole_moment (void) |
Using this backend, calculate the dipole moment in C-m (1 D = 3.33564e-30 C-m) More... | |
virtual CoolPropDbl | calc_second_partial_deriv (parameters Of1, parameters Wrt1, parameters Constant1, parameters Wrt2, parameters Constant2) |
Calculate the second partial derivative using the given backend. More... | |
virtual CoolPropDbl | calc_reduced_density (void) |
Using this backend, calculate the reduced density (rho/rhoc) More... | |
virtual CoolPropDbl | calc_reciprocal_reduced_temperature (void) |
Using this backend, calculate the reciprocal reduced temperature (Tc/T) More... | |
virtual CoolPropDbl | calc_Bvirial (void) |
Using this backend, calculate the second virial coefficient. More... | |
virtual CoolPropDbl | calc_Cvirial (void) |
Using this backend, calculate the third virial coefficient. More... | |
virtual CoolPropDbl | calc_dBvirial_dT (void) |
Using this backend, calculate the derivative dB/dT. More... | |
virtual CoolPropDbl | calc_dCvirial_dT (void) |
Using this backend, calculate the derivative dC/dT. More... | |
virtual CoolPropDbl | calc_compressibility_factor (void) |
Using this backend, calculate the compressibility factor Z \( Z = p/(\rho R T) \). More... | |
virtual CoolPropDbl | calc_Ttriple (void) |
Using this backend, get the triple point temperature in K. More... | |
virtual CoolPropDbl | calc_p_triple (void) |
Using this backend, get the triple point pressure in Pa. More... | |
virtual CoolPropDbl | calc_T_critical (void) |
Using this backend, get the critical point temperature in K. More... | |
virtual CoolPropDbl | calc_T_reducing (void) |
Using this backend, get the reducing point temperature in K. More... | |
virtual CoolPropDbl | calc_p_critical (void) |
Using this backend, get the critical point pressure in Pa. More... | |
virtual CoolPropDbl | calc_p_reducing (void) |
Using this backend, get the reducing point pressure in Pa. More... | |
virtual CoolPropDbl | calc_rhomolar_critical (void) |
Using this backend, get the critical point molar density in mol/m^3. More... | |
virtual CoolPropDbl | calc_rhomass_critical (void) |
Using this backend, get the critical point mass density in kg/m^3 - Added for IF97Backend which is mass based. More... | |
virtual CoolPropDbl | calc_rhomolar_reducing (void) |
Using this backend, get the reducing point molar density in mol/m^3. More... | |
virtual void | calc_phase_envelope (const std::string &type) |
Using this backend, construct the phase envelope, the variable type describes the type of phase envelope to be built. More... | |
virtual CoolPropDbl | calc_hmass_excess (void) |
virtual CoolPropDbl | calc_smass_excess (void) |
virtual CoolPropDbl | calc_cp0mass (void) |
virtual CoolPropDbl | calc_umass_excess (void) |
virtual CoolPropDbl | calc_gibbsmass (void) |
virtual CoolPropDbl | calc_gibbsmass_excess (void) |
virtual CoolPropDbl | calc_helmholtzmass (void) |
virtual CoolPropDbl | calc_helmholtzmass_excess (void) |
virtual CoolPropDbl | calc_volumemass_excess (void) |
virtual void | update_states (void) |
Update the states after having changed the reference state for enthalpy and entropy. More... | |
virtual CoolPropDbl | calc_saturation_ancillary (parameters param, int Q, parameters given, double value) |
virtual phases | calc_phase (void) |
Using this backend, calculate the phase. More... | |
virtual void | calc_specify_phase (phases phase) |
Using this backend, specify the phase to be used for all further calculations. More... | |
virtual void | calc_unspecify_phase (void) |
Using this backend, unspecify the phase. More... | |
virtual std::vector< std::string > | calc_fluid_names (void) |
Using this backend, get a vector of fluid names. More... | |
virtual const CoolProp::SimpleState & | calc_state (const std::string &state) |
virtual const CoolProp::PhaseEnvelopeData & | calc_phase_envelope_data (void) |
virtual std::vector< CoolPropDbl > | calc_mole_fractions_liquid (void) |
virtual std::vector< CoolPropDbl > | calc_mole_fractions_vapor (void) |
virtual const std::vector< CoolPropDbl > | calc_mass_fractions (void) |
virtual CoolPropDbl | calc_first_saturation_deriv (parameters Of1, parameters Wrt1) |
virtual CoolPropDbl | calc_second_saturation_deriv (parameters Of1, parameters Wrt1, parameters Wrt2) |
virtual CoolPropDbl | calc_first_two_phase_deriv (parameters Of, parameters Wrt, parameters Constant) |
virtual CoolPropDbl | calc_second_two_phase_deriv (parameters Of, parameters Wrt, parameters Constant, parameters Wrt2, parameters Constant2) |
virtual CoolPropDbl | calc_first_two_phase_deriv_splined (parameters Of, parameters Wrt, parameters Constant, CoolPropDbl x_end) |
virtual CoolPropDbl | calc_saturated_liquid_keyed_output (parameters key) |
virtual CoolPropDbl | calc_saturated_vapor_keyed_output (parameters key) |
virtual void | calc_ideal_curve (const std::string &type, std::vector< double > &T, std::vector< double > &p) |
virtual CoolPropDbl | calc_T (void) |
Using this backend, get the temperature. More... | |
virtual CoolPropDbl | calc_rhomolar (void) |
Using this backend, get the molar density in mol/m^3. More... | |
virtual double | calc_tangent_plane_distance (const double T, const double p, const std::vector< double > &w, const double rhomolar_guess) |
Using this backend, calculate the tangent plane distance for a given trial composition. More... | |
virtual void | calc_true_critical_point (double &T, double &rho) |
Using this backend, return true critical point where dp/drho|T = 0 and d2p/drho^2|T = 0. More... | |
virtual void | calc_conformal_state (const std::string &reference_fluid, CoolPropDbl &T, CoolPropDbl &rhomolar) |
virtual void | calc_viscosity_contributions (CoolPropDbl &dilute, CoolPropDbl &initial_density, CoolPropDbl &residual, CoolPropDbl &critical) |
virtual void | calc_conductivity_contributions (CoolPropDbl &dilute, CoolPropDbl &initial_density, CoolPropDbl &residual, CoolPropDbl &critical) |
virtual std::vector< CriticalState > | calc_all_critical_points (void) |
virtual void | calc_build_spinodal () |
virtual SpinodalData | calc_get_spinodal_data () |
virtual void | calc_criticality_contour_values (double &L1star, double &M1star) |
virtual void | mass_to_molar_inputs (CoolProp::input_pairs &input_pair, CoolPropDbl &value1, CoolPropDbl &value2) |
Convert mass-based input pair to molar-based input pair; If molar-based, do nothing. More... | |
virtual void | calc_change_EOS (const std::size_t i, const std::string &EOS_name) |
Change the equation of state for a given component to a specified EOS. More... | |
Additional Inherited Members | |
Static Public Member Functions inherited from CoolProp::AbstractState | |
static AbstractState * | factory (const std::string &backend, const std::string &fluid_names) |
A factory function to return a pointer to a new-allocated instance of one of the backends. More... | |
static AbstractState * | factory (const std::string &backend, const std::vector< std::string > &fluid_names) |
A factory function to return a pointer to a new-allocated instance of one of the backends. More... | |
CoolProp::IncompressibleBackend::IncompressibleBackend | ( | ) |
Definition at line 25 of file IncompressibleBackend.cpp.
|
inlinevirtual |
Definition at line 40 of file IncompressibleBackend.h.
CoolProp::IncompressibleBackend::IncompressibleBackend | ( | IncompressibleFluid * | fluid | ) |
The instantiator
fluid | object, mostly for testing purposes |
Definition at line 30 of file IncompressibleBackend.cpp.
CoolProp::IncompressibleBackend::IncompressibleBackend | ( | const std::string & | fluid_name | ) |
The instantiator
fluid_name | the string with the fluid name |
Definition at line 45 of file IncompressibleBackend.cpp.
CoolProp::IncompressibleBackend::IncompressibleBackend | ( | const std::vector< std::string > & | component_names | ) |
The instantiator
component_names | The vector of strings of the fluid components, without file ending |
Definition at line 55 of file IncompressibleBackend.cpp.
|
inlinevirtual |
Get a string representation of the backend - for instance "HelmholtzEOSMixtureBackend" for the core mixture model in CoolProp
Must be overloaded by the backend to provide the backend's name
Implements CoolProp::AbstractState.
Definition at line 41 of file IncompressibleBackend.h.
|
inline |
Definition at line 192 of file IncompressibleBackend.h.
|
inlinevirtual |
Using this backend, calculate the thermal conductivity in W/m/K.
Reimplemented from CoolProp::AbstractState.
Definition at line 204 of file IncompressibleBackend.h.
|
inlinevirtual |
Reimplemented from CoolProp::AbstractState.
Definition at line 195 of file IncompressibleBackend.h.
|
inlinevirtual |
Reimplemented from CoolProp::AbstractState.
Definition at line 198 of file IncompressibleBackend.h.
|
inlinevirtual |
Using this backend, get the description of the fluid.
Reimplemented from CoolProp::AbstractState.
Definition at line 283 of file IncompressibleBackend.h.
|
protected |
Partial derivative of enthalpy with respect to pressure at constant temperature and composition.
Partial derivative of enthalpy with respect to pressure at constant temperature and composition
\[ \left( \frac{\partial h}{\partial p} \right)_T = v - T \left( \frac{\partial v}{\partial T} \right)_p = \rho^{-1} \left( 1 + T \rho^{-1} \left( \frac{\partial \rho}{\partial T} \right)_p \right) \]
Definition at line 564 of file IncompressibleBackend.cpp.
|
inlineprotected |
Partial derivative of enthalpy with respect to temperature at constant pressure and composition.
Definition at line 241 of file IncompressibleBackend.h.
|
inlineprotected |
Partial derivative of enthalpy with respect to temperature at constant pressure and composition integrated in temperature
Definition at line 253 of file IncompressibleBackend.h.
|
inlineprotected |
Partial derivative of density with respect to temperature at constant pressure and composition.
Definition at line 233 of file IncompressibleBackend.h.
|
protected |
Partial derivative of entropy with respect to pressure at constant temperature and composition.
Partial derivative of entropy with respect to pressure at constant temperature and composition
\[ \left( \frac{\partial s}{\partial p} \right)_T = - \left( \frac{\partial v}{\partial T} \right)_p = \rho^{-2} \left( \frac{\partial \rho}{\partial T} \right)_p \right) \]
Definition at line 559 of file IncompressibleBackend.cpp.
|
inlineprotected |
Partial derivative of entropy with respect to temperature at constant pressure and composition.
Definition at line 237 of file IncompressibleBackend.h.
|
inlineprotected |
Partial derivative of entropy with respect to temperature at constant pressure and composition integrated in temperature
Definition at line 247 of file IncompressibleBackend.h.
|
protectedvirtual |
Calculate the first partial derivative for the desired derivative.
Reimplemented from CoolProp::AbstractState.
Definition at line 529 of file IncompressibleBackend.cpp.
|
inlinevirtual |
Get the maximum fraction (mole, mass, volume) for incompressible fluid.
Reimplemented from CoolProp::AbstractState.
Definition at line 277 of file IncompressibleBackend.h.
|
inlinevirtual |
Get the minimum fraction (mole, mass, volume) for incompressible fluid.
Reimplemented from CoolProp::AbstractState.
Definition at line 274 of file IncompressibleBackend.h.
|
virtual |
... and continue with the ones that depend on reference conditions.
Reimplemented from CoolProp::AbstractState.
Definition at line 513 of file IncompressibleBackend.cpp.
|
virtual |
We start with the functions that do not need a reference state.
Calculate T given pressure and internal energy **
Reimplemented from CoolProp::AbstractState.
Definition at line 501 of file IncompressibleBackend.cpp.
|
inlinevirtual |
Using this backend, get the name of the fluid.
Reimplemented from CoolProp::AbstractState.
Definition at line 280 of file IncompressibleBackend.h.
|
inlinevirtual |
We start with the functions that do not need a reference state.
Reimplemented from CoolProp::AbstractState.
Definition at line 189 of file IncompressibleBackend.h.
|
virtual |
Reimplemented from CoolProp::AbstractState.
Definition at line 516 of file IncompressibleBackend.cpp.
|
inlinevirtual |
Reimplemented from CoolProp::AbstractState.
Definition at line 207 of file IncompressibleBackend.h.
|
inlinevirtual |
Constants from the fluid object.
Reimplemented from CoolProp::AbstractState.
Definition at line 268 of file IncompressibleBackend.h.
|
inlinevirtual |
Using this backend, calculate the minimum temperature in K.
Reimplemented from CoolProp::AbstractState.
Definition at line 271 of file IncompressibleBackend.h.
|
virtual |
Reimplemented from CoolProp::AbstractState.
Definition at line 508 of file IncompressibleBackend.cpp.
|
inlinevirtual |
Using this backend, calculate the viscosity in Pa-s.
Reimplemented from CoolProp::AbstractState.
Definition at line 201 of file IncompressibleBackend.h.
void CoolProp::IncompressibleBackend::check_status | ( | ) |
Check if the mole fractions have been set, etc.
Definition at line 288 of file IncompressibleBackend.cpp.
|
virtual |
Clear all the cached values.
Additional cached elements used for the partial derivatives
Reimplemented from CoolProp::AbstractState.
Definition at line 155 of file IncompressibleBackend.cpp.
double CoolProp::IncompressibleBackend::cmass | ( | void | ) |
Return the mass constant pressure specific heat in J/kg/K.
Definition at line 318 of file IncompressibleBackend.cpp.
double CoolProp::IncompressibleBackend::dhdpatTx | ( | void | ) |
Definition at line 347 of file IncompressibleBackend.cpp.
double CoolProp::IncompressibleBackend::dhdTatPx | ( | void | ) |
Definition at line 331 of file IncompressibleBackend.cpp.
double CoolProp::IncompressibleBackend::dhdTatPxdT | ( | void | ) |
Definition at line 339 of file IncompressibleBackend.cpp.
CoolPropDbl CoolProp::IncompressibleBackend::DmassP_flash | ( | CoolPropDbl | rhomass, |
CoolPropDbl | p | ||
) |
Calculate T given pressure and density.
These functions should be protected, but that requires new tests. I'll leave that as a TODO item for now.
rhomass | The mass density in kg/m^3 |
p | The pressure in Pa |
rhomass | The mass density in kg/m^3 |
p | The pressure in Pa |
Definition at line 395 of file IncompressibleBackend.cpp.
double CoolProp::IncompressibleBackend::drhodTatPx | ( | void | ) |
Definition at line 323 of file IncompressibleBackend.cpp.
double CoolProp::IncompressibleBackend::dsdpatTx | ( | void | ) |
Definition at line 343 of file IncompressibleBackend.cpp.
double CoolProp::IncompressibleBackend::dsdTatPx | ( | void | ) |
Definition at line 327 of file IncompressibleBackend.cpp.
double CoolProp::IncompressibleBackend::dsdTatPxdT | ( | void | ) |
Definition at line 335 of file IncompressibleBackend.cpp.
|
inlinevirtual |
Return a string from the backend for the mixture/fluid - backend dependent - could be CAS #, name, etc.
Reimplemented from CoolProp::AbstractState.
Definition at line 77 of file IncompressibleBackend.h.
|
inlinevirtual |
Get the mole fractions of the fluid.
Implements CoolProp::AbstractState.
Definition at line 96 of file IncompressibleBackend.h.
double CoolProp::IncompressibleBackend::h_ref | ( | void | ) |
Return the mass enthalpy in J/kg.
Definition at line 368 of file IncompressibleBackend.cpp.
double CoolProp::IncompressibleBackend::hmass | ( | void | ) |
Return the mass enthalpy in J/kg.
Definition at line 303 of file IncompressibleBackend.cpp.
double CoolProp::IncompressibleBackend::hmass_ref | ( | void | ) |
Return the mass enthalpy in J/kg.
Definition at line 379 of file IncompressibleBackend.cpp.
CoolPropDbl CoolProp::IncompressibleBackend::HmassP_flash | ( | CoolPropDbl | hmass, |
CoolPropDbl | p | ||
) |
Calculate T given pressure and enthalpy.
hmass | The mass enthalpy in J/kg |
p | The pressure in Pa |
Definition at line 404 of file IncompressibleBackend.cpp.
double CoolProp::IncompressibleBackend::p_ref | ( | void | ) |
Return the pressure in Pa.
Definition at line 358 of file IncompressibleBackend.cpp.
CoolPropDbl CoolProp::IncompressibleBackend::PSmass_flash | ( | CoolPropDbl | p, |
CoolPropDbl | smass | ||
) |
Calculate T given pressure and entropy.
smass | The mass entropy in J/kg/K |
p | The pressure in Pa |
Definition at line 436 of file IncompressibleBackend.cpp.
CoolPropDbl CoolProp::IncompressibleBackend::raw_calc_hmass | ( | double | T, |
double | p, | ||
double | x | ||
) |
Functions that can be used with the solver, they miss the reference values!
Definition at line 521 of file IncompressibleBackend.cpp.
CoolPropDbl CoolProp::IncompressibleBackend::raw_calc_smass | ( | double | T, |
double | p, | ||
double | x | ||
) |
Definition at line 524 of file IncompressibleBackend.cpp.
double CoolProp::IncompressibleBackend::rhomass | ( | void | ) |
Return the mass density in kg/m^3.
We have to override some of the functions from the AbstractState. The incompressibles are only mass-based and do not support conversion from molar to specific quantities. We also have a few new chaced variables that we need.
Definition at line 298 of file IncompressibleBackend.cpp.
double CoolProp::IncompressibleBackend::s_ref | ( | void | ) |
Return the molar entropy in J/mol/K.
Definition at line 373 of file IncompressibleBackend.cpp.
|
protected |
Set the fractions.
fractions | The vector of fractions of the components converted to the correct input |
Definition at line 193 of file IncompressibleBackend.cpp.
|
virtual |
Set the mass fractions.
mass_fractions | The vector of mass fractions of the components |
Implements CoolProp::AbstractState.
Definition at line 238 of file IncompressibleBackend.cpp.
|
virtual |
Set the mole fractions.
mole_fractions | The vector of mole fractions of the components |
Implements CoolProp::AbstractState.
Definition at line 212 of file IncompressibleBackend.cpp.
void CoolProp::IncompressibleBackend::set_reference_state | ( | double | T0 = 20 + 273.15 , |
double | p0 = 101325 , |
||
double | x0 = 0.0 , |
||
double | h0 = 0.0 , |
||
double | s0 = 0.0 |
||
) |
Update the reference values and clear the state.
Reference values, no need to calculate them each time
Definition at line 176 of file IncompressibleBackend.cpp.
|
virtual |
Set the volume fractions.
volu_fractions | The vector of volume fractions of the components |
Reimplemented from CoolProp::AbstractState.
Definition at line 264 of file IncompressibleBackend.cpp.
double CoolProp::IncompressibleBackend::smass | ( | void | ) |
Return the molar entropy in J/mol/K.
Definition at line 308 of file IncompressibleBackend.cpp.
double CoolProp::IncompressibleBackend::smass_ref | ( | void | ) |
Return the molar entropy in J/mol/K.
Definition at line 384 of file IncompressibleBackend.cpp.
double CoolProp::IncompressibleBackend::T_ref | ( | void | ) |
Return the temperature in K.
Definition at line 353 of file IncompressibleBackend.cpp.
double CoolProp::IncompressibleBackend::umass | ( | void | ) |
Return the molar internal energy in J/mol.
Definition at line 313 of file IncompressibleBackend.cpp.
|
virtual |
Updating function for incompressible fluid.
In this function we take a pair of thermodynamic states, those defined in the input_pairs enumeration and update all the internal variables that we can.
input_pair | Integer key from CoolProp::input_pairs to the two inputs that will be passed to the function |
value1 | First input value |
value2 | Second input value |
Implements CoolProp::AbstractState.
Definition at line 59 of file IncompressibleBackend.cpp.
|
inlinevirtual |
Implements CoolProp::AbstractState.
Definition at line 59 of file IncompressibleBackend.h.
|
inlinevirtual |
Implements CoolProp::AbstractState.
Definition at line 56 of file IncompressibleBackend.h.
|
inlinevirtual |
Implements CoolProp::AbstractState.
Definition at line 62 of file IncompressibleBackend.h.
double CoolProp::IncompressibleBackend::x_ref | ( | void | ) |
Return the composition.
Definition at line 363 of file IncompressibleBackend.cpp.
|
protected |
Additional cached elements used for the partial derivatives.
Definition at line 27 of file IncompressibleBackend.h.
|
protected |
Definition at line 28 of file IncompressibleBackend.h.
|
protected |
Definition at line 28 of file IncompressibleBackend.h.
|
protected |
Definition at line 28 of file IncompressibleBackend.h.
|
protected |
Definition at line 28 of file IncompressibleBackend.h.
|
protected |
Definition at line 28 of file IncompressibleBackend.h.
|
protected |
Definition at line 28 of file IncompressibleBackend.h.
|
protected |
Definition at line 28 of file IncompressibleBackend.h.
|
protected |
Bulk values, state variables.
Definition at line 20 of file IncompressibleBackend.h.
|
protected |
Definition at line 23 of file IncompressibleBackend.h.
|
protected |
Definition at line 27 of file IncompressibleBackend.h.
|
protected |
Definition at line 24 of file IncompressibleBackend.h.
|
protected |
Definition at line 23 of file IncompressibleBackend.h.
|
protected |
Definition at line 27 of file IncompressibleBackend.h.
|
protected |
Definition at line 23 of file IncompressibleBackend.h.
|
protected |
Definition at line 27 of file IncompressibleBackend.h.
|
protected |
Definition at line 24 of file IncompressibleBackend.h.
|
protected |
Reference values, no need to calculate them each time.
Definition at line 23 of file IncompressibleBackend.h.
|
protected |
Definition at line 27 of file IncompressibleBackend.h.
|
protected |
Definition at line 23 of file IncompressibleBackend.h.
|
protected |
Definition at line 30 of file IncompressibleBackend.h.