Welcome to seawater’s documentation!

Python Seawater

Introduction:

This python package contains a python translation for two MatlabTM Toolboxes.

  1. The CSIRO seawater toolbox (SEAWATER-3.3) for calculating the properties of sea water. Uses the formulas from Unesco’s joint panel on oceanographic tables and standards, UNESCO 1981 and UNESCO 1983 (EOS-80).

The EOS-80 library is considered now obsolete; it is provided here for compatibility with old scripts, and to allow a smooth transition to the new TEOS-10.

  1. The Gibbs Sea Water (GSW v2.0).

A oceanographic toolbox of the International Thermodynamic Equation of Seawater - 2010 or TEOS-10.

Contains the functions for evaluating the thermodynamic properties of pure water (using IAPWS-09) and seawater (using IAPWS-08 for the saline part).

The author has no intention to do things in a “pythonic-way”, it is just a “work around” from someone that couldn’t afford MatlabTM anymore.

GSW toolbox was rewritten in OO approach, i.e., there is a Gibbs class that contains all (SA,t,p) functions as methods.

Ex.:
>>> from seawater.gibbs import Gibbs
>>> SA = [34.5075, 34.7165, 34.8083, 34.8465, 34.8636, 34.8707, 34.8702]
>>> t  = [27.9620, 4.4726, 2.1178, 1.6031, 1.4601, 1.4753, 1.5998]
>>> p  = [0., 1010., 2025., 3045., 4069., 5098., 6131.]
>>> STP = Gibbs(SA, t, p)
>>> STP.beta_const_pt()
For more information:
http://pypi.python.org/pypi/seawater/

Modules:

gibbs Seawater Documentation

This page contains the Csiro Module documentation.

The gibbs module

seawater.gibbs.CT_from_pt(SA, pt)

Calculates Conservative Temperature of seawater from potential temperature (whose reference sea pressure is zero dbar).

Parameters :

SA : array_like

Absolute salinity [g kg -1]

pt : array_like

potential temperature referenced to a sea pressure of zero dbar [^\circ C (ITS-90)]

Returns :

CT : array_like

Conservative Temperature [^\circ C (TEOS-10)]

See also

TODO

Notes

TODO

References

[R1]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. See section 3.3.

Modifications: 2010-08-05. David Jackett, Trevor McDougall and Paul Barker. 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
>>> pt = [28.7832, 28.4209, 22.7850, 10.2305, 6.8292, 4.3245]
>>> gsw.CT_from_pt(SA, pt)
array([ 28.80992302,  28.43914426,  22.78624661,  10.22616561,
         6.82718342,   4.32356518])
class seawater.gibbs.Dict2Struc(adict)

Bases: object

Open variables from a dictionary in a “matlab-like-structure”

class seawater.gibbs.Gibbs(SA=None, t=None, p=None)

Class that aggregate all SA, t, p functions.

Parameters :

SA : array_like

Absolute salinity [g kg -1]

t : array_like

in situ temperature [^\circ C (ITS-90)]

p : array_like

pressure [dbar]

Methods

adiabatic_lapse_rate()

See adiabatic_lapse_rate docstring

alpha_wrt_CT()

See alpha_wrt_CT docstring

alpha_wrt_pt()

See alpha_wrt_pt docstring

alpha_wrt_t()

See alpha_wrt_t docstring

beta_const_CT()

See beta_const_CT docstring

beta_const_pt()

See beta_const_pt docstring

beta_const_t()

See beta_const_t docstring

chem_potential_relative()

See chem_potential_relative docstring

chem_potential_salt()

See chem_potential_salt docstring

chem_potential_water()

See chem_potential_water docstring

conservative_t()

See conservative_t docstring

cp()

See cp docstring

enthalpy()

See enthalpy docstring

entropy()

See entropy docstring

helmholtz_energy()

See helmholtz_energy docstring

internal_energy()

See internal_energy docstring

ionic_strength()

See ionic_strength docstring

isochoric_heat_cap()

See isochoric_heat_cap docstring

kappa()

See kappa docstring

kappa_const_t()

See kappa_const_t docstring

molality()

See molality docstring

osmotic_coefficient()

See osmotic_coefficient docstring

pot_rho(pr=0)

See pot_rho docstring

potential_t(pr=0)

See potential_t docstring

rho()

See rho docstring

sound_speed()

See sound_speed docstring

specvol()

See specvol docstring

specvol_anom()

See specvol_anom docstring

seawater.gibbs.SA_Sstar_from_SP(SP, p, lon, lat)

Calculates Absolute Salinity and Preformed Salinity from Practical Salinity.

Parameters :

SP : array_like

salinity [psu (PSS-78)], unitless

p : array_like

pressure [dbar]

lon : array_like

decimal degrees east [0..+360] or [-180..+180]

lat : array_like

decimal degrees (+ve N, -ve S) [-90..+90]

Returns :

SA : array_like

Absolute salinity [g kg -1]

Sstar : array_like

Preformed Salinity [g kg -1]

in_ocean : False, if [lon, lat] are a long way from the ocean

True, [lon, lat] are in the ocean.

See also

_delta_SA, _SA_from_SP_Baltic

Notes

The in_ocean flag is only set when the observation is well and truly on dry land; often the warning flag is not set until one is several hundred kilometers inland from the coast.

Since SP is non-negative by definition, this function changes any negative input values of SP to be zero.

References

[R2]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. See section 2.5 and appendices A.4 and A.5.
[R3]McDougall, T.J., D.R. Jackett and F.J. Millero, 2010: An algorithm for estimating Absolute Salinity in the global ocean. Submitted to Ocean Science. A preliminary version is available at Ocean Sci. Discuss., 6, 215-242.

http://www.ocean-sci-discuss.net/6/215/2009/osd-6-215-2009-print.pdf

Modifications: 2010-07-23. David Jackett, Trevor McDougall & Paul Barker. 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> SP = [34.5487, 34.7275, 34.8605, 34.6810, 34.5680, 34.5600]
>>> p = [10, 50, 125, 250, 600, 1000]
>>> lon, lat =  188, 4
>>> gsw.SA_Sstar_from_SP(SP, p, lon, lat)
(array([ 34.71177971,  34.89152372,  35.02554774,  34.84723008,
        34.7366296 ,  34.73236186]), array([ 34.7115532 ,  34.89116101,  35.02464926,  34.84359277,
        34.7290336 ,  34.71967638]), array([ True,  True,  True,  True,  True,  True], dtype=bool))
seawater.gibbs.SA_from_SP(SP, p, lon, lat)

Calculates Absolute Salinity from Practical Salinity.

Parameters :

SP : array_like

salinity [psu (PSS-78)], unitless

p : array_like

pressure [dbar]

lon : array_like

decimal degrees east [0..+360] or [-180..+180]

lat : array_like

decimal degrees (+ve N, -ve S) [-90..+90]

Returns :

SA : array_like

Absolute salinity [g kg -1]

in_ocean : False, if [lon, lat] are a long way from the ocean

True, [lon, lat] are in the ocean.

See also

_delta_SA, _SA_from_SP_Baltic

Notes

Since SP is non-negative by definition, this function changes any negative input values of SP to be zero.

References

[R4]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. See section 2.5 and appendices A.4 and A.5.
[R5]McDougall, T.J., D.R. Jackett and F.J. Millero, 2010: An algorithm for estimating Absolute Salinity in the global ocean. Submitted to Ocean Science. A preliminary version is available at Ocean Sci. Discuss., 6, 215-242.

http://www.ocean-sci-discuss.net/6/215/2009/osd-6-215-2009-print.pdf

Modifications: 2010-07-23. David Jackett, Trevor McDougall & Paul Barker. 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> SP = [34.5487, 34.7275, 34.8605, 34.6810, 34.5680, 34.5600]
>>> p = [10, 50, 125, 250, 600, 1000]
>>> lon, lat = 188, 4
>>> gsw.SA_from_SP(SP, p, lon, lat)
(array([ 34.71177971,  34.89152372,  35.02554774,  34.84723008,
        34.7366296 ,  34.73236186]), array([ True,  True,  True,  True,  True,  True], dtype=bool))
seawater.gibbs.SA_from_Sstar(Sstar, p, lon, lat)

Calculates Absolute Salinity from Preformed Salinity.

Parameters :

Sstar : array_like

Preformed Salinity [g kg -1]

p : array_like

pressure [dbar]

lon : array_like

decimal degrees east [0..+360] or [-180..+180]

lat : array_like

decimal degrees (+ve N, -ve S) [-90..+90]

Returns :

SA : array_like

Absolute salinity [g kg -1]

in_ocean : False, if [lon, lat] are a long way from the ocean

True, [lon, lat] are in the ocean.

See also

_delta_SA

Notes

The in_ocean flag is only set when the observation is well and truly on dry land; often the warning flag is not set until one is several hundred kilometers inland from the coast.

References

[R6]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.

Modifications: 2010-07-23. David Jackett, Trevor McDougall and Paul Barker. 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> Sstar = [34.7115, 34.8912, 35.0247, 34.8436, 34.7291, 34.7197]
>>> p = [10, 50, 125, 250, 600, 1000]
>>> lon, lat = 188, 4
>>> gsw.SA_from_Sstar(Sstar, p, lon, lat)
(array([ 34.71172651,  34.89156271,  35.02559848,  34.84723731,
        34.736696  ,  34.73238548]), array([ True,  True,  True,  True,  True,  True], dtype=bool))
seawater.gibbs.SA_from_rho(rho, t, p)

Calculates the Absolute Salinity of a seawater sample, for given values of its density, in situ temperature and sea pressure (in dbar).

One use for this function is in the laboratory where a measured value of the in situ density \rho of a seawater sample may have been made at the laboratory temperature t and at atmospheric pressure p. The present function will return the Absolute Salinity SA of this seawater sample.

Parameters :

rho : array_like

in situ density [kg m -3]

t : array_like

in situ temperature [^\circ C (ITS-90)]

p : array_like

pressure [dbar]

Returns :

SA : array_like

Absolute salinity [g kg -1]

See also

TODO

Notes

This is expressed on the Reference-Composition Salinity Scale of Millero et al. (2008).

After two iterations of a modified Newton-Raphson iteration, the error in SA is typically no larger than 2 ^\times 10 -13 [g kg -1]

References

[R7]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. See section 2.5.
[R8]Millero, F. J., R. Feistel, D. G. Wright, and T. J. McDougall, 2008: The composition of Standard Seawater and the definition of the Reference-Composition Salinity Scale, Deep-Sea Res. I, 55, 50-72.

Modifications: 2010-08-23. Trevor McDougall & Paul Barker. 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> rho = [1021.839, 1022.262, 1024.426, 1027.792, 1029.839, 1032.002]
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
>>> p = [10, 50, 125, 250, 600, 1000]
>>> gsw.SA_from_rho(rho, t, p)
array([ 34.71022966,  34.89057683,  35.02332421,  34.84952096,
        34.73824809,  34.73188384])
seawater.gibbs.SP_from_SA(SA, p, lon, lat)

Calculates Practical Salinity from Absolute Salinity.

Parameters :

SA : array_like

Absolute salinity [g kg -1]

p : array_like

pressure [dbar]

lon : array_like

decimal degrees east [0..+360] or [-180..+180]

lat : array_like

decimal degrees (+ve N, -ve S) [-90..+90]

Returns :

SP : array_like

salinity [psu (PSS-78)], unitless

in_ocean : False, if [lon, lat] are a long way from the ocean

True, [lon, lat] are in the ocean.

See also

_delta_SA, _SP_from_SA_Baltic

Notes

The in_ocean flag is only set when the observation is well and truly on dry land; often the warning flag is not set until one is several hundred kilometers inland from the coast.

References

[R9]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.

Modifications: 2010-07-23. David Jackett, Trevor McDougall and Paul Barker. 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
>>> p = [10, 50, 125, 250, 600, 1000]
>>> lon, lat =  188, 4
>>> gsw.SP_from_SA(SA, p, lon, lat)
(array([ 34.54872019,  34.72747639,  34.86055202,  34.68097006,
        34.56797054,  34.56003796]), array([ True,  True,  True,  True,  True,  True], dtype=bool))
seawater.gibbs.SP_from_Sstar(Sstar, p, lon, lat)

Calculates Practical Salinity from Preformed Salinity.

Parameters :

Sstar : array_like

Preformed Salinity [g kg -1]

p : array_like

pressure [dbar]

lon : array_like

decimal degrees east [0..+360] or [-180..+180]

lat : array_like

decimal degrees (+ve N, -ve S) [-90..+90]

Returns :

SP : array_like

salinity [psu (PSS-78)], unitless

in_ocean : False, if [lon, lat] are a long way from the ocean

True, [lon, lat] are in the ocean.

See also

_delta_SA, _SP_from_SA_Baltic

Notes

The in_ocean flag is only set when the observation is well and truly on dry land; often the warning flag is not set until one is several hundred kilometers inland from the coast.

References

[R10]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. See section 2.5 and appendices A.4 and A.5.

Modifications: 2010-07-23. David Jackett, Trevor McDougall and Paul Barker. 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> Sstar = [34.7115, 34.8912, 35.0247, 34.8436, 34.7291, 34.7197]
>>> p = [10, 50, 125, 250, 600, 1000]
>>> lon, lat =  188, 4
>>> gsw.SP_from_Sstar(Sstar, p, lon, lat)
(array([ 34.54864705,  34.72753881,  34.8605505 ,  34.68100719,
        34.56806609,  34.56002351]), array([ True,  True,  True,  True,  True,  True], dtype=bool))
