Fluid properties#

class barotropy.properties.fluid_properties.Fluid(name, backend='HEOS', exceptions=True, identifier=None)[source]#

Bases: object

Represents a fluid with various thermodynamic properties computed via CoolProp.

This class provides a convenient interface to CoolProp for various fluid property calculations.

Critical and triple point properties are computed upon initialization and stored internally for convenience.

Attributes:
namestr

Name of the fluid.

backendstr

Backend used for CoolProp, default is ‘HEOS’.

exceptionsbool

Determines if exceptions should be raised during state calculations. Default is True.

converged_flagbool

Flag indicating whether properties calculations converged.

propertiesdict

Dictionary of various fluid properties. Accessible directly as attributes (e.g., fluid.p for pressure).

critical_pointFluidState

Properties at the fluid’s critical point.

triple_point_liquidFluidState

Properties at the fluid’s triple point in the liquid state.

triple_point_vaporFluidState

Properties at the fluid’s triple point in the vapor state.

Methods

get_state(input_type, prop_1, prop_2):

Set the thermodynamic state of the fluid using specified property inputs.

Examples

Calculating properties with Fluid.get_state()

>>> fluid = bpy.Fluid(name="Water", backend="HEOS")
>>> state = fluid.get_state(bpy.PT_INPUTS, 101325, 300)
>>> print(f"Water density is {state.rho:0.2f} kg/m3 at p={state.p:0.2f} Pa and T={state.T:0.2f} K")
Water density is 996.56 kg/m3 at p=101325.00 Pa and T=300.00 K
>>> fluid = bpy.Fluid(name="Air", backend="HEOS")
>>> state = fluid.get_state(bpy.PT_INPUTS, 101325, 300)
>>> print(f"Air heat capacity ratio is {state.gamma:0.2f} at p={state.p:0.2f} Pa and T={state.T:0.2f} K")
Air heat capacity ratio is 1.40 at p=101325.00 Pa and T=300.00 K

Accessing critical point properties:

>>> fluid.critical_point.p  # Retrieves critical pressure
>>> fluid.critical_point['T']  # Retrieves critical temperature

Accessing triple point properties:

>>> fluid.triple_point_liquid.h  # Retrieves liquid enthalpy at the triple point
>>> fluid.triple_point_vapor.s  # Retrieves vapor entropy at the triple point
compute_properties_meanline(input_type, prop_1, prop_2)[source]#

Extract fluid properties for meanline model

get_property(propname)[source]#

Get the value of a single property

get_state(input_type, prop_1, prop_2, generalize_quality=False, supersaturation=False)[source]#

Set the thermodynamic state of the fluid using the CoolProp low level interface.

This method updates the thermodynamic state of the fluid in the CoolProp abstractstate object using the given input properties. It then calculates either single-phase or two-phase properties based on the current phase of the fluid.

If the calculation of properties fails, converged_flag is set to False, indicating an issue with the property calculation. Otherwise, it’s set to True.

Parameters:
input_typeint

The variable pair used to define the thermodynamic state. This should be one of the predefined input pairs in CoolProp, such as PT_INPUTS for pressure and temperature.

prop_1float

The first property value corresponding to the input type.

prop_2float

The second property value corresponding to the input type.

Returns:
barotropy.State

A State object containing the fluid properties

Raises:
Exception

If throw_exceptions attribute is set to True and an error occurs during property calculation, the original exception is re-raised.

get_state_blending(prop_1, prop_1_value, prop_2, prop_2_value, rhoT_guess_equilibrium, rhoT_guess_metastable, blending_variable, blending_onset, blending_width, phase_change, supersaturation=True, generalize_quality=True, solver_algorithm='hybr', solver_tolerance=1e-06, solver_max_iterations=100, print_convergence=False)[source]#

Calculate fluid properties by blending equilibrium and metastable properties

Note

For a detailed description of input arguments and calculation methods, see the documentation of the function compute_properties.

get_state_equilibrium(prop_1, prop_1_value, prop_2, prop_2_value, rhoT_guess=None, supersaturation=True, generalize_quality=True, solver_algorithm='hybr', solver_tolerance=1e-06, solver_max_iterations=100, print_convergence=False)[source]#

