Contents
USAGE:
[N2, N2_p, N2_specvol, N2_alpha, N2_beta, dSA, dCT, dp] = ...
gsw_Nsquared_min(SA,CT,p,{lat})
DESCRIPTION:
Calculates the minimum buoyancy frequency squared (N2)(i.e. the
Brunt-Vaisala frequency squared) at the mid pressure from the equation,
( beta x d(SA) - alpha x d(CT) )
N2 = g2 x ---------------------------------
specvol_local x dP
Note. This routine uses specvol from "gsw_specvol", which is the
computationally efficient 75-term expression for specific volume in
terms of SA, CT and p (Roquet et al., 2015).
Note also that the pressure increment, dP, in the above formula is in
Pa, so that it is 104 times the pressure increment dp in dbar.
Note that the 75-term equation has been fitted in a restricted range of
parameter space, and is most accurate inside the "oceanographic funnel"
described in McDougall et al. (2003). The GSW library function
"gsw_infunnel(SA,CT,p)" is avaialble to be used if one wants to test if
some of one's data lies outside this "funnel".
INPUT:
SA = Absolute Salinity [ g/kg ]
CT = Conservative Temperature [ deg C ]
p = sea pressure [ dbar ]
( i.e. absolute pressure - 10.1325 dbar )
OPTIONAL:
lat = latitude in decimal degrees north [ -90 ... +90 ]
Note. If lat is not supplied, a default gravitational acceleration
of 9.7963 m/s2 (Griffies, 2004) will be applied.
SA & CT need to have the same dimensions.
p & lat may have dimensions 1x1 or Mx1 or 1xN or MxN, where SA & CT
are MxN.
OUTPUT:
N2 = Brunt-Vaisala Frequency squared (M-1xN) [rad2 s-2 ]
N2_p = pressure of minimum N2 [ dbar ]
N2_specvol = specific volume at the minimum N2 [ kg m-3 ]
N2_alpha = thermal expansion coefficient with respect [ K-1 ]
to Conservative Temperature at the minimum N2
N2_beta = saline contraction coefficient at constant [ kg g-1 ]
Conservative Temperature at the minimum N2
dSA = difference in salinity between bottles [ g kg-1 ]
dCT = difference in Conservative Temperature between [ deg C ]
bottles
dp = difference in pressure between bottles [ dbar ]
The units of N2 are radians2 s-2 however in may textbooks this is
abreviated to s-2 as radians does not have a unit. To convert the
frequency to hertz, cycles sec-1, divide the frequency by 2π, ie N/(2π).
EXAMPLE:
SA = [34.7118; 34.8915; 35.0256; 34.8472; 34.7366; 34.7324;]
CT = [28.8099; 28.4392; 22.7862; 10.2262; 6.8272; 4.3236;]
p = [ 10; 50; 125; 250; 600; 1000;]
lat = 4;
[N2, N2_p, N2_specvol, N2_alpha, N2_beta, dSA, dCT, dp] = ...
gsw_Nsquared_min(SA,CT,p,lat)
N2 =
1.0e-003 *
0.0608
0.2204
0.1606
0.0116
0.0079
N2_p =
50
125
250
600
1000
N2_specvol =
1.0e-03 *
0.9782
0.9762
0.9730
0.9710
0.9690
N2_alpha =
1.0e-03 *
0.3227
0.2811
0.1732
0.1463
0.1294
N2_beta =
1.0e-03 *
0.7176
0.7262
0.7505
0.7551
0.7571
dSA =
0.1797
0.1341
-0.1784
-0.1106
-0.0042
dCT =
-0.3707
-5.6530
-12.5600
-3.3990
-2.5036
dp =
40
75
125
350
400
AUTHOR:
Trevor McDougall and Paul Barker. [ help@teos-10.org ]
VERSION NUMBER:
3.06.13 (4th August, 2021)
REFERENCES:
Griffies, S. M., 2004: Fundamentals of Ocean Climate Models. Princeton,
NJ: Princeton University Press, 518 pp + xxxiv.
IOC, SCOR and IAPSO, 2010: The international thermodynamic equation of
seawater - 2010: Calculation and use of thermodynamic properties.
Intergovernmental Oceanographic Commission, Manuals and Guides No. 56,
UNESCO (English), 196 pp. Available from the TEOS-10 web site.
See section 3.10 and Eqn. (3.10.2) of this TEOS-10 Manual.
McDougall, T.J., D.R. Jackett, D.G. Wright and R. Feistel, 2003:
Accurate and computationally efficient algorithms for potential
temperature and density of seawater. J. Atmosph. Ocean. Tech., 20,
pp. 730-741.
Roquet, F., G. Madec, T.J. McDougall and P.M. Barker, 2015: Accurate
polynomial expressions for the density and specific volume of seawater
using the TEOS-10 standard. Ocean Modelling, 90, pp. 29-43.
http://dx.doi.org/10.1016/j.ocemod.2015.04.002
The software is available from http://www.TEOS-10.org