seawater.gibbs.Sstar_from_SA(SA, p, lon, lat)

Converts Preformed Salinity from Absolute Salinity.

Parameters :

SA : array_like

Absolute salinity [g kg -1]

p : array_like

pressure [dbar]

lon : array_like

decimal degrees east [0..+360] or [-180..+180]

lat : array_like

decimal degrees (+ve N, -ve S) [-90..+90]

Returns :

Sstar : array_like

Preformed Salinity [g kg -1]

in_ocean : False, if [lon, lat] are a long way from the ocean

True, [lon, lat] are in the ocean.

See also

_delta_SA

Notes

The in_ocean flag is only set when the observation is well and truly on dry land; often the warning flag is not set until one is several hundred kilometers inland from the coast.

References

[R11]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. See section 2.5 and appendices A.4 and A.5.

Modifications: 2010-07-23. David Jackett, Trevor McDougall and Paul Barker. 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
>>> p = [10, 50, 125, 250, 600, 1000]
>>> lon, lat =  188, 4
>>> gsw.Sstar_from_SA(SA, p, lon, lat)
(array([ 34.71157349,  34.89113729,  35.02470152,  34.84356269,
        34.729004  ,  34.71971452]), array([ True,  True,  True,  True,  True,  True], dtype=bool))
seawater.gibbs.Sstar_from_SP(SP, p, lon, lat)

Calculates Preformed Salinity from Absolute Salinity.

Parameters :

SP : array_like

salinity [psu (PSS-78)], unitless

p : array_like

pressure [dbar]

lon : array_like

decimal degrees east [0..+360] or [-180..+180]

lat : array_like

decimal degrees (+ve N, -ve S) [-90..+90]

Returns :

Sstar : array_like

Preformed Salinity [g kg -1]

in_ocean : False, if [lon, lat] are a long way from the ocean

True, [lon, lat] are in the ocean.

See also

_delta_SA, _SA_from_SP_Baltic

Notes

The in_ocean flag is only set when the observation is well and truly on dry land; often the warning flag is not set until one is several hundred kilometers inland from the coast.

Since SP is non-negative by definition, this function changes any negative input values of SP to be zero.

References

[R12]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. See section 2.5 and appendices A.4 and A.5.

,, [2] McDougall, T.J., D.R. Jackett and F.J. Millero, 2010: An algorithm for estimating Absolute Salinity in the global ocean. Submitted to Ocean Science. A preliminary version is available at Ocean Sci. Discuss., 6, 215-242.

Modifications: 2010-07-23. David Jackett, Trevor McDougall and Paul Barker. 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> SP = [34.5487, 34.7275, 34.8605, 34.6810, 34.5680, 34.5600]
>>> p = [10, 50, 125, 250, 600, 1000]
>>> lon, lat =  188, 4
>>> gsw.Sstar_from_SP(SP, p, lon, lat)
(array([ 34.7115532 ,  34.89116101,  35.02464926,  34.84359277,
        34.7290336 ,  34.71967638]), array([ True,  True,  True,  True,  True,  True], dtype=bool))
seawater.gibbs.adiabatic_lapse_rate(SA, t, p)

Calculates the adiabatic lapse rate of sea water.

Parameters :

SA : array_like

Absolute salinity [g kg -1]

t : array_like

in situ temperature [^\circ C (ITS-90)]

p : array_like

pressure [dbar]

Returns :

adiabatic_lapse_rate : array_like

Adiabatic lapse rate [K Pa -1]

See also

TODO

Notes

The output is in unit of degrees Celsius per Pa, (or equivalently K/Pa) not in units of K/dbar

References

[R14]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. See Eqn. (2.22.1).

Modifications: 2010-08-26. Trevor McDougall and Paul Barker 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
>>> p = [10, 50, 125, 250, 600, 1000]
>>> gsw.adiabatic_lapse_rate(SA, t, p)
array([  2.40350282e-08,   2.38496700e-08,   2.03479880e-08,
         1.19586543e-08,   9.96170718e-09,   8.71747270e-09])
seawater.gibbs.alpha_wrt_CT(SA, t, p)

Calculates the thermal expansion coefficient of seawater with respect to Conservative Temperature.

Parameters :

SA : array_like

Absolute salinity [g kg -1]

t : array_like

in situ temperature [^\circ C (ITS-90)]

p : array_like

pressure [dbar]

Returns :

alpha_wrt_CT : array_like

thermal expansion coefficient [K -1]

See also

TODO

Notes

TODO

References

[R15]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. See Eqn. (2.18.3).

Modifications: 2010-07-23. David Jackett, Trevor McDougall and Paul Barker 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
>>> p = [10, 50, 125, 250, 600, 1000]
>>> gsw.alpha_wrt_CT(SA, t, p)
array([ 0.00032471,  0.00032272,  0.00028118,  0.00017314,  0.00014627,
        0.00012943])
seawater.gibbs.alpha_wrt_pt(SA, t, p)

Calculates the thermal expansion coefficient of seawater with respect to potential temperature, with a reference pressure of zero.

Parameters :

SA : array_like

Absolute salinity [g kg -1]

t : array_like

in situ temperature [^\circ C (ITS-90)]

p : array_like

pressure [dbar]

Returns :

alpha_wrt_pt : array_like

thermal expansion coefficient [K -1]

See also

TODO

Notes

TODO

References

[R16]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. See Eqn. (2.18.2).

Modifications: 2010-07-23. David Jackett, Trevor McDougall and Paul Barker 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
>>> p = [10, 50, 125, 250, 600, 1000]
>>> gsw.alpha_wrt_pt(SA, t, p)
array([ 0.00032562,  0.00032355,  0.00028164,  0.00017314,  0.00014623,
        0.00012936])
seawater.gibbs.alpha_wrt_t(SA, t, p)

Calculates the thermal expansion coefficient of seawater with respect to in situ temperature.

Parameters :

SA : array_like

Absolute salinity [g kg -1]

t : array_like

in situ temperature [^\circ C (ITS-90)]

p : array_like

pressure [dbar]

Returns :

alpha_wrt_t : array_like

thermal expansion coefficient [K -1]

See also

TODO

Notes

TODO

References

[R17]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. See Eqn. (2.18.1)
[R18]McDougall, T.J., D.R. Jackett and F.J. Millero, 2010: An algorithm for estimating Absolute Salinity in the global ocean. Submitted to Ocean Science. A preliminary version is available at Ocean Sci. Discuss., 6, 215-242.

Modifications: 2010-07-23. David Jackett, Trevor McDougall and Paul Barker 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
>>> p = [10, 50, 125, 250, 600, 1000]
>>> gsw.alpha_wrt_t(SA, t, p)
array([ 0.0003256 ,  0.00032345,  0.00028141,  0.00017283,  0.00014557,
        0.00012836])
seawater.gibbs.beta_const_CT(SA, t, p)

Calculates the saline (i.e. haline) contraction coefficient of seawater at constant Conservative Temperature.

Parameters :

SA : array_like

Absolute salinity [g kg -1]

t : array_like

in situ temperature [^\circ C (ITS-90)]

p : array_like

pressure [dbar]

Returns :

beta_const_CT : array_like

saline contraction coefficient [kg g -1]

See also

TODO

Notes

TODO

References

[R19]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. See Eqn. (2.19.3)

Modifications: 2010-07-23. David Jackett, Trevor McDougall and Paul Barker 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
>>> p = [10, 50, 125, 250, 600, 1000]
>>> gsw.beta_const_CT(SA, t, p)
array([ 0.00071749,  0.00071765,  0.00072622,  0.00075051,  0.00075506,
        0.00075707])
seawater.gibbs.beta_const_pt(SA, t, p)

Calculates the saline (i.e. haline) contraction coefficient of seawater at constant potential temperature with a reference pressure of 0 dbar.

Parameters :

SA : array_like

Absolute salinity [g kg -1]

t : array_like

in situ temperature [^\circ C (ITS-90)]

p : array_like

pressure [dbar]

Returns :

beta_const_pt : array_like

saline contraction coefficient [kg g -1]

See also

TODO

Notes

TODO

References

[R20]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. See Eqn. (2.19.2)

Modifications: 2010-07-23. David Jackett, Trevor McDougall and Paul Barker 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
>>> p = [10, 50, 125, 250, 600, 1000]
>>> gsw.beta_const_pt(SA, t, p)
array([ 0.00073112,  0.00073106,  0.00073599,  0.00075375,  0.00075712,
        0.00075843])
seawater.gibbs.beta_const_t(SA, t, p)

Calculates the saline (i.e. haline) contraction coefficient of seawater at constant in situ temperature.

Parameters :

SA : array_like

Absolute salinity [g kg -1]

t : array_like

in situ temperature [^\circ C (ITS-90)]

p : array_like

pressure [dbar]

Returns :

beta_const_t : array_like

saline contraction coefficient [kg g -1]

See also

TODO

Notes

TODO

References

[R21]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. See Eqn. (2.19.1)

Modifications: 2010-07-23. David Jackett, Trevor McDougall and Paul Barker 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
>>> p = [10, 50, 125, 250, 600, 1000]
>>> gsw.beta_const_t(SA, t, p)
array([ 0.00073112,  0.00073107,  0.00073602,  0.00075381,  0.00075726,
        0.00075865])
seawater.gibbs.chem_potential_relative(SA, t, p)

Calculates the adiabatic lapse rate of sea water.

Parameters :

SA : array_like

Absolute salinity [g kg -1]

t : array_like

in situ temperature [^\circ C (ITS-90)]

p : array_like

pressure [dbar]

Returns :

chem_potential_relative : array_like

relative chemical potential [J kg -1]

See also

TODO

Notes

TODO

References

[R22]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.

Modifications: 2010-08-26. Trevor McDougall and Paul Barker 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
>>> p = [10, 50, 125, 250, 600, 1000]
>>> gsw.chem_potential_relative(SA, t, p)
array([ 79.4254481 ,  79.25989214,  74.69154859,  65.64063719,
        61.22685656,  57.21298557])
seawater.gibbs.chem_potential_salt(SA, t, p)

Calculates the chemical potential of salt in seawater.

Parameters :

SA : array_like

Absolute salinity [g kg -1]

t : array_like

in situ temperature [^\circ C (ITS-90)]

p : array_like

pressure [dbar]

Returns :

chem_potential_salt : array_like

chemical potential of salt in seawater [J kg -1]

See also

TODO

Notes

TODO

References

[R23]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. See section 2.9.

Modifications: 2010-09-28. Trevor McDougall and Paul Barker 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
>>> p = [10, 50, 125, 250, 600, 1000]
>>> gsw.chem_potential_salt(SA, t, p)
array([-8466.13569818, -7928.8256562 , -5029.28859129,  -568.42714556,
        3396.79366004,  7612.64743154])
seawater.gibbs.chem_potential_water(SA, t, p)

Calculates the chemical potential of water in seawater.

Parameters :

SA : array_like

Absolute salinity [g kg -1]

t : array_like

in situ temperature [^\circ C (ITS-90)]

p : array_like

pressure [dbar]

Returns :

chem_potential_water : array_like

chemical potential of water in seawater [J kg -1]

See also

TODO

Notes

TODO

References

[R24]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.

Modifications: 2010-09-28. Trevor McDougall and Paul Barker 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
>>> p = [10, 50, 125, 250, 600, 1000]
>>> gsw.chem_potential_water(SA, t, p)
array([-8545.56114628, -8008.08554834, -5103.98013987,  -634.06778275,
        3335.56680347,  7555.43444597])
seawater.gibbs.conservative_t(SA, t, p)

Calculates Conservative Temperature of seawater from in situ temperature.

Parameters :

SA : array_like

Absolute salinity [g kg -1]

t : array_like

in situ temperature [^\circ C (ITS-90)]

p : array_like

pressure [dbar]

Returns :

CT : array_like

Conservative Temperature [^\circ C (TEOS-10)]

See also

TODO

Notes

TODO

References

[R25]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. See section 3.3.

Modifications: 2010-08-26. David Jackett, Trevor McDougall and Paul Barker 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
>>> p = [10, 50, 125, 250, 600, 1000]
>>> gsw.conservative_t(SA, t, p)
array([ 28.80991983,  28.43922782,  22.78617689,  10.22618927,
         6.82721363,   4.32357575])
seawater.gibbs.cp(SA, t, p)

Calculates the isobaric heat capacity of seawater.

SA : array_like
Absolute salinity [g kg -1]
t : array_like
in situ temperature [^\circ C (ITS-90)]
p : array_like
pressure [dbar]
Returns :

cp : array_like

heat capacity of seawater [J kg -1 K -1]

See also

TODO

Notes

TODO

References

[R26]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.

Modifications: 2010-07-23. David Jackett, Trevor McDougall and Paul Barker 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
>>> p = [10, 50, 125, 250, 600, 1000]
>>> gsw.cp(SA, t, p)
array([ 4002.88800396,  4000.98028393,  3995.54646889,  3985.07676902,
        3973.59384348,  3960.18408479])
seawater.gibbs.datadir

Section A: Library functions

seawater.gibbs.distance(lon, lat, p=None)