Calculate fluid properties according to thermodynamic equilibrium.

Note

For a detailed description of input arguments and calculation methods, see the documentation of the function compute_properties.

get_state_metastable(prop_1, prop_1_value, prop_2, prop_2_value, rhoT_guess=None, supersaturation=True, generalize_quality=True, solver_algorithm='hybr', solver_tolerance=1e-06, solver_max_iterations=100, print_convergence=False)[source]#

Calculate fluid properties assuming phase metastability

Note

For a detailed description of input arguments and calculation methods, see the documentation of the function compute_properties.

plot_phase_diagram(x_prop='s', y_prop='T', axes=None, N=100, plot_saturation_line=True, plot_critical_point=True, plot_triple_point_liquid=False, plot_triple_point_vapor=False, plot_spinodal_line=False, spinodal_line_color=array([0., 0., 0.]), spinodal_line_width=1.25, spinodal_line_method='bfgs', spinodal_line_use_previous=False, plot_quality_isolines=False, plot_pseudocritical_line=False, quality_levels=array([0.1, 0.2, 0.3, 0.4, 0.5, 0.6, 0.7, 0.8, 0.9, 1.]), quality_labels=False, show_in_legend=False, x_scale='linear', y_scale='linear')[source]#
class barotropy.properties.fluid_properties.FluidState(properties, fluid_name)[source]#

Bases: object

An immutable class representing the thermodynamic state of a fluid.

This class is designed to provide a read-only representation of a fluid’s state, with properties accessible through both attribute-style and dictionary-style access. The state of a fluid is defined at initialization and cannot be altered thereafter, ensuring the integrity of the data.

Instances of this class store fluid properties and the fluid’s name, providing methods to access these properties but not to modify them.

Parameters:
propertiesdict

A dictionary where keys are property names (as strings) and values are the corresponding properties of the fluid. This dictionary is converted to an immutable internal representation.

fluid_namestr

The name of the fluid.

Attributes:
propertiesdict

An internal dictionary storing the properties of the fluid state. This attribute is immutable and cannot be modified after the object’s initialization.

fluid_namestr

The name of the fluid. Immutable after initialization.

Methods

to_dict():

Returns a copy of the fluid properties as a dictionary.

keys():

Returns the keys of the fluid properties.

items():

Returns the items (key-value pairs) of the fluid properties.

fluid_name#
get(key, default=None)[source]#
items()[source]#
keys()[source]#
to_dict()[source]#
values()[source]#
barotropy.properties.fluid_properties.compute_property_grid(fluid, input_pair, range_1, range_2, generalize_quality=False, supersaturation=False)[source]#

Compute fluid properties over a specified range and store them in a dictionary.

This function creates a meshgrid of property values based on the specified ranges and input pair, computes the properties of the fluid at each point on the grid, and stores the results in a dictionary where each key corresponds to a fluid property.

Parameters:
fluidFluid object

An instance of the Fluid class.

input_pairtuple

The input pair specifying the property type (e.g., PT_INPUTS for pressure-temperature).

range1tuple

The range linspace(min, max, n) for the first property of the input pair.

range2tuple

The range linspace(min, max, n) for the second property of the input pair.

Returns:
properties_dictdict

A dictionary where keys are property names and values are 2D numpy arrays of computed properties.

grid1, grid2numpy.ndarray

The meshgrid arrays for the first and second properties.

barotropy.properties.fluid_properties.compute_property_grid_rhoT(fluid, rho_array, T_array)[source]#
barotropy.properties.fluid_properties.compute_pseudocritical_line(fluid, N_points=100)[source]#
barotropy.properties.fluid_properties.compute_quality_grid(fluid, num_points, quality_levels)[source]#
barotropy.properties.fluid_properties.compute_saturation_line(fluid, N=100)[source]#

Compute the saturation line for a given fluid.

Parameters:
fluidobject

The fluid object containing thermodynamic properties and methods.

Nint, optional

Number of points to compute along the saturation line. Default is 100.

