CoolProp.AbstractState_all_critical_pointsMethod
abstractState_all_critical_points(handle::Clong, length::Integer, T::Array{Float64}, p::Array{Float64}, rhomolar::Array{Float64}, stable::Array{Clong})

Calculate all the critical points for a given composition.

Arguments

  • handle: The integer handle for the state class stored in memory
  • length: The length of the buffers passed to this function
  • T: The pointer to the array of temperature (K)
  • p: The pointer to the array of pressure (Pa)
  • rhomolar: The pointer to the array of molar density (m^3/mol)
  • stable: The pointer to the array of boolean flags for whether the critical point is stable (1) or unstable (0)

Note

If there is an error in an update call for one of the inputs, no change in the output array will be made

Ref

CoolProp::AbstractStateallcriticalpoints(const long handle, const long length, double* T, double* p, double* rhomolar, long* stable, long* errcode, char* messagebuffer, const long buffer_length);

CoolProp.AbstractState_build_phase_envelopeMethod
AbstractState_build_phase_envelope(handle::Clong, level::AbstractString)

Build the phase envelope.

Arguments

  • handle: The integer handle for the state class stored in memory
  • level: How much refining of the phase envelope ("none" to skip refining (recommended) or "veryfine")

Note

If there is an error in an update call for one of the inputs, no change in the output array will be made

Ref

CoolPRop::AbstractStatebuildphaseenvelope(const long handle, const char* level, long* errcode, char* messagebuffer, const long buffer_length);

CoolProp.AbstractState_build_spinodalMethod
AbstractState_build_spinodal(handle::Clong)

Build the spinodal.

Arguments

  • handle: The integer handle for the state class stored in memory

Ref

CoolProp::AbstractStatebuildspinodal(const long handle, long* errcode, char* messagebuffer, const long bufferlength);

CoolProp.AbstractState_factoryMethod
AbstractState_factory(backend::AbstractString, fluids::AbstractString)

Generate an AbstractState instance return an integer handle to the state class generated to be used in the other low-level accessor functions.

Arguments

  • backend: The backend you will use could be: ["HEOS", "REFPROP", "INCOMP", "IF97", "TREND", "HEOS&TTSE", "HEOS&BICUBIC", "SRK", "PR", "VTPR"] etc.
  • fluids: '&' delimited list of fluids. To get a list of possible values call get_global_param_string(key) with key one of the following: ["FluidsList", "incompressible_list_pure", "incompressible_list_solution", "mixture_binary_pairs_list", "predefined_mixtures"], also there is a list in CoolProp online documentation List of Fluids, or simply type ?CoolProp_fluids

Example

julia> HEOS = AbstractState_factory("HEOS", "R245fa");
julia> TTSE = AbstractState_factory("HEOS&TTSE", "R245fa");
julia> BICU = AbstractState_factory("HEOS&BICUBIC", "R245fa");
julia> SRK = AbstractState_factory("SRK", "R245fa");
julia> PR = AbstractState_factory("PR", "R245fa");
CoolProp.AbstractState_first_partial_derivMethod
AbstractState_first_partial_deriv(handle::Clong, of::Clong, wrt::Clong, constant::Clong)

Calculate the first partial derivative in homogeneous phases from the AbstractState using integer values for the desired parameters.

Arguments

  • handle: The integer handle for the state class stored in memory
  • of: The parameter of which the derivative is being taken
  • Wrt: The derivative with with respect to this parameter
  • Constant: The parameter that is not affected by the derivative

Example

julia> as = AbstractState_factory("HEOS", "Water");
julia> AbstractState_update(as, "PQ_INPUTS", 15e5, 0);
julia> AbstractState_first_partial_deriv(as, get_param_index("Hmolar"), get_param_index("P"), get_param_index("S"))
2.07872526058326e-5
julia> AbstractState_first_partial_deriv(as, get_param_index("Hmolar"), get_param_index("P"), get_param_index("D"))
5.900781297636475e-5

Ref

double CoolProp::AbstractStatefirstpartialderiv(const long handle, const long Of, const long Wrt, const long Constant, long* errcode, char* messagebuffer, const long buffer_length);

CoolProp.AbstractState_first_saturation_derivMethod
AbstractState_first_saturation_deriv(handle::Clong, of::Clong, wrt::Clong)

Calculate a saturation derivative from the AbstractState using integer values for the desired parameters.

Arguments

  • handle: The integer handle for the state class stored in memory
  • of: The parameter of which the derivative is being taken
  • wrt: The derivative with with respect to this parameter

Example

julia> as = AbstractState_factory("HEOS", "Water");
julia> AbstractState_update(as, "PQ_INPUTS", 15e5, 0);
julia> AbstractState_first_saturation_deriv(as, get_param_index("Hmolar"), get_param_index("P"))
0.0025636362140578207

Ref

double CoolProp::AbstractStatefirstsaturationderiv(const long handle, const long Of, const long Wrt, long* errcode, char* messagebuffer, const long buffer_length);

CoolProp.AbstractState_freeMethod
AbstractState_free(handle::Clong)

Release a state class generated by the low-level interface wrapper.

Arguments

  • handle: The integer handle for the state class stored in memory
CoolProp.AbstractState_get_fugacityMethod
AbstractState_get_fugacity(handle::Clong, i::Integer)

Return the fugacity of species i

Arguments

  • handle: The integer handle for the state class stored in memory
  • i: The index of the species

Example