Calculates the distance in met res between successive points in the vectors lon and lat, computed using the Haversine formula on a spherical earth of radius 6,371 km, being the radius of a sphere having the same volume as Earth. For a spherical Earth of radius 6,371,000 m, one nautical mile is 1,853.2488 m, thus one degree of latitude is 111,194.93 m.

Haversine formula:
R = earth’s radius (mean radius = 6,371 km)

a = \sin^2(\delta \text{lat}/2) + \cos(\text{lat}_1) \cos(\text{lat}_2) \sin^2(\delta \text{lon}/2)

c = 2 \times \text{atan2}(\sqrt{a}, \sqrt{(1-a)})

d = R \times c

Parameters :

lon : array_like

decimal degrees east [0..+360] or [-180 ... +180]

lat : array_like

latitude in decimal degrees north [-90..+90]

p : number or array_like. Default p = 0

pressure [dbar]

Returns :

dist: array_like :

distance between points on a spherical Earth at pressure (p) [m]

See also

TODO

Notes

Distances are probably good to better than 1% of the “true” distance on the ellipsoidal earth. The check value below differ from the original online docs at “http://www.teos-10.org/pubs/gsw/html/gsw_distance.html” but agree with the result.

References

[R27]http://www.eos.ubc.ca/~rich/map.html

Modifications: 2000-11-06. Rich Pawlowicz 2010-07-28. Paul Barker and Trevor McDougall 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> lon = [159, 220]
>>> lat = [-35, 35]
>>> gsw.distance(lon, lat)
array([[ 10030974.652916]])
>>> p = [200, 1000]
>>> gsw.distance(lon, lat, p)
array([[ 10030661.63878009]])
seawater.gibbs.enthalpy(SA, t, p)

Calculates the specific enthalpy of seawater.

The specific enthalpy of seawater h is given by:

h(SA, t, p) = g + (T_0 + t)\eta = g - (T_0 + t) \frac{\partial g}{\partial T}\Big|_{SA,p}

Parameters :

SA : array_like

Absolute salinity [g kg -1]

t : array_like

in situ temperature [^\circ C (ITS-90)]

p : array_like

pressure [dbar]

Returns :

enthalpy : array_like

specific enthalpy [J kg -1]

See also

TODO

Notes

TODO

References

[R28]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. See appendix A.11.

Modifications: 2010-08-26. David Jackett, Trevor McDougall and Paul Barker. 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
>>> p = [10, 50, 125, 250, 600, 1000]
>>> gsw.enthalpy(SA, t, p)
array([ 115103.26047838,  114014.8036012 ,   92179.9209311 ,
         43255.32838089,   33087.21597002,   26970.5880448 ])
seawater.gibbs.entropy(SA, t, p)

Calculates specific entropy of seawater.

The specific entropy of seawater \eta is given by:

\eta(SA, t, p) = -g_T = \frac{\partial g}{\partial T}\Big|_{SA,p}

When taking derivatives with respect to in situ temperature, the symbol T will be used for temperature in order that these derivatives not be confused with time derivatives.

Parameters :

SA : array_like

Absolute salinity [g kg -1]

t : array_like

in situ temperature [^\circ C (ITS-90)]

p : array_like

pressure [dbar]

Returns :

entropy : array_like

specific entropy [J kg -1 K -1]

See also

TODO

Notes

TODO

References

[R29]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.

Modifications: 2010-07-23. David Jackett, Trevor McDougall and Paul Barker. 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
>>> p = [10, 50, 125, 250, 600, 1000]
>>> gsw.entropy(SA, t, p)
array([ 400.38942528,  395.43817843,  319.8664982 ,  146.79088159,
         98.64734087,   62.79150873])
seawater.gibbs.entropy_from_t(SA, t, t_type='pt')

Calculates specific entropy of seawater.

Parameters :

SA : array_like

Absolute salinity [g kg -1]

t : array_like

temperature [^\circ C]

t_type : str, optional

‘pt’ for potential temperature [^\circ C (ITS-90)] ‘CT’ for Conservative Temperature [^\circ C (TEOS-10)]

Returns :

entropy : array_like

specific entropy [J kg -1 K -1]

See also

TODO

Notes

TODO

References

[R30]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. See appendix A.10.

Modifications: 2010-10-13. Trevor McDougall & Paul Barker 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
>>> pt = [28.7832, 28.4210, 22.7850, 10.2305, 6.8292, 4.3245]
>>> gsw.entropy_from_t(SA, pt)
array([ 400.38946744,  395.43839949,  319.86743859,  146.79054828,
         98.64691006,   62.79135672])
>>> CT = [28.8099, 28.4392, 22.7862, 10.2262, 6.8272, 4.3236]
>>> gsw.entropy_from_t(SA, CT, 'CT')
array([ 400.38916315,  395.43781023,  319.86680989,  146.79103279,
         98.64714648,   62.79185763])
seawater.gibbs.grav(lat, p=0)

Calculates acceleration due to gravity as a function of latitude and as a function of pressure in the ocean.

Parameters :

lat : array_like

latitude in decimal degrees north [-90..+90]

p : number or array_like. Default p = 0

pressure [dbar]

Returns :

g : array_like

gravity [m s 2]

See also

TODO

Notes

In the ocean z is negative.

References

[R31]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.
[R32]Moritz (2000) Goedetic reference system 1980. J. Geodesy, 74, 128-133.
[R33]Saunders, P.M., and N.P. Fofonoff (1976) Conversion of pressure to depth in the ocean. Deep-Sea Res.,pp. 109 - 111.

Modifications: 2010-07-23. Trevor McDougall & Paul Barker 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> lat = [-90, -60, -30, 0]
>>> p = 0
>>> gsw.grav(lat, p)
array([ 9.83218621,  9.81917886,  9.79324926,  9.780327  ])
>>> gsw.grav(45)
9.8061998770458008
seawater.gibbs.helmholtz_energy(SA, t, p)

Calculates the Helmholtz energy of seawater.

The specific Helmholtz energy of seawater f is given by:

f(SA, t, p) = g - (p + P_0) \nu = g - (p + P_0) \frac{\partial g}{\partial P}\Big|_{SA,T}

Parameters :

SA : array_like

Absolute salinity [g kg -1]

t : array_like

in situ temperature [^\circ C (ITS-90)]

p : array_like

pressure [dbar]

Returns :

Helmholtz_energy : array_like

Helmholtz energy [J kg -1]

See also

TODO

Notes

TODO

References

[R34]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. See section 2.13.

Modifications: 2010-08-26. Trevor McDougall 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
>>> p = [10, 50, 125, 250, 600, 1000]
>>> gsw.helmholtz_energy(SA, t, p)
array([-5985.58288209, -5830.81845224, -3806.96617841,  -877.66369421,
        -462.17033905,  -245.50407205])
seawater.gibbs.internal_energy(SA, t, p)

Calculates the Helmholtz energy of seawater.

The specific internal energy of seawater u is given by:

u(SA, t, p) = g + (T_0 + t)\eta - (p + P_0)\nu = g - (T_0 + t)\frac{\partial g}{\partial T}\Big|_{SA,p} - (p + P_0)\frac{\partial g}{\partial P}\Big|_{SA,T}

where T_0 is the Celsius zero point, 273.15 K and P_0 = 101 325 Pa is the standard atmosphere pressure.

Parameters :

SA : array_like

Absolute salinity [g kg -1]

t : array_like

in situ temperature [^\circ C (ITS-90)]

p : array_like

pressure [dbar]

Returns :

internal_energy (u) : array_like

specific internal energy [J kg -1]

See also

TODO

Notes

TODO

References

[R35]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. See Eqn. (2.11.1)

Modifications: 2010-08-22. Trevor McDougall 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
>>> p = [10, 50, 125, 250, 600, 1000]
>>> gsw.internal_energy(SA, t, p)
array([ 114906.23847309,  113426.57417062,   90860.81858842,
         40724.34005719,   27162.66600185,   17182.50522667])
seawater.gibbs.ionic_strength(SA)

Calculates the ionic strength of seawater.

Parameters :

SA : array_like

Absolute salinity [g kg -1]

Returns :

ionic_strength : array_like

ionic strength of seawater [mol kg -1]

See also

TODO

Notes

TODO

References

[R36]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. See Table L.1.
[R37]Millero, F. J., R. Feistel, D. G. Wright, and T. J. McDougall, 2008: The composition of Standard Seawater and the definition of the Reference-Composition Salinity Scale, Deep-Sea Res. I, 55, 50-72. see Eqns. 5.9 and 5.12.

Modifications: 2010-09-28. Trevor McDougall & Paul Barker 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
>>> gsw.ionic_strength(SA)
array([ 0.71298118,  0.71680567,  0.71966059,  0.71586272,  0.71350891,
        0.71341953])
seawater.gibbs.isochoric_heat_cap(SA, t, p)

Calculates the isochoric heat capacity of seawater.

Parameters :

SA : array_like

Absolute salinity [g kg -1]

t : array_like

in situ temperature [^\circ C (ITS-90)]

p : array_like

pressure [dbar]

Returns :

isochoric_heat_cap : array_like

isochoric heat capacity [J kg -1 K -1]

See also

TODO

Notes

TODO

References

[R38]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. See section 2.21.

Modifications: 2010-08-26. Trevor McDougall 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
>>> p = [10, 50, 125, 250, 600, 1000]
>>> gsw.isochoric_heat_cap(SA, t, p)
array([ 3928.13708702,  3927.27381633,  3941.36418525,  3966.26126146,
        3960.50903222,  3950.13901342])
seawater.gibbs.kappa(SA, t, p)

Calculates the isentropic compressibility of seawater.

When the entropy and Absolute Salinity are held constant while the pressure is changed, the isentropic and isohaline compressibility kappa is obtained:

\kappa(SA, t, p) = \rho^{-1}\frac{\partial \rho}{\partial P}\Big|_{SA,\eta} = -\nu^{-1}\frac{\partial \nu}{\partial P}\Big|_{SA,\eta} = \rho^{-1}\frac{\partial \rho}{\partial P}\Big|_{SA,\theta} = -\nu^{-1}\frac{\partial \nu}{\partial P}\Big|_{SA,\theta} =
-\frac{ (g_{TP}^2 - g_{TT} g_{PP} ) }{g_P g_{TT}}

The isentropic and isohaline compressibility is sometimes called simply the isentropic compressibility (or sometimes the “adiabatic compressibility”), on the unstated understanding that there is also no transfer of salt during the isentropic or adiabatic change in pressure.

Parameters :

SA : array_like

Absolute salinity [g kg -1]

t : array_like

in situ temperature [^\circ C (ITS-90)]

p : array_like

pressure [dbar]

Returns :

kappa : array_like

Isentropic compressibility [Pa -1]

See also

TODO

Notes

The output is Pascal and not dbar.

References

[R39]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. See Eqns. (2.16.1) and the row for kappa in Table P.1 of appendix P

Modifications: 2010-07-23. David Jackett, Trevor McDougall and Paul Barker 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
>>> p = [10, 50, 125, 250, 600, 1000]
>>> gsw.kappa(SA, t, p)
array([  4.11245799e-10,   4.11029072e-10,   4.16539558e-10,
         4.35668338e-10,   4.38923693e-10,   4.40037576e-10])
seawater.gibbs.kappa_const_t(SA, t, p)

Calculates isothermal compressibility of seawater at constant in situ temperature.

\kappa^t(SA, t, p) = \rho^{-1}\frac{\partial \rho}{\partial P}\Big|_{SA,T} = -\nu^{-1}\frac{\partial \nu}{\partial P}\Big|_{SA,T} = -\frac{g_{PP}}{g_P}

Parameters :

SA : array_like

Absolute salinity [g kg -1]

t : array_like

in situ temperature [^\circ C (ITS-90)]

p : array_like

pressure [dbar]

Returns :

kappa : array_like

Isothermal compressibility [Pa -1]

See also

TODO

Notes

This is the compressibility of seawater at constant in situ temperature.

References

[R40]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. See Eqn. (2.15.1)

Modifications: 2010-07-23. David Jackett, Trevor McDougall and Paul Barker 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
>>> p = [10, 50, 125, 250, 600, 1000]
>>> gsw.kappa_const_t(SA, t, p)
array([  4.19071646e-10,   4.18743202e-10,   4.22265764e-10,
         4.37735100e-10,   4.40373818e-10,   4.41156577e-10])
seawater.gibbs.molality(SA)

Calculates the molality of seawater.

Parameters :

SA : array_like

Absolute salinity [g kg -1]

Returns :

molal : array_like

seawater molality [mol kg -1]

See also

TODO

Notes

TODO

References

[R41]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.

Modifications: 2010-09-28. Trevor McDougall & Paul Barker 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
>>> gsw.molality(SA)
array([ 1.14508476,  1.15122708,  1.15581223,  1.14971265,  1.14593231,
        1.14578877])
seawater.gibbs.osmotic_coefficient(SA, t, p)

Calculates the osmotic coefficient of seawater.

Parameters :

SA : array_like

Absolute salinity [g kg -1]

t : array_like

in situ temperature [^\circ C (ITS-90)]

p : array_like

pressure [dbar]

Returns :

osmotic_coefficient : array_like

osmotic coefficient of seawater [unitless]

See also

TODO

Notes

TODO

References

[R42]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.