Returns:
saturation_liqdict of lists

Dictionary containing the liquid saturation properties.

saturation_vapdict of lists

Dictionary containing the vapor saturation properties.

barotropy.properties.fluid_properties.compute_spinodal_line(fluid, N=50, method='bfgs', use_previous_as_initial_guess=False, supersaturation=False)[source]#

Compute the spinodal line for a given fluid.

Parameters:
fluidobject

The fluid object containing thermodynamic properties and methods.

Nint, optional

Number of points to compute along the spinodal line. Default is 50.

methodstr, optional

The optimization method to solve the spinodal point problem (‘bfgs’ or ‘slsqp’). Default is ‘bfgs’.

use_previous_as_initial_guessbool, optional

Whether to use the previous point as the initial guess for the next point. Default is False.

supersaturationbool, optional

Whether to compute supersaturation properties. Default is False.

Returns:
spinodal_liqdict of lists

Dictionary containing the liquid spinodal properties.

spinodal_vapdict of lists

Dictionary containing the vapor spinodal properties.

barotropy.properties.fluid_properties.compute_spinodal_point(temperature, fluid, branch, rho_guess=None, N_trial=100, method='bfgs', tolerance=1e-06, supersaturation=False, print_convergence=False)[source]#

Compute the vapor or liquid spinodal point of a fluid at a given temperature.

Parameters:
temperaturefloat

Temperature of the fluid (K).

fluidbarotropy.Fluid

The fluid for which the spinodal point is to be calculated.

branchstr

Branch of the spinodal line used to determine the density initial guess. Options: ‘liquid’ or ‘vapor’.

rho_guessfloat, optional

Initial guess for the density. If provided, this value will be used directly. If not provided, the density initial guess will be generated based on a number of trial points.

N_trialint, optional

Number of trial points to generate the density initial guess. Default is 50.

methodstr, optional

The optimization method to solve the problem (‘bfgs’ or ‘slsqp’).

tolerancefloat, optional

Tolerance for the solver termination. Defaults to 1e-6.

print_convergencebool, optional

If True, displays the convergence progress. Defaults to False.

Returns:
barotropy.State

A State object containing the fluid properties

Raises:
ValueError

If the input temperature is higher than the critical temperature or lower than the triple temperature.

Notes

When a single-phase fluid undergoes a thermodynamic process and enters the two-phase region it can exist in a single-phase state that is different from the equilibrium two-phase state. Such states are know as metastable states and they are only possible in the thermodynamic region between the saturation lines and the spinodal lines. If the thermodynamic process continues and crosses the spinodal lines metastable states become unstable and the transition to two-distinct phase occurs rapidly. Therefore, the spinodal line represents the locus of points that separates the region where a mixture is thermodynamically unstable and prone to phase separation from the region where metastable states are physically possible.

In mathematical terms, the spinodal line is defined as the loci of thermodynamic states in which the isothermal bulk modulus of the fluid is zero:

\[K_T = \rho \left( \frac{\partial p}{\partial \rho} \right)_T = 0\]

More precisely, a vapor spinodal point is the first local maximum of a isotherm line in a pressure-density diagram as the density increases. Conversely, a liquid spinodal point is the first local minimum of a isotherm line in a pressure-density diagram as the density decreases. The spinodal lines and phase envelope of carbon dioxide according to the HEOS developed by [Span and Wagner, 1996] are illustrated in the figure below

Pressure-density diagram and spinodal points for carbon dioxide.

Some equations of state are not well-posed and do not satisfy the condition \(K_T=0\) within the two-phase region. This is exemplified by the nitrogen HEOS developed by [Span et al., 2000].

Pressure-density diagram and "pseudo" spinodal points for carbon dioxide.

As seen in the figure, this HEOS is not well-posed because there are isotherms that do not have a local minimum/maximum corresponding to a state with zero isothermal bulk modulus. In such cases, this function returns the inflection point of the isotherms (corresponding to the point closest to zero bulk modulus) as the spinodal point.