julia> handle = AbstractState_factory("HEOS", "Water&Ethanol");
julia> pq_inputs = get_input_pair_index("PQ_INPUTS");
julia> t = get_param_index("T");
julia> AbstractState_set_fractions(handle, [0.4, 0.6]);
julia> AbstractState_update(handle, pq_inputs, 101325, 0);
julia> AbstractState_get_fugacity(handle, 0)
30227.119385400914
julia> AbstractState_free(handle);
CoolProp.AbstractState_get_fugacity_coefficientMethod
AbstractState_get_fugacity_coefficient(handle::Clong, i::Real)

Return the fugacity coefficient of species i

Arguments

  • handle: The integer handle for the state class stored in memory
  • i: The index of the species

Example

julia> handle = AbstractState_factory("HEOS", "Water&Ethanol");
julia> pq_inputs = get_input_pair_index("PQ_INPUTS");
julia> t = get_param_index("T");
julia> AbstractState_set_fractions(handle, [0.4, 0.6]);
julia> AbstractState_update(handle, pq_inputs, 101325, 0);
julia> AbstractState_get_fugacity_coefficient(handle, 0)
0.7457961851803392
julia> AbstractState_free(handle);
CoolProp.AbstractState_get_mole_fractionsMethod
AbstractState_get_mole_fractions(handle::Clong, fractions::Array{Float64})

Get the molar fractions for the AbstractState.

Arguments

  • handle: The integer handle for the state class stored in memory
  • fractions: The array of fractions

Example

julia> handle = AbstractState_factory("HEOS", "Water&Ethanol");
julia> pq_inputs = get_input_pair_index("PQ_INPUTS");
julia> t = get_param_index("T");
julia> AbstractState_set_fractions(handle, [0.4, 0.6]);
julia> AbstractState_update(handle, pq_inputs, 101325, 0);
julia> mole_frac = [0.0, 0.0]
julia> AbstractState_get_mole_fractions(handle, mole_frac)
julia> @show mole_frac
mole_frac = [0.4, 0.6]
julia> AbstractState_free(handle);
CoolProp.AbstractState_get_mole_fractions_satStateMethod
AbstractState_get_mole_fractions_satState(handle::Clong, saturated_state::AbstractString,
                                          fractions::Array{Float64})

Get the molar fractions for the AbstractState and the desired saturated State.

Arguments

  • handle: The integer handle for the state class stored in memory
  • saturated_state: The string specifying the state (liquid or gas)
  • fractions: The array of fractions

Example

julia> handle = AbstractState_factory("HEOS", "Water&Ethanol");
julia> pq_inputs = get_input_pair_index("PQ_INPUTS");
julia> t = get_param_index("T");
julia> AbstractState_set_fractions(handle, [0.4, 0.6]);
julia> AbstractState_update(handle, pq_inputs, 101325, 0);
julia> gas_frac = [0.0, 0.0]
julia> AbstractState_get_mole_fractions_satState(handle, "gas", gas_frac)
julia> @show gas_frac
gas_frac = [0.30304421184833846, 0.6969557881516616]
julia> liq_frac = [0.0, 0.0]
julia> AbstractState_get_mole_fractions_satState(handle, "liquid", liq_frac)
julia> @show liq_frac
liq_frac = [0.4, 0.6]
julia> AbstractState_free(handle);
CoolProp.AbstractState_get_phase_envelope_dataMethod
AbstractState_get_phase_envelope_data(handle::Clong, length::Integer, T::Array{Float64}, p::Array{Float64}, rhomolar_vap::Array{Float64}, rhomolar_liq::Array{Float64}, x::Array{Float64}, y::Array{Float64})

Get data from the phase envelope for the given mixture composition.

Arguments

  • handle: The integer handle for the state class stored in memory
  • length: The number of elements stored in the arrays (both inputs and outputs MUST be the same length)
  • T: The pointer to the array of temperature (K)
  • p: The pointer to the array of pressure (Pa)
  • rhomolar_vap: The pointer to the array of molar density for vapor phase (m^3/mol)
  • rhomolar_liq: The pointer to the array of molar density for liquid phase (m^3/mol)
  • x: The compositions of the "liquid" phase (WARNING: buffer should be Ncomp*Npoints in length, at a minimum, but there is no way to check buffer length at runtime)
  • y: The compositions of the "vapor" phase (WARNING: buffer should be Ncomp*Npoints in length, at a minimum, but there is no way to check buffer length at runtime)

Example

julia> HEOS=AbstractState_factory("HEOS","Methane&Ethane");
julia> length=200;
julia> t=zeros(length);p=zeros(length);x=zeros(2*length);y=zeros(2*length);rhomolar_vap=zeros(length);rhomolar_liq=zeros(length);
julia> AbstractState_set_fractions(HEOS, [0.2, 1 - 0.2])
julia> AbstractState_build_phase_envelope(HEOS, "none")
julia> AbstractState_get_phase_envelope_data(HEOS, length, t, p, rhomolar_vap, rhomolar_liq, x, y)
julia> AbstractState_free(HEOS)

Note

If there is an error in an update call for one of the inputs, no change in the output array will be made

Ref

CoolProp::AbstractStategetphaseenvelopedata(const long handle, const long length, double* T, double* p, double* rhomolarvap, double* rhomolarliq, double* x, double* y, long* errcode, char* messagebuffer, const long bufferlength);

CoolProp.AbstractState_get_spinodal_dataMethod
AbstractState_get_spinodal_data(handle::Clong, length::Integer, tau::Array{Float64}, dalta::Array{Float64}, m1::Array{Float64})