Modifications: 2010-09-28. Trevor McDougall and Paul Barker 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
>>> p = [10, 50, 125, 250, 600, 1000]
>>> gsw.osmotic_coefficient(SA,t , p)
array([ 0.90284718,  0.90298624,  0.90238866,  0.89880927,  0.89801054,
        0.89767912])
seawater.gibbs.p_from_z(z, lat)

Calculates sea pressure from height using computationally-efficient 25-term expression for density, in terms of SA, CT and p.

Parameters :

lat : array_like

latitude in decimal degrees north [-90..+90]

z : array_like

height [m]

Returns :

p : array_like

pressure [dbar]

See also

_specvol_SSO_0_CT25, _enthalpy_SSO_0_CT25

Notes

Height (z) is NEGATIVE in the ocean. Depth is -z. Depth is not used in the gibbs library.

References

[R43]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.

,, [2] McDougall T. J., D. R. Jackett, P. M. Barker, C. Roberts-Thomson, R. Feistel and R. W. Hallberg, 2010: A computationally efficient 25-term expression for the density of seawater in terms of Conservative Temperature, and related properties of seawater.

[R45]Moritz (2000) Goedetic reference system 1980. J. Geodesy, 74, 128-133.
[R46]Saunders, P. M., 1981: Practical conversion of pressure to depth. Journal of Physical Oceanography, 11, 573-574.

Modifications: 2010-08-26. Trevor McDougall, Claire Roberts-Thomson and Paul Barker. 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> z = [10, 50, 125, 250, 600, 1000]
>>> lat = 4.
>>> gsw.p_from_z(z, lat)
array([  -10.05521794,   -50.2711751 ,  -125.6548857 ,  -251.23284504,
        -602.44050752, -1003.07609807])
>>> z = [-9.94460074, -49.71817465, -124.2728275, -248.47044828, -595.82618014, -992.0931748]
>>> gsw.p_from_z(z, lat)
array([   10.,    50.,   125.,   250.,   600.,  1000.])
seawater.gibbs.pot_rho(SA, t, p, pr=0)

Calculates potential density of seawater.

Parameters :

SA : array_like

Absolute salinity [g kg -1]

t : array_like

in situ temperature [^\circ C (ITS-90)]

p : array_like

pressure [dbar]

pr : int, float, optional

reference pressure, default = 0

Returns :

pot_rho : array_like

potential density [kg m -3]

See also

TODO

Notes

TODO

References

[R47]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. See section 3.4.

Modifications: 2010-07-23. David Jackett, Trevor McDougall and Paul Barker 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
>>> p = [10, 50, 125, 250, 600, 1000]
>>> gsw.pot_rho(SA, t, p)
array([ 1021.79814581,  1022.05248442,  1023.89358365,  1026.66762112,
        1027.10723087,  1027.40963126])
>>> gsw.pot_rho(SA, t, p, pr=1000)
array([ 1025.95554512,  1026.21306986,  1028.12563226,  1031.1204547 ,
        1031.63768355,  1032.00240412])
seawater.gibbs.potential_t(SA, t, p, pr=0)

Calculates potential temperature with the general reference pressure, pr, from in situ temperature.

Parameters :

SA : array_like

Absolute salinity [g kg -1]

t : array_like

in situ temperature [^\circ C (ITS-90)]

p : array_like

pressure [dbar]

pr : int, float, optional

reference pressure, default = 0

Returns :

pt : array_like

potential temperature [^\circ C (ITS-90)]

See also

_entropy_part

Notes

This function calls “entropy_part” which evaluates entropy except for the parts which are a function of Absolute Salinity alone. A faster routine exists pt0_from_t(SA,t,p) if pr is indeed zero dbar.

References

[R48]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. See section 3.1.
[R49]McDougall T. J., D. R. Jackett, P. M. Barker, C. Roberts-Thomson, R. Feistel and R. W. Hallberg, 2010: A computationally efficient 25-term expression for the density of seawater in terms of Conservative Temperature, and related properties of seawater. To be submitted to Ocean Science Discussions.

Modifications: 2010-08-26. Trevor McDougall, David Jackett, Claire Roberts-Thomson and Paul Barker. 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
>>> p = [10, 50, 125, 250, 600, 1000]
>>> gsw.potential_t(SA, t, p)
array([ 28.78319682,  28.42098334,  22.7849304 ,  10.23052366,
         6.82923022,   4.32451057])
>>> gsw.potential_t(SA, t, p, pr = 1000)
array([ 29.02665528,  28.662375  ,  22.99149634,  10.35341725,
         6.92732954,   4.4036    ])
seawater.gibbs.pt0_from_t(SA, t, p)

Calculates potential temperature with reference pressure, pr = 0 dbar. The present routine is computationally faster than the more general function “potential_t(SA, t, p, pr)” which can be used for any reference pressure value.

Parameters :

SA : array_like

Absolute salinity [g kg -1]

t : array_like

in situ temperature [^\circ C (ITS-90)]

p : array_like

pressure [dbar]

Returns :

pt0 : array_like

potential temperature relative to 0 dbar [^\circ C (ITS-90)]

See also

_entropy_part, _gibbs_pt0_pt0, _entropy_part_zerop

Notes

TODO

References

[R50]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. See section 3.1.
[R51]McDougall T. J., D. R. Jackett, P. M. Barker, C. Roberts-Thomson, R. Feistel and R. W. Hallberg, 2010: A computationally efficient 25-term expression for the density of seawater in terms of Conservative Temperature, and related properties of seawater. To be submitted to Ocean Science Discussions.

Modifications: 2010-08-26. Trevor McDougall, David Jackett, Claire Roberts-Thomson and Paul Barker. 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
>>> p = [10, 50, 125, 250, 600, 1000]
>>> gsw.pt0_from_t(SA, t, p)
array([ 28.78319682,  28.42098334,  22.7849304 ,  10.23052366,
         6.82923022,   4.32451057])
seawater.gibbs.pt_from_CT(SA, CT)

Calculates potential temperature (with a reference sea pressure of zero dbar) from Conservative Temperature.

Parameters :

SA : array_like

Absolute salinity [g kg -1]

CT : array_like

Conservative Temperature [^\circ C (TEOS-10)]

Returns :

pt : array_like

potential temperature referenced to a sea pressure of zero dbar [^\circ C (ITS-90)]

See also

TODO

Notes

This function uses 1.5 iterations through a modified Newton-Raphson (N-R) iterative solution procedure, starting from a rational-function-based initial condition for both pt and dCT_dpt.

References

[R52]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. See sections 3.1 and 3.3.
[R53]McDougall T. J., D. R. Jackett, P. M. Barker, C. Roberts-Thomson, R. Feistel and R. W. Hallberg, 2010: A computationally efficient 25-term expression for the density of seawater in terms of Conservative Temperature, and related properties of seawater. To be submitted to Ocean Science Discussions.

Modifications: 2010-08-26. Trevor McDougall, David Jackett, Claire Roberts-Thomson and Paul Barker. 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> 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]
>>> gsw.pt_from_CT(SA, CT)
array([ 28.78317705,  28.4209556 ,  22.78495347,  10.23053439,
         6.82921659,   4.32453484])
seawater.gibbs.rho(SA, t, p)

Calculates in situ density of seawater from Absolute Salinity and in situ temperature.

Parameters :

SA : array_like

Absolute salinity [g kg -1]

t : array_like

in situ temperature [^\circ C (ITS-90)]

p : array_like

pressure [dbar]

Returns :

rho : array_like

in situ density [kg m -3]

See also

TODO

Notes

TODO

References

[R54]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. See section 2.8.

Modifications: 2010-07-23. David Jackett, Trevor McDougall & Paul Barker 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
>>> p = [10, 50, 125, 250, 600, 1000]
>>> gsw.rho(SA, t, p)
array([ 1021.84017319,  1022.26268993,  1024.42771594,  1027.79020181,
        1029.83771473,  1032.00240412])
seawater.gibbs.sound_speed(SA, t, p)

Calculates the speed of sound in seawater.

The speed of sound in seawater c is given by:

c(SA, t, p) = \sqrt{ \partial P  / \partial \rho |_{SA,\eta}} = \sqrt{(\rho\kappa)^{-1}} = g_P \sqrt{g_{TT}/(g^2_{TP} - g_{TT}g_{PP})}

Note that in these expressions, since sound speed is in m s :sup`-1` and density has units of kg m -3 it follows that the pressure of the partial derivatives must be in Pa and the isentropic compressibility kappa must have units of Pa -1. The sound speed c produced by both the SIA and the GSW software libraries (appendices M and N) has units of m s -1.

Parameters :

SA : array_like

Absolute salinity [g kg -1]

t : array_like

in situ temperature [^\circ C (ITS-90)]

p : array_like

pressure [dbar]

Returns :

sound_speed : array_like

speed of sound in seawater [m s -1]

See also

TODO

Notes

TODO

References

[R55]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. See Eqn. (2.17.1)

Modifications: 2010-07-23. David Jackett, Paul Barker and Trevor McDougall. 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
>>> p = [10, 50, 125, 250, 600, 1000]
>>> gsw.sound_speed(SA, t, p)
array([ 1542.61580359,  1542.70353407,  1530.84497914,  1494.40999692,
        1487.37710252,  1483.93460908])
seawater.gibbs.specvol(SA, t, p)

Calculates the specific volume of seawater.

Parameters :

SA : array_like

Absolute salinity [g kg -1]

t : array_like

in situ temperature [^\circ C (ITS-90)]

p : array_like

pressure [dbar]

Returns :

specvol : array_like

specific volume [m 3 kg -1]

See also

TODO

Notes

TODO

References

[R56]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. See section 2.7.

Modifications: 2010-08-26. David Jackett & Paul Barker. 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
>>> p = [10, 50, 125, 250, 600, 1000]
>>> gsw.specvol(SA, t, p)
array([ 0.00097863,  0.00097822,  0.00097615,  0.00097296,  0.00097103,
        0.00096899])
seawater.gibbs.specvol_anom(SA, t, p)

Calculates specific volume anomaly from Absolute Salinity, in situ temperature and pressure, using the full TEOS-10 Gibbs function.

The reference value of Absolute Salinity is SSO and the reference value of Conservative Temperature is equal to 0 ^\circ C.

Parameters :

SA : array_like

Absolute salinity [g kg -1]

t : array_like

in situ temperature [^\circ C (ITS-90)]

p : array_like

pressure [dbar]

Returns :

specvol_anom : array_like

specific volume anomaly [m 3 kg -1]

See also

TODO

Notes

TODO

References

[R57]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. See Eqn. (3.7.3)

Modifications: 2010-08-26. Trevor McDougall & Paul Barker 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
>>> t = [28.7856, 28.4329, 22.8103, 10.2600, 6.8863, 4.4036]
>>> p = [10, 50, 125, 250, 600, 1000]
>>> gsw.specvol_anom(SA, t, p)
array([  6.01044463e-06,   5.78602432e-06,   4.05564999e-06,
         1.42198662e-06,   1.04351837e-06,   7.63964850e-07])
seawater.gibbs.t_from_CT(SA, CT, p)

Calculates in situ temperature from Conservative Temperature of seawater.

Parameters :

SA : array_like

Absolute salinity [g kg -1]

CT : array_like

Conservative Temperature [^\circ C (TEOS-10)]

p : array_like

pressure [dbar]

Returns :

t : array_like

in situ temperature [^\circ C (ITS-90)]

See also

TODO

Notes

TODO

References

[R58]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. See sections 3.1 and 3.3.

Modifications: 2010-08-26. Trevor McDougall & Paul Barker 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> 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]
>>> gsw.t_from_CT(SA, CT, p)
array([ 28.78558023,  28.43287225,  22.81032309,  10.26001075,
         6.8862863 ,   4.40362445])
seawater.gibbs.t_from_entropy(SA, entropy, t_type='pt')

Calculates potential temperature with reference pressure pr = 0 dbar or Conservative temperature from entropy.

Parameters :

SA : array_like

Absolute salinity [g kg -1]

entropy : array_like

specific entropy [J kg -1 K -1]

t_type : str, optional

‘pt’ for potential temperature [^\circ C (ITS-90)] ‘CT’ for Conservative Temperature [^\circ C (TEOS-10)]

Returns :

t : array_like

potential temperature [^\circ C (ITS-90)] (Default) or Conservative Temperature [^\circ C (TEOS-10)]

See also

_gibbs_pt0_pt0

Notes

TODO

References

[R59]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. See appendix A.10.

Modifications: 2010-10-13. Trevor McDougall and Paul Barker. 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
>>> entropy = [400.3892, 395.4378, 319.8668, 146.7910, 98.6471, 62.7919]
>>> gsw.t_from_entropy(SA, entropy)
array([ 28.78317983,  28.42095483,  22.78495274,  10.23053207,
         6.82921333,   4.32453778])
>>> gsw.t_from_entropy(SA, entropy, t_type='CT')
array([ 28.80990279,  28.43919923,  22.78619927,  10.22619767,
         6.82719674,   4.32360295])
seawater.gibbs.z_from_p(p, lat)

Calculates height from sea pressure using the computationally-efficient 25-term expression for density in terms of SA, CT and p.

Parameters :

lat : array_like

latitude in decimal degrees north [-90..+90]

p : array_like

pressure [dbar]

Returns :

z : array_like

height [m]

See also

_enthalpy_SSO_0_CT25

Notes

At sea level z = 0, and since z (HEIGHT) is defined to be positive upwards, it follows that while z is positive in the atmosphere, it is NEGATIVE in the ocean.