barotropy.properties.fluid_properties.compute_spinodal_point_general(prop_type, prop_value, fluid, branch, rho_guess=None, N_trial=100, method='bfgs', tolerance=1e-06, print_convergence=False, supersaturation=False)[source]#

General function to compute the spinodal point for a given property name and value.

This function uses the underlying compute_spinodal_point function and iterates on temperature until the specified property at the spinodal point matches the given value.

Parameters:
prop_typestr

The type of property to match (e.g., ‘rho’, ‘p’).

prop_valuefloat

The value of the property to match at the spinodal point.

fluidobject

The fluid object containing thermodynamic properties and methods.

branchstr

The branch of the spinodal line. Options: ‘liquid’ or ‘vapor’.

rho_guessfloat, optional

Initial guess for the density. If provided, this value will be used directly. If not provided, the density initial guess will be generated based on a number of trial points.

N_trialint, optional

Number of trial points to generate the density initial guess. Default is 100.

methodstr, optional

The optimization method to solve the problem (‘bfgs’ or ‘slsqp’). Default is ‘bfgs’.

tolerancefloat, optional

Tolerance for the solver termination. Defaults to 1e-6.

print_convergencebool, optional

If True, displays the convergence progress. Defaults to False.

Returns:
barotropy.State

A State object containing the fluid properties

Raises:
ValueError

If the scalar root to determine the spinodal point fails to converge.

Notes

This function uses a root-finding algorithm to iterate on temperature until the property specified by prop_type at the spinodal point matches prop_value. The compute_spinodal_point function is called iteratively to evaluate the spinodal properties at different temperatures.

Examples

>>> state = compute_spinodal_point_general(
...     prop_type='density',
...     prop_value=500,
...     T_guess=300,
...     fluid=my_fluid,
...     branch='liquid',
...     rho_guess=10,
...     N_trial=150,
...     method='slsqp',
...     tolerance=1e-7,
...     print_convergence=True,
... )
barotropy.properties.fluid_properties.states_to_dict(states)[source]#

Convert a list of state objects into a dictionary.

Each key is a field name of the state objects, and each value is a Numpy array of all the values for that field.

Parameters:
states_gridlist of FluidState

A 1D grid where each element is a state object with the same keys.

Returns:
dict

A dictionary where keys are field names and values are 1D arrays of field values.

barotropy.properties.fluid_properties.states_to_dict_2d(states)[source]#

Convert a 2D list (grid) of state objects into a dictionary.

Each key is a field name of the state objects, and each value is a 2D Numpy array of all the values for that field.

Parameters:
states_gridlist of lists of FluidState

A 2D grid where each element is a state object with the same keys.

Returns:
dict

A dictionary where keys are field names and values are 2D arrays of field values.

Core calculations#

barotropy.properties.core_calculations.blend_properties(phase_change, props_equilibrium, props_metastable, blending_variable, blending_onset, blending_width)[source]#

Calculate blending between equilibrum and metastable fluid properties

Parameters:
phase_changestr

The type of phase change (e.g., ‘condensation’, ‘evaporation’, ‘flashing’, ‘cavitation’). Cavitation, flashing, and evaporation do the same calculations, they are aliases added for convenience.

props_equilibriumdict

The equilibrium properties.

props_metastabledict

The metastable properties.

blending_variablestr

The variable used for blending.

blending_onsetfloat

The onset value for blending.

blending_widthfloat

The width value for blending.

Returns:
dict

Blended thermodynamic properties.

barotropy.properties.core_calculations.calculate_generalized_quality(abstract_state, alpha=10)[source]#

Calculate the generalized quality of a fluid, extending the quality calculation beyond the two-phase region if necessary.

Below the critical temperature, the quality is calculated from the specific volume of the saturated liquid and vapor states. Above the critical temperature, an artificial two-phase region is defined around the critical density line to allow for a finite-width quality computation.

The quality, \(Q\), is given by:

\[Q = \frac{v - v(T, Q=0)}{v(T, Q=1) - v(T, Q=0)}\]

where \(v=1/\rho\) is the specific volume and \(T\) is the temperature.