Get data for the spinodal curve.

Arguments

  • handle: The integer handle for the state class stored in memory
  • length: The number of elements stored in the arrays (all outputs MUST be the same length)
  • tau: The pointer to the array of reciprocal reduced temperature
  • delta: The pointer to the array of reduced density
  • m1: The pointer to the array of M1 values (when L1=M1=0, critical point)

Note

If there is an error, no change in the output arrays will be made

Example

julia> HEOS=AbstractStatefactory("HEOS","Methane&Ethane"); julia> AbstractStatesetfractions(HEOS, [0.1, 0.9]); julia> AbstractStatebuildspinodal(HEOS); julia> tau, delta, m1 = AbstractStategetspinodaldata(HEOS, 127);

Ref

CoolProp::AbstractStategetspinodaldata(const long handle, const long length, double* tau, double* delta, double* M1, long* errcode, char* messagebuffer, const long buffer_length);

CoolProp.AbstractState_keyed_outputMethod
AbstractState_keyed_output(handle::Clong, param::Clong)

Get an output value from the AbstractState using an integer value for the desired output value.

Arguments

  • handle: The integer handle for the state class stored in memory
  • param::Clong: param The integer value for the parameter you want

Note

See AbstractState_output

CoolProp.AbstractState_outputMethod
AbstractState_output(handle::Clong, param::AbstractString)

Get an output value from the AbstractState using an integer value for the desired output value. It is a convenience function that call AbstractState_keyed_output

Arguments

  • handle: The integer handle for the state class stored in memory
  • param::AbstractString: The name for the parameter you want
CoolProp.AbstractState_set_binary_interaction_doubleMethod
AbstractState_set_binary_interaction_double(handle::Clong, i::Int, j::Int, parameter::AbstractString, value::Float64)

Set binary interraction parameter for different mixtures model e.g.: "linear", "Lorentz-Berthelot"

Arguments

  • handle: The integer handle for the state class stored in memory
  • i: indice of the first fluid of the binary pair
  • j: indice of the second fluid of the binary pair
  • parameter: string wit the name of the parameter, e.g.: "betaT", "gammaT", "betaV", "gammaV"
  • value: the value of the binary interaction parameter

Example

julia> handle = AbstractState_factory("HEOS", "Water&Ethanol");
julia> AbstractState_set_binary_interaction_double(handle, 0, 1, "betaT", 0.987);
julia> pq_inputs = get_input_pair_index("PQ_INPUTS");
julia> t = get_param_index("T");
julia> AbstractState_set_fractions(handle, [0.4, 0.6]);
julia> AbstractState_update(handle, pq_inputs, 101325, 0);
julia> AbstractState_keyed_output(handle, t)
349.32634425309755
julia> AbstractState_free(handle);
CoolProp.AbstractState_set_cubic_alpha_CMethod
AbstractState_set_cubic_alpha_C(handle::Clong, i::Integer, parameter::AbstractString, c1::Real, c2::Real, c3::Real)

Set cubic's alpha function parameters.

Arguments

  • handle: The integer handle for the state class stored in memory
  • i: indice of the fluid the parameter should be applied too (for mixtures)
  • parameter: the string specifying the alpha function to use, e.g. "TWU" for the Twu or "MC" for Mathias-Copeman alpha function.
  • c1: the first parameter for the alpha function
  • c2: the second parameter for the alpha function
  • c3: the third parameter for the alpha function
CoolProp.AbstractState_set_fluid_parameter_doubleMethod
AbstractState_set_fluid_parameter_double(handle::Clong, i::Integer, parameter::AbstractString, value::Real)

Set some fluid parameter (ie volume translation for cubic). Currently applied to the whole fluid not to components.

Arguments

  • handle: The integer handle for the state class stored in memory
  • i: indice of the fluid the parameter should be applied to (for mixtures)
  • parameter: string wit the name of the parameter, e.g. "c", "cm", "c_m" for volume translation parameter.
  • value: the value of the parameter
CoolProp.AbstractState_set_fractionsMethod
AbstractState_set_fractions(handle::Clong, fractions::Array{Float64})

Set the fractions (mole, mass, volume) for the AbstractState.

Arguments

  • handle: The integer handle for the state class stored in memory
  • fractions: The array of fractions

Example

julia> handle = AbstractState_factory("HEOS", "Water&Ethanol");
julia> pq_inputs = get_input_pair_index("PQ_INPUTS");
julia> t = get_param_index("T");
julia> AbstractState_set_fractions(handle, [0.4, 0.6]);
julia> AbstractState_update(handle, pq_inputs, 101325, 0);
julia> AbstractState_keyed_output(handle, t)
352.3522212991724
julia> AbstractState_free(handle);
CoolProp.AbstractState_specify_phaseMethod
AbstractState_specify_phase(handle::Clong, phase::AbstractString)

Specify the phase to be used for all further calculations.

Arguments

  • handle: The integer handle for the state class stored in memory
  • phase: The string with the phase to use. Possible candidates are listed below:
Phase nameCondition
phase_liquid
phase_gas
phase_twophase
phase_supercritical
phasesupercriticalgasp < pc, T > Tc
phasesupercriticalliquidp > pc, T < Tc
phasecriticalpointp = pc, T = Tc
phase_unknown
phasenotimposed

Example

julia> heos = AbstractState_factory("HEOS", "Water");
# Do a flash call that is a very low density state point, definitely vapor
julia> @time AbstractState_update(heos, "DmolarT_INPUTS", 1e-6, 300);
  0.025233 seconds (5.23 k allocations: 142.283 KB)