References

[R60]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.

,, [2] McDougall T. J., D. R. Jackett, P. M. Barker, C. Roberts-Thomson, R. Feistel and R. W. Hallberg, 2010: A computationally efficient 25-term expression for the density of seawater in terms of Conservative Temperature, and related properties of seawater.

[R62]Moritz (2000) Goedetic reference system 1980. J. Geodesy, 74, 128-133.

Modifications: 2010-08-26. Trevor McDougall, Claire Roberts-Thomson & Paul Barker. 2010-12-09. Filipe Fernandes, Python translation from gsw toolbox.

Examples

>>> import seawater.gibbs as gsw
>>> p = [10, 50, 125, 250, 600, 1000]
>>> lat = 4
>>> gsw.z_from_p(p, lat)
array([  -9.94460074,  -49.71817465, -124.2728275 , -248.47044828,
       -595.82618014, -992.0931748 ])

csiro Seawater Documentation

This page contains the Csiro Module documentation.

The csiro module

seawater.csiro.T68conv(T90)

Convert ITS-90 temperature to IPTS-68

T68  = T90 * 1.00024

Parameters :

t : array_like

temperature [^\circ C (ITS-90)]

Returns :

t : array_like

temperature [^\circ C (IPTS-68)]

See also

TODO

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 ^\circ C for conversion between IPTS-68 and ITS-90 over the oceanographic temperature range.

References

[R66]Saunders, P. M., 1991: The International Temperature Scale of 1990, ITS-90. WOCE Newsletter, No. 10, WOCE International Project Office, Southampton, United Kingdom, 10.
Modifications: Filipe Fernandes, 2010
10-11-24. Filipe Fernandes, first version.

Examples

>>> import seawater.csiro as sw
>>> sw.T68conv(19.995201151723585)
20.0
seawater.csiro.T90conv(t, t_type='T68')

Convert IPTS-68 or IPTS-48 to temperature to ITS-90

..math:

T90 = T68 / 1.00024

T90 = T48 - (4.4e-6) * T48 * (100-T48) ) / 1.00024

Parameters :

t : array_like

temperature [^\circ C (IPTS-68) or (IPTS-48)]

t_type : string, optional

‘T68’ (default) or ‘T48’

Returns :

T90 : array_like

temperature [^\circ C (ITS-90)]

See also

TODO

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 ^\circ C for conversion between IPTS-68 and ITS-90 over the oceanographic temperature range.

References

[R67]Saunders, P. M., 1991: The International Temperature Scale of 1990, ITS-90. WOCE Newsletter, No. 10, WOCE International Project Office, Southampton, United Kingdom, 10.
[R68]International Temperature Scales of 1948, 1968 and 1990, an ICES note, available from http://www.ices.dk/ocean/procedures/its.htm
Modifications: Filipe Fernandes, 2010
2010-11-24. Filipe Fernandes, first version. 2010-12-24. Filipe Fernandes, added T48.

Examples

>>> import seawater.csiro as sw
>>> sw.T90conv(20.004799999999999)
20.0
>>> sw.T90conv(20., t_type='T48')
19.988162840918179
seawater.csiro.adtg(s, t, p)

Calculates adiabatic temperature gradient as per UNESCO 1983 routines.

Parameters :

s(p) : array_like

salinity [psu (PSS-78)]

t(p) : array_like

temperature [^\circ C (ITS-90)]

p : array_like

pressure [db]. The shape can be “broadcasted”

Returns :

adtg : array_like

adiabatic temperature gradient [^\circ C db -1]

See also

ptmp

Notes

TODO: Pressure broadcast feature need to be tested.

References

[R69]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. http://unesdoc.unesco.org/images/0005/000598/059832eb.pdf
[R70]Bryden, H. 1973. New Polynomials for thermal expansion, adiabatic temperature gradient and potential temperature of sea water. Deep-Sea Res. Vol20,401-408. doi:10.1016/0011-7471(73)90063-6
Modifications: Phil Morgan, Lindsay Pender (Lindsay.Pender@csiro.au)
99-06-25. Lindsay Pender, Fixed transpose of row vectors. 03-12-12. Lindsay Pender, Converted to ITS-90. 10-01-14. Filipe Fernandes, Python translation. 10-08-16. Filipe Fernandes, Reformulated docstring.

Examples

Data from UNESCO 1983 p45

>>> import seawater.csiro as sw
>>> t = sw.T90conv([[ 0,  0,  0,  0,  0,  0], [10, 10, 10, 10, 10, 10], [20, 20, 20, 20, 20, 20], [30, 30, 30, 30, 30, 30], [40, 40, 40, 40, 40, 40]])
>>> s = [[25, 25, 25, 35, 35, 35], [25, 25, 25, 35, 35, 35], [25, 25, 25, 35, 35, 35], [25, 25, 25, 35, 35, 35], [25, 25, 25, 35, 35, 35]]
>>> p = [0, 5000, 10000, 0, 5000, 10000]
>>> sw.adtg(s, t, p)
array([[  1.68710000e-05,   1.04700000e-04,   1.69426000e-04,
          3.58030000e-05,   1.17956500e-04,   1.77007000e-04],
       [  1.00194580e-04,   1.60959050e-04,   2.06874170e-04,
          1.14887280e-04,   1.71364200e-04,   2.12991770e-04],
       [  1.73819840e-04,   2.13534000e-04,   2.44483760e-04,
          1.84273240e-04,   2.21087800e-04,   2.49137960e-04],
       [  2.41720460e-04,   2.64764100e-04,   2.82959590e-04,
          2.47934560e-04,   2.69466550e-04,   2.86150390e-04],
       [  3.07870120e-04,   3.16988600e-04,   3.23006480e-04,
          3.09844920e-04,   3.18839700e-04,   3.24733880e-04]])
seawater.csiro.alpha(s, t, p, pt=False)

Calculate the thermal expansion coefficient.

Parameters :

s(p) : array_like

salinity [psu (PSS-78)]

t(p) : array_like

temperature or potential temperature [^\circ C (ITS-90)]

p : array_like

pressure [db]. The shape can be “broadcasted”

pt : bool

True if temperature is potential, default is False

Returns :

alpha : array_like

thermal expansion coeff (\alpha) [^\circ C -1]

See also

beta, aonb

Notes

TODO: Pressure broadcast feature need to be tested.

References

[R71]McDougall, Trevor J., 1987: Neutral Surfaces. J. Phys. Oceanogr., 17, 1950-1964 doi: 10.1175/1520-0485(1987)017<1950:NS>2.0.CO;2
Modifications: N.L. Bindoff 1993, Lindsay Pender (Lindsay.Pender@csiro.au)
93-04-22. Phil Morgan, Help display modified to suit library 93-04-23. Phil Morgan, Input argument checking 94-10-15. Phil Morgan, Pass S,T,P and keyword for ‘ptmp’ 99-06-25. Lindsay Pender, Fixed transpose of row vectors. 03-12-12. Lindsay Pender, Converted to ITS-90. 10-01-14. Filipe Fernandes, Python translation. 10-08-16. Filipe Fernandes, Reformulated docstring.

Examples

Data from McDougall 1987

>>> import seawater.csiro as sw
>>> s, t, p = 40, 10, 4000
>>> sw.alpha(s, t, p, pt=True)
0.00025061316481624323
seawater.csiro.aonb(s, t, p, pt=False)

Calculate \alpha/\beta.

Parameters :

s(p) : array_like

salinity [psu (PSS-78)]

t(p) : array_like

temperature or potential temperature [^\circ C (ITS-90)]

p : array_like

pressure [db]. The shape can be “broadcasted”

pt : bool

True if temperature is potential, default is False

Returns :

aonb : array_like

\alpha/\beta [psu ^\circ C -1]

See also

alpha, beta

Notes

TODO: Pressure broadcast feature need to be tested.

TODO: Test pt=False

References

[R72]McDougall, Trevor J., 1987: Neutral Surfaces. J. Phys. Oceanogr., 17, 1950-1964 doi: 10.1175/1520-0485(1987)017<1950:NS>2.0.CO;2
Modifications: N.L. Bindoff 1993, Lindsay Pender (Lindsay.Pender@csiro.au)
93-04-22. Phil Morgan, Help display modified to suit library 93-04-23. Phil Morgan, Input argument checking 94-10-15. Phil Morgan, Pass S,T,P and keyword for ‘ptmp’ 99-06-25. Lindsay Pender, Fixed transpose of row vectors. 03-12-12. Lindsay Pender, Converted to ITS-90. 10-01-14. Filipe Fernandes, Python translation. 10-08-16. Filipe Fernandes, Reformulated docstring.

Examples

Data from McDouogall 1987

>>> import seawater.csiro as sw
>>> s, t, p = 40, 10, 4000
>>> sw.aonb(s, t, p, pt=True)
0.347650567047807
seawater.csiro.beta(s, t, p, pt=False)

Calculate the saline contraction coefficient \beta as defined by T.J. McDougall.

Parameters :

s(p) : array_like

salinity [psu (PSS-78)]

t(p) : array_like

temperature or potential temperature [^\circ C (ITS-90)]

p : array_like

pressure [db]. The shape can be “broadcasted”

pt : bool

True if temperature is potential, default is False

Returns :

beta : array_like

saline Contraction Coefficient [psu -1]

Notes

TODO: Pressure broadcast feature need to be tested.

TODO: Test pt=False for alpha, beta and aonb.

References

[R73]McDougall, Trevor J., 1987: Neutral Surfaces. J. Phys. Oceanogr., 17, 1950-1964 doi: 10.1175/1520-0485(1987)017<1950:NS>2.0.CO;2
Modifications: N.L. Bindoff 1993, Lindsay Pender (Lindsay.pender@csiro.au)
93-04-22. Phil Morgan, Help display modified to suit library 93-04-23. Phil Morgan, Input argument checking 94-10-15. Phil Morgan, Pass S,T,P and keyword for ‘ptmp’ 99-06-25. Lindsay Pender, Fixed transpose of row vectors. 03-12-12. Lindsay Pender, Converted to ITS-90. 10-01-14. Filipe Fernandes, Python translation. 10-08-16. Filipe Fernandes, Reformulated docstring.

Examples

Data from McDouogall 1987

>>> import seawater.csiro as sw
>>> s, t, p = 40, 10, 4000
>>> sw.beta(s, t, p, pt=True)
0.00072087661741618932
seawater.csiro.bfrq(s, t, p, lat=None)

Calculates Brünt-Väisälä Frequency squared (N 2) at the mid depths from the equation:

N^{2} = \frac{-g}{\sigma_{\theta}} \frac{d\sigma_{\theta}}{dz}

Also calculates Potential Vorticity from:

q=f \frac{N^2}{g}

Parameters :

s(p) : array_like

salinity [psu (PSS-78)]

t(p) : array_like

temperature or potential temperature [^\circ C (ITS-90)]

p : array_like

pressure [db]. The shape can be “broadcasted”

lat : number or array_like, optional

latitude in decimal degrees north [-90..+90]. Will grav instead of the default g = 9.8 m 2 s -1) and d(z) instead of d(p)

Returns :

n2 : array_like

Brünt-Väisälä Frequency squared (M-1xN) [rad s -2]

q : array_like

planetary potential vorticity (M-1xN) [ m s -1]

p_ave : array_like

mid pressure between P grid (M-1xN) [db]

See also

pden, dens

Notes

The value of gravity is a global constant if lat is not provided.

TODO: Pressure broadcast feature need to be tested.

References

[R74]A.E. Gill 1982. p.54 eqn 3.7.15 “Atmosphere-Ocean Dynamics” Academic Press: New York. ISBN: 0-12-283522-0
[R75]Jackett, David R., Trevor J. Mcdougall, 1995: Minimal Adjustment of Hydrographic Profiles to Achieve Static Stability. J. Atmos. Oceanic Technol., 12, 381-389. doi: 10.1175/1520-0426(1995)012<0381:MAOHPT>2.0.CO;2
Modifications: Phil Morgan 93-06-24, Lindsay Pender (Lindsay.Pender@csiro.au)
Greg Johnson (gjohnson@pmel.noaa.gov) added potential vorticity calculation 03-12-12. Lindsay Pender, Converted to ITS-90. 06-04-19. Lindsay Pender, Corrected sign of PV. 10-01-14. Filipe Fernandes, Python translation. 10-08-17. Filipe Fernandes, Reformulated docstring.

Examples

>>> import numpy as np
>>> import seawater.csiro as sw
>>> s = np.array([[0, 0, 0], [15, 15, 15], [30, 30, 30],[35,35,35]])
>>> t = np.repeat(15, s.size).reshape(s.shape)
>>> p = [0, 250, 500, 1000]
>>> lat = [30,32,35]
>>> sw.bfrq(s, t, p, lat)[0]
array([[  4.51543648e-04,   4.51690708e-04,   4.51920753e-04],
       [  4.45598092e-04,   4.45743207e-04,   4.45970207e-04],
       [  7.40996788e-05,   7.41238078e-05,   7.41615525e-05]])
seawater.csiro.cndr(s, t, p)

Calculates conductivity ratio.

Parameters :

s(p) : array_like

salinity [psu (PSS-78)]

t(p) : array_like

temperature [^\circ C (ITS-90)]

p : array_like

pressure [db]. The shape can be “broadcasted”

Returns :

cndr : array_like