Additionally, this function applies smoothing limiters to ensure the quality is bounded between -1 and 2. These limiters prevent the quality value from taking arbitrarily large values, enhancing stability and robustness of downstream calculation. The limiters use the logsumexp method for smooth transitions, with a specified alpha value controlling the smoothness.

Parameters:
abstract_stateCoolProp.AbstractState

CoolProp abstract state of the fluid.

alphafloat

Smoothing factor of the quality-calculation limiters.

Returns:
float

The calculated quality of the fluid.

barotropy.properties.core_calculations.calculate_mixture_properties(props_1, props_2, y_1, y_2)[source]#

Calculate the thermodynamic properties of the mixture.

TODO: add model equations and explanation

Parameters:
state_1dict

Thermodynamic properties of fluid 1.

state_2dict

Thermodynamic properties of fluid 2.

y_1float

Mass fraction of fluid 1.

y_2float

Mass fraction of fluid 2.

Returns:
mixture_propertiesdict

A dictionary containing the mixture’s properties.

barotropy.properties.core_calculations.calculate_subcooling(abstract_state)[source]#

Calculate the degree of subcooling for a given CoolProp abstract state.

This function computes the subcooling of a fluid using the CoolProp library’s abstract state class. It handles both subcritical and supercritical conditions to provide a measure of subcooling for any thermodynamic state. This results in a continuous variation of subcooling in the two-phase region, which is necessary to achieve reliable convergence of systems of equations and optimization problems involving the degree of subcooling.

Calculation cases:
  • Physically meaningful subcooling for subcritical states in the liquid region:

\[\Delta T = T(p, Q=0) - T \quad \text{for} \quad h < h(p, Q=0)\]
  • Artificial subcooling for subcritical states in the vapor and two-phase regions:

\[\Delta T = \frac{h(p, Q=0) - h}{c_p(p, Q=0)}\]
  • Artificial subcooling for supercritical states defined using the critical density line:

\[\Delta T = T(p, \rho_{\text{crit}}) - T\]
Parameters:
abstract_stateCoolProp.AbstractState

The abstract state of the fluid for which the subcooling is to be calculated.

Returns:
float

The degree of subcooling in Kelvin.

Examples

>>> import CoolProp as CP
>>> abstract_state = CP.AbstractState("HEOS", "water")
>>> abstract_state.update(CP.PT_INPUTS, 101325, 25+273.15)
>>> subcooling = bpy.calculate_subcooling(abstract_state)
>>> print(f"Subcooling is {subcooling:+0.3f} K")
Subcooling is +74.974 K
>>> abstract_state = CP.AbstractState("HEOS", "water")
>>> abstract_state.update(CP.PQ_INPUTS, 101325, 0.05)
>>> subcooling = bpy.calculate_subcooling(abstract_state)
>>> print(f"Subcooling is {subcooling:+0.3f} K")
Subcooling is -26.763 K
barotropy.properties.core_calculations.calculate_superheating(abstract_state)[source]#

Calculate the degree of superheating for a given CoolProp abstract state.

This function computes the superheating of a fluid using the CoolProp library’s abstract state class. It handles both subcritical and supercritical conditions to provide a measure of superheating for any thermodynamic state. This results in a continuous variation of superheating in the two-phase region, which is necessary to achieve in reliable convergence of systems of equations and optimization problems involving the degree of superheating.

Calculation cases:
  • Physically meaningful superheating for subcritical states in the vapor region:

\[\Delta T = T - T(p, Q=1) \quad \text{for} \quad h > h(p, Q=1)\]
  • Artificial superheating for subcritical states in the liquid and two-phase regions:

\[\Delta T = \frac{h - h(p, Q=1)}{c_p(p, Q=1)}\]
  • Artificial superheating for supercritical states defined using the critical density line:

\[\Delta T = T - T(p, \rho_{\text{crit}})\]
Parameters:
abstract_stateCoolProp.AbstractState

The abstract state of the fluid for which the superheating is to be calculated.

Returns:
float

The degree of superheating in Kelvin.

Examples

