1
# -*- coding: utf-8 -*-
3
from __future__ import division
7
from library import gibbs
8
from absolute_salinity_sstar_ct import CT_from_t
9
from gsw.utilities import match_args_return, strip_mask
10
from conversions import pt_from_CT, pt_from_t, pt0_from_t
11
from constants import Kelvin, db2Pascal, P0, SSO, cp0, R, sfac, M_S
17
'alpha_wrt_CT_t_exact',
18
'alpha_wrt_pt_t_exact',
20
'beta_const_CT_t_exact',
21
'beta_const_pt_t_exact',
24
'specvol_anom_t_exact',
25
'sound_speed_t_exact',
27
'kappa_const_t_exact',
28
'internal_energy_t_exact',
30
'dynamic_enthalpy_t_exact',
31
'SA_from_rho_t_exact',
36
'isochoric_heat_cap_t_exact',
37
'chem_potential_relative_t_exact',
38
'chem_potential_water_t_exact',
39
'chem_potential_salt_t_exact',
40
'Helmholtz_energy_t_exact',
41
'adiabatic_lapse_rate_t_exact',
42
'osmotic_coefficient_t_exact',
43
'osmotic_pressure_t_exact'
50
def Helmholtz_energy_t_exact(SA, t, p):
51
r"""Calculates the Helmholtz energy of seawater.
53
The specific Helmholtz energy of seawater :math:`f` is given by:
56
f(SA, t, p) = g - (p + P_0) \nu =
57
g - (p + P_0) \frac{\partial g}{\partial P}\Big|_{SA,T}
62
Absolute salinity [g kg :sup:`-1`]
64
in situ temperature [:math:`^\circ` C (ITS-90)]
70
Helmholtz_energy : array_like
71
Helmholtz energy [J kg :sup:`-1`]
84
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
85
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
86
>>> p = [10, 50, 125, 250, 600, 1000]
87
>>> gsw.Helmholtz_energy_t_exact(SA, t, p)
88
array([-5985.58288209, -5830.81845224, -3806.96617841, -877.66369421,
89
-462.17033905, -245.50407205])
93
.. [1] IOC, SCOR and IAPSO, 2010: The international thermodynamic equation
94
of seawater - 2010: Calculation and use of thermodynamic properties.
95
Intergovernmental Oceanographic Commission, Manuals and Guides No. 56,
96
UNESCO (English), 196 pp. See section 2.13.
99
2011-03-29. Trevor McDougall
102
return (gibbs(n0, n0, n0, SA, t, p) -
103
(db2Pascal * p + P0) * gibbs(n0, n0, n1, SA, t, p))
107
def rho_t_exact(SA, t, p):
108
r"""Calculates in situ density of seawater from Absolute Salinity and in
114
Absolute salinity [g kg :sup:`-1`]
116
in situ temperature [:math:`^\circ` C (ITS-90)]
122
rho_t_exact : array_like
123
in situ density [kg m :sup:`-3`]
136
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
137
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
138
>>> p = [10, 50, 125, 250, 600, 1000]
139
>>> gsw.rho(SA, t, p)
140
array([ 1021.84017319, 1022.26268993, 1024.42771594, 1027.79020181,
141
1029.83771473, 1032.00240412])
145
.. [1] IOC, SCOR and IAPSO, 2010: The international thermodynamic equation
146
of seawater - 2010: Calculation and use of thermodynamic properties.
147
Intergovernmental Oceanographic Commission, Manuals and Guides No. 56,
148
UNESCO (English), 196 pp. See section 2.8.
151
2011-03-29. Paul Barker, David Jackett and Trevor McDougal
154
return 1. / gibbs(n0, n0, n1, SA, t, p)
158
def sigma0_pt0_exact(SA, pt0):
159
r"""Calculates potential density anomaly with reference sea pressure of
160
zero (0) dbar. The temperature input to this function is potential
161
temperature referenced to zero dbar.
166
Absolute salinity [g kg :sup:`-1`]
168
potential temperature [:math:`^\circ` C (ITS-90)]
169
with respect to a reference sea pressure of 0 dbar
173
sigma0_pt0_exact : array_like
174
potential density anomaly [kg m :sup:`-3`]
175
respect to a reference pressure of 0 dbar
188
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
189
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
190
>>> p = [10, 50, 125, 250, 600, 1000]
191
>>> gsw.rho(SA, t, p)
192
array([ 1021.84017319, 1022.26268993, 1024.42771594, 1027.79020181,
193
1029.83771473, 1032.00240412])
197
.. [1] IOC, SCOR and IAPSO, 2010: The international thermodynamic equation
198
of seawater - 2010: Calculation and use of thermodynamic properties.
199
Intergovernmental Oceanographic Commission, Manuals and Guides No. 56,
200
UNESCO (English), 196 pp. See Eqn. (3.6.1).
203
2011-03-29. Trevor McDougal and Paul Barker.
205
SA = np.maximum(SA, 0) # Ensure that SA is non-negative.
211
g03 = (100015.695367145 +
212
y * (-270.983805184062 +
213
y * (1455.0364540468 +
214
y * (-672.50778314507 +
215
y * (397.968445406972 +
216
y * (-194.618310617595 +
217
y * (63.5113936641785 -
218
y * 9.63108119393062)))))))
220
g08 = x2 * (-3310.49154044839 +
221
x * (199.459603073901 +
222
x * (-54.7919133532887 +
223
x * 36.0284195611086 -
224
y * 22.6683558512829) +
225
y * (-175.292041186547 +
226
y * (383.058066002476 +
227
y * (-460.319931801257 +
228
y * 234.565187611355)))) +
229
y * (729.116529735046 +
230
y * (-860.764303783977 +
231
y * (694.244814133268 +
232
y * (-297.728741987187)))))
234
"""The above code is exactly the same as the following two lines of code.
235
sigma0_pt_exact = rho_t_exact(SA, pt0, 0.) - 1000
238
return 100000000. / (g03 + g08) - 1000.0
242
def enthalpy_t_exact(SA, t, p):
243
r"""Calculates the specific enthalpy of seawater.
245
The specific enthalpy of seawater :math:`h` is given by:
248
h(SA, t, p) = g + (T_0 + t)\eta =
249
g - (T_0 + t) \frac{\partial g}{\partial T}\Big|_{SA,p}
254
Absolute salinity [g kg :sup:`-1`]
256
in situ temperature [:math:`^\circ` C (ITS-90)]
262
enthalpy : array_like
263
specific enthalpy [J kg :sup:`-1`]
276
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
277
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
278
>>> p = [10, 50, 125, 250, 600, 1000]
279
>>> gsw.enthalpy(SA, t, p)
280
array([ 115103.26047838, 114014.8036012 , 92179.9209311 ,
281
43255.32838089, 33087.21597002, 26970.5880448 ])
285
.. [1] IOC, SCOR and IAPSO, 2010: The international thermodynamic equation
286
of seawater - 2010: Calculation and use of thermodynamic properties.
287
Intergovernmental Oceanographic Commission, Manuals and Guides No. 56,
288
UNESCO (English), 196 pp. See appendix A.11.
291
2011-03-29. David Jackett, Trevor McDougall and Paul Barker.
294
return (gibbs(n0, n0, n0, SA, t, p) -
295
(t + Kelvin) * gibbs(n0, n1, n0, SA, t, p))
299
def specvol_t_exact(SA, t, p):
300
r"""Calculates the specific volume of seawater.
305
Absolute salinity [g kg :sup:`-1`]
307
in situ temperature [:math:`^\circ` C (ITS-90)]
314
specific volume [m :sup:`3` kg :sup:`-1`]
327
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
328
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
329
>>> p = [10, 50, 125, 250, 600, 1000]
330
>>> gsw.specvol(SA, t, p)
331
array([ 0.00097863, 0.00097822, 0.00097615, 0.00097296, 0.00097103,
336
.. [1] IOC, SCOR and IAPSO, 2010: The international thermodynamic equation
337
of seawater - 2010: Calculation and use of thermodynamic properties.
338
Intergovernmental Oceanographic Commission, Manuals and Guides No. 56,
339
UNESCO (English), 196 pp. See section 2.7.
342
2011-03-23. David Jackett and Paul Barker.
345
return gibbs(n0, n0, n1, SA, t, p)
349
def entropy_t_exact(SA, t, p):
350
r"""Calculates specific entropy of seawater.
352
The specific entropy of seawater :math:`\eta` is given by:
355
\eta(SA, t, p) = -g_T = \frac{\partial g}{\partial T}\Big|_{SA,p}
357
When taking derivatives with respect to *in situ* temperature, the symbol
358
:math:`T` will be used for temperature in order that these derivatives not
359
be confused with time derivatives.
364
Absolute salinity [g kg :sup:`-1`]
366
in situ temperature [:math:`^\circ` C (ITS-90)]
373
specific entropy [J kg :sup:`-1` K :sup:`-1`]
386
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
387
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
388
>>> p = [10, 50, 125, 250, 600, 1000]
389
>>> gsw.entropy_t_exact(SA, t, p)
390
array([ 400.38942528, 395.43817843, 319.8664982 , 146.79088159,
391
98.64734087, 62.79150873])
395
.. [1] IOC, SCOR and IAPSO, 2010: The international thermodynamic equation
396
of seawater - 2010: Calculation and use of thermodynamic properties.
397
Intergovernmental Oceanographic Commission, Manuals and Guides No. 56,
398
UNESCO (English), 196 pp.
401
2011-03-29. David Jackett, Trevor McDougall and Paul Barker.
404
return -gibbs(n0, n1, n0, SA, t, p)
408
def cp_t_exact(SA, t, p):
409
r"""Calculates the isobaric heat capacity of seawater.
414
Absolute salinity [g kg :sup:`-1`]
416
in situ temperature [:math:`^\circ` C (ITS-90)]
422
cp_t_exact : array_like
423
heat capacity of seawater [J kg :sup:`-1` K :sup:`-1`]
436
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
437
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
438
>>> p = [10, 50, 125, 250, 600, 1000]
439
>>> gsw.cp_t_exact(SA, t, p)
440
array([ 4002.88800396, 4000.98028393, 3995.54646889, 3985.07676902,
441
3973.59384348, 3960.18408479])
443
.. [1] IOC, SCOR and IAPSO, 2010: The international thermodynamic equation
444
of seawater - 2010: Calculation and use of thermodynamic properties.
445
Intergovernmental Oceanographic Commission, Manuals and Guides No. 56,
446
UNESCO (English), 196 pp.
449
2011-03-29. David Jackett, Trevor McDougall and Paul Barker
452
return -(t + Kelvin) * gibbs(n0, n2, n0, SA, t, p)
456
def sound_speed_t_exact(SA, t, p):
457
r"""Calculates the speed of sound in seawater.
459
The speed of sound in seawater :math:`c` is given by:
462
c(SA, t, p) = \sqrt{ \partial P / \partial \rho |_{SA,\eta}} =
463
\sqrt{(\rho\kappa)^{-1}} =
464
g_P \sqrt{g_{TT}/(g^2_{TP} - g_{TT}g_{PP})}
466
Note that in these expressions, since sound speed is in m s :sup`-1` and
467
density has units of kg m :sup:`-3` it follows that the pressure of the
468
partial derivatives must be in Pa and the isentropic compressibility
469
:math:`kappa` must have units of Pa :sup:`-1`. The sound speed c produced
470
by both the SIA and the GSW software libraries (appendices M and N) has
471
units of m s :sup:`-1`.
476
Absolute salinity [g kg :sup:`-1`]
478
in situ temperature [:math:`^\circ` C (ITS-90)]
484
sound_speed : array_like
485
speed of sound in seawater [m s :sup:`-1`]
498
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
499
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
500
>>> p = [10, 50, 125, 250, 600, 1000]
501
>>> gsw.sound_speed_t_exact(SA, t, p)
502
array([ 1542.61580359, 1542.70353407, 1530.84497914, 1494.40999692,
503
1487.37710252, 1483.93460908])
507
.. [1] IOC, SCOR and IAPSO, 2010: The international thermodynamic equation
508
of seawater - 2010: Calculation and use of thermodynamic properties.
509
Intergovernmental Oceanographic Commission, Manuals and Guides No. 56,
510
UNESCO (English), 196 pp. See Eqn. (2.17.1)
513
2011-03-29. David Jackett, Paul Barker and Trevor McDougall.
516
return (gibbs(n0, n0, n1, SA, t, p) * np.sqrt(gibbs(n0, n2, n0, SA, t, p) /
517
(gibbs(n0, n1, n1, SA, t, p) ** 2 - gibbs(n0, n2, n0, SA, t, p) *
518
gibbs(n0, n0, n2, SA, t, p))))
522
def specvol_anom_t_exact(SA, t, p):
523
r"""Calculates specific volume anomaly from Absolute Salinity, in situ
524
temperature and pressure, using the full TEOS-10 Gibbs function.
526
The reference value of Absolute Salinity is SSO and the reference value of
527
Conservative Temperature is equal to 0 :math:`^\circ` C.
532
Absolute salinity [g kg :sup:`-1`]
534
in situ temperature [:math:`^\circ` C (ITS-90)]
540
specvol_anom_t_exact : array_like
541
specific volume anomaly [m :sup:`3` kg :sup:`-1`]
554
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
555
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
556
>>> p = [10, 50, 125, 250, 600, 1000]
557
>>> gsw.specvol_anom_t_exact(SA, t, p)
558
array([ 6.01044463e-06, 5.78602432e-06, 4.05564999e-06,
559
1.42198662e-06, 1.04351837e-06, 7.63964850e-07])
563
.. [1] IOC, SCOR and IAPSO, 2010: The international thermodynamic equation
564
of seawater - 2010: Calculation and use of thermodynamic properties.
565
Intergovernmental Oceanographic Commission, Manuals and Guides No. 56,
566
UNESCO (English), 196 pp. See Eqn. (3.7.3)
569
2011-03-23. Trevor McDougall and Paul Barker
572
pt_zero = pt_from_CT(SSO, 0)
573
t_zero = pt_from_t(SSO, pt_zero, 0, p)
574
return (gibbs(n0, n0, n1, SA, t, p) -
575
gibbs(n0, n0, n1, SSO, t_zero, p))
579
def chem_potential_relative_t_exact(SA, t, p):
580
r"""Calculates the adiabatic lapse rate of seawater.
585
Absolute salinity [g kg :sup:`-1`]
587
in situ temperature [:math:`^\circ` C (ITS-90)]
593
chem_potential_relative : array_like
594
relative chemical potential [J kg :sup:`-1`]
607
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
608
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
609
>>> p = [10, 50, 125, 250, 600, 1000]
610
>>> gsw.chem_potential_relative_t_exact(SA, t, p)
611
array([ 79.4254481 , 79.25989214, 74.69154859, 65.64063719,
612
61.22685656, 57.21298557])
616
.. [1] IOC, SCOR and IAPSO, 2010: The international thermodynamic equation
617
of seawater - 2010: Calculation and use of thermodynamic properties.
618
Intergovernmental Oceanographic Commission, Manuals and Guides No. 56,
619
UNESCO (English), 196 pp.
622
2011-03-29. Trevor McDougall and Paul Barker
625
return gibbs(n1, n0, n0, SA, t, p)
629
def internal_energy_t_exact(SA, t, p):
630
r"""Calculates the Helmholtz energy of seawater.
632
The specific internal energy of seawater :math:`u` is given by:
635
u(SA, t, p) = g + (T_0 + t)\eta - (p + P_0)\nu =
636
g - (T_0 + t)\frac{\partial g}{\partial T}\Big|_{SA,p} -
637
(p + P_0)\frac{\partial g}{\partial P}\Big|_{SA,T}
639
where :math:`T_0` is the Celsius zero point, 273.15 K and
640
:math:`P_0` = 101 325 Pa is the standard atmosphere pressure.
645
Absolute salinity [g kg :sup:`-1`]
647
in situ temperature [:math:`^\circ` C (ITS-90)]
653
internal_energy (u) : array_like
654
specific internal energy [J kg :sup:`-1`]
667
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
668
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
669
>>> p = [10, 50, 125, 250, 600, 1000]
670
>>> gsw.internal_energy_t_exact(SA, t, p)
671
array([ 114906.23847309, 113426.57417062, 90860.81858842,
672
40724.34005719, 27162.66600185, 17182.50522667])
676
.. [1] IOC, SCOR and IAPSO, 2010: The international thermodynamic equation
677
of seawater - 2010: Calculation and use of thermodynamic properties.
678
Intergovernmental Oceanographic Commission, Manuals and Guides No. 56,
679
UNESCO (English), 196 pp. See Eqn. (2.11.1)
682
2011-03-29. Trevor McDougall
685
return (gibbs(n0, n0, n0, SA, t, p) -
686
(Kelvin + t) * gibbs(n0, n1, n0, SA, t, p) -
687
(db2Pascal * p + P0) * gibbs(n0, n0, n1, SA, t, p))
691
def kappa_const_t_exact(SA, t, p):
692
r"""Calculates isothermal compressibility of seawater at constant in situ
697
\rho^{-1}\frac{\partial \rho}{\partial P}\Big|_{SA,T} =
698
-\nu^{-1}\frac{\partial \nu}{\partial P}\Big|_{SA,T} =
704
Absolute salinity [g kg :sup:`-1`]
706
in situ temperature [:math:`^\circ` C (ITS-90)]
713
Isothermal compressibility [Pa :sup:`-1`]
721
This is the compressibility of seawater at constant in situ temperature.
726
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
727
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
728
>>> p = [10, 50, 125, 250, 600, 1000]
729
>>> gsw.kappa_const_t_exact(SA, t, p)
730
array([ 4.19071646e-10, 4.18743202e-10, 4.22265764e-10,
731
4.37735100e-10, 4.40373818e-10, 4.41156577e-10])
735
.. [1] IOC, SCOR and IAPSO, 2010: The international thermodynamic equation
736
of seawater - 2010: Calculation and use of thermodynamic properties.
737
Intergovernmental Oceanographic Commission, Manuals and Guides No. 56,
738
UNESCO (English), 196 pp. See Eqn. (2.15.1)
741
2011-03-29. David Jackett, Trevor McDougall and Paul Barker
744
return -gibbs(n0, n0, n2, SA, t, p) / gibbs(n0, n0, n1, SA, t, p)
748
def alpha_wrt_t_exact(SA, t, p):
749
r"""Calculates the thermal expansion coefficient of seawater with respect
750
to in situ temperature.
755
Absolute salinity [g kg :sup:`-1`]
757
in situ temperature [:math:`^\circ` C (ITS-90)]
763
alpha_wrt_t : array_like
764
thermal expansion coefficient [K :sup:`-1`]
777
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
778
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
779
>>> p = [10, 50, 125, 250, 600, 1000]
780
>>> gsw.alpha_wrt_t_exact(SA, t, p)
781
array([ 0.0003256 , 0.00032345, 0.00028141, 0.00017283, 0.00014557,
786
.. [1] IOC, SCOR and IAPSO, 2010: The international thermodynamic equation
787
of seawater - 2010: Calculation and use of thermodynamic properties.
788
Intergovernmental Oceanographic Commission, Manuals and Guides No. 56,
789
UNESCO (English), 196 pp. See Eqn. (2.18.1)
791
.. [2] McDougall, T.J., D.R. Jackett and F.J. Millero, 2010: An algorithm
792
for estimating Absolute Salinity in the global ocean. Submitted to Ocean
793
Science. A preliminary version is available at Ocean Sci. Discuss.,
797
2011-03-29. David Jackett, Trevor McDougall and Paul Barker
800
return gibbs(n0, n1, n1, SA, t, p) / gibbs(n0, n0, n1, SA, t, p)
804
def isochoric_heat_cap_t_exact(SA, t, p):
805
r"""Calculates the isochoric heat capacity of seawater.
810
Absolute salinity [g kg :sup:`-1`]
812
in situ temperature [:math:`^\circ` C (ITS-90)]
818
isochoric_heat_cap : array_like
819
isochoric heat capacity [J kg :sup:`-1` K :sup:`-1`]
832
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
833
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
834
>>> p = [10, 50, 125, 250, 600, 1000]
835
>>> gsw.isochoric_heat_cap_t_exact(SA, t, p)
836
array([ 3928.13708702, 3927.27381633, 3941.36418525, 3966.26126146,
837
3960.50903222, 3950.13901342])
841
.. [1] IOC, SCOR and IAPSO, 2010: The international thermodynamic equation
842
of seawater - 2010: Calculation and use of thermodynamic properties.
843
Intergovernmental Oceanographic Commission, Manuals and Guides No. 56,
844
UNESCO (English), 196 pp. See section 2.21.
847
2011-03-29. Trevor McDougall
850
return (-(Kelvin + t) * (gibbs(n0, n2, n0, SA, t, p) -
851
gibbs(n0, n1, n1, SA, t, p) ** 2 / gibbs(n0, n0, n2, SA, t, p)))
855
def kappa_t_exact(SA, t, p):
856
r"""Calculates the isentropic compressibility of seawater.
858
When the entropy and Absolute Salinity are held constant while the pressure
859
is changed, the isentropic and isohaline compressibility
860
:math:`kappa` is obtained:
864
\rho^{-1}\frac{\partial \rho}{\partial P}\Big|_{SA,\eta} =
865
-\nu^{-1}\frac{\partial \nu}{\partial P}\Big|_{SA,\eta} =
866
\rho^{-1}\frac{\partial \rho}{\partial P}\Big|_{SA,\theta} =
867
-\nu^{-1}\frac{\partial \nu}{\partial P}\Big|_{SA,\theta} =
868
-\frac{ (g_{TP}^2 - g_{TT} g_{PP} ) }{g_P g_{TT}}
870
The isentropic and isohaline compressibility is sometimes called simply the
871
isentropic compressibility (or sometimes the "adiabatic compressibility"),
872
on the unstated understanding that there is also no transfer of salt during
873
the isentropic or adiabatic change in pressure.
878
Absolute salinity [g kg :sup:`-1`]
880
in situ temperature [:math:`^\circ` C (ITS-90)]
887
Isentropic compressibility [Pa :sup:`-1`]
895
The output is Pascal and not dbar.
900
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
901
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
902
>>> p = [10, 50, 125, 250, 600, 1000]
903
>>> gsw.kappa_t_exact(SA, t, p)
904
array([ 4.11245799e-10, 4.11029072e-10, 4.16539558e-10,
905
4.35668338e-10, 4.38923693e-10, 4.40037576e-10])
909
.. [1] IOC, SCOR and IAPSO, 2010: The international thermodynamic equation
910
of seawater - 2010: Calculation and use of thermodynamic properties.
911
Intergovernmental Oceanographic Commission, Manuals and Guides No. 56,
912
UNESCO (English), 196 pp. See Eqns. (2.16.1) and the row for kappa in
913
Table P.1 of appendix P
916
2011-03-23. David Jackett, Trevor McDougall and Paul Barker
919
return ((gibbs(n0, n1, n1, SA, t, p) ** 2 - gibbs(n0, n2, n0, SA, t, p) *
920
gibbs(n0, n0, n2, SA, t, p)) / (gibbs(n0, n0, n1, SA, t, p) *
921
gibbs(n0, n2, n0, SA, t, p)))
925
def SA_from_rho_t_exact(rho, t, p):
926
r"""Calculates the Absolute Salinity of a seawater sample, for given values
927
of its density, in situ temperature and sea pressure (in dbar).
929
One use for this function is in the laboratory where a measured value of
930
the in situ density :math:`\rho` of a seawater sample may have been made at
931
the laboratory temperature :math:`t` and at atmospheric pressure :math:`p`.
932
The present function will return the Absolute Salinity SA of this seawater
938
in situ density [kg m :sup:`-3`]
940
in situ temperature [:math:`^\circ` C (ITS-90)]
947
Absolute salinity [g kg :sup:`-1`]
955
This is expressed on the Reference-Composition Salinity Scale of
956
Millero et al. (2008).
958
After two iterations of a modified Newton-Raphson iteration,
959
the error in SA is typically no larger than
960
2 :math:`^\times` 10 :sup:`-13` [g kg :sup:`-1`]
965
>>> rho = [1021.839, 1022.262, 1024.426, 1027.792, 1029.839, 1032.002]
966
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
967
>>> p = [10, 50, 125, 250, 600, 1000]
968
>>> gsw.SA_from_rho_t_exact(rho, t, p)
969
array([ 34.71022966, 34.89057683, 35.02332421, 34.84952096,
970
34.73824809, 34.73188384])
974
.. [1] IOC, SCOR and IAPSO, 2010: The international thermodynamic equation
975
of seawater - 2010: Calculation and use of thermodynamic properties.
976
Intergovernmental Oceanographic Commission, Manuals and Guides No. 56,
977
UNESCO (English), 196 pp. See section 2.5.
979
.. [2] Millero, F. J., R. Feistel, D. G. Wright, and T. J. McDougall, 2008:
980
The composition of Standard Seawater and the definition of the
981
Reference-Composition Salinity Scale, Deep-Sea Res. I, 55, 50-72.
984
2011-03-28. Trevor McDougall and Paul Barker.
987
v_lab = np.ones_like(rho) / rho
988
v_0 = gibbs(n0, n0, n1, 0, t, p)
989
v_120 = gibbs(n0, n0, n1, 120, t, p)
991
# Initial estimate of SA.
992
SA = 120 * (v_lab - v_0) / (v_120 - v_0)
993
Ior = np.logical_or(SA < 0, SA > 120)
995
# Initial estimate of v_SA, SA derivative of v
996
v_SA = (v_120 - v_0) / 120
998
for k in range(0, 2):
1000
delta_v = gibbs(n0, n0, n1, SA_old, t, p) - v_lab
1001
# Half way the mod. N-R method (McDougall and Wotherspoon, 2012)
1002
SA = SA_old - delta_v / v_SA
1003
SA_mean = 0.5 * (SA + SA_old)
1004
v_SA = gibbs(n1, n0, n1, SA_mean, t, p)
1005
SA = SA_old - delta_v / v_SA
1007
SA[Ior] = np.ma.masked
1013
def t_from_rho_exact(rho, SA, p):
1014
r"""Calculates the in-situ temperature of a seawater sample, for given
1015
values of its density, Absolute Salinity and sea pressure (in dbar).
1021
in situ density [kg m :sup:`-3`]
1023
Absolute salinity [g kg :sup:`-1`]
1030
in situ temperature [:math:`^\circ` C (ITS-90)]
1031
t_multiple : array_like
1032
in situ temperature [:math:`^\circ` C (ITS-90)]
1040
At low salinities, in brackish water, there are two possible temperatures
1041
for a single density. This program will output both valid solutions
1042
(t, t_multiple), if there is only one possible solution the second variable
1052
.. [1] IOC, SCOR and IAPSO, 2010: The international thermodynamic equation
1053
of seawater - 2010: Calculation and use of thermodynamic properties.
1054
Intergovernmental Oceanographic Commission, Manuals and Guides No. 56,
1055
UNESCO (English), 196 pp.
1058
2011-04-21. Trevor McDougall and Paul Barker.
1061
"""alpha_limit is the positive value of the thermal expansion coefficient
1062
which is used at the freezing temperature to distinguish between I_salty
1066
"""rec_half_rho_TT is a constant representing the reciprocal of half the
1067
second derivative of density with respect to temperature near the
1068
temperature of maximum density."""
1069
rec_half_rho_TT = -110.0
1071
t = np.zeros_like(SA) + np.NaN
1072
t_multiple = np.zeros_like(SA) + np.NaN
1074
I_SA = np.logical_or(SA < 0, SA > 42)
1075
I_p = np.logical_or(p < -1.5, p > 12000)
1076
I_SA_p = np.logical_or(I_SA, I_p)
1078
SA[I_SA_p] = np.ma.masked
1080
rho_40 = rho_t_exact(SA, 40 * np.ones_like(SA), p)
1082
I_rho_light = (rho - rho_40) < 0
1084
SA[I_rho_light] = np.ma.masked
1086
t_max_rho = t_maxdensity_exact(SA, p)
1087
rho_max = rho_t_exact(SA, t_max_rho, p)
1088
rho_extreme = rho_max
1089
t_freezing = t_freezing(SA, p) # Assumes seawater is saturated with air.
1090
rho_freezing = rho_t_exact(SA, t_freezing, p)
1092
I_fr_gr_max = (t_freezing - t_max_rho) > 0
1093
rho_extreme[I_fr_gr_max] = rho_freezing[I_fr_gr_max]
1095
I_rho_dense = rho > rho_extreme
1096
SA[I_rho_dense] = np.ma.masked
1098
# FIXME: Is this needed?
1099
I_bad = np.isnan(SA * p * rho)
1100
SA[I_bad] = np.ma.masked
1102
alpha_freezing = alpha_wrt_t_exact(SA, t_freezing, p)
1104
I_salty = alpha_freezing > alpha_limit
1106
t_diff = 40. * np.ones_like(I_salty) - t_freezing(I_salty)
1108
top = (rho_40[I_salty] - rho_freezing[I_salty] +
1109
rho_freezing[I_salty] * alpha_freezing[I_salty] * t_diff)
1111
a = top / (t_diff ** 2)
1112
b = -rho_freezing[I_salty] * alpha_freezing[I_salty]
1113
c = rho_freezing[I_salty] - rho[I_salty]
1114
sqrt_disc = np.sqrt(b ** 2 - 4 * a * c)
1115
# The value of t[I_salty] is the initial guess `t` in the range of I_salty.
1116
t[I_salty] = t_freezing[I_salty] + 0.5 * (-b - sqrt_disc) / a
1118
I_fresh = alpha_freezing <= alpha_limit
1119
t_diff = 40 * np.ones_like[I_fresh] - t_max_rho[I_fresh]
1120
factor = ((rho_max[I_fresh] - rho[I_fresh]) /
1121
(rho_max[I_fresh] - rho_40[I_fresh]))
1122
delta_t = t_diff * np.sqrt(factor)
1124
I_fresh_NR = delta_t > 5
1125
t[I_fresh[I_fresh_NR]] = (t_max_rho[I_fresh[I_fresh_NR]] +
1126
delta_t[I_fresh_NR])
1128
I_quad = delta_t <= 5
1129
t_a = np.zeros_like(SA) + np.NaN
1130
# Set the initial value of the quadratic solution roots.
1131
t_a[I_fresh[I_quad]] = (t_max_rho[I_fresh[I_quad]] +
1132
np.sqrt(rec_half_rho_TT * (rho[I_fresh[I_quad]] -
1133
rho_max[I_fresh[I_quad]])))
1135
for Number_of_iterations in range(0, 5):
1137
rho_old = rho_t_exact(SA, t_old, p)
1138
factorqa = (rho_max - rho) / (rho_max - rho_old)
1139
t_a = t_max_rho + (t_old - t_max_rho) * np.sqrt(factorqa)
1141
t_a[t_freezing - t_a < 0] = np.ma.masked
1143
t_b = np.zeros_like(SA) + np.NaN
1144
# Set the initial value of the quadratic solution routes.
1145
t_b[I_fresh[I_quad]] = (t_max_rho[I_fresh[I_quad]] -
1146
np.sqrt(rec_half_rho_TT * (rho[I_fresh[I_quad]] -
1147
rho_max[I_fresh[I_quad]])))
1148
for Number_of_iterations in range(0, 6):
1150
rho_old = rho_t_exact(SA, t_old, p)
1151
factorqb = (rho_max - rho) / (rho_max - rho_old)
1152
t_b = t_max_rho + (t_old - t_max_rho) * np.sqrt(factorqb)
1154
# After seven iterations of this quadratic iterative procedure,
1155
# the error in rho is no larger than 4.6x10^-13 kg/m^3.
1156
t_b[t_freezing - t_b < 0] = np.ma.masked
1158
# Begin the modified Newton-Raphson iterative method, which will only
1159
# operate on non-masked data.
1161
v_lab = np.ones_like(rho) / rho
1162
v_t = gibbs(0, 1, 1, SA, t, p)
1163
for Number_of_iterations in range(0, 3):
1165
delta_v = gibbs(0, 0, 1, SA, t_old, p) - v_lab
1166
t = t_old - delta_v / v_t # Half way through the modified N-R method.
1167
t_mean = 0.5 * (t + t_old)
1168
v_t = gibbs(0, 1, 1, SA, t_mean, p)
1169
t = t_old - delta_v / v_t
1171
I_quad = ~np.isnan(t_a)
1172
t[I_quad] = t_a[I_quad]
1174
I_quad = ~np.isnan(t_b)
1175
t_multiple[I_quad] = t_b[I_quad]
1177
# After three iterations of this modified Newton-Raphson iteration,
1178
# the error in rho is no larger than 4.6x10^-13 kg/m^3.
1180
return t, t_multiple
1184
def pot_rho_t_exact(SA, t, p, p_ref=0):
1185
r"""Calculates potential density of seawater.
1190
Absolute salinity [g kg :sup:`-1`]
1192
in situ temperature [:math:`^\circ` C (ITS-90)]
1195
p_ref : int, float, optional
1196
reference pressure, default = 0
1200
pot_rho : array_like
1201
potential density [kg m :sup:`-3`]
1214
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
1215
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
1216
>>> p = [10, 50, 125, 250, 600, 1000]
1217
>>> gsw.pot_rho_t_exact(SA, t, p)
1218
array([ 1021.79814581, 1022.05248442, 1023.89358365, 1026.66762112,
1219
1027.10723087, 1027.40963126])
1220
>>> gsw.pot_rho(SA, t, p, p_ref=1000)
1221
array([ 1025.95554512, 1026.21306986, 1028.12563226, 1031.1204547 ,
1222
1031.63768355, 1032.00240412])
1226
.. [1] IOC, SCOR and IAPSO, 2010: The international thermodynamic equation
1227
of seawater - 2010: Calculation and use of thermodynamic properties.
1228
Intergovernmental Oceanographic Commission, Manuals and Guides No. 56,
1229
UNESCO (English), 196 pp. See section 3.4.
1232
2011-03-29. David Jackett, Trevor McDougall and Paul Barker
1235
pt = pt_from_t(SA, t, p, p_ref=p_ref)
1237
return rho_t_exact(SA, pt, p_ref)
1241
def alpha_wrt_CT_t_exact(SA, t, p):
1242
r"""Calculates the thermal expansion coefficient of seawater with respect
1243
to Conservative Temperature.
1248
Absolute salinity [g kg :sup:`-1`]
1250
in situ temperature [:math:`^\circ` C (ITS-90)]
1256
alpha_wrt_CT : array_like
1257
thermal expansion coefficient [K :sup:`-1`]
1270
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
1271
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
1272
>>> p = [10, 50, 125, 250, 600, 1000]
1273
>>> gsw.alpha_wrt_CT_t_exact(SA, t, p)
1274
array([ 0.00032471, 0.00032272, 0.00028118, 0.00017314, 0.00014627,
1279
.. [1] IOC, SCOR and IAPSO, 2010: The international thermodynamic equation
1280
of seawater - 2010: Calculation and use of thermodynamic properties.
1281
Intergovernmental Oceanographic Commission, Manuals and Guides No. 56,
1282
UNESCO (English), 196 pp. See Eqn. (2.18.3).
1285
2011-03-29. Trevor McDougall and Paul Barker
1288
pt0 = pt0_from_t(SA, t, p)
1289
factor = -cp0 / ((Kelvin + pt0) * gibbs(n0, n2, n0, SA, t, p))
1290
return factor * (gibbs(n0, n1, n1, SA, t, p) / gibbs(n0, n0, n1, SA, t, p))
1294
def alpha_wrt_pt_t_exact(SA, t, p):
1295
r"""Calculates the thermal expansion coefficient of seawater with respect
1296
to potential temperature, with a reference pressure of zero.
1301
Absolute salinity [g kg :sup:`-1`]
1303
in situ temperature [:math:`^\circ` C (ITS-90)]
1309
alpha_wrt_pt : array_like
1310
thermal expansion coefficient [K :sup:`-1`]
1323
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
1324
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
1325
>>> p = [10, 50, 125, 250, 600, 1000]
1326
>>> gsw.alpha_wrt_pt_t_exact(SA, t, p)
1327
array([ 0.00032562, 0.00032355, 0.00028164, 0.00017314, 0.00014623,
1332
.. [1] IOC, SCOR and IAPSO, 2010: The international thermodynamic equation
1333
of seawater - 2010: Calculation and use of thermodynamic properties.
1334
Intergovernmental Oceanographic Commission, Manuals and Guides No. 56,
1335
UNESCO (English), 196 pp. See Eqn. (2.18.2).
1338
2011-03-29. David Jackett, Trevor McDougall and Paul Barker
1341
pt0 = pt0_from_t(SA, t, p)
1342
factor = gibbs(n0, n2, n0, SA, pt0, 0) / gibbs(n0, n2, n0, SA, t, p)
1343
return factor * (gibbs(n0, n1, n1, SA, t, p) / gibbs(n0, n0, n1, SA, t, p))
1347
def beta_const_CT_t_exact(SA, t, p):
1348
r"""Calculates the saline (i.e. haline) contraction coefficient of seawater
1349
at constant Conservative Temperature.
1354
Absolute salinity [g kg :sup:`-1`]
1356
in situ temperature [:math:`^\circ` C (ITS-90)]
1362
beta_const_CT : array_like
1363
saline contraction coefficient [kg g :sup:`-1`]
1376
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
1377
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
1378
>>> p = [10, 50, 125, 250, 600, 1000]
1379
>>> gsw.beta_const_CT_t_exact(SA, t, p)
1380
array([ 0.00071749, 0.00071765, 0.00072622, 0.00075051, 0.00075506,
1385
.. [1] IOC, SCOR and IAPSO, 2010: The international thermodynamic equation
1386
of seawater - 2010: Calculation and use of thermodynamic properties.
1387
Intergovernmental Oceanographic Commission, Manuals and Guides No. 56,
1388
UNESCO (English), 196 pp. See Eqn. (2.19.3)
1391
2010-07-23. David Jackett, Trevor McDougall and Paul Barker
1394
# TODO: Original GSW-V3 re-implements gibbs, check what to do here!
1396
pt0 = pt0_from_t(SA, t, p)
1398
factora = (gibbs(n1, n1, n0, SA, t, p) - gibbs(n1, n0, n0, SA, pt0, 0) /
1400
factor = (factora / (gibbs(n0, n0, n1, SA, t, p) *
1401
gibbs(n0, n2, n0, SA, t, p)))
1403
return (gibbs(n0, n1, n1, SA, t, p) * factor -
1404
gibbs(n1, n0, n1, SA, t, p) / gibbs(n0, n0, n1, SA, t, p))
1408
def beta_const_pt_t_exact(SA, t, p):
1409
r"""Calculates the saline (i.e. haline) contraction coefficient of seawater
1410
at constant potential temperature with a reference pressure of 0 dbar.
1415
Absolute salinity [g kg :sup:`-1`]
1417
in situ temperature [:math:`^\circ` C (ITS-90)]
1423
beta_const_pt : array_like
1424
saline contraction coefficient [kg g :sup:`-1`]
1437
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
1438
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
1439
>>> p = [10, 50, 125, 250, 600, 1000]
1440
>>> gsw.beta_const_pt_t_exact(SA, t, p)
1441
array([ 0.00073112, 0.00073106, 0.00073599, 0.00075375, 0.00075712,
1446
.. [1] IOC, SCOR and IAPSO, 2010: The international thermodynamic equation
1447
of seawater - 2010: Calculation and use of thermodynamic properties.
1448
Intergovernmental Oceanographic Commission, Manuals and Guides No. 56,
1449
UNESCO (English), 196 pp. See Eqn. (2.19.2)
1452
2011-04-10. Trevor McDougall and Paul Barker
1454
# NOTE: The original Matlab toolbox re-implement some code here. Why?
1456
pt0 = pt0_from_t(SA, t, p)
1458
factora = gibbs(n1, n1, n0, SA, t, p) - gibbs(n1, n1, n0, SA, pt0, 0)
1460
factor = (factora / (gibbs(n0, n0, n1, SA, t, p) *
1461
gibbs(n0, n2, n0, SA, t, p)))
1463
return (gibbs(n0, n1, n1, SA, t, p) * factor -
1464
gibbs(n1, n0, n1, SA, t, p) / gibbs(n0, n0, n1, SA, t, p))
1468
def beta_const_t_exact(SA, t, p):
1469
r"""Calculates the saline (i.e. haline) contraction coefficient of seawater
1470
at constant in situ temperature.
1475
Absolute salinity [g kg :sup:`-1`]
1477
in situ temperature [:math:`^\circ` C (ITS-90)]
1483
beta_const_t : array_like
1484
saline contraction coefficient [kg g :sup:`-1`]
1497
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
1498
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
1499
>>> p = [10, 50, 125, 250, 600, 1000]
1500
>>> gsw.beta_const_t_exact(SA, t, p)
1501
array([ 0.00073112, 0.00073107, 0.00073602, 0.00075381, 0.00075726,
1506
.. [1] IOC, SCOR and IAPSO, 2010: The international thermodynamic equation
1507
of seawater - 2010: Calculation and use of thermodynamic properties.
1508
Intergovernmental Oceanographic Commission, Manuals and Guides No. 56,
1509
UNESCO (English), 196 pp. See Eqn. (2.19.1)
1512
2011-03-29. David Jackett, Trevor McDougall and Paul Barker
1515
return -gibbs(n1, n0, n1, SA, t, p) / gibbs(n0, n0, n1, SA, t, p)
1519
def chem_potential_water_t_exact(SA, t, p):
1520
r"""Calculates the chemical potential of water in seawater.
1525
Absolute salinity [g kg :sup:`-1`]
1527
in situ temperature [:math:`^\circ` C (ITS-90)]
1533
chem_potential_water : array_like
1534
chemical potential of water in seawater
1548
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
1549
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
1550
>>> p = [10, 50, 125, 250, 600, 1000]
1551
>>> gsw.chem_potential_water_t_exact(SA, t, p)
1552
array([-8545.56114628, -8008.08554834, -5103.98013987, -634.06778275,
1553
3335.56680347, 7555.43444597])
1557
.. [1] IOC, SCOR and IAPSO, 2010: The international thermodynamic equation
1558
of seawater - 2010: Calculation and use of thermodynamic properties.
1559
Intergovernmental Oceanographic Commission, Manuals and Guides No. 56,
1560
UNESCO (English), 196 pp.
1563
2011-03-29. Trevor McDougall and Paul Barker
1565
SA, t, p, mask = strip_mask(SA, t, p)
1567
# FIXME: Ugly copy from gibbs, why?
1572
z = p * 1e-4 # Pressure (p) is sea pressure in units of dbar.
1574
g03_g = (101.342743139674 + z * (100015.695367145 +
1575
z * (-2544.5765420363 + z * (284.517778446287 +
1576
z * (-33.3146754253611 + (4.20263108803084 -
1577
0.546428511471039 * z) * z)))) +
1578
y * (5.90578347909402 + z * (-270.983805184062 +
1579
z * (776.153611613101 + z * (-196.51255088122 +
1580
(28.9796526294175 - 2.13290083518327 * z) * z))) +
1581
y * (-12357.785933039 + z * (1455.0364540468 +
1582
z * (-756.558385769359 + z * (273.479662323528 +
1583
z * (-55.5604063817218 + 4.34420671917197 * z)))) +
1584
y * (736.741204151612 + z * (-672.50778314507 +
1585
z * (499.360390819152 + z * (-239.545330654412 +
1586
(48.8012518593872 - 1.66307106208905 * z) * z))) +
1587
y * (-148.185936433658 + z * (397.968445406972 +
1588
z * (-301.815380621876 + (152.196371733841 -
1589
26.3748377232802 * z) * z)) +
1590
y * (58.0259125842571 + z * (-194.618310617595 +
1591
z * (120.520654902025 + z * (-55.2723052340152 +
1592
6.48190668077221 * z))) +
1593
y * (-18.9843846514172 + y * (3.05081646487967 -
1594
9.63108119393062 * z) +
1595
z * (63.5113936641785 + z * (-22.2897317140459 +
1596
8.17060541818112 * z)))))))))
1598
g08_g = x2 * (1416.27648484197 +
1599
x * (-2432.14662381794 + x * (2025.80115603697 +
1600
y * (543.835333000098 + y * (-68.5572509204491 +
1601
y * (49.3667694856254 + y * (-17.1397577419788 +
1602
2.49697009569508 * y))) - 22.6683558512829 * z) +
1603
x * (-1091.66841042967 - 196.028306689776 * y +
1604
x * (374.60123787784 - 48.5891069025409 * x +
1605
36.7571622995805 * y) + 36.0284195611086 * z) +
1606
z * (-54.7919133532887 + (-4.08193978912261 -
1607
30.1755111971161 * z) * z)) +
1608
z * (199.459603073901 + z * (-52.2940909281335 +
1609
(68.0444942726459 - 3.41251932441282 * z) * z)) +
1610
y * (-493.407510141682 + z * (-175.292041186547 +
1611
(83.1923927801819 - 29.483064349429 * z) * z) +
1612
y * (-43.0664675978042 + z * (383.058066002476 +
1613
z * (-54.1917262517112 + 25.6398487389914 * z)) +
1614
y * (-10.0227370861875 - 460.319931801257 * z + y *
1615
(0.875600661808945 + 234.565187611355 * z))))) +
1616
y * (168.072408311545))
1618
g_SA_part = (8645.36753595126 +
1619
x * (-7296.43987145382 + x * (8103.20462414788 +
1620
y * (2175.341332000392 + y * (-274.2290036817964 +
1621
y * (197.4670779425016 + y * (-68.5590309679152 +
1622
9.98788038278032 * y))) - 90.6734234051316 * z) +
1623
x * (-5458.34205214835 - 980.14153344888 * y +
1624
x * (2247.60742726704 - 340.1237483177863 * x +
1625
220.542973797483 * y) + 180.142097805543 * z) +
1626
z * (-219.1676534131548 + (-16.32775915649044 -
1627
120.7020447884644 * z) * z)) +
1628
z * (598.378809221703 + z * (-156.8822727844005 +
1629
(204.1334828179377 - 10.23755797323846 * z) * z)) +
1630
y * (-1480.222530425046 + z * (-525.876123559641 +
1631
(249.57717834054571 - 88.449193048287 * z) * z) +
1632
y * (-129.1994027934126 + z * (1149.174198007428 +
1633
z * (-162.5751787551336 + 76.9195462169742 * z)) +
1634
y * (-30.0682112585625 - 1380.9597954037708 * z + y *
1635
(2.626801985426835 + 703.695562834065 * z))))) +
1636
y * (1187.3715515697959))
1638
chem_potential_water = g03_g + g08_g - 0.5 * sfac * SA * g_SA_part
1640
return np.ma.array(chem_potential_water, mask=mask, copy=False)
1644
def chem_potential_salt_t_exact(SA, t, p):
1645
r"""Calculates the chemical potential of salt in seawater.
1650
Absolute salinity [g kg :sup:`-1`]
1652
in situ temperature [:math:`^\circ` C (ITS-90)]
1658
chem_potential_salt : array_like
1659
chemical potential of salt in seawater [J kg :sup:`-1`]
1672
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
1673
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
1674
>>> p = [10, 50, 125, 250, 600, 1000]
1675
>>> gsw.chem_potential_salt_t_exact(SA, t, p)
1676
array([-8466.13569818, -7928.8256562 , -5029.28859129, -568.42714556,
1677
3396.79366004, 7612.64743154])
1681
.. [1] IOC, SCOR and IAPSO, 2010: The international thermodynamic equation
1682
of seawater - 2010: Calculation and use of thermodynamic properties.
1683
Intergovernmental Oceanographic Commission, Manuals and Guides No. 56,
1684
UNESCO (English), 196 pp. See section 2.9.
1687
2010-03-29. Trevor McDougall and Paul Barker
1690
return (chem_potential_relative_t_exact(SA, t, p) +
1691
chem_potential_water_t_exact(SA, t, p))
1695
def adiabatic_lapse_rate_t_exact(SA, t, p):
1696
r"""Calculates the adiabatic lapse rate of seawater.
1701
Absolute salinity [g kg :sup:`-1`]
1703
in situ temperature [:math:`^\circ` C (ITS-90)]
1709
adiabatic_lapse_rate : array_like
1710
Adiabatic lapse rate [K Pa :sup:`-1`]
1718
The output is in unit of degrees Celsius per Pa, (or equivalently K/Pa) not
1724
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
1725
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
1726
>>> p = [10, 50, 125, 250, 600, 1000]
1727
>>> gsw.adiabatic_lapse_rate_t_exact(SA, t, p)
1728
array([ 2.40350282e-08, 2.38496700e-08, 2.03479880e-08,
1729
1.19586543e-08, 9.96170718e-09, 8.71747270e-09])
1733
.. [1] IOC, SCOR and IAPSO, 2010: The international thermodynamic equation
1734
of seawater - 2010: Calculation and use of thermodynamic properties.
1735
Intergovernmental Oceanographic Commission, Manuals and Guides No. 56,
1736
UNESCO (English), 196 pp. See Eqn. (2.22.1).
1739
2011-03-29. Trevor McDougall and Paul Barker
1742
return -gibbs(n0, n1, n1, SA, t, p) / gibbs(n0, n2, n0, SA, t, p)
1746
def osmotic_coefficient_t_exact(SA, t, p):
1747
r"""Calculates the osmotic coefficient of seawater.
1752
Absolute salinity [g kg :sup:`-1`]
1754
in situ temperature [:math:`^\circ` C (ITS-90)]
1760
osmotic_coefficient : array_like
1761
osmotic coefficient of seawater [unitless]
1774
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
1775
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
1776
>>> p = [10, 50, 125, 250, 600, 1000]
1777
>>> gsw.osmotic_coefficient_t_exact(SA,t , p)
1778
array([ 0.90284718, 0.90298624, 0.90238866, 0.89880927, 0.89801054,
1783
.. [1] IOC, SCOR and IAPSO, 2010: The international thermodynamic equation
1784
of seawater - 2010: Calculation and use of thermodynamic properties.
1785
Intergovernmental Oceanographic Commission, Manuals and Guides No. 56,
1786
UNESCO (English), 196 pp.
1789
2011-04-01. Trevor McDougall and Paul Barker.
1790
2012-11-15. Trevor McDougall and Paul Barker.
1793
SA = np.maximum(SA, 0)
1795
part = k * (1000 - SA) / (Kelvin + t)
1800
# Note that the input pressure (p) is sea pressure in units of dbar.
1803
oc = (7.231916621570606e1, 1.059039593127674e1, -3.025914794694813e1,
1804
5.040733670521486e1, -4.074543321119333e1, 1.864215613820487e1,
1805
-3.022566485046178, -6.138647522851840, 1.353207379758663e1,
1806
-7.316560781114737, 1.829232499785750, -5.358042980767074e-1,
1807
-1.705887283375562, -1.246962174707332e-1, 1.228376913546017,
1808
1.089364009088042e-2, -4.264828939262248e-1, 6.213127679460041e-2,
1809
2.481543497315280, -1.363368964861909, -5.640491627443773e-1,
1810
1.344724779893754, -2.180866793244492, 4.765753255963401,
1811
-5.726993916772165, 2.918303792060746, -6.506082399183509e-1,
1812
-1.015695507663942e-1, 1.035024326471108, -6.742173543702397e-1,
1813
8.465642650849419e-1, -7.508472135244717e-1, -3.668086444057845e-1,
1814
3.189939162107803e-1, -4.245629194309487e-2)
1816
tl = (oc[0] + oc[1] * y + x * (oc[2] + x * (oc[3] + x * (oc[4] + x *
1817
(oc[5] + oc[6] * x))) + y * (oc[7] + x * (oc[8] + x *
1818
(oc[9] + oc[10] * x)) + y * (oc[11] + oc[12] * x + y * (oc[13] +
1819
oc[14] * x + y * (oc[15] + x * (oc[16] + oc[17] * y))))) + z *
1820
(oc[18] + x * (oc[19] + oc[20] * y + oc[21] * x) + y * (oc[22] + y *
1821
(oc[23] + y * (oc[24] + oc[25] * y))) + z * (oc[26] + oc[27] * x + y *
1822
(oc[28] + oc[29] * y) + z * (oc[30] + oc[31] * x + y * (oc[32] +
1823
oc[33] * y) + oc[34] * z)))))
1829
def dynamic_enthalpy_t_exact(SA, t, p):
1830
r"""Calculates the dynamic enthalpy of seawater from Absolute Salinity, in
1831
situ temperature and pressure. Dynamic enthalpy was defined by Young
1832
(2010) as the difference between enthalpy and potential enthalpy. Note that
1833
this function uses the full TEOS-10 Gibbs function (i.e. the sum of the
1834
IAPWS-09 and IAPWS-08 Gibbs functions, see the TEOS-10 Manual, IOC et al.
1840
Absolute salinity [g kg :sup:`-1`]
1842
in situ temperature [:math:`^\circ` C (ITS-90)]
1848
dynamic_enthalpy_t_exact : array_like
1849
dynamic enthalpy [J :sup:`-1`]
1865
.. [1] IOC, SCOR and IAPSO, 2010: The international thermodynamic equation
1866
of seawater - 2010: Calculation and use of thermodynamic properties.
1867
Intergovernmental Oceanographic Commission, Manuals and Guides No. 56,
1868
UNESCO (English), 196 pp.
1870
.. [2] Young, W.R., 2010: Dynamic enthalpy, Conservative Temperature, and
1871
the seawater. Boussinesq approximation. Journal of Physical Oceanography,
1875
2011-04-11. Trevor McDougall and Paul Barker
1878
CT = CT_from_t(SA, t, p)
1880
return enthalpy_t_exact(SA, t, p) - cp0 * CT
1884
def t_maxdensity_exact(SA, p):
1885
r"""Calculates the in-situ temperature of maximum density of seawater.
1886
This function returns the in-situ temperature at which the density of
1887
seawater is a maximum, at given Absolute Salinity, SA, and sea pressure, p
1893
Absolute salinity [g kg :sup:`-1`]
1899
t_maxdensity_exact : array_like
1900
max in-situ temperature [:math:`^\circ` C]
1916
.. [1] IOC, SCOR and IAPSO, 2010: The international thermodynamic equation
1917
of seawater - 2010: Calculation and use of thermodynamic properties.
1918
Intergovernmental Oceanographic Commission, Manuals and Guides No. 56,
1919
UNESCO (English), 196 pp. See section 3.42.
1922
2011-04-03. Trevor McDougall and Paul Barker
1925
# The temperature increment for calculating the gibbs_PTT derivative.
1927
t = 3.978 - 0.22072 * SA # The initial guess of t_maxden.
1928
gibbs_PTT = 1.1e-8 # The initial guess for g_PTT.
1930
for Number_of_iterations in range(0, 3):
1932
gibbs_PT = gibbs(n0, n1, n1, SA, t_old, p)
1933
# Half way through the mod. method (McDougall and Wotherspoon, 2012)
1934
t = t_old - gibbs_PT / gibbs_PTT
1935
t_mean = 0.5 * (t + t_old)
1936
gibbs_PTT = (gibbs(n0, n1, n1, SA, t_mean + dt, p) -
1937
gibbs(n0, n1, n1, SA, t_mean - dt, p)) / (dt + dt)
1938
t = t_old - gibbs_PT / gibbs_PTT
1940
# After three iterations of this modified Newton-Raphson iteration, the
1941
# error in t_maxdensity_exact is typically no larger than 1x10^-15 deg C.
1947
def osmotic_pressure_t_exact(SA, t, pw):
1948
r"""Calculates the osmotic pressure of seawater.
1953
Absolute salinity [g kg :sup:`-1`]
1955
in situ temperature [:math:`^\circ` C (ITS-90)]
1957
sea pressure of the pure water side [dbar]
1961
osmotic_pressure_t_exact : array_like
1962
dynamic osmotic pressure of seawater [dbar]
1978
.. [1] IOC, SCOR and IAPSO, 2010: The international thermodynamic equation
1979
of seawater - 2010: Calculation and use of thermodynamic properties.
1980
Intergovernmental Oceanographic Commission, Manuals and Guides No. 56,
1981
UNESCO (English), 196 pp. See section 3.41.
1984
2011-05-26. Trevor McDougall and Paul Barker
1986
SA = np.maximum(SA, 0)
1987
gibbs_pure_water = gibbs(0, 0, 0, 0, t, pw)
1989
# Initial guess of p, in dbar.
1992
# Initial guess of df/dp.
1993
df_dp = -db2Pascal * (gibbs(n0, n0, n1, SA, t, p) -
1994
SA * gibbs(n1, n0, n1, SA, t, p))
1996
for Number_of_iterations in range(0, 2):
1998
f = gibbs_pure_water - chem_potential_water_t_exact(SA, t, p_old)
1999
# This is half way through the modified N-R method.
2000
p = p_old - f / df_dp
2001
p_mean = 0.5 * (p + p_old)
2002
df_dp = -db2Pascal * (gibbs(0, 0, 1, SA, t, p_mean) -
2003
SA * gibbs(1, 0, 1, SA, t, p_mean))
2004
p = p_old - f / df_dp
2006
# After two iterations though the modified Newton-Raphson technique the
2007
# maximum error is 6x10^-12 dbar.
2009
# Osmotic pressure of seawater in dbar.
2013
if __name__ == '__main__':