# Specify the phase - for some inputs (especially density-temperature), this will result in a
# more direct evaluation of the equation of state without checking the saturation boundary
julia> AbstractState_specify_phase(heos, "phase_gas");
# We try it again - a bit faster
julia> @time AbstractState_update(heos, "DmolarT_INPUTS", 1e-6, 300);
  0.000050 seconds (5 allocations: 156 bytes)
julia> AbstractState_free(heos);
CoolProp.AbstractState_unspecify_phaseMethod
AbstractState_unspecify_phase(handle::Clong)

Unspecify the phase to be used for all further calculations.

Arguments

  • handle: The integer handle for the state class stored in memory
CoolProp.AbstractState_updateMethod
AbstractState_update(handle::Clong, input_pair::Clong, value1::Real, value2::Real)
AbstractState_update(handle::Clong, input_pair::AbstractString, value1::Real, value2::Real)

Update the state of the AbstractState.

Arguments

  • handle: The integer handle for the state class stored in memory
  • input_pair::Clong: The integer value for the input pair obtained from getinputpair_index(param::AbstractString)
  • input_pair::AbstractString: The name of an input pair
  • value1: The first input value
  • value2: The second input value

Example

julia> handle = AbstractState_factory("HEOS", "Water&Ethanol");
julia> pq_inputs = get_input_pair_index("PQ_INPUTS");
julia> t = get_param_index("T");
julia> AbstractState_set_fractions(handle, [0.4, 0.6]);
julia> AbstractState_update(handle, pq_inputs, 101325, 0);
julia> AbstractState_keyed_output(handle, t)
352.3522212991724
julia> AbstractState_free(handle);
CoolProp.AbstractState_update_and_1_outMethod
AbstractState_update_and_1_out(handle::Clong, input_pair::Clong, value1::Array{Float64}, value2::Array{Float64}, length::Integer, output::Clong, out::Array{Float64})
AbstractState_update_and_1_out(handle::Clong, input_pair::AbstractString, value1::Array{Float64}, value2::Array{Float64}, length::Integer, output::AbstractString, out::Array{Float64})

Update the state of the AbstractState and get one output value (temperature, pressure, molar density, molar enthalpy and molar entropy) from the AbstractState using pointers as inputs and output to allow array computation.

Arguments

  • handle: The integer handle for the state class stored in memory
  • input_pair::Clong: The integer value for the input pair obtained from getinputpair_index
  • input_pair::AbstractString:
  • value1: The pointer to the array of the first input parameters
  • value2: The pointer to the array of the second input parameters
  • length: The number of elements stored in the arrays (both inputs and outputs MUST be the same length)
  • output: The indice for the output desired
  • out: The array for output
CoolProp.AbstractState_update_and_5_outMethod
AbstractState_update_and_5_out(handle::Clong, input_pair::Clong, value1::Array{Float64}, value2::Array{Float64}, length::Integer, outputs::Array{Clong}, out1::Array{Float64}, out2::Array{Float64}, out3::Array{Float64}, out4::Array{Float64}, out5::Array{Float64})
AbstractState_update_and_5_out{S<:AbstractString}(handle::Clong, input_pair::AbstractString, value1::Array{Float64}, value2::Array{Float64}, length::Integer, outputs::Array{S}, out1::Array{Float64}, out2::Array{Float64}, out3::Array{Float64}, out4::Array{Float64}, out5::Array{Float64})

Update the state of the AbstractState and get an output value five common outputs (temperature, pressure, molar density, molar enthalpy and molar entropy) from the AbstractState using pointers as inputs and output to allow array computation.

Arguments

  • handle: The integer handle for the state class stored in memory
  • input_pair::Clong: The integer value for the input pair obtained from getinputpair_index
  • input_pair::AbstractString:
  • value1: The pointer to the array of the first input parameters
  • value2: The pointer to the array of the second input parameters
  • length: The number of elements stored in the arrays (both inputs and outputs MUST be the same length)
  • outputs: The 5-element vector of indices for the outputs desired
  • out1: The array for the first output
  • out2: The array for the second output
  • out3: The array for the third output
  • out4: The array for the fourth output
  • out5: The array for the fifth output
CoolProp.AbstractState_update_and_common_outMethod
AbstractState_update_and_common_out(handle::Clong, input_pair::Clong, value1::Array{Float64}, value2::Array{Float64}, length::Integer, T::Array{Float64}, p::Array{Float64}, rhomolar::Array{Float64}, hmolar::Array{Float64}, smolar::Array{Float64})
AbstractState_update_and_common_out(handle::Clong, input_pair::AbstractString, value1::Array{Float64}, value2::Array{Float64}, length::Integer, T::Array{Float64}, p::Array{Float64}, rhomolar::Array{Float64}, hmolar::Array{Float64}, smolar::Array{Float64})

Update the state of the AbstractState and get an output value five common outputs (temperature, pressure, molar density, molar enthalpy and molar entropy) from the AbstractState using pointers as inputs and output to allow array computation.

Arguments

  • handle: The integer handle for the state class stored in memory
  • input_pair::Clong: The integer value for the input pair obtained from getinputpair_index
  • input_pair::AbstractString:
  • value1: The pointer to the array of the first input parameters
  • value2: The pointer to the array of the second input parameters
  • length: The number of elements stored in the arrays (both inputs and outputs MUST be the same length)
  • T: The pointer to the array of temperature
  • p: The pointer to the array of pressure
  • rhomolar: Array of molar density
  • hmolar: The array of molar enthalpy
  • smolar: Trray of molar entropy