>>> import CoolProp as CP
>>> abstract_state = CP.AbstractState("HEOS", "water")
>>> abstract_state.update(CP.PT_INPUTS, 101325, 120 + 273.15)
>>> superheating = calculate_superheating(abstract_state)
>>> print(f"Superheating is {superheating:+0.3f} K")
Superheating is +20.026 K
>>> abstract_state = CP.AbstractState("HEOS", "water")
>>> abstract_state.update(CP.PQ_INPUTS, 101325, 0.95)
>>> superheating = calculate_superheating(abstract_state)
>>> print(f"Superheating is {superheating:+0.3f} K")
Superheating is -54.244 K
barotropy.properties.core_calculations.calculate_supersaturation(abstract_state, props)[source]#

Evaluate degree of supersaturation and supersaturation ratio.

The supersaturation degree is defined as the actual temperature minus the saturation temperature at the corresponding pressure:

\[\Delta T = T - T_{\text{sat}}(p)\]

The supersaturation ratio is defined as the actual pressure divided by the saturation pressure at the corresponding temperature:

\[S = \frac{p}{p_{\text{sat}}(T)}\]

The metastable liquid and metastable vapor regions are illustrated in the pressure-temperature diagram. In the metastable liquid region, \(\Delta T > 0\) and \(S < 1\), indicating that the liquid temperature is higher than the equilibrium temperature and will tend to evaporate. Conversely, in the metastable vapor region, \(\Delta T < 0\) and \(S > 1\), indicating that the vapor temperature is lower than the equilibrium temperature and will tend to condense.

Pressure-density diagram and spinodal points for carbon dioxide.

Note

Supersaturation properties are only computed for subcritical pressures and temperatures. If the fluid is in the supercritical region, the function will return NaN for the supersaturation properties.

Parameters:
abstract_stateCoolProp.AbstractState

The abstract state of the fluid for which the properties are to be calculated.

propsdict

Dictionary to store the computed properties.

Returns:
dict

Thermodynamic properties including supersaturation properties

barotropy.properties.core_calculations.compute_properties(abstract_state, prop_1, prop_1_value, prop_2, prop_2_value, calculation_type, rhoT_guess_metastable=None, rhoT_guess_equilibrium=None, phase_change=None, blending_variable=None, blending_onset=None, blending_width=None, supersaturation=False, generalize_quality=False, solver_algorithm='hybr', solver_tolerance=1e-06, solver_max_iterations=100, print_convergence=False)[source]#

Determine the thermodynamic state for the given fluid property pair by iterating on the density-temperature native inputs,

This function uses a non-linear root finder to determine the solution of the nonlinear system given by:

\[\begin{split}R_1(\rho,\, T) = y_1 - y_1(\rho,\, T) = 0 \\ R_2(\rho,\, T) = y_2 - y_2(\rho,\, T) = 0\end{split}\]

where \((y_1,\, y_2)\) is the given thermodynamic property pair (e.g., enthalpy and pressure), while density and temperature \((\rho,\, T)\) are the independent variables that the solver iterates until the residual of the problem is driven to zero.

TODO The calculations \(y_1(\rho,\, T)\) and \(y_1(\rho,\, T)\) are performed by [dependns on input]

equilibrium calculations (coolprop) evaluating the Helmholtz energy equation of state. blending of both

Parameters:
prop_1str

The the type of the first thermodynamic property.

prop_1_valuefloat

The the numerical value of the first thermodynamic property.

prop_2str

The the type of the second thermodynamic property.

prop_2_valuefloat

The the numerical value of the second thermodynamic property.

rho_guessfloat

Initial guess for density

T_guessfloat

Initial guess for temperature

calculation_typestr

The type of calculation to perform. Valid options are ‘equilibrium’, ‘metastable’, and ‘blending’.

supersaturationbool, optional

If True, calculates supersaturation variables. Defaults to False.

generalize_qualitybool, optional

If True, generalizes quality outside two-phase region. Defaults to False.

blending_variablestr, optional

The variable used for blending properties. Required if calculation_type is ‘blending’.

blending_onsetfloat, optional

The onset value for blending. Required if calculation_type is ‘blending’.

blending_widthfloat, optional

The width value for blending. Required if calculation_type is ‘blending’.

