Properties of Frazil ice i.t.o. potential enthalpy
Source:R/gsw.R
gsw_frazil_properties_potential.Rd
Calculation of Absolute Salinity, Conservative Temperature, and ice mass fraction based on bulk Absolute Salinity, bulk potential enthalpy, and pressure
Implementation Note
This R function uses a wrapper to a C function contained within the GSW-C system as updated 2022-10-11 at https://github.com/TEOS-10/GSW-C with git commit `657216dd4f5ea079b5f0e021a4163e2d26893371`.
The C function uses data from the library/gsw_data_v3_0.mat
file provided in the GSW-Matlab source code, version 3.06-11.
Unfortunately, this version of the mat file is no longer displayed on the
TEOS-10.org website. Therefore, in the interests of making GSW-R be
self-contained, a copy was downloaded from
http://www.teos-10.org/software/gsw_matlab_v3_06_11.zip on 2022-05-25,
the .mat file was stored in the developer/create_data directory of
https://github.com/TEOS-10/GSW-R, and then the dataset used in GSW-R
was created based on that .mat file.
Please consult http://www.teos-10.org to learn more about the various TEOS-10 software systems.
See also
Other things related to enthalpy:
gsw_CT_from_enthalpy()
,
gsw_dynamic_enthalpy()
,
gsw_enthalpy()
,
gsw_enthalpy_CT_exact()
,
gsw_enthalpy_diff()
,
gsw_enthalpy_first_derivatives()
,
gsw_enthalpy_first_derivatives_CT_exact()
,
gsw_enthalpy_ice()
,
gsw_enthalpy_t_exact()
,
gsw_frazil_properties_potential_poly()
,
gsw_pot_enthalpy_from_pt_ice()
,
gsw_pot_enthalpy_from_pt_ice_poly()
,
gsw_pot_enthalpy_ice_freezing()
,
gsw_pot_enthalpy_ice_freezing_poly()
,
gsw_pt_from_pot_enthalpy_ice()
,
gsw_pt_from_pot_enthalpy_ice_poly()
,
gsw_specvol_first_derivatives()
,
gsw_specvol_first_derivatives_wrt_enthalpy()
Examples
SA_bulk <- c( 34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324)
h_pot_bulk <- c(-4.5544e4, -4.6033e4, -4.5830e4, -4.5589e4, -4.4948e4, -4.4027e4)
p <- c( 10, 50, 125, 250, 600, 1000)
r <- gsw_frazil_properties_potential(SA_bulk, h_pot_bulk, p)
stopifnot(all.equal(r$SA_final, c(39.098258701462051, 39.343217598625756, 39.434254585716296,
39.159536295126657, 38.820511558004590, 38.542322667924459)))
stopifnot(all.equal(r$CT_final, c(-2.155553336670014, -2.200844802695826, -2.264077329325076,
-2.344567015865174, -2.598559540430464, -2.900814843304696)))
stopifnot(all.equal(r$w_Ih_final, c(0.112190640891586, 0.113150826758543, 0.111797588975174,
0.110122251260246, 0.105199838799201, 0.098850365110330)))