Example

julia> handle = AbstractState_factory("HEOS", "Water&Ethanol");
julia> pq_inputs = get_input_pair_index("PQ_INPUTS");
julia> AbstractState_set_fractions(handle, [0.4, 0.6]);
julia> T = [0.0]; p = [0.0]; rhomolar = [0.0]; hmolar = [0.0]; smolar = [0.0];
julia> AbstractState_update_and_common_out(handle, pq_inputs, [101325.0], [0.0], 1, T, p, rhomolar, hmolar, smolar);
julia> AbstractState_free(handle);
CoolProp.F2KMethod
F2K(tf::Real)

Convert from degrees Fahrenheit to Kelvin (useful primarily for testing).

CoolProp.HAPropsSIMethod
HAPropsSI(output::AbstractString, name1::AbstractString, value1::Real, name2::AbstractString, value2::Real, name3::AbstractString, value3::Real)

DLL wrapper of the HAPropsSI function.

Arguments

  • output: Output name for desired property, accepetd values are listed in the following table
  • name1, name2, name3: Input name of given state values, one must be "P"

Note

Here, all outputs calculated as functions of these three: "Y"(Water mole fraction), "T" and "P", as "P" is mandatory so more performance achieved when "Y" and "T" is given (or at least one of them).

Parameter(s) nameDescriptionUnitFormula
Omega HumRat WHumidity ratio
psi_w YWater mole fractionmol_w/mol
Tdp T_dp DewPoint DDew point temperatureK
Twb T_wb WetBulb BWet bulb temperatureK
Enthalpy H HdaEnthalpyJ/kg_da
HhaEnthalpy per kg of humid airJ/kg_ha
InternalEnergy U UdaInternal energyJ/kg_da
UhaInternal energy per kg of humid airJ/kg_ha
Entropy S SdaEntropyJ/kg_da/K
ShaEntropy per kg of humid airJ/kg_ha/K
RH RelHum RRelative humidity
Tdb T_db TTemperatureK
PPressurePa
V VdaSpecific volumem^3/kg_da$MolarVolume*(1+HumidityRatio)/M_ha$
VhaSpecific volume per kg of humid airm^3/kg_ha$MolarVolume/M_ha$
mu Visc MViscosity
k Conductivity KConductivity
C cpHeat cap. const. press.J/kg_da/K
Cha cp_haHeat cap. const. press. per kg of humid airJ/kg_ha/K
CVHeat cap. const. vol.J/kg_da/K
CVha cv_haHeat cap. const. vol. per kg of humid airJ/kg_ha/K
P_wPartial pressure of water
isentropic_exponentIsentropic exponent
speedofsoundSpeed of sound$sqrt(1/M_ha*cp/cv*dpdrho__constT)$
ZCompressibility factor$P*MolarVolume/(R*T)$

Example

# Enthalpy (J per kg dry air) as a function of temperature, pressure,
# and relative humidity at STP
julia> h = HAPropsSI("H", "T", 298.15, "P", 101325, "R", 0.5)
50423.45039247604
# Temperature of saturated air at the previous enthalpy
julia> T = HAPropsSI("T", "P", 101325, "H", h, "R", 1.0)
290.9620891952412
# Temperature of saturated air - order of inputs doesn't matter
julia> T = HAPropsSI("T", "H", h, "R", 1.0, "P", 101325)
290.9620891952412

Ref

HumidAir::HAPropsSI(const char* OutputName, const char* Input1Name, double Input1, const char* Input2Name, double Input2, const char* Input3Name, double Input3);

CoolProp.K2FMethod
K2F(tk::Real)

Convert from Kelvin to degrees Fahrenheit (useful primarily for testing).

CoolProp.PhaseSIMethod
PhaseSI(name1::AbstractString, value1::Real, name2::AbstractString, value2::Real, fluid::AbstractString)

Return a string representation of the phase. Valid states are: "liquid", "supercritical", "supercriticalgas", "supercriticalliquid", "critical_point", "gas", "twophase"

Arguments

  • fluid::AbstractString: The name of the fluid that is part of CoolProp, for instance "n-Propane", to get a list of different passible fulid types call get_global_param_string(key) with key one of the following: ["FluidsList", "incompressible_list_pure", "incompressible_list_solution", "mixture_binary_pairs_list", "predefined_mixtures"], also there is a list in CoolProp online documentation List of Fluids
  • name1::AbstractString: The name of parameter for first state point
  • value1::Real: Value of the first state point
  • name2::AbstractString: The name of parameter for second state point
  • value2::Real: Value of the second state point

Example

julia> PhaseSI("T", PropsSI("TCRIT", "Water"), "P", PropsSI("PCRIT", "Water"), "Water")
"critical_point"
julia> PhaseSI("T", 300, "Q", 1, "Water")
"twophase"
julia> PhaseSI("P", "T", 300, "Q", 1, "Water")
3536.806750422325
julia> PhaseSI("T", 300, "P", 3531, "Water")
"gas"
julia> PhaseSI("T", 300, "P", 3541, "Water")
"liquid"

Ref

CoolProp::PhaseSI(const std::string &, double, const std::string &, double, const std::string&)

CoolProp.PropsSIMethod
PropsSI(output::AbstractString, name1::AbstractString, value1::Real, name2::AbstractString, value2::Real, fluid::AbstractString)