solver_algorithmstr, optional

Method to be used for solving the nonlinear system. Defaults to ‘hybr’.

solver_tolerancefloat, optional

Tolerance for the solver termination. Defaults to 1e-6.

solver_max_iterationsinteger, optional

Maximum number of iterations of the solver. Defaults to 100.

print_convergencebool, optional

If True, displays the convergence progress. Defaults to False.

Returns:
dict

Thermodynamic properties corresponding to the given input pair.

barotropy.properties.core_calculations.compute_properties_1phase(abstract_state, generalize_quality=False, supersaturation=False)[source]#

Extract single-phase properties from CoolProp abstract state

barotropy.properties.core_calculations.compute_properties_2phase(abstract_state, supersaturation=False)[source]#

Compute two-phase fluid properties from CoolProp abstract state

Get two-phase properties from mixing rules and single-phase CoolProp properties

Homogeneous equilibrium model

State formulas for T=T, p=p, mfrac/vfrac(rho), h-s-g-u-cp-cv, mu-k, a

barotropy.properties.core_calculations.compute_properties_coolprop(abstract_state, input_type, prop_1, prop_2, generalize_quality=False, supersaturation=False)[source]#

Set the thermodynamic state of the fluid based on input properties.

This method updates the thermodynamic state of the fluid in the CoolProp abstractstate object using the given input properties. It then calculates either single-phase or two-phase properties based on the current phase of the fluid.

If the calculation of properties fails, converged_flag is set to False, indicating an issue with the property calculation. Otherwise, it’s set to True.

Aliases of the properties are also added to the Fluid.properties dictionary for convenience.

Parameters:
input_typeint

The variable pair used to define the thermodynamic state. This should be one of the predefined input pairs in CoolProp, such as PT_INPUTS for pressure and temperature.

The valid options for the argument ‘input_type’ are summarized below.

Input pair name

Input pair mapping

QT_INPUTS

1

PQ_INPUTS

2

QSmolar_INPUTS

3

QSmass_INPUTS

4

HmolarQ_INPUTS

5

HmassQ_INPUTS

6

DmolarQ_INPUTS

7

DmassQ_INPUTS

8

PT_INPUTS

9

DmassT_INPUTS

10

DmolarT_INPUTS

11

HmolarT_INPUTS

12

HmassT_INPUTS

13

SmolarT_INPUTS

14

SmassT_INPUTS

15

TUmolar_INPUTS

16

TUmass_INPUTS

17

DmassP_INPUTS

18

DmolarP_INPUTS

19

HmassP_INPUTS

20

HmolarP_INPUTS

21

PSmass_INPUTS

22

PSmolar_INPUTS

23

PUmass_INPUTS

24

PUmolar_INPUTS

25

HmassSmass_INPUTS

26

HmolarSmolar_INPUTS

27

SmassUmass_INPUTS

28

SmolarUmolar_INPUTS

29

DmassHmass_INPUTS

30

DmolarHmolar_INPUTS

31

DmassSmass_INPUTS

32

DmolarSmolar_INPUTS

33

DmassUmass_INPUTS

34

DmolarUmolar_INPUTS

35

prop_1float

The first property value corresponding to the input type.

prop_2float

The second property value corresponding to the input type.

Returns:
dict

A dictionary object containing the fluid properties

barotropy.properties.core_calculations.compute_properties_metastable_rhoT(abstract_state, rho, T, supersaturation=False, generalize_quality=False)[source]#

Compute the thermodynamic properties of a fluid using temperature-density calls to the Helmholtz energy equation of state (HEOS).

Parameters:
abstract_stateCoolProp.AbstractState

The abstract state of the fluid for which the properties are to be calculated.

rhofloat

Density of the fluid (kg/m³).

Tfloat

Temperature of the fluid (K).

supersaturationbool, optional

Whether to evaluate supersaturation properties. Default is False.

Returns:
dict

Thermodynamic properties computed at the given density and temperature.

Notes

The Helmholtz energy equation of state (HEOS) expresses the Helmholtz energy as an explicit function of temperature and density:

