Skip to contents

First Derivatives of Enthalpy wrt CT

Usage

gsw_enthalpy_first_derivatives_CT_exact(SA, CT, p)

Arguments

SA

Absolute Salinity [ g/kg ]. The valid range for most `gsw` functions is 0 to 42 g/kg.

CT

Conservative Temperature [ degC ].

p

sea pressure [dbar], i.e. absolute pressure [dbar] minus 10.1325 dbar

Value

a list containing h_SA [ (J/kg)/(g/kg) ], the derivative of enthalpy wrt Absolute Salinity, and h_CT [ (J/kg)/degC ], the derivative of enthalpy wrt Conservative Temperature.

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.

Bugs

The HTML documentation suggests that this function returns 3 values, but there are only 2 returned values in the C code used here (and the matlab code on which that is based). Also, the d/dSA check values given the HTML are not reproduced by the present function. This was reported on Mar 18, 2017 as https://github.com/TEOS-10/GSW-Matlab/issues/7. See https://github.com/TEOS-10/GSW-R/issues/34

Examples

SA <- c(34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324)
CT <- c(28.7856, 28.4329, 22.8103, 10.2600,  6.8863,  4.4036)
p <-  c(     10,      50,     125,     250,     600,    1000)
d <- gsw_enthalpy_first_derivatives_CT_exact(SA, CT, p)
stopifnot(all.equal(d$h_SA, c(-0.070224183838619, -0.351159869043798, -0.887036550157504,
                              -1.829626251448858, -4.423522691827955, -7.405211691293971)))
stopifnot(all.equal(d$h_CT/1e3, c(3.991899712269790, 3.992025674159605, 3.992210402650973,
                                  3.992283991748418, 3.992685275917238, 3.993014370250710)))