Return a value that depends on the thermodynamic state.

For pure and pseudo-pure fluids, two state points are required to fix the state. The equations of state are based on T and ρ as state variables, so T, ρ will always be the fastest inputs. P, T will be a bit slower (3-10 times), and then comes inputs where neither T nor ρ are given, like p, h. They will be much slower. If speed is an issue, you can look into table-based interpolation methods using TTSE or bicubic interpolation.

Arguments

  • fluid::AbstractString: The name of the fluid that is part of CoolProp, for instance "n-Propane", to get a list of different passible fulid types call get_global_param_string(key) with key one of the following: ["FluidsList", "incompressible_list_pure", "incompressible_list_solution", "mixture_binary_pairs_list", "predefined_mixtures"], also there is a list in CoolProp online documentation List of Fluids
  • output::AbstractString: The name of parameter to evaluate. to see a list type ?CoolProp_parameters
  • name1::AbstractString: The name of parameter for first state point
  • value1::Real: Value of the first state point
  • name2::AbstractString: The name of parameter for second state point
  • value2::Real: Value of the second state point

Example

julia> PropsSI("D", "T", 300, "P", 101325, "n-Butane")
2.4325863624814326
julia> PropsSI("D", "T", 300, "P", 101325, "INCOMP::DEB") # incompressible pure
857.1454
julia> PropsSI("D", "T", 300, "P", 101325, "INCOMP::LiBr[0.23]") # incompressible mass-based binary mixture
1187.5438243617214
julia> PropsSI("D", "T", 300, "P", 101325, "INCOMP::ZM[0.23]") # incompressible volume-based binary mixtures
1028.7273860290911
julia> PropsSI("Dmass", "T", 300, "P", 101325, "R125[0.5]&R32[0.5]")
3.5413381483914512
julia> split(get_global_param_string("mixture_binary_pairs_list"), ', ')[1] # a random binary pair
"100-41-4&106-42-3"
julia> PropsSI("Dmass", "T", 300, "P", 101325, "100-41-4[0.5]&106-42-3[0.5]") # ethylbenzene[0.5]&p-Xylene[0.5]
857.7381127561846

Ref

CoolProp::PropsSI(const std::string &, const std::string &, double, const std::string &, double, const std::string&)

CoolProp.PropsSIMethod
PropsSI(fluid::AbstractString, output::AbstractString)

Return a value that does not depend on the thermodynamic state - this is a convenience function that does the call PropsSI(output, "", 0, "", 0, fluid).

Arguments

  • fluid::AbstractString: The name of the fluid that is part of CoolProp, for instance "n-Propane", to get a list of possible values types call get_global_param_string(key) with key one of the following: ["FluidsList", "incompressible_list_pure", "incompressible_list_solution", "mixture_binary_pairs_list", "predefined_mixtures"], also there is a list in CoolProp online documentation List of Fluids, or simply type ?CoolProp_fluids
  • output::AbstractString: The name of parameter to evaluate. to see a list type ?CoolProp_parameters

Example

julia> PropsSI("n-Butane", "rhomolar_critical")
3922.769612987809

Ref

CoolProp::Props1SI(std::string, std::string)

CoolProp.cair_satMethod
cair_sat(t::Real)

Humid air saturation specific in [kJ/kg-K] heat at 1 atmosphere, based on a correlation from EES.

Arguments

  • t: T [K] good from 250K to 300K, no error bound checking is carried out.

Ref

HumidAir::cair_sat(double);

Note

Equals partial derivative of enthalpy with respect to temperature at constant relative humidity of 100 percent and pressure of 1 atmosphere.

CoolProp.get_debug_levelMethod
get_debug_level()

Get the debug level.

Return value

Level The level of the verbosity for the debugging output (0-10) 0: no debgging output

CoolProp.get_fluid_param_stringMethod
get_fluid_param_string(fluid::AbstractString, param::AbstractString)

Get a string for a value from a fluid (numerical values for the fluid can be obtained from Props1SI function).

Arguments

  • fluid: The name of the fluid that is part of CoolProp, for instance "n-Propane"
  • param: A string, can be in one of the terms described in the following table
ParamNameDescription
"aliases"A comma separated list of aliases for the fluid
"CAS", "CAS_number"The CAS number
"ASHRAE34"The ASHRAE standard 34 safety rating
"REFPROPName", "REFPROP_name"The name of the fluid used in REFPROP
"Bibtex-XXX"A BibTeX key, where XXX is one of the bibtex keys used in get_BibTeXKey
"pure""true" if the fluid is pure, "false" otherwise
"formula"The chemical formula of the fluid in LaTeX form if available, "" otherwise

Note

A tabular output for this function is available with ?CoolProp_fluids

CoolProp.get_global_param_stringMethod
get_global_param_string(key::AbstractString)

Get a globally-defined string.

Ref

ref CoolProp::getglobalparam_string

Arguments

  • key: A string represents parameter name, could be one of ["version", "gitrevision", "errstring", "warnstring", "FluidsList", "incompressiblelistpure", "incompressiblelistsolution", "mixturebinarypairslist", "parameterlist", "predefinedmixtures", "HOME", "cubicfluids_schema"]
CoolProp.get_input_pair_indexMethod
get_input_pair_index(pair::AbstractString)

Get the index for an input pair "PTINPUTS", "HmolarQINPUTS", etc, for abstractstate_update() function.

Arguments

  • pair::AbstractString: The name of an input pair, described in the following table