\[\Phi = \Phi(\rho, T)\]

In dimensionless form, the Helmholtz energy is given by:

\[\phi(\delta, \tau) = \frac{\Phi(\delta, \tau)}{R T}\]

where:

  • \(\phi\) is the dimensionless Helmholtz energy

  • \(R\) is the gas constant of the fluid

  • \(\delta = \rho / \rho_c\) is the reduced density

  • \(\tau = T_c / T\) is the inverse of the reduced temperature

Thermodynamic properties can be derived from the Helmholtz energy and its partial derivatives. The following table summarizes the expressions for various properties:

Helmholtz energy thermodynamic relations#

Property name

Mathematical relation

Pressure

\[Z = \frac{p}{\rho R T} = \delta \cdot \phi_{\delta}\]

Entropy

\[\frac{s}{R} = \tau \cdot \phi_{\tau} - \phi\]

Internal energy

\[\frac{u}{R T} = \tau \cdot \phi_{\tau}\]

Enthalpy

\[\frac{h}{R T} = \tau \cdot \phi_{\tau} + \delta \cdot \phi_{\delta}\]

Isochoric heat capacity

\[\frac{c_v}{R} = -\tau^2 \cdot \phi_{\tau \tau}\]

Isobaric heat capacity

\[\frac{c_p}{R} = -\tau^2 \cdot \phi_{\tau \tau} + \frac{(\delta \cdot \phi_{\delta} - \tau \cdot \delta \cdot \phi_{\tau \delta})^2}{2 \cdot \delta \cdot \phi_{\delta} + \delta^2 \cdot \phi_{\delta \delta}}\]

Isobaric expansivity

\[\alpha \cdot T = \frac{\delta \cdot \phi_{\delta} - \tau \cdot \delta \cdot \phi_{\tau \delta}}{2 \cdot \delta \cdot \phi_{\delta} + \delta^2 \cdot \phi_{\delta \delta}}\]

Isothermal compressibility

\[\rho R T \beta = \left(2 \cdot \delta \cdot \phi_{\delta} + \delta^2 \cdot \phi_{\delta \delta} \right)^{-1}\]

Isothermal bulk modulus

\[\frac{K_T}{\rho R T} = 2 \cdot \delta \cdot \phi_{\delta} + \delta^2 \cdot \phi_{\delta \delta}\]

Isentropic bulk modulus

\[\frac{K_s}{\rho R T} = 2 \cdot \delta \cdot \phi_{\delta} + \delta^2 \ \cdot \phi_{\delta \delta} - \frac{(\delta \cdot \phi_{\delta} - \tau \cdot \delta \cdot \phi_{\tau \delta})^2}{\tau^2 \cdot \phi_{\tau \tau}}\]

Joule-Thompson coefficient

\[\rho R \mu_{\mathrm{JT}} = - \frac{\delta \cdot \phi_{\delta} + \tau \cdot \delta \cdot \phi_{\tau \delta} + \delta^2 \cdot \phi_{\delta \delta}}{(\delta \cdot \phi_{\delta} - \tau \cdot \delta \cdot \phi_{\tau \delta})^2 - \tau^2 \cdot \phi_{\tau \tau} \cdot (2 \cdot \delta \cdot \phi_{\delta} + \delta^2 \cdot \phi_{\delta \delta})}\]

Where the following abbreviations are used:

  • \(\phi_{\delta} = \left(\frac{\partial \phi}{\partial \delta}\right)_{\tau}\)

  • \(\phi_{\tau} = \left(\frac{\partial \phi}{\partial \tau}\right)_{\delta}\)

  • \(\phi_{\delta \delta} = \left(\frac{\partial^2 \phi}{\partial \delta^2}\right)_{\tau, \tau}\)

  • \(\phi_{\tau \tau} = \left(\frac{\partial^2 \phi}{\partial \tau^2}\right)_{\delta, \delta}\)

  • \(\phi_{\tau \delta} = \left(\frac{\partial^2 \phi}{\partial \tau \delta}\right)_{\delta, \tau}\)

This function can be used to estimate metastable properties using the equation of state beyond the saturation lines.