library
¶
- seawater.library.T68conv(T90)¶
Convert ITS-90 temperature to IPTS-68.
\(T68 = T90 * 1.00024\)
- Parameters:
- tarray_like
temperature [℃ (ITS-90)]
- Returns:
- tarray_like
temperature [℃ (IPTS-68)]
Notes
The International Practical Temperature Scale of 1968 (IPTS-68) need to be correct to the ITS-90. This linear transformation is accurate within 0.5 ℃ for conversion between IPTS-68 and ITS-90 over the oceanographic temperature range.
References
[1]Saunders, P. M., 1991: The International Temperature Scale of 1990, ITS-90. WOCE Newsletter, No. 10, WOCE International Project Office, Southampton, United Kingdom, 10.
Examples
>>> import seawater as sw >>> T68conv(19.995201151723585) 20.0
- seawater.library.T90conv(t, t_type='T68')¶
Convert IPTS-68 or IPTS-48 to temperature to ITS-90.
T48 apply to all data collected prior to 31/12/1967. T68 apply to all data collected between 01/10/1968 and 31/12/1989.
\[ \begin{align}\begin{aligned}T90 = T68 / 1.00024\\T90 = T48 - (4.4e-6) * T48 * (100-T48) ) / 1.00024\end{aligned}\end{align} \]- Parameters:
- tarray_like
temperature [℃ (IPTS-68) or (IPTS-48)]
- t_typestring, optional
‘T68’ (default) or ‘T48’
- Returns:
- T90array_like
temperature [℃ (ITS-90)]
Notes
The International Practical Temperature Scale of 1968 (IPTS-68) need to be correct to the ITS-90. This linear transformation is accurate within 0.5 ℃ for conversion between IPTS-68 and ITS-90 over the oceanographic temperature range.
References
[1]Saunders, P. M., 1991: The International Temperature Scale of 1990, ITS-90. WOCE Newsletter, No. 10, WOCE International Project Office, Southampton, United Kingdom, 10.
Examples
>>> T90conv(20.004799999999999) 20.0 >>> T90conv(20.0, t_type="T48") 19.98816284091818
- seawater.library.cndr(s, t, p)¶
Calculates conductivity ratio.
- Parameters:
- s(p)array_like
salinity [psu (PSS-78)]
- t(p)array_like
temperature [℃ (ITS-90)]
- parray_like
pressure [db]
- Returns:
- cndrarray_like
conductivity ratio. R = C(s,t,p) / C(35,15(IPTS-68),0) [no units]
References
[1]Fofonoff, P. and Millard, R.C. Jr UNESCO 1983. Algorithms for computation of fundamental properties of seawater. UNESCO Tech. Pap. in Mar. Sci., No. 44, 53 pp. Eqn.(31) p.39. https://unesdoc.unesco.org/ark:/48223/pf0000059832_eng
Examples
>>> # Data from UNESCO 1983 p9. >>> import seawater as sw >>> t = T90conv([0, 10, 0, 10, 10, 30]) >>> p = [0, 0, 1000, 1000, 0, 0] >>> s = [25, 25, 25, 25, 40, 40] >>> sw.cndr(s, t, p) array([0.49800825, 0.65499015, 0.50624434, 0.66297496, 1.00007311, 1.52996697])
- seawater.library.salds(rtx, delt)¶
Calculates Salinity differential (\(\\frac{dS}{d(\\sqrt{Rt})}\)) at constant temperature.
- Parameters:
- rtxarray_like
\(\\sqrt{rt}\)
- deltarray_like
t-15 [℃ (IPTS-68)]
- Returns:
- dsarray_like
\(\\frac{dS}{d rtx}\)
References
[1]Fofonoff, P. and Millard, R.C. Jr UNESCO 1983. Algorithms for computation of fundamental properties of seawater. UNESCO Tech. Pap. in Mar. Sci., No. 44, 53 pp. Eqn.(31) p.39. https://unesdoc.unesco.org/ark:/48223/pf0000059832_eng
Examples
>>> # Data from UNESCO 1983 p9. >>> import numpy as np >>> import seawater as sw >>> delt = T90conv([15, 20, 5]) - 15 >>> rtx = np.array([1, 1.0568875, 0.81705885]) ** 0.5 >>> sw.salds(rtx, delt) array([78.31921607, 81.5689307 , 68.19023687])
- seawater.library.salrp(r, t, p)¶
Equation for Rp used in calculating salinity. UNESCO 1983 polynomial.
\[\begin{split}Rp(S,T,P) = \\frac{C(S,T,P)}{C(S,T,0)}\end{split}\]- Parameters:
- rarray_like
conductivity ratio \(R = \\frac{C(S,T,P)}{C(35,15(IPTS-68),0)}\)
- tarray_like
temperature [℃ (ITS-90)]
- parray_like
pressure [db]
- Returns:
- rparray_like
conductivity ratio \(Rp(S,T,P) = \\frac{C(S,T,P)}{C(S,T,0)}\)
References
[1]Fofonoff, P. and Millard, R.C. Jr UNESCO 1983. Algorithms for computation of fundamental properties of seawater. UNESCO Tech. Pap. in Mar. Sci., No. 44, 53 pp. Eqn.(31) p.39. https://unesdoc.unesco.org/ark:/48223/pf0000059832_eng
Examples
>>> import seawater as sw >>> r = [1, 1.2, 0.65] >>> t = T90conv([15, 20, 5]) >>> p = [0, 2000, 1500] >>> sw.salrp(r, t, p) array([1. , 1.01694294, 1.02048638])
- seawater.library.salrt(t)¶
Equation for rt used in calculating salinity. UNESCO 1983 polynomial.
\[\begin{split}rt(t) = \\frac{C(35,t,0)}{C(35,15(\\textrm{IPTS-68}), 0)}\end{split}\]- Parameters:
- tarray_like
temperature [℃ (ITS-90)]
- Returns:
- rtarray_like
- conductivity ratio [no units]
References
[1]Fofonoff, P. and Millard, R.C. Jr UNESCO 1983. Algorithms for computation of fundamental properties of seawater. UNESCO Tech. Pap. in Mar. Sci., No. 44, 53 pp. Eqn.(31) p.39. https://unesdoc.unesco.org/ark:/48223/pf0000059832_eng
Examples
>>> # Data from UNESCO 1983 p9. >>> import seawater as sw >>> t = T90conv([15, 20, 5]) >>> sw.salrt(t) array([1. , 1.11649272, 0.77956585])
- seawater.library.sals(rt, t)¶
Salinity of sea water as a function of Rt and T. UNESCO 1983 polynomial.
- Parameters:
- rtarray_like
\(rt(s,t) = \\frac{C(s,t,0)}{C(35, t(\\textrm{IPTS-68}), 0)}\)
- tarray_like
temperature [℃ (ITS-90)]
- Returns:
- sarray_like
salinity [psu (PSS-78)]
References
[1]Fofonoff, P. and Millard, R.C. Jr UNESCO 1983. Algorithms for computation of fundamental properties of seawater. UNESCO Tech. Pap. in Mar. Sci., No. 44, 53 pp. Eqn.(31) p.39. https://unesdoc.unesco.org/ark:/48223/pf0000059832_eng
Examples
>>> # Data from UNESCO 1983 p9. >>> import seawater as sw >>> t = T90conv([15, 20, 5]) >>> rt = [1, 1.0568875, 0.81705885] >>> sw.sals(rt, t) array([35. , 37.24562718, 27.99534701])
- seawater.library.seck(s, t, p=0)¶
Secant Bulk Modulus (K) of Sea Water using Equation of state 1980. UNESCO polynomial implementation.
- Parameters:
- s(p)array_like
salinity [psu (PSS-78)]
- t(p)array_like
temperature [℃ (ITS-90)]
- parray_like
pressure [db].
- Returns:
- karray_like
secant bulk modulus [bars]
References
[1]Fofonoff, P. and Millard, R.C. Jr UNESCO 1983. Algorithms for computation of fundamental properties of seawater. UNESCO Tech. Pap. in Mar. Sci., No. 44, 53 pp. Eqn.(31) p.39. https://unesdoc.unesco.org/ark:/48223/pf0000059832_eng
[2]Millero, F.J. and Poisson, A. International one-atmosphere equation of state of seawater. Deep-Sea Res. 1981. Vol28A(6) pp625-629. doi:10.1016/0198-0149(81)90122-9
Examples
>>> # Data from Unesco Tech. Paper in Marine Sci. No. 44, p22. >>> import seawater as sw >>> s = [0, 0, 0, 0, 35, 35, 35, 35] >>> t = T90conv([0, 0, 30, 30, 0, 0, 30, 30]) >>> p = [0, 10000, 0, 10000, 0, 10000, 0, 10000] >>> sw.seck(s, t, p) array([19652.21 , 22977.2115 , 22336.0044572 , 25656.8196222 , 21582.27006823, 24991.99729129, 23924.21823158, 27318.32472464])
- seawater.library.smow(t)¶
Density of Standard Mean Ocean Water (Pure Water) using EOS 1980.
- Parameters:
- tarray_like
temperature [℃ (ITS-90)]
- Returns:
- dens(t)array_like
density [kg m 3]
References
[1]Fofonoff, P. and Millard, R.C. Jr UNESCO 1983. Algorithms for computation of fundamental properties of seawater. UNESCO Tech. Pap. in Mar. Sci., No. 44, 53 pp. Eqn.(31) p.39. https://unesdoc.unesco.org/ark:/48223/pf0000059832_eng
[2]Millero, F.J. and Poisson, A. International one-atmosphere equation of state of seawater. Deep-Sea Res. 1981. Vol28A(6) pp625-629. doi:10.1016/0198-0149(81)90122-9
Examples
>>> # Data from UNESCO Tech. Paper in Marine Sci. No. 44, p22. >>> import seawater as sw >>> t = T90conv([0, 0, 30, 30, 0, 0, 30, 30]) >>> sw.smow(t) array([999.842594 , 999.842594 , 995.65113374, 995.65113374, 999.842594 , 999.842594 , 995.65113374, 995.65113374])