Input PairDescription
QT_INPUTSMolar quality, Temperature in K
QSmolar_INPUTSMolar quality, Entropy in J/mol/K
QSmass_INPUTSMolar quality, Entropy in J/kg/K
HmolarQ_INPUTSEnthalpy in J/mol, Molar quality
HmassQ_INPUTSEnthalpy in J/kg, Molar quality
DmassQ_INPUTSMolar density kg/m^3, Molar quality
DmolarQ_INPUTSMolar density in mol/m^3, Molar quality
PQ_INPUTSPressure in Pa, Molar quality
PT_INPUTSPressure in Pa, Temperature in K
DmassT_INPUTSMass density in kg/m^3, Temperature in K
DmolarT_INPUTSMolar density in mol/m^3, Temperature in K
HmassT_INPUTSEnthalpy in J/kg, Temperature in K
HmolarT_INPUTSEnthalpy in J/mol, Temperature in K
SmassT_INPUTSEntropy in J/kg/K, Temperature in K
SmolarT_INPUTSEntropy in J/mol/K, Temperature in K
TUmass_INPUTSTemperature in K, Internal energy in J/kg
TUmolar_INPUTSTemperature in K, Internal energy in J/mol
DmassP_INPUTSMass density in kg/m^3, Pressure in Pa
DmolarP_INPUTSMolar density in mol/m^3, Pressure in Pa
HmassP_INPUTSEnthalpy in J/kg, Pressure in Pa
HmolarP_INPUTSEnthalpy in J/mol, Pressure in Pa
PSmass_INPUTSPressure in Pa, Entropy in J/kg/K
PSmolar_INPUTSPressure in Pa, Entropy in J/mol/K
PUmass_INPUTSPressure in Pa, Internal energy in J/kg
PUmolar_INPUTSPressure in Pa, Internal energy in J/mol
DmassHmass_INPUTSMass density in kg/m^3, Enthalpy in J/kg
DmolarHmolar_INPUTSMolar density in mol/m^3, Enthalpy in J/mol
DmassSmass_INPUTSMass density in kg/m^3, Entropy in J/kg/K
DmolarSmolar_INPUTSMolar density in mol/m^3, Entropy in J/mol/K
DmassUmass_INPUTSMass density in kg/m^3, Internal energy in J/kg
DmolarUmolar_INPUTSMolar density in mol/m^3, Internal energy in J/mol
HmassSmass_INPUTSEnthalpy in J/kg, Entropy in J/kg/K
HmolarSmolar_INPUTSEnthalpy in J/mol, Entropy in J/mol/K
SmassUmass_INPUTSEntropy in J/kg/K, Internal energy in J/kg
SmolarUmolar_INPUTSEntropy in J/mol/K, Internal energy in J/mol

Example

julia> get_input_pair_index("PT_INPUTS")
9
CoolProp.get_param_indexMethod
get_param_index(param::AbstractString)

Get the index as a long for a parameter "T", "P", etc, for abstractstate_keyed_output() function.

Arguments

  • param: A string represents parameter name, to see full list check "Table of string inputs to PropsSI function": http://www.coolprop.org/coolprop/HighLevelAPI.html#parameter-table, or simply type get_global_param_string("parameter_list")
CoolProp.get_parameter_information_stringMethod
get_parameter_information_string(key::AbstractString, outtype::AbstractString)
get_parameter_information_string(key::AbstractString)

Get information for a parameter.

Arguments

  • key: A string represents parameter name, to see full list check "Table of string inputs to PropsSI function": http://www.coolprop.org/coolprop/HighLevelAPI.html#parameter-table, or simply type get_global_param_string("parameter_list")
  • outtype="long": Output type, could be one of the ["IO", "short", "long", "units"], with a default value of "long"

Example

julia> get_parameter_information_string("HMOLAR")
"Molar specific enthalpy"
julia> get_parameter_information_string("HMOLAR", "units")
"J/mol"

Note

A tabular output for this function is available with ?CoolProp_parameters

CoolProp.saturation_ancillaryMethod
saturation_ancillary(fluid_name::AbstractString, output::AbstractString, quality::Integer, input::AbstractString, value::Real)

Extract a value from the saturation ancillary.

Arguments

  • fluid_name: The name of the fluid to be used - HelmholtzEOS backend only
  • output: The desired output variable ("P" for instance for pressure)
  • quality: The quality, 0 or 1
  • input: The input variable ("T")
  • value: The input value

Example

julia> saturation_ancillary("R410A","I",1,"T", 300) 0.004877519938463293

Ref

double saturationancillary(const char* fluidname, const char* output, int Q, const char* input, double value);

CoolProp.set_configMethod
set_config(key::AbstractString, val::AbstractString)

Set configuration string.

Arguments

  • key::AbstractString: The key to configure, following table shows possible key values, its default setting and its usage.
  • val::AbstractString: The value to set to the key
KeyDefaultDescription
ALTERNATIVETABLESDIRECTORY""If provided, this path will be the root directory for the tabular data. Otherwise, {HOME}/.CoolProp/Tables is used
ALTERNATIVEREFPROPPATH""An alternative path to be provided to the directory that contains REFPROP's fluids and mixtures directories. If provided, the SETPATH function will be called with this directory prior to calling any REFPROP functions.
ALTERNATIVEREFPROPHMXBNCPATH""An alternative path to the HMX.BNC file. If provided, it will be passed into REFPROP's SETUP or SETMIX routines
VTPRUNIFACPATH""The path to the directory containing the UNIFAC JSON files. Should be slash terminated
CoolProp.set_configMethod
set_config(key::AbstractString, val::Bool)