conductivity ratio. R = C(s,t,p)/C(35,15(IPTS-68),0) [no units]

See also

salds, sals, salrt

Notes

TODO: Pressure broadcast feature need to be tested.

References

[R76]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. http://unesdoc.unesco.org/images/0005/000598/059832eb.pdf
Modifications: Phil Morgan 93-04-21, Lindsay Pender (Lindsay.Pender@csiro.au)
99-06-25. Lindsay Pender, Fixed transpose of row vectors. 03-12-12. Lindsay Pender, Converted to ITS-90. 10-01-14. Filipe Fernandes, Python translation. 10-08-19. Filipe Fernandes, Reformulated docstring.

Examples

Data from UNESCO 1983 p9

>>> import seawater.csiro as sw
>>> t = sw.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.csiro.cor(lat)

Calculates the Coriolis factor f defined by:

f = 2 \Omega \sin(lat)

where:

\Omega = \frac{2 \pi}{\textrm{sidereal day}} = 7.2921150e^{-5} \textrm{ radians sec}^{-1}

Parameters :

lat : array_like

latitude in decimal degrees north [-90..+90].

Returns :

f : array_like

Coriolis factor [s -1]

See also

inertial_period

Notes

TODO

References

[R77]
  1. Pond & G.Pickard 2nd Edition 1986 Introductory Dynamical Oceanography Pergamon Press Sydney. ISBN 0-08-028728-X
[R78]A.E. Gill 1982. p.54 eqn 3.7.15 “Atmosphere-Ocean Dynamics” Academic Press: New York. ISBN: 0-12-283522-0
[R79]Groten, E., 2004: Fundamental Parameters and Current (2004) Best Estimates of the Parameters of Common Relevance to Astronomy, Geodesy, and Geodynamics. Journal of Geodesy, 77, pp. 724-797.
Modifications: Phil Morgan 93-04-20 (morgan@ml.csiro.au)
10-01-14. Filipe Fernandes, Python translation. 10-08-19. Filipe Fernandes, Reformulated docstring.

Examples

>>> import seawater.csiro as sw
>>> sw.cor(45)
0.00010312607931384281
seawater.csiro.cp(s, t, p)

Heat Capacity of Sea Water using UNESCO 1983 polynomial.

Parameters :

s(p) : array_like

salinity [psu (PSS-78)]

t(p) : array_like

temperature [^\circ C (ITS-90)]

p : array_like

pressure [db]. The shape can be “broadcasted”

Returns :

cp : array_like

specific heat capacity [J kg -1 C -1]

See also

TODO

Notes

TODO

References

[R80]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. http://unesdoc.unesco.org/images/0005/000598/059832eb.pdf
Modifications: Phil Morgan, Lindsay Pender (Lindsay.Pender@csiro.au)
99-06-25. Lindsay Pender, Fixed transpose of row vectors. 03-12-12. Lindsay Pender, Converted to ITS-90. 10-01-14. Filipe Fernandes, Python translation. 10-08-25. Filipe Fernandes, Reformulated docstring.

Examples

Data from Pond and Pickard Intro. Dynamical Oceanography 2nd ed. 1986

>>> import seawater.csiro as sw
>>> t = T90conv([[0, 0, 0, 0, 0, 0], [10, 10, 10, 10, 10, 10], [20, 20, 20, 20, 20, 20], [30, 30, 30, 30, 30, 30], [40, 40, 40, 40, 40, 40]])
>>> s = [[25, 25, 25, 35, 35, 35], [25, 25, 25, 35, 35, 35], [25, 25, 25, 35, 35, 35], [25, 25, 25, 35, 35, 35], [25, 25, 25, 35, 35, 35]]
>>> p = [0, 5000, 10000, 0, 5000, 10000]
>>> sw.cp(s, t, p)
array([[ 4048.4405375 ,  3896.25585   ,  3807.7330375 ,  3986.53309476,
         3849.26094605,  3769.11791286],
       [ 4041.8276691 ,  3919.5550066 ,  3842.3111366 ,  3986.34061786,
         3874.72665865,  3804.415624  ],
       [ 4044.8438591 ,  3938.5978466 ,  3866.7400391 ,  3993.85441786,
         3894.99294519,  3828.29059113],
       [ 4049.0984351 ,  3952.0375476 ,  3882.9855526 ,  4000.68382238,
         3909.24271128,  3844.32151784],
       [ 4051.2244911 ,  3966.1132036 ,  3905.9162711 ,  4003.46192541,
         3923.89463092,  3868.28959814]])
seawater.csiro.dens(s, t, p)

Density of Sea Water using UNESCO 1983 (EOS 80) polynomial.

Parameters :

s(p) : array_like

salinity [psu (PSS-78)]

t(p) : array_like

temperature [^\circ C (ITS-90)]

p : array_like

pressure [db]. The shape can be “broadcasted”

Returns :

dens : array_like

density [kg m 3]

See also

dens0, seck

Notes

TODO

References

[R81]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. http://unesdoc.unesco.org/images/0005/000598/059832eb.pdf
[R82]Millero, F.J., Chen, C.T., Bradshaw, A., and Schleicher, K. A new high pressure equation of state for seawater. Deap-Sea Research., 1980, Vol27A, pp255-264. doi:10.1016/0198-0149(80)90016-3
Modifications: Phil Morgan 92-11-05, Lindsay Pender (Lindsay.Pender@csiro.au)
99-06-25. Lindsay Pender, Fixed transpose of row vectors. 03-12-12. Lindsay Pender, Converted to ITS-90. 10-01-14. Filipe Fernandes, Python translation. 10-08-25. Filipe Fernandes, Reformulated docstring.

Examples

Data from Unesco Tech. Paper in Marine Sci. No. 44, p22

>>> import seawater.csiro 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.dens(s, t, p)
array([  999.842594  ,  1045.33710972,   995.65113374,  1036.03148891,
        1028.10633141,  1070.95838408,  1021.72863949,  1060.55058771])
seawater.csiro.dens0(s, t)

Density of Sea Water at atmospheric pressure.

Parameters :

s(p=0) : array_like

salinity [psu (PSS-78)]

t(p=0) : array_like

temperature [^\circ C (ITS-90)]

Returns :

dens0(s, t) : array_like

density [kg m 3] of salt water with properties (s, t, p=0) 0 db gauge pressure

See also

svan, dens, smow, seck, pden

Notes

TODO

References

[R83]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. http://unesdoc.unesco.org/images/0005/000598/059832eb.pdf
[R84]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
Modifications: Phil Morgan 92-11-05, Lindsay Pender (Lindsay.Pender@csiro.au)
03-12-12. Lindsay Pender, Converted to ITS-90. 10-01-14. Filipe Fernandes, Python translation. 10-08-25. Filipe Fernandes, Reformulated docstring.

Examples

Data from UNESCO Tech. Paper in Marine Sci. No. 44, p22

>>> import seawater.csiro as sw
>>> s = [0, 0, 0, 0, 35, 35, 35, 35]
>>> t = T90conv([0, 0, 30, 30, 0, 0, 30, 30])
>>> sw.dens0(s, t)
array([  999.842594  ,   999.842594  ,   995.65113374,   995.65113374,
        1028.10633141,  1028.10633141,  1021.72863949,  1021.72863949])
seawater.csiro.depth(p, lat)

Calculates depth in meters from pressure in dbars.

Parameters :

p : array_like

pressure [db].

lat : number or array_like

latitude in decimal degrees north [-90..+90]. The shape can be “broadcasted”

Returns :

z : array_like

depth [meters]

Notes

Original matlab seawater name is dpth and not depth.

References

[R85]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. http://unesdoc.unesco.org/images/0005/000598/059832eb.pdf
Modifications: Phil Morgan 92-04-06 (morgan@ml.csiro.au)
99-06-25. Lindsay Pender, Fixed transpose of row vectors. 10-01-14. Filipe Fernandes, Python translation. 10-08-19. Filipe Fernandes, Reformulated docstring.

Examples

UNESCO 1983 data p30

>>> import seawater.csiro as sw
>>> lat = [0, 30, 45, 90]
>>> p   = [[  500,   500,   500,  500], [ 5000,  5000,  5000, 5000], [10000, 10000, 10000, 10000]]
>>> sw.depth(p, lat)
array([[  496.65299239,   495.99772917,   495.3427354 ,   494.03357499],
       [ 4915.04099112,  4908.55954332,  4902.08075214,  4889.13132561],
       [ 9725.47087508,  9712.6530721 ,  9699.84050403,  9674.23144056]])
seawater.csiro.dist(lon, lat, units='km')

Calculate distance between two positions on globe using the “Plane Sailing” method. Also uses simple geometry to calculate the bearing of the path between position pairs.

Parameters :

lon : array_like

decimal degrees (+ve E, -ve W) [-180..+180]

lat : array_like

decimal degrees (+ve N, -ve S) [- 90.. +90]

units : string, optional

default kilometers

Returns :

dist : array_like

distance between positions in units

phaseangle : array_like

angle of line between stations with x axis (East). Range of values are -180..+180. (E=0, N=90, S=-90)

See also

TODO

Notes

Usually used to create a distance vector to plot hydrographic data. However, pay attention to the phaseangle to avoid apples and oranges!

Also not that the input order for the matlab version is lat,lon (alphabetic order), while this version is lon,lat (geometric order).

References

[R86]The PLANE SAILING method as described in “CELESTIAL NAVIGATION” 1989 by Dr. P. Gormley. The Australian Antarctic Division.
Modifications: Phil Morgan and Steve Rintoul 92-02-10
99-06-25. Lindsay Pender, Function name change from distance to sw_dist. 99-06-25. Lindsay Pender, Fixed transpose of row vectors. 10-01-14. Filipe Fernandes, Python translation. 10-08-25. Filipe Fernandes, Reformulated docstring.

Examples

>>> import seawater.csiro as sw
>>> lon = [35, 35]
>>> lat = [41, 40]
>>> sw.dist(lon, lat)
(array([ 111.12]), array([-90.]))

Create a distance vector

>>> lon = np.arange(30,40,1)
>>> lat = 35
>>> np.cumsum(np.append(0, sw.dist(lon, lat, units='km')[0]))
array([   0.        ,   91.02417516,  182.04835032,  273.07252548,
        364.09670065,  455.12087581,  546.14505097,  637.16922613,
        728.19340129,  819.21757645])
seawater.csiro.fp(s, p)

Freezing point of Sea Water using UNESCO 1983 polynomial.

Parameters :

s : array_like

salinity [psu (PSS-78)]

p : array_like

pressure [db]

Returns :

fp : array_like

freezing point temperature [^\circ C (ITS-90)]

See also

TODO

Notes

TODO

References

[R87]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. http://unesdoc.unesco.org/images/0005/000598/059832eb.pdf
Modifications: Phil Morgan 93-04-20, Lindsay Pender (Lindsay.Pender@csiro.au)
99-06-25. Lindsay Pender, Fixed transpose of row vectors. 03-12-12. Lindsay Pender, Converted to ITS-90. 10-01-14. Filipe Fernandes, Python translation. 10-08-19. Filipe Fernandes, Reformulated docstring.

Examples

UNESCO DATA p.30

>>> import seawater.csiro as sw
>>> s = [[5, 10, 15, 20, 25, 30, 35, 40], [5, 10, 15, 20, 25, 30, 35, 40]]
>>> p = [[ 0, 0, 0, 0, 0, 0, 0, 0], [500, 500, 500, 500, 500, 500, 500, 500]]
>>> sw.fp(s, p)
array([[-0.27369757, -0.54232831, -0.81142026, -1.0829461 , -1.35804594,
        -1.63748903, -1.9218401 , -2.2115367 ],
       [-0.65010724, -0.91873798, -1.18782992, -1.45935577, -1.73445561,
        -2.01389869, -2.29824976, -2.58794636]])
seawater.csiro.gpan(s, t, p)

Geopotential Anomaly calculated as the integral of svan from the the sea surface to the bottom. THUS RELATIVE TO SEA SURFACE.

Parameters :

s(p) : array_like

salinity [psu (PSS-78)]

t(p) : array_like

temperature [^\circ C (ITS-90)]

p : array_like

pressure [db]. The shape can be “broadcasted”

axis : int, optional

Axis along which the means are computed. The default is to compute using axis=1 or stations. Usually arrays are 2D MxN where M is depth and M stations.

Returns :

gpan : array_like

geopotential anomaly [m 3 kg -1 Pa == m 2 s -2 == J kg -1]

See also

svan, gvel

Notes

Adapted method from Pond and Pickard (p76) to calculate gpan relative to sea surface whereas P&P calculated relative to the deepest common depth. Note that older literature may use units of “dynamic decimeter” for above.

TODO: example with values that make some sense

TODO: pass axis as argument

References

[R88]
  1. Pond & G.Pickard 2nd Edition 1986 Introductory Dynamical Oceanography Pergamon Press Sydney. ISBN 0-08-028728-X
Modifications: Phil Morgan 92-11-05, Lindsay Pender (Lindsay.Pender@csiro.au)
03-12-12. Lindsay Pender, Converted to ITS-90. 10-01-14. Filipe Fernandes, Python translation. 10-08-25. Filipe Fernandes, Reformulated docstring.

Examples

Data from Unesco Tech. Paper in Marine Sci. No. 44, p22

>>> import numpy as np
>>> import seawater.csiro as sw
>>> s = np.array([[0, 0, 0], [15, 15, 15], [30, 30, 30],[35,35,35]])
>>> t = np.repeat(15, s.size).reshape(s.shape)
>>> p = [0, 250, 500, 1000]
>>> sw.gpan(s, t, p)
array([[   0.        ,    0.        ,    0.        ],
       [  56.35465209,   56.35465209,   56.35465209],
       [  84.67266947,   84.67266947,   84.67266947],
       [ 104.95799186,  104.95799186,  104.95799186]])
seawater.csiro.grav(lat, z=0)

Calculates acceleration due to gravity as function of latitude.

Parameters :

lat : array_like

latitude in decimal degrees north [-90..+90].

z : number or array_like. Default z = 0

height in meters (+ve above sea surface, -ve below).

Returns :

g : array_like

gravity [m s 2]

See also

bfrq

Notes

Original matlab name is g and not grav.

References

[R89]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. http://unesdoc.unesco.org/images/0005/000598/059832eb.pdf
[R90]A.E. Gill 1982. p.54 eqn 3.7.15 “Atmosphere-Ocean Dynamics” Academic Press: New York. ISBN: 0-12-283522-0
Modifications: Phil Morgan 93-04-20 (morgan@ml.csiro.au)
10-01-14. Filipe Fernandes, Python translation. 10-08-19. Filipe Fernandes, Reformulated docstring.

Examples

>>> import seawater.csiro as sw
>>> sw.grav(45, z=0)
9.8061898752053995
seawater.csiro.gvel(ga, distm, lat)

Calculates geostrophic velocity given the geopotential anomaly and position of each station.

Parameters :

ga : array_like

geopotential anomaly relative to the sea surface.

dist : array_like

distance between stations [meters]

lat : array_like

lat to use for constant f determination

Returns :

vel : array_like

geostrophic velocity relative to the sea surface. dimension will be MxN-1 (N: stations)

See also

svan, gpan

Notes

The original matlab version had gvel and gvel2 only, here the logic is “gvel2”, where one must compute the distance first.

TODO: dim(m, nstations-1) or pass axis?

TODO: add example with a reference level.

TODO: example with values that make some sense.

References

[R91]
  1. Pond & G.Pickard 2nd Edition 1986 Introductory Dynamical Oceanography Pergamon Press Sydney. ISBN 0-08-028728-X
Modifications: Phil Morgan 1992/03/26 (morgan@ml.csiro.au)
10-01-14. Filipe Fernandes, Python translation. 10-08-25. Filipe Fernandes, Reformulated docstring.

Examples

>>> import numpy as np
>>> import seawater.csiro as sw
>>> lon = [30, 30, 30]
>>> lat = [30, 32, 35]
>>> s = np.array([[0, 1, 2], [15, 16, 17], [30, 31, 32],[35,35,35]])
>>> t = np.repeat(15, s.size).reshape(s.shape)
>>> p = [0, 250, 500, 1000]
>>> ga = sw.gpan(s,t,p)
>>> distm = 1000.0 * sw.dist(lon, lat, units='km')[0]
>>> sw.gvel(ga, distm, lat)
array([[-0.        , -0.        ],
       [ 0.11385677,  0.07154215],
       [ 0.22436555,  0.14112761],
       [ 0.33366412,  0.20996272]])
seawater.csiro.pden(s, t, p, pr=0)

Calculates potential density of water mass relative to the specified reference pressure by pden = dens(S, ptmp, PR).

Parameters :

s(p) : array_like

salinity [psu (PSS-78)]

t(p) : array_like

temperature [^\circ C (ITS-90)]

p : array_like

pressure [db]. The shape can be “broadcasted”

pr : number

reference pressure [db], default = 0

Returns :

pden : array_like

potential density relative to the ref. pressure [kg m :sup:3]

See also

ptmp, dens

Notes

The reference pressure is in “oceanographic” standards, so 0 db means at surface or 1 atm.

References

[R92]A.E. Gill 1982. p.54 eqn 3.7.15 “Atmosphere-Ocean Dynamics” Academic Press: New York. ISBN: 0-12-283522-0
Modifications: Phil Morgan 1992/04/06, Lindsay Pender (Lindsay.Pender@csiro.au)
03-12-12. Lindsay Pender, Converted to ITS-90. 10-01-14. Filipe Fernandes, Python translation. 10-08-19. Filipe Fernandes, Reformulated docstring.

Examples

Data from Unesco Tech. Paper in Marine Sci. No. 44, p22

>>> import seawater.csiro 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.pden(s, t, p)
array([  999.842594  ,   999.79523994,   995.65113374,   996.36115932,
        1028.10633141,  1028.15738545,  1021.72863949,  1022.59634627])

\sigma_{4} (at 4000 db)

>>> sw.pden(s, t, p, 4000) - 1000
array([ 19.2895493 ,  19.33422519,  12.43271053,  13.27563816,
        46.30976432,  46.48818851,  37.76150878,  38.74500757])
seawater.csiro.pres(depth, lat)

Calculates pressure in dbars from depth in meters.

Parameters :

depth : array_like

depth [meters]

lat : array_like

latitude in decimal degrees north [-90..+90] The shape can be “broadcasted”

Returns :

p : array_like

pressure [db]

See also

pressure

Notes

TODO: lat broadcast feature need to be tested.

References

[R93]Saunders, Peter M., 1981: Practical Conversion of Pressure to Depth. J. Phys. Oceanogr., 11, 573-574. doi: 10.1175/1520-0485(1981)011<0573:PCOPTD>2.0.CO;2
Modifications: Phil Morgan 93-06-25 (morgan@ml.csiro.au)
99-06-25. Lindsay Pender, Fixed transpose of row vectors. 10-01-14. Filipe Fernandes, Python translation. 10-08-25. Filipe Fernandes, Reformulated docstring.

Examples

>>> import seawater.csiro as sw
>>> depth, lat = 7321.45, 30
>>> sw.pres(depth,lat)
7500.0065130118019
seawater.csiro.ptmp(s, t, p, pr=0)

Calculates potential temperature as per UNESCO 1983 report.

Parameters :

s(p) : array_like

salinity [psu (PSS-78)]

t(p) : array_like

temperature [^\circ C (ITS-90)]

p : array_like

pressure [db]. The shape can be “broadcasted”

pr : array_like

reference pressure [db], default = 0

Returns :

pt : array_like

potential temperature relative to PR [^\circ C (ITS-90)]

See also

adtg, pden

Notes

TODO

References

[R94]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. http://unesdoc.unesco.org/images/0005/000598/059832eb.pdf
[R95]Bryden, H. 1973. New Polynomials for thermal expansion, adiabatic temperature gradient and potential temperature of sea water. Deep-Sea Res. Vol20,401-408. doi:10.1016/0011-7471(73)90063-6
Modifications: Phil Morgan 92-04-06, Lindsay Pender (Lindsay.Pender@csiro.au)
99-06-25. Lindsay Pender, Fixed transpose of row vectors. 03-12-12. Lindsay Pender, Converted to ITS-90. 10-01-14. Filipe Fernandes, Python translation. 10-08-25. Filipe Fernandes, Reformulated docstring.

Examples

>>> import seawater.csiro as sw
>>> t = T90conv([[0, 0, 0, 0, 0, 0], [10, 10, 10, 10, 10, 10], [20, 20, 20, 20, 20, 20], [30, 30, 30, 30, 30, 30], [40, 40, 40, 40, 40, 40]])
>>> s = [[25, 25, 25, 35, 35, 35], [25, 25, 25, 35, 35, 35],  [25, 25, 25, 35, 35, 35], [25, 25, 25, 35, 35, 35], [25, 25, 25, 35, 35, 35]]
>>> p = [0, 5000, 10000, 0, 5000, 10000]
>>> sw.T68conv(sw.ptmp(s, t, p, pr=0))
array([[  0.        ,  -0.30614418,  -0.96669485,   0.        ,
         -0.3855565 ,  -1.09741136],
       [ 10.        ,   9.35306331,   8.46840949,  10.        ,
          9.29063461,   8.36425752],
       [ 20.        ,  19.04376281,  17.94265   ,  20.        ,
         18.99845171,  17.86536441],
       [ 30.        ,  28.75124632,  27.43529911,  30.        ,
         28.72313484,  27.38506197],
       [ 40.        ,  38.46068173,  36.92544552,  40.        ,
         38.44979906,  36.90231661]])
seawater.csiro.salds(rtx, delt)

Calculates Salinity differential (\frac{dS}{d(\sqrt{Rt})}) at constant temperature.

Parameters :

rtx : array_like

\sqrt{rt}

delt : array_like

t-15 [^\circ C (IPTS-68)]

Returns :

ds : array_like

\frac{dS}{d rtx}

See also

cndr, salt

Notes

TODO

References

[R96]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. http://unesdoc.unesco.org/images/0005/000598/059832eb.pdf
Modifications: Phil Morgan 93-04-21, Lindsay Pender (Lindsay.Pender@csiro.au)
10-01-14. Filipe Fernandes, Python translation. 10-08-19. Filipe Fernandes, Reformulated docstring.

Examples

Data from UNESCO 1983 p9

>>> import numpy as np
>>> import seawater.csiro 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.csiro.salrp(r, t, p)

Equation for Rp used in calculating salinity. UNESCO 1983 polynomial.

Rp(S,T,P) = \frac{C(S,T,P)}{C(S,T,0)}

Parameters :

r : array_like

conductivity ratio R = \frac{C(S,T,P)}{C(35,15(IPTS-68),0)} [no units]

t : array_like

temperature [^\circ C (ITS-90)]

p : array_like

pressure [db]

Returns :

rp : array_like

conductivity ratio Rp(S,T,P) = \frac{C(S,T,P)}{C(S,T,0)} [no units]

See also

salt

Notes

TODO

References

[R97]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. http://unesdoc.unesco.org/images/0005/000598/059832eb.pdf
Modifications: Phil Morgan 93-04-17, Lindsay Pender (Lindsay.Pender@csiro.au)
03-12-12. Lindsay Pender, Converted to ITS-90. 10-01-14. Filipe Fernandes, Python translation. 10-08-19. Filipe Fernandes, Reformulated docstring.

Examples

>>> import seawater.csiro 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.csiro.salrt(t)

Equation for rt used in calculating salinity. UNESCO 1983 polynomial.

rt(t) = \frac{C(35,t,0)}{C(35,15(\textrm{IPTS-68}), 0)}

Parameters :

t : array_like

temperature [^\circ C (ITS-90)]

Returns :

rt : array_like

conductivity ratio [no units] :

See also

salt

Notes

TODO

References

[R98]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. http://unesdoc.unesco.org/images/0005/000598/059832eb.pdf
Modifications: Phil Morgan 93-04-17, Lindsay Pender (Lindsay.Pender@csiro.au)
03-12-12. Lindsay Pender, Converted to ITS-90. 10-01-14. Filipe Fernandes, Python translation. 10-08-19. Filipe Fernandes, Reformulated docstring.

Examples

Data from UNESCO 1983 p9

>>> import seawater.csiro as sw
>>> t = T90conv([15, 20, 5])
>>> sw.salrt(t)
array([ 1.        ,  1.11649272,  0.77956585])
seawater.csiro.sals(rt, t)

Salinity of sea water as a function of Rt and T. UNESCO 1983 polynomial.

Parameters :

rt : array_like

rt(s,t) = \frac{C(s,t,0)}{C(35, t(\textrm{IPTS-68}), 0)}

t : array_like

temperature [^\circ C (ITS-90)]

Returns :

s : array_like

salinity [psu (PSS-78)]

See also

salt

Notes

TODO

References

[R99]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. http://unesdoc.unesco.org/images/0005/000598/059832eb.pdf
Modifications: Phil Morgan 93-04-17, Lindsay Pender (Lindsay.Pender@csiro.au)
03-12-12. Lindsay Pender, Converted to ITS-90. 10-01-14. Filipe Fernandes, Python translation. 10-08-19. Filipe Fernandes, Reformulated docstring.

Examples

Data from UNESCO 1983 p9

>>> import seawater.csiro as sw
>>> t = T90conv([15, 20, 5])
>>> rt   = [  1, 1.0568875, 0.81705885]
>>> sw.sals(rt, t)
array([ 35.        ,  37.24562718,  27.99534701])
seawater.csiro.salt(r, t, p)

Calculates Salinity from conductivity ratio. UNESCO 1983 polynomial.

Parameters :

r : array_like

conductivity ratio R = \frac{C(S,T,P)}{C(35,15(IPTS-68),0)} [no units]

t : array_like

temperature [^\circ C (ITS-90)]

p : array_like

pressure [db]

Returns :

s : array_like

salinity [psu (PSS-78)]

See also

sals, salrt, salrp

Notes

TODO

References

[R100]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. http://unesdoc.unesco.org/images/0005/000598/059832eb.pdf
Modifications: Phil Morgan 93-04-17, Lindsay Pender (Lindsay.Pender@csiro.au)
03-12-12. Lindsay Pender, Converted to ITS-90. 10-01-14. Filipe Fernandes, Python translation. 10-08-19. Filipe Fernandes, Reformulated docstring.

Examples

Data from UNESCO 1983 p9

>>> import seawater.csiro as sw
>>> r = [1, 1.2, 0.65]
>>> t = sw.T90conv([15, 20, 5])
>>> p = [0, 2000, 1500]
>>> sw.salt(r, t, p)
array([ 34.99999992,  37.24562765,  27.99534693])
seawater.csiro.satAr(s, t)

Solubility (saturation) of Argon (Ar) in sea water.

Parameters :

s : array_like

salinity [psu (PSS-78)]

t : array_like

temperature [^\circ C (ITS-90)]

Returns :

satAr : array_like

solubility of Ar [ml l -1]

See also

satO2, satN2

Notes

TODO

References

[R101]Weiss, R. F. 1970. The Solubility of Nitrogen, Oxygen and Argon in Water and Seawater Deep-Sea Research Vol. 17, p. 721-735 doi:10.1016/0011-7471(70)90037-9
Modifications: Phil Morgan 97-11-05, Lindsay Pender (Lindsay.Pender@csiro.au)
99-06-25. Lindsay Pender, Fixed transpose of row vectors. 03-12-12. Lindsay Pender, Converted to ITS-90. 10-01-14. Filipe Fernandes, Python translation. 10-08-25. Filipe Fernandes, Reformulated docstring.

Examples

Data from Weiss 1970

>>> import seawater.csiro as sw
>>> t = T90conv([[ -1, -1], [ 10, 10], [ 20, 20], [ 40, 40]])
>>> s = [[ 20, 40], [ 20, 40], [ 20, 40], [ 20, 40]]
>>> sw.satAr(s, t)
array([[ 0.4455784 ,  0.38766011],
       [ 0.33970659,  0.29887756],
       [ 0.27660227,  0.24566428],
       [ 0.19861429,  0.17937698]])
seawater.csiro.satN2(s, t)

Solubility (saturation) of Nitrogen (N2) in sea water.

Parameters :

s : array_like

salinity [psu (PSS-78)]

t : array_like

temperature [^\circ C (ITS-90)]

Returns :

satN2 : array_like

solubility of N2 [ml l -1]

See also

satO2, satAr

Notes

TODO

References

[R102]Weiss, R. F. 1970. The Solubility of Nitrogen, Oxygen and Argon in Water and Seawater Deep-Sea Research Vol. 17, p. 721-735 doi:10.1016/0011-7471(70)90037-9
Modifications: Phil Morgan 97-11-05, Lindsay Pender (Lindsay.Pender@csiro.au)
99-06-25. Lindsay Pender, Fixed transpose of row vectors. 03-12-12. Lindsay Pender, Converted to ITS-90. 10-01-14. Filipe Fernandes, Python translation. 10-08-25. Filipe Fernandes, Reformulated docstring.

Examples

Data from Weiss 1970

>>> import seawater.csiro as sw
>>> t = T90conv([[ -1, -1], [ 10, 10], [ 20, 20], [ 40, 40]])
>>> s = [[ 20, 40], [ 20, 40], [ 20, 40], [ 20, 40]]
>>> sw.satN2(s, t)
array([[ 16.27952432,  14.00784526],
       [ 12.64036196,  11.01277257],
       [ 10.46892822,   9.21126859],
       [  7.78163876,   6.95395099]])
seawater.csiro.satO2(s, t)

Solubility (saturation) of Oxygen (O2) in sea water.

Parameters :

s : array_like

salinity [psu (PSS-78)]

t : array_like

temperature [^\circ C (ITS-68)]

Returns :

satO2 : array_like

solubility of O2 [ml l -1 ]

Notes

TODO

References

[R103]Weiss, R. F. 1970. The Solubility of Nitrogen, Oxygen and Argon in Water and Seawater Deep-Sea Research Vol. 17, p. 721-735 doi:10.1016/0011-7471(70)90037-9
Modifications: Phil Morgan 97-11-05, Lindsay Pender (Lindsay.Pender@csiro.au)
99-06-25. Lindsay Pender, Fixed transpose of row vectors. 03-12-12. Lindsay Pender, Converted to ITS-90. 10-01-14. Filipe Fernandes, Python translation. 10-08-25. Filipe Fernandes, Reformulated docstring.

Examples

Data from Weiss 1970

>>> import seawater.csiro as sw
>>> t = T90conv([[ -1, -1], [ 10, 10], [ 20, 20], [ 40, 40]])
>>> s = [[ 20, 40], [ 20, 40], [ 20, 40], [ 20, 40]]
>>> sw.satO2(s, t)
array([[ 9.162056  ,  7.98404249],
       [ 6.95007741,  6.12101928],
       [ 5.64401453,  5.01531004],
       [ 4.0495115 ,  3.65575811]])
seawater.csiro.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 [^\circ C (ITS-90)]

p : array_like

pressure [db]. The shape can be “broadcasted”

Returns :

k : array_like

secant bulk modulus [bars]

See also

dens

Notes

TODO

References

[R104]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. http://unesdoc.unesco.org/images/0005/000598/059832eb.pdf
[R105]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
Modifications: Phil Morgan 92-11-05, Lindsay Pender (Lindsay.Pender@csiro.au)
99-06-25. Lindsay Pender, Fixed transpose of row vectors. 03-12-12. Lindsay Pender, Converted to ITS-90. 10-01-14. Filipe Fernandes, Python translation. 10-08-19. Filipe Fernandes, Reformulated docstring.

Examples

Data from Unesco Tech. Paper in Marine Sci. No. 44, p22

>>> import seawater.csiro 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.csiro.smow(t)

Density of Standard Mean Ocean Water (Pure Water) using EOS 1980.

Parameters :

t : array_like

temperature [^\circ C (ITS-90)]

Returns :

dens(t) : array_like

density [kg m 3]

See also

svan, dens, dens0, seck, pden

Notes

Standard Mean Ocean Water (SMOW) is the water collected in the deep ocean used as a reference.

References

[R106]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. http://unesdoc.unesco.org/images/0005/000598/059832eb.pdf
[R107]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
Modifications: Phil Morgan 92-11-05, Lindsay Pender (Lindsay.Pender@csiro.au)
99-06-25. Lindsay Pender, Fixed transpose of row vectors. 03-12-12. Lindsay Pender, Converted to ITS-90. 10-01-14. Filipe Fernandes, Python translation. 10-08-25. Filipe Fernandes, Reformulated docstring.

Examples

Data from UNESCO Tech. Paper in Marine Sci. No. 44, p22

>>> import seawater.csiro 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])
seawater.csiro.svan(s, t, p=0)

Specific Volume Anomaly calculated as svan = 1/dens(s, t, p) - 1/dens(35, 0, p).

Note that it is often quoted in literature as 1e8*units.

Parameters :

s(p) : array_like

salinity [psu (PSS-78)]

t(p) : array_like

temperature [^\circ C (ITS-90)]

p : array_like

pressure [db]. The shape can be “broadcasted”

Returns :

svan : array_like

specific volume anomaly [m 3 kg -1]

See also

dens

Notes

TODO

References

[R108]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. http://unesdoc.unesco.org/images/0005/000598/059832eb.pdf
[R109]
  1. Pond & G.Pickard 2nd Edition 1986 Introductory Dynamical Oceanography Pergamon Press Sydney. ISBN 0-08-028728-X
Modifications: Phil Morgan 92-11-05, Lindsay Pender (Lindsay.Pender@csiro.au)
99-06-25. Lindsay Pender, Fixed transpose of row vectors. 03-12-12. Lindsay Pender, Converted to ITS-90. 10-01-14. Filipe Fernandes, Python translation. 10-08-25. Filipe Fernandes, Reformulated docstring.

Examples

Data from Unesco Tech. Paper in Marine Sci. No. 44, p22

>>> import seawater.csiro 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.svan(s, t, p)
array([  2.74953924e-05,   2.28860986e-05,   3.17058231e-05,
         3.14785290e-05,   0.00000000e+00,   0.00000000e+00,
         6.07141523e-06,   9.16336113e-06])
seawater.csiro.svel(s, t, p)

Sound Velocity in sea water using UNESCO 1983 polynomial.

Parameters :

s(p) : array_like

salinity [psu (PSS-78)]

t(p) : array_like

temperature [^\circ C (ITS-90)]

p : array_like

pressure [db]. The shape can be “broadcasted”

Returns :

svel : array_like

sound velocity [m/s]

See also

TODO

Notes

TODO: Pressure broadcast feature need to be tested.

TODO: Add equation to docstring.

References

[R110]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. http://unesdoc.unesco.org/images/0005/000598/059832eb.pdf
Modifications: Phil Morgan 93-04-20, Lindsay Pender (Lindsay.Pender@csiro.au)
99-06-25. Lindsay Pender, Fixed transpose of row vectors. 03-12-12. Lindsay Pender, Converted to ITS-90. 10-01-14. Filipe Fernandes, Python translation. 10-08-25. Filipe Fernandes, Reformulated docstring.

Examples

Data from Pond and Pickard Intro. Dynamical Oceanography 2nd ed. 1986

>>> import seawater.csiro as sw
>>> t = T90conv([[  0,  0,  0,  0,  0,  0], [ 10, 10, 10, 10, 10, 10], [ 20, 20, 20, 20, 20, 20], [ 30, 30, 30, 30, 30, 30], [ 40, 40, 40, 40, 40, 40]])
>>> s = [[ 25, 25, 25, 35, 35, 35], [ 25, 25, 25, 35, 35, 35], [ 25, 25, 25, 35, 35, 35], [ 25, 25, 25, 35, 35, 35], [ 25, 25, 25, 35, 35, 35]]
>>> p = [ 0, 5000, 10000, 0, 5000, 10000]
>>> sw.svel(s, t, p)
array([[ 1435.789875  ,  1520.358725  ,  1610.4074    ,  1449.13882813,
         1533.96863705,  1623.15007097],
       [ 1477.68316464,  1561.30635914,  1647.39267114,  1489.82233602,
         1573.40946928,  1658.99115504],
       [ 1510.31388348,  1593.59671798,  1676.80967748,  1521.4619731 ,
         1604.4762822 ,  1687.18305631],
       [ 1535.21434752,  1618.95631952,  1700.60547902,  1545.59485539,
         1628.97322783,  1710.06294277],
       [ 1553.44506636,  1638.02522336,  1719.15088536,  1563.20925247,
         1647.29949576,  1727.83176404]])
seawater.csiro.swvel(length, depth)

Calculates surface wave velocity.

length : array_like
wave length
depth : array_like
water depth [meters]
Returns :

speed : array_like

surface wave speed [m s -1]

See also

TODO

Notes

TODO: add my wave function to extras

Examples

>>> import seawater.csiro as sw
>>> sw.swvel(10, 100)
3.9493270848342941
Modifications: Lindsay Pender 2005
10-01-14. Filipe Fernandes, Python translation. 10-08-25. Filipe Fernandes, Reformulated docstring.
seawater.csiro.temp(s, pt, p, pr=0)

Calculates temperature from potential temperature at the reference pressure PR and in situ pressure P.

Parameters :

s(p) : array_like

salinity [psu (PSS-78)]

pt(p) : array_like

potential temperature [^\circ C (ITS-90)]

p : array_like

pressure [db]. The shape can be “broadcasted”

pr : array_like

reference pressure [db]

Returns :

temp : array_like

temperature [^\circ C (ITS-90)]

See also

ptmp

Notes

TODO

References

[R111]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. http://unesdoc.unesco.org/images/0005/000598/059832eb.pdf
[R112]Bryden, H. 1973. New Polynomials for thermal expansion, adiabatic temperature gradient and potential temperature of sea water. Deep-Sea Res. Vol20,401-408. doi:10.1016/0011-7471(73)90063-6
Modifications: Phil Morgan 92-04-06, Lindsay Pender (Lindsay.Pender@csiro.au)
03-12-12. Lindsay Pender, Converted to ITS-90. 10-01-14. Filipe Fernandes, Python translation. 10-08-19. Filipe Fernandes, Reformulated docstring.

Examples

>>> import seawater.csiro as sw
>>> s, t, p = 35, 15, 100
>>> sw.temp(s, sw.ptmp(s, t, p), p)
15.0
seawater.csiro.test(fileout)

Execute test routines to test and verify SEAWATER Library routines for your platform. Prints output to file.

Notes

This is only to reproduce sw_test.m from the original. A better more complete test is performed via doctest.

Modifications: Phil Morgan, Lindsay Pender (Lindsay.Pender@csiro.au)
03-12-12. Lindsay Pender, Converted to ITS-90. 10-01-14. Filipe Fernandes, Python translation.

extras Documentation

This page contains the Extras Package documentation.

Waves Documentation

This page contains the Waves Package documentation.

The waves Package

class seawater.extras.waves.Waves(h, T=None, L=None, thetao=None, Ho=None, lat=None)

Solves the wave dispersion relationship.

\Omega = ...

Parameters :

h : array_like, str

Water depth [m] or ‘deep’, ‘shallow’ as keywords

T : array_like

Wave period [s]

L : array_like

Wave length [m]

thetao : array_like

TODO

Ho : array_like

TODO

Returns :

omega : array_like

Wave frequency

TODO: hoLo, hoL, Lo, L, k, T, Co, C, Cg, G, Ks, Kr, theta, H :

Notes

Compare values with: http://www.coastal.udel.edu/faculty/rad/wavetheory.html

References

TODO

Modifications: Filipe Fernandes, 2010
10-01-26. Filipe Fernandes, first version.

Examples

>>> from seawater.extras.waves import Waves
>>> Waves(h=10, T=5, L=None)

Indices and tables