jaxprop.coolprop.fluid_properties module#

class jaxprop.coolprop.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
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 specified input properties. It supports both scalar and array inputs for prop_1 and prop_2, automatically broadcasting them to a common shape and evaluating the fluid properties element-wise. The resulting thermodynamic properties are returned as arrays matching the broadcasted input dimensions.

Depending on the specified inputs and fluid conditions, the calculation may involve either single-phase or two-phase states.

Parameters:
input_typeint

The variable pair used to define the thermodynamic state. This should be one of the predefined CoolProp input pairs (e.g., PT_INPUTS for pressure and temperature).

prop_1float or array_like

The first property value corresponding to the input type.

prop_2float or array_like

The second property value corresponding to the input type.

generalize_qualitybool, optional

If True, extends quality-based inputs beyond the strict two-phase limits.

supersaturationbool, optional

If True, allows evaluation of metastable (supersaturated) states.

Returns:
dict

A dictionary containing arrays of canonical thermodynamic properties, each with the same shape as the broadcasted inputs. Metadata fields such as fluid_name and identifier are included as scalars.

Raises:
Exception

If throw_exceptions is enabled 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_onset, blending_width, phase_change, blending_variable='quality_mass', 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.5, 0.5, 0.5]), spinodal_line_width=1.25, spinodal_line_method='slsqp', spinodal_line_use_previous=False, plot_quality_isolines=False, plot_two_phase_patch=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', dT_crit=1.0)[source]#
jaxprop.coolprop.fluid_properties.build_alias_map()[source]#

Build a mapping from alias names to canonical CoolProp fluid names.

Returns:
dict

Dictionary where keys are aliases and values are canonical fluid names.

jaxprop.coolprop.fluid_properties.compute_property_grid_rhoT(fluid, rho_array, T_array)[source]#
jaxprop.coolprop.fluid_properties.compute_pseudocritical_line(fluid, N_points=100)[source]#

Compute pseudocritical line (approximate, defined here at critical density).

jaxprop.coolprop.fluid_properties.compute_quality_grid(fluid, num_points, quality_levels, dT_crit=1.0)[source]#
jaxprop.coolprop.fluid_properties.compute_saturation_line(fluid, N=100, dT_crit=0.5)[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_liqFluidState

Batched FluidState for the liquid saturation line.

saturation_vapFluidState

Batched FluidState for the vapor saturation line.

jaxprop.coolprop.fluid_properties.compute_spinodal_line(fluid, N=50, method='slsqp', use_previous_as_initial_guess=False, supersaturation=False, dT_crit=0.25, tolerance=1e-08)[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_liqFluidState

Batched FluidState for the liquid spinodal line.

spinodal_vapFluidState

Batched FluidState for the vapor spinodal line.

jaxprop.coolprop.fluid_properties.compute_spinodal_point(temperature, fluid, branch, rho_guess=None, N_trial=100, method='slsqp', tolerance=1e-08, 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.

jaxprop.coolprop.fluid_properties.compute_spinodal_point_general(prop_type, prop_value, fluid, branch, rho_guess=None, N_trial=100, method='slsqp', tolerance=1e-08, 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,
... )
jaxprop.coolprop.fluid_properties.is_pure_substance(AS) bool[source]#
jaxprop.coolprop.fluid_properties.print_fluid_names()[source]#

Print all canonical fluids and their aliases in a structured format.

jaxprop.coolprop.fluid_properties.safe_alias_split(raw_string)[source]#

Split a CoolProp alias string safely, preserving multi-part names that include commas (e.g., chemical names like ‘1,2-dichloroethane’).

Parameters:
raw_stringstr

Comma-separated alias string from CoolProp.

Returns:
list of str

Cleaned list of alias names, with spurious splits recombined.

jaxprop.coolprop.fluid_properties.try_initialize_fluid(name: str, backend: str = 'HEOS')[source]#

Attempt to initialize a CoolProp AbstractState for the given fluid name. If initialization fails, raise a ValueError with close name suggestions.

Parameters:
namestr

Fluid name or alias to initialize.

backendstr, optional

CoolProp backend to use (default is ‘HEOS’).

Returns:
CP.AbstractState

Initialized CoolProp AbstractState object.

Raises:
ValueError

If the fluid name is not recognized by CoolProp, with close suggestions.