Set configuration value as a boolean.

Arguments

  • key::AbstractString: The key to configure, following table shows possible key values, its default setting and its usage.
  • val::Bool: The value to set to the key
KeyDefaultDescription
NORMALIZEGASCONSTANTStrueIf true, for mixtures, the molar gas constant (R) will be set to the CODATA value
CRITICALWITHIN1UKtrueIf true, any temperature within 1 uK of the critical temperature will be considered to be AT the critical point
CRITICALSPLINESENABLEDtrueIf true, the critical splines will be used in the near-vicinity of the critical point
SAVERAWTABLESfalseIf true, the raw, uncompressed tables will also be written to file
REFPROPDONTESTIMATEINTERACTIONPARAMETERSfalseIf true, if the binary interaction parameters in REFPROP are estimated, throw an error rather than silently continuing
REFPROPIGNOREERRORESTIMATEDINTERACTION_PARAMETERSfalseIf true, if the binary interaction parameters in REFPROP are unable to be estimated, silently continue rather than failing
REFPROPUSEGERGfalseIf true, rather than using the highly-accurate pure fluid equations of state, use the pure-fluid EOS from GERG-2008
REFPROPUSEPENGROBINSONfalseIf true, rather than using the highly-accurate pure fluid equations of state, use the Peng-Robinson EOS
DONTCHECKPROPERTY_LIMITSfalseIf true, when possible, CoolProp will skip checking whether values are inside the property limits
HENRYSLAWTOGENERATEVLE_GUESSESfalseIf true, when doing water-based mixture dewpoint calculations, use Henry's Law to generate guesses for liquid-phase composition
OVERWRITE_FLUIDSfalseIf true, and a fluid is added to the fluids library that is already there, rather than not adding the fluid (and probably throwing an exception), overwrite it
OVERWRITEDEPARTUREFUNCTIONfalseIf true, and a departure function to be added is already there, rather than not adding the departure function (and probably throwing an exception), overwrite it
OVERWRITEBINARYINTERACTIONfalseIf true, and a pair of binary interaction pairs to be added is already there, rather than not adding the binary interaction pair (and probably throwing an exception), overwrite it
CoolProp.set_configMethod
set_config(key::AbstractString, val::Real)

Set configuration numerical value as double.

Arguments

  • key::AbstractString: The key to configure, following table shows possible key values, its default setting and its usage.
  • val::Real: The value to set to the key
KeyDefaultDescription
MAXIMUMTABLEDIRECTORYSIZEIN_GB1.0The maximum allowed size of the directory that is used to store tabular data
PHASEENVELOPESTARTINGPRESSUREPA100.0Starting pressure [Pa] for phase envelope construction
RUCODATA8.3144598The value for the ideal gas constant in J/mol/K according to CODATA 2014. This value is used to harmonize all the ideal gas constants. This is especially important in the critical region.
SPINODALMINIMUMDELTA0.5The minimal delta to be used in tracing out the spinodal; make sure that the EOS has a spinodal at this value of delta=rho/rho_r
CoolProp.set_debug_levelMethod
set_debug_level(level::Integer)

Set the debug level.

Arguments

  • level::Integer: The level of the verbosity for the debugging output (0-10) 0: no debgging output
CoolProp.set_reference_stateMethod
set_reference_state(ref::AbstractString, reference_state::AbstractString)

Set the reference state based on a string representation.

#Arguments

  • fluid::AbstractString The name of the fluid (Backend can be provided like "REFPROP::Water", or if no backend is provided, "HEOS" is the assumed backend)
  • reference_state::AbstractString The reference state to use, one of:
Reference StateDescription
"IIR"h = 200 kJ/kg, s=1 kJ/kg/K at 0C saturated liquid
"ASHRAE"h = 0, s = 0 @ -40C saturated liquid
"NBP"h = 0, s = 0 @ 1.0 bar saturated liquid
"DEF"Reset to the default reference state for the fluid
"RESET"Remove the offset

#Example

julia> h0=-15870000.0; # J/kg
julia> s0= 3887.0; #J/kg
julia> rho0=997.1;
julia> T0=298.15;
julia> M = PropsSI("molemass", "Water");
julia> set_reference_stateD("Water", T0, rho0/M, h0*M, s0*M);
julia> PropsSI("H", "T", T0, "P", 101325, "Water")
-1.5870107493843542e7
julia> set_reference_stateS("Water", "DEF");
julia> PropsSI("H", "T", T0, "P", 101325, "Water")
104920.1198093371

#Ref setreferencestateS(const std::string& FluidName, const std::string& reference_state)

#Note The changing of the reference state should be part of the initialization of your program, and it is not recommended to change the reference state during the course of making calculations

CoolProp.set_reference_stateMethod
set_reference_state(fluid::AbstractString, T::Real, rhomolar::Real, hmolar0::Real, smolar0::Real)

Set the reference state based on a thermodynamic state point specified by temperature and molar density.

#Arguments

  • fluid::AbstractString The name of the fluid
  • T::Real Temperature at reference state [K]
  • rhomolar::Real Molar density at reference state [mol/m^3]
  • hmolar0::Real Molar enthalpy at reference state [J/mol]
  • smolar0::Real Molar entropy at reference state [J/mol/K]

#Ref setreferencestateD(const char* Ref, double T, double rhomolar, double hmolar0, double smolar0)