gsw.gibbs package

Submodules

gsw.gibbs.constants module

Constants

gsw.gibbs.constants.C3515 = 42.914

Conductivity of 42.914 [mmho cm -1 == mS cm -1] at Salinity 35 psu, Temperature 15 \(^\circ\) C [ITPS 68] and Pressure 0 db.

References

[R7]Culkin and Smith, 1980: Determination of the Concentration of Potassium Chloride Solution Having the Same Electrical Conductivity, at 15C and Infinite Frequency, as Standard Seawater of Salinity 35.0000 (Chlorinity 19.37394), IEEE J. Oceanic Eng, 5, 22-23.
[R8]Unesco, 1983: Algorithms for computation of fundamental properties of seawater. Unesco Technical Papers in Marine Science, 44, 53 pp.
gsw.gibbs.constants.Kelvin = 273.15

The Celsius zero point; 273.15 K. That is T = t + T0 where T is the Absolute Temperature (in degrees K) and t is temperature in degrees C.

gsw.gibbs.constants.M_S = 0.0314038218

Mole-weighted average atomic weight of the elements of Reference-Composition sea salt, in units of kg mol -1. Strictly speaking, the formula below applies only to seawater of Reference Composition. If molality is required to an accuracy of better than 0.1% we suggest you contact the authors for further guidance.

gsw.gibbs.constants.OMEGA = 7.292115e-05
\(\Omega = \frac{2\pi}{\textrm{sidereal day}}\) =
7.292e-5.radians sec -1

1 sidereal day = 23.9344696 hours

Changed to a more precise value at Groten 2004

References

[R9]A.E. Gill 1982. p.54 eqn 3.7.15 “Atmosphere-Ocean Dynamics” Academic Press: New York. ISBN: 0-12-283522-0. page: 597
[R10]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.
gsw.gibbs.constants.P0 = 101325

Absolute Pressure of one standard atmosphere in Pa, 101325 Pa.

gsw.gibbs.constants.R = 8.314472

The molar gas constant = 8.314472 m 2 kg s:sup:-21 K :sup:-1` mol -1.

gsw.gibbs.constants.SSO = 35.16504

SSO is the Standard Ocean Reference Salinity (35.16504 g/kg.)

SSO is the best estimate of the Absolute Salinity of Standard Seawater when the seawater sample has a Practical Salinity, SP, of 35 (Millero et al., 2008), and this number is a fundamental part of the TEOS-10 definition of seawater.

gsw.gibbs.constants.SonCl = 1.80655

The ratio of Practical Salinity, SP, to Chlorinity, 1.80655 kg/g for Reference Seawater (Millero et al., 2008). This is the ratio that was used by the JPOTS committee in their construction of the 1978 Practical Salinity Scale (PSS-78) to convert between the laboratory measurements of seawater samples (which were measured in Chlorinity) to Practical Salinity.

gsw.gibbs.constants.T0 = 273.15

The Celsius zero point; 273.15 K. That is T = t + T0 where T is the Absolute Temperature (in degrees K) and t is temperature in degrees C.

gsw.gibbs.constants.cp0 = 3991.86795711963

The “specific heat” for use with Conservative Temperature. cp0 is the ratio of potential enthalpy to Conservative Temperature. See Eqn. (3.3.3) and Table D.5 from IOC et al. (2010).

gsw.gibbs.constants.db2Pascal = 10000.0

Decibar to pascal.

gsw.gibbs.constants.earth_radius = 6371000.0

Mean radius of earth A.E. Gill.

gsw.gibbs.constants.gamma = 2.26e-07

Gamma (A.E. Gill).

gsw.gibbs.constants.r1 = 0.35

TODO

gsw.gibbs.constants.sfac = 0.0248826675584615

sfac = 1 / (40 * uPS) = 1 / (40. * (SSO / 35.))

gsw.gibbs.constants.soffset = 0.5971840214030754

soffset = deltaS * sfac.

gsw.gibbs.constants.uPS = 1.0047154285714286

The unit conversion factor for salinities (35.16504/35) g/kg (Millero et al., 2008). Reference Salinity SR is uPS times Practical Salinity SP.

Ratio, unit conversion factor for salinities [g kg -1]

References

[R11]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 section 6, Eqn.(6.1).
gsw.gibbs.constants.valence_factor = 1.2452898

This function returns the valence factor of sea salt of Reference Composition, 1.2452898. This valence factor is exact, and follows from the definition of the Reference-Composition Salinity Scale 2008 of Millero et al. (2008). The valence factor is the mole-weighted square of the charges, Z, of the ions comprising Reference Composition sea salt.

gsw.gibbs.conversions module

gsw.gibbs.conversions.CT_from_t(SA, t, p)[source]

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:

TODO

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 3.3.

Examples

TODO

gsw.gibbs.conversions.SA_from_SP(SP, p, lon, lat)[source]

Calculates Absolute Salinity from Practical Salinity.

Parameters:

SP : array_like

salinity (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 : masked array

Absolute salinity [g kg -1]

Notes

The mask 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

[R13]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.
[R14]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

Examples

>>> import gsw
>>> SP = [34.5487, 34.7275, 34.8605, 34.6810, 34.5680, 34.5600]
>>> p = [10, 50, 125, 250, 600, 1000]
>>> lon = 188
>>> lat = 4
>>> gsw.SA_from_SP(SP, p, lon, lat)
array([ 34.71177834,  34.89152262,  35.02554486,  34.84722903,
        34.73662847,  34.73236307])
gsw.gibbs.conversions.Sstar_from_SP(SP, p, lon, lat)[source]

Calculates Preformed Salinity from Absolute Salinity.

Parameters:

SP : array_like

salinity (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 : masked array

Preformed Salinity [g kg -1]

Notes

The mask 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

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

Examples

>>> import gsw
>>> SP = [34.5487, 34.7275, 34.8605, 34.6810, 34.5680, 34.5600]
>>> p = [10, 50, 125, 250, 600, 1000]
>>> lon = 188
>>> lat = 4
>>> gsw.Sstar_from_SP(SP, p, lon, lat)
array([ 34.71155368,  34.8911614 ,  35.02465027,  34.84359314,
        34.729034  ,  34.71967596])
gsw.gibbs.conversions.Abs_Pressure_from_p(p)[source]

Calculates Absolute Pressure from sea pressure. Note that Absolute Pressure is in Pa NOT dbar.

p : array_like
sea pressure [dbar]
Returns:

Absolute_Pressure : array_like

Absolute Pressure [Pa]

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.2.1).

Examples

>>> import gsw
>>> gsw.Abs_Pressure_from_p([0, 500, 1000])
array([   101325.,   5101325.,  10101325.])
gsw.gibbs.conversions.CT_from_entropy(SA, entropy)[source]

Calculates Conservative Temperature with entropy as an input variable.

Parameters:

SA : array_like

Absolute salinity [g kg -1]

entropy : array_like

specific entropy [J kg -1 K -1]

Returns:

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

References

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

Examples

>>> import 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.CT_from_entropy(SA, entropy)
array([ 28.80990279,  28.43919923,  22.78619927,  10.22619767,
         6.82719674,   4.32360295])
gsw.gibbs.conversions.CT_from_pt(SA, pt)[source]

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 (ITS-90)]

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 section 3.3.

Examples

>>> import 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])
gsw.gibbs.conversions.SA_Sstar_from_SP(SP, p, lon, lat)[source]

TODO: docstring

gsw.gibbs.conversions.SA_from_Sstar(Sstar, p, lon, lat)[source]

TODO: docstring

gsw.gibbs.conversions.SP_from_SA(SA, p, lon, lat)[source]

TODO: docstring

gsw.gibbs.conversions.SP_from_SR(SR)[source]

Calculates Practical Salinity from Reference Salinity.

SR : array_like
Reference Salinity [g kg -1]
Returns:

SP : array_like

Practical Salinity (PSS-78) [unitless]

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.

Examples

TODO

gsw.gibbs.conversions.SR_from_SP(SP)[source]

Calculates Reference Salinity from Practical Salinity.

SP : array_like
Practical Salinity (PSS-78) [unitless]
Returns:

SR : array_like

Reference Salinity [g kg -1]

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.

Examples

TODO

gsw.gibbs.conversions.Sstar_from_SA(SA, p, lon, lat)[source]

TODO: docstring

gsw.gibbs.conversions.deltaSA_from_SP(SP, p, lon, lat)[source]
gsw_deltaSA_from_SP Absolute Salinity Anomaly
from Practical Salinity
USAGE:
deltaSA = gsw_deltaSA_from_SP(SP,p,long,lat)
DESCRIPTION:
Calculates Absolute Salinity Anomaly from Practical Salinity. Since SP is non-negative by definition, this function changes any negative input values of SP to be zero.
INPUT:

SP = Practical Salinity (PSS-78) [ unitless ] p = sea pressure [ dbar ]

( i.e. absolute pressure - 10.1325 dbar )
long = longitude in decimal degrees [ 0 ... +360 ]
or [ -180 ... +180 ]

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

p, lat & long may have dimensions 1x1 or Mx1 or 1xN or MxN, where SP is MxN.

OUTPUT:
deltaSA = Absolute Salinity Anomaly [ g/kg ]
AUTHOR:
Trevor McDougall & Paul Barker [ help@teos-10.org ]

VERSION NUMBER: 3.03 (29th April, 2013)

REFERENCES:
IOC, SCOR and IAPSO, 2010: The international thermodynamic equation of

seawater - 2010: Calculation and use of thermodynamic properties. Intergovernmental Oceanographic Commission, Manuals and Guides No. 56, UNESCO (English), 196 pp. Available from http://www.TEOS-10.org

See section 2.5 and appendices A.4 and A.5 of this TEOS-10 Manual.
McDougall, T.J., D.R. Jackett, F.J. Millero, R. Pawlowicz and
P.M. Barker, 2012: A global algorithm for estimating Absolute Salinity. Ocean Science, 8, 1117-1128. http://www.ocean-sci.net/8/1117/2012/os-8-1117-2012.pdf
gsw.gibbs.conversions.depth_from_z(z)[source]

Calculates depth from height, z. Note that in general height is negative in the ocean.

z : array_like
height [m]
Returns:

depth : array_like

depth [m]

gsw.gibbs.conversions.entropy_from_CT(SA, CT)[source]

Calculates specific entropy of seawater.

Parameters:

SA : array_like

Absolute salinity [g kg -1]

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

Returns:

entropy : array_like

specific entropy [J kg -1 K -1]

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

Examples

>>> import 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.entropy_from_CT(SA, CT)
array([ 400.38916315,  395.43781023,  319.86680989,  146.79103279,
         98.64714648,   62.79185763])
gsw.gibbs.conversions.entropy_from_pt(SA, pt)[source]

Calculates specific entropy of seawater.

Parameters:

SA : array_like

Absolute salinity [g kg -1]

pt : array_like

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

Returns:

entropy : array_like

specific entropy [J kg -1 K -1]

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 appendix A.10.

Examples

>>> import 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_pt(SA, pt)
array([ 400.38946744,  395.43839949,  319.86743859,  146.79054828,
         98.64691006,   62.79135672])
gsw.gibbs.conversions.entropy_from_t(SA, t, p)[source]
gsw.gibbs.conversions.ionic_strength_from_SA(SA)[source]

Calculates the ionic strength of seawater from Absolute Salinity.

Parameters:

SA : array_like

Absolute salinity [g kg -1]

Returns:

ionic_strength : array_like

ionic strength of seawater [mol kg -1]

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. See Table L.1.
[R25]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.

Examples

>>> import gsw
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
>>> gsw.ionic_strength_from_SA(SA)
array([ 0.71298118,  0.71680567,  0.71966059,  0.71586272,  0.71350891,
        0.71341953])
gsw.gibbs.conversions.molality_from_SA(SA)[source]

Calculates the molality of seawater from Absolute Salinity.

Parameters:

SA : array_like

Absolute salinity [g kg -1]

Returns:

molality : array_like

seawater molality [mol kg -1]

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.

Examples

>>> import gsw
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
>>> gsw.molality_from_SA(SA)
array([ 1.14508476,  1.15122708,  1.15581223,  1.14971265,  1.14593231,
        1.14578877])
gsw.gibbs.conversions.p_from_Abs_Pressure(Absolute_Pressure)[source]

Calculates sea pressure from Absolute Pressure. Note that Absolute Pressure is in Pa NOT dbar.

Absolute_Pressure : array_like
Absolute Pressure [Pa]
Returns:

p : array_like

sea pressure [dbar]

References

[R27]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.2.1).

Examples

TODO

gsw.gibbs.conversions.p_from_z(z, lat, geo_strf_dyn_height=0)[source]

Calculates sea pressure from height using computationally-efficient 48-term expression for density, in terms of SA, CT and p (McDougall et al., 2011). Dynamic height anomaly, geo_strf_dyn_height, if provided, must be computed with its pr=0 (the surface.)

Parameters:

z : array_like

height [m]

lat : array_like

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

geo_strf_dyn_height : float, optional

dynamic height anomaly [ m 2 s -2 ] The reference pressure (p_ref) of geo_strf_dyn_height must be zero (0) dbar.

Returns:

p : array_like

pressure [dbar]

Notes

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

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.
[R29]McDougall T.J., P.M. Barker, R. Feistel and D.R. Jackett, 2011: A computationally efficient 48-term expression for the density of seawater in terms of Conservative Temperature, and related properties of seawater.
[R30]Moritz (2000) Goedetic reference system 1980. J. Geodesy, 74, 128-133.
[R31]Saunders, P. M., 1981: Practical conversion of pressure to depth. Journal of Physical Oceanography, 11, 573-574.

Examples

>>> import gsw
>>> z = [-10., -50., -125., -250., -600., -1000.]
>>> lat = 4.
>>> gsw.p_from_z(z, lat)
array([   10.05572704,    50.28354425,   125.73185732,   251.54028663,
         604.2099135 ,  1007.9900587 ])
>>> -gsw.z_from_p(gsw.p_from_z(z, lat), lat)
array([  10.,   50.,  125.,  250.,  600., 1000.])
gsw.gibbs.conversions.pot_enthalpy_from_pt(SA, pt)[source]

Calculates the potential enthalpy 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:

pot_enthalpy : array_like

potential enthalpy [J kg -1]

Notes

TODO

References

[R32]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.2.

Examples

>>> import 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.pot_enthalpy_from_pt(SA, pt)
array([ 115005.40853458,  113525.30870246,   90959.68769935,
         40821.50280454,   27253.21472227,   17259.10131183])
gsw.gibbs.conversions.pt0_from_t(SA, t, p)[source]

Calculates potential temperature with reference pressure, pr = 0 dbar. The present routine is computationally faster than the more general function “pt_from_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

pt_from_t has the same result (only slower)

References

[R33]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.
[R34]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.

Examples

>>> import 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])
gsw.gibbs.conversions.pt_from_CT(SA, CT)[source]

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 (ITS-90)]

Returns:

pt : array_like

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

See also

specvol_anom

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

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

Examples

>>> import 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])
gsw.gibbs.conversions.pt_from_entropy(SA, entropy)[source]

Calculates potential temperature with reference pressure p_ref = 0 dbar and with entropy as an input variable.

Parameters:

SA : array_like

Absolute salinity [g kg -1]

entropy : array_like

specific entropy [J kg -1 K -1]

Returns:

pt : array_like

potential temperature [\(^\circ\) C (ITS-90)] with reference sea pressure (p_ref) = 0 dbar.

References

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

Examples

>>> import 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.pt_from_entropy(SA, entropy)
array([ 28.78317983,  28.42095483,  22.78495274,  10.23053207,
         6.82921333,   4.32453778])
gsw.gibbs.conversions.pt_from_t(SA, t, p, p_ref=0)[source]

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]

p_ref : int, float, optional

reference pressure, default = 0

Returns:

pt : array_like

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

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 p_ref is indeed zero dbar.

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 3.1.
[R39]McDougall T.J., P.M. Barker, R. Feistel and D.R. Jackett, 2011: A computationally efficient 48-term expression for the density of seawater in terms of Conservative Temperature, and related properties of seawater.

Examples

>>> import 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.pt_from_t(SA, t, p)
array([ 28.78319682,  28.42098334,  22.7849304 ,  10.23052366,
         6.82923022,   4.32451057])
>>> gsw.pt_from_t(SA, t, p, p_ref=1000)
array([ 29.02665528,  28.662375  ,  22.99149634,  10.35341725,
         6.92732954,   4.4036    ])
gsw.gibbs.conversions.t90_from_t48(t48)[source]

Converts IPTS-48 temperature to International Temperature Scale 1990 (ITS-90) temperature. This conversion should be applied to all in-situ data collected prior to 31/12/1967.

t48 : array_like
in-situ temperature [\(^\circ\) C (ITPS-48)]
Returns:

t90 : array_like

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

References

[R40]International Temperature Scales of 1948, 1968 and 1990, an ICES note, available from http://www.ices.dk/ocean/procedures/its.htm
gsw.gibbs.conversions.t90_from_t68(t68)[source]

Converts IPTS-68 temperature to International Temperature Scale 1990 (ITS-90) temperature. This conversion should be applied to all in-situ data collected between 1/1/1968 and 31/12/1989.

t68 : array_like
in-situ temperature [\(^\circ\) C (ITPS-68)]
Returns:

t90 : array_like

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

References

[R41]International Temperature Scales of 1948, 1968 and 1990, an ICES note, available from http://www.ices.dk/ocean/procedures/its.htm
gsw.gibbs.conversions.t_from_CT(SA, CT, p)[source]

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 (ITS-90)]

p : array_like

pressure [dbar]

Returns:

t : array_like

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

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

Examples

>>> import 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])
gsw.gibbs.conversions.t_from_entropy(SA, entropy, p)[source]
gsw_t_from_entropy in-situ temperature
as a function of entropy
USAGE:
t = gsw_t_from_entropy(SA,entropy,p)
DESCRIPTION:
Calculates in-situ temperature with entropy as an input variable.
INPUT:

SA = Absolute Salinity [ g/kg ] entropy = specific entropy [ J/(kg*K) ] p = sea pressure [ dbar ]

( i.e. absolute pressure - 10.1325 dbar )

SA & entropy need to have the same dimensions. p may have dimensions 1x1 or Mx1 or 1xN or MxN, where SA & entropy are MxN.

OUTPUT:
t = in-situ temperature (ITS-90) [ deg C ]
AUTHOR:
Trevor McDougall and Paul Barker [ help@teos-10.org ]

VERSION NUMBER: 3.03 (29th April, 2013)

REFERENCES:
IOC, SCOR and IAPSO, 2010: The international thermodynamic equation of
seawater - 2010: Calculation and use of thermodynamic properties. Intergovernmental Oceanographic Commission, Manuals and Guides No. 56, UNESCO (English), 196 pp. Available from http://www.TEOS-10.org See appendix A.10 of this TEOS-10 Manual.
gsw.gibbs.conversions.z_from_depth(depth)[source]

Calculates height, z, from depth. Note that in general height is negative in the ocean.

depth : array_like
depth [m]
Returns:

z : array_like

height [m]

gsw.gibbs.conversions.z_from_p(p, lat, geo_strf_dyn_height=0)[source]

Calculates height from sea pressure using the computationally-efficient 48-term expression for density in terms of SA, CT and p (McDougall et al., 2011). Dynamic height anomaly, geo_strf_dyn_height, if provided, must be computed with its pr=0 (the surface).

Parameters:

p : array_like

pressure [dbar]

lat : array_like

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

geo_strf_dyn_height : float, optional

dynamic height anomaly [ m 2 s -2 ]

Returns:

z : array_like

height [m]

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

[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.
[R44]McDougall T.J., P.M. Barker, R. Feistel and D.R. Jackett, 2011: A computationally efficient 48-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.

Examples

>>> import gsw
>>> p = [10, 50, 125, 250, 600, 1000]
>>> lat = 4
>>> gsw.z_from_p(p, lat)
array([  -9.94458313,  -49.71808883, -124.27262301, -248.47007032,
       -595.82544461, -992.09217796])

gsw.gibbs.density_enthalpy_48 module

gsw.gibbs.density_enthalpy_48.alpha(SA, CT, p)[source]

Calculates the thermal expansion coefficient of seawater with respect to Conservative Temperature using the computationally-efficient 48-term expression for density in terms of SA, CT and p (McDougall et al., 2011)

Parameters:

SA : array_like

Absolute Salinity [g/kg]

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

p : array_like

sea pressure [dbar]

Returns:

alpha : array_like

thermal expansion coefficient [K \(-1\)] with respect to Conservative Temperature

Notes

The 48-term equation has been fitted in a restricted range of parameter space, and is most accurate inside the “oceanographic funnel” described in McDougall et al. (2011). The GSW library function “infunnel(SA, CT, p)” is available to be used if one wants to test if some of one’s data lies outside this “funnel”.

References

[R46]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).
[R47]McDougall T.J., P.M. Barker, R. Feistel and D.R. Jackett, 2011: A computationally efficient 48-term expression for the density of seawater in terms of Conservative Temperature, and related properties of seawater.

Examples

TODO

gsw.gibbs.density_enthalpy_48.alpha_on_beta(SA, CT, p)[source]
gsw.gibbs.density_enthalpy_48.alpha_on_beta_CT(SA, CT, p)
gsw.gibbs.density_enthalpy_48.beta(SA, CT, p)[source]

Calculates the saline (i.e. haline) contraction coefficient of seawater at constant Conservative Temperature using the computationally-efficient 48-term expression for density in terms of SA, CT and p (McDougall et al., 2011).

Parameters:

SA : array_like

Absolute Salinity [g/kg]

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

p : array_like

sea pressure [dbar]

Returns:

beta : array_like

saline contraction coefficient [kg g \(-1\)] at constant Conservative Temperature.

Notes

The 48-term equation has been fitted in a restricted range of parameter space, and is most accurate inside the “oceanographic funnel” described in McDougall et al. (2011). The GSW library function “infunnel(SA, CT, p)” is available to be used if one wants to test if some of one’s data lies outside this “funnel”.

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 Eqn. (2.19.3).
[R49]McDougall T.J., P.M. Barker, R. Feistel and D.R. Jackett, 2011: A computationally efficient 48-term expression for the density of seawater in terms of Conservative Temperature, and related properties of seawater.

Examples

TODO

gsw.gibbs.density_enthalpy_48.dynamic_enthalpy(SA, CT, p)[source]

Calculates dynamic enthalpy of seawater using the computationally- efficient 48-term expression for density in terms of SA, CT and p (McDougall et al., 2011). Dynamic enthalpy is defined as enthalpy minus potential enthalpy (Young, 2010).

Parameters:

SA : array_like

Absolute Salinity [g/kg]

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

p : array_like

sea pressure [dbar]

Returns:

dynamic_enthalpy : array_like

dynamic enthalpy [J/kg]

Notes

The 48-term equation has been fitted in a restricted range of parameter space, and is most accurate inside the “oceanographic funnel” described in McDougall et al. (2011). The GSW library function “infunnel(SA, CT, p)” is available to be used if one wants to test if some of one’s data lies outside this “funnel”.

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.2
[R51]McDougall T.J., P.M. Barker, R. Feistel and D.R. Jackett, 2011: A computationally efficient 48-term expression for the density of seawater in terms of Conservative Temperature, and related properties of seawater.
[R52]Young, W.R., 2010: Dynamic enthalpy, Conservative Temperature, and the seawater Boussinesq approximation. Journal of Physical Oceanography, 40, 394-400.

Examples

TODO

gsw.gibbs.density_enthalpy_48.enthalpy(SA, CT, p)[source]

Calculates specific enthalpy of seawater using the computationally- efficient 48-term expression for density in terms of SA, CT and p (McDougall et al., 2011)

Parameters:

SA : array_like

Absolute Salinity [g/kg]

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

p : array_like

sea pressure [dbar]

Returns:

enthalpy : array_like

specific enthalpy [J/kg]

Notes

The 48-term equation has been fitted in a restricted range of parameter space, and is most accurate inside the “oceanographic funnel” described in McDougall et al. (2011). The GSW library function “infunnel(SA, CT, p)” is available to be used if one wants to test if some of one’s data lies outside this “funnel”.

References

[R53]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. (A.30.6).
[R54]McDougall T.J., P.M. Barker, R. Feistel and D.R. Jackett, 2011: A computationally efficient 48-term expression for the density of seawater in terms of Conservative Temperature, and related properties of seawater.

Examples

TODO

gsw.gibbs.density_enthalpy_48.enthalpy_diff(SA, CT, p_shallow, p_deep)[source]

Calculates the difference of the specific enthalpy of seawater between two different pressures, p_deep (the deeper pressure) and p_shallow (the shallower pressure), at the same values of SA and CT. This function uses the computationally-efficient 48-term expression for density in terms of SA, CT and p (McDougall et al., 2011). The output (enthalpy_diff_CT) is the specific enthalpy evaluated at (SA, CT, p_deep) minus the specific enthalpy at (SA, CT, p_shallow).

Parameters:

SA : array_like

Absolute Salinity [g/kg]

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

p_shallow : array_like

lower sea pressure [dbar]

p_deep : array_like

upper sea pressure [dbar]

Returns:

enthalpy_diff : array_like

difference of specific enthalpy [J/kg] (deep minus shallow)

Notes

The 48-term equation has been fitted in a restricted range of parameter space, and is most accurate inside the “oceanographic funnel” described in McDougall et al. (2011). The GSW library function “infunnel(SA, CT, p)” is available to be used if one wants to test if some of one’s data lies outside this “funnel”.

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 Eqns. (3.32.2) and (A.30.6).
[R56]McDougall T.J., P.M. Barker, R. Feistel and D.R. Jackett, 2011: A computationally efficient 48-term expression for the density of seawater in terms of Conservative Temperature, and related properties of seawater.

Examples

TODO

gsw.gibbs.density_enthalpy_48.internal_energy(SA, CT, p)[source]

Calculates specific internal energy of seawater using the computationally-efficient 48-term expression for density in terms of SA, CT and p (McDougall et al., 2011).

Parameters:

SA : array_like

Absolute Salinity [g/kg]

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

p : array_like

sea pressure [dbar]

Returns:

internal_energy : array_like

specific internal energy [J/kg]

Notes

The 48-term equation has been fitted in a restricted range of parameter space, and is most accurate inside the “oceanographic funnel” described in McDougall et al. (2011). The GSW library function “infunnel(SA, CT, p)” is available to be used if one wants to test if some of one’s data lies outside this “funnel”.

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.
[R58]McDougall T.J., P.M. Barker, R. Feistel and D.R. Jackett, 2011: A computationally efficient 48-term expression for the density of seawater in terms of Conservative Temperature, and related properties of seawater.

Examples

TODO

gsw.gibbs.density_enthalpy_48.rho(SA, CT, p)[source]

Calculates in-situ density from Absolute Salinity and Conservative Temperature, using the computationally-efficient 48-term expression for density in terms of SA, CT and p (McDougall et al., 2011).

Parameters:

SA : array_like

Absolute Salinity [g/kg]

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

p : array_like

sea pressure [dbar]

Returns:

rho : array_like

in-situ density [kg/m**3]

Notes

The potential density with respect to reference pressure, pr, is obtained by calling this function with the pressure argument being pr (i.e. “rho(SA,CT,pr)”).

The 48-term equation has been fitted in a restricted range of parameter space, and is most accurate inside the “oceanographic funnel” described in McDougall et al. (2011). The GSW library function “infunnel(SA,CT,p)” is available to be used if one wants to test if some of one’s data lies outside this “funnel”.

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.20 and appendix K.
[R60]McDougall T.J., P.M. Barker, R. Feistel and D.R. Jackett, 2011: A computationally efficient 48-term expression for the density of seawater in terms of Conservative Temperature, and related properties of seawater.

Examples

TODO

gsw.gibbs.density_enthalpy_48.rho_alpha_beta(SA, CT, p)[source]

Calculates in-situ density, the appropriate thermal expansion coefficient and the appropriate saline contraction coefficient of seawater from Absolute Salinity and Conservative Temperature. This function uses the computationally-efficient 48-term expression for density in terms of SA, CT and p (McDougall et al., 2011).

The potential density (pot_rho) with respect to reference pressure p_ref is obtained by calling this function with the pressure argument being p_ref as in pot_rho, _, _] = rho_alpha_beta(SA, CT, p_ref).

Parameters:

SA : array_like

Absolute Salinity [g/kg]

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

p : array_like

sea pressure [dbar]

Returns:

rho : array_like

in-situ density [kg/m**3]

alpha : array_like

thermal expansion coefficient [K \(-1\)] with respect to Conservative Temperature

beta : array_like

saline contraction coefficient [kg g \(-1\)] at constant Conservative Temperature.

Notes

The 48-term equation has been fitted in a restricted range of parameter space, and is most accurate inside the “oceanographic funnel” described in McDougall et al. (2011). The GSW library function “infunnel(SA, CT, p)” is available to be used if one wants to test if some of one’s data lies outside this “funnel”.

References

[R61]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.20 and appendix K.
[R62]McDougall T.J., P.M. Barker, R. Feistel and D.R. Jackett, 2011: A computationally efficient 48-term expression for the density of seawater in terms of Conservative Temperature, and related properties of seawater.

Examples

TODO

gsw.gibbs.density_enthalpy_48.SA_from_rho(rho, CT, p)[source]

Calculates the Absolute Salinity of a seawater sample, for given values of its density, Conservative Temperature and sea pressure (in dbar). This function uses the computationally-efficient 48-term expression for density in terms of SA, CT and p (McDougall et al., 2011).

Parameters:

rho : array_like

density of a seawater sample [kg/m**3] This input has not had 1000 kg/m^3 subtracted from it (e.g. 1026 kg m**-3), that is, it is density, NOT density anomaly.

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

p : array_like

sea pressure [dbar]

Returns:

SA : array_like

Absolute Salinity [g/kg]

Notes

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

The 48-term equation has been fitted in a restricted range of parameter space, and is most accurate inside the “oceanographic funnel” described in McDougall et al. (2011). The GSW library function “infunnel(SA, CT, p)” is available to be used if one wants to test if some of one’s data lies outside this “funnel”.

References

[R63]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
[R64]McDougall T.J., P.M. Barker, R. Feistel and D.R. Jackett, 2011: A computationally efficient 48-term expression for the density of seawater in terms of Conservative Temperature, and related properties of seawater.
[R65]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.

Examples

TODO

gsw.gibbs.density_enthalpy_48.sigma0(SA, CT)[source]

Calculates potential density anomaly with reference pressure of 0 dbar, this being this particular potential density minus 1000 kg/m^3. This function has inputs of Absolute Salinity and Conservative Temperature. This function uses the computationally-efficient 48-term expression for density in terms of SA, CT and p (McDougall et al., 2011).

Parameters:

SA : array_like

Absolute Salinity [g/kg]

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

Returns:

sigma0 : array_like

potential density anomaly with [kg/m**3] respect to a reference pressure of 0 dbar

gsw.gibbs.density_enthalpy_48.sigma1(SA, CT)[source]

Calculates potential density anomaly with reference pressure of 1000 dbar, this being this particular potential density minus 1000 kg/m^3. This function has inputs of Absolute Salinity and Conservative Temperature.

Parameters:

SA : array_like

Absolute Salinity [g/kg]

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

Returns:

sigma1 : array_like

potential density anomaly with [kg/m**3] respect to a reference pressure of 1000 dbar

gsw.gibbs.density_enthalpy_48.sigma2(SA, CT)[source]

Calculates potential density anomaly with reference pressure of 2000 dbar, this being this particular potential density minus 1000 kg/m^3. This function has inputs of Absolute Salinity and Conservative Temperature.

Parameters:

SA : array_like

Absolute Salinity [g/kg]

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

Returns:

sigma2 : array_like

potential density anomaly with [kg/m**3] respect to a reference pressure of 2000 dbar

gsw.gibbs.density_enthalpy_48.sigma3(SA, CT)[source]

Calculates potential density anomaly with reference pressure of 3000 dbar, this being this particular potential density minus 1000 kg/m^3. This function has inputs of Absolute Salinity and Conservative Temperature.

Parameters:

SA : array_like

Absolute Salinity [g/kg]

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

Returns:

sigma1 : array_like

potential density anomaly with [kg/m**3] respect to a reference pressure of 3000 dbar

gsw.gibbs.density_enthalpy_48.sigma4(SA, CT)[source]

Calculates potential density anomaly with reference pressure of 4000 dbar, this being this particular potential density minus 1000 kg/m^3. This function has inputs of Absolute Salinity and Conservative Temperature.

Parameters:

SA : array_like

Absolute Salinity [g/kg]

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

Returns:

sigma1 : array_like

potential density anomaly with [kg/m**3] respect to a reference pressure of 4000 dbar

gsw.gibbs.density_enthalpy_48.sound_speed(SA, CT, p)[source]

Calculates the speed of sound in seawater. This function has inputs of Absolute Salinity and Conservative Temperature. This function uses the computationally-efficient 48-term expression for density in terms of SA, CT and p (McDougall et al., 2011).

Parameters:

SA : array_like

Absolute Salinity [g/kg]

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

p : array_like

sea pressure [dbar]

Returns:

sound_speed : array_like

speed of sound in seawater [m/s]

Notes

Approximate with a r.m.s. of 6.7 cm s^-1.

The 48-term equation has been fitted in a restricted range of parameter space, and is most accurate inside the “oceanographic funnel” described in McDougall et al. (2011). The GSW library function “infunnel(SA, CT, p)” is available to be used if one wants to test if some of one’s data lies outside this “funnel”.

References

[R66]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).
[R67]McDougall T.J., P.M. Barker, R. Feistel and D.R. Jackett, 2011: A computationally efficient 48-term expression for the density of seawater in terms of Conservative Temperature, and related properties of seawater.

Examples

TODO

gsw.gibbs.density_enthalpy_48.specvol(SA, CT, p)[source]

Calculates specific volume from Absolute Salinity, Conservative Temperature and pressure, using the computationally-efficient 48-term expression for density (McDougall et al., 2011).

Parameters:

SA : array_like

Absolute Salinity [g/kg]

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

p : array_like

sea pressure [dbar]

Returns:

specvol : array_like

specific volume [m**3/kg]

Notes

The 48-term equation has been fitted in a restricted range of parameter space, and is most accurate inside the “oceanographic funnel” described in McDougall et al. (2011). The GSW library function “infunnel(SA,CT,p)” is available to be used if one wants to test if some of one’s data lies outside this “funnel”.

References

[R68]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.7.2).
[R69]McDougall T.J., P.M. Barker, R. Feistel and D.R. Jackett, 2011: A computationally efficient 48-term expression for the density of seawater in terms of Conservative Temperature, and related properties of seawater.

Examples

TODO

gsw.gibbs.density_enthalpy_48.specvol_anom(SA, CT, p)[source]

Calculates specific volume anomaly from Absolute Salinity, Conservative Temperature and pressure. It uses the computationally-efficient 48-term expression for density as a function of SA, CT and p (McDougall et al., 2011). The reference value of Absolute Salinity is SSO and the reference value of Conservative Temperature is equal to 0 degrees C.

Parameters:

SA : array_like

Absolute Salinity [g/kg]

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

p : array_like

sea pressure [dbar]

Returns:

specvol_anom : array_like

specific volume anomaly [m**3/kg]

Notes

The 48-term equation has been fitted in a restricted range of parameter space, and is most accurate inside the “oceanographic funnel” described in McDougall et al. (2011). The GSW library function “infunnel(SA,CT,p)” is available to be used if one wants to test if some of one’s data lies outside this “funnel”.

References

[R70]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).
[R71]McDougall T.J., P.M. Barker, R. Feistel and D.R. Jackett, 2011: A computationally efficient 48-term expression for the density of seawater in terms of Conservative Temperature, and related properties of seawater.

Examples

TODO

gsw.gibbs.density_enthalpy_48_ct module

gsw.gibbs.density_enthalpy_ct_exact module

gsw.gibbs.density_enthalpy_ct_exact.alpha_CT_exact(SA, CT, p)[source]

Calculates the thermal expansion coefficient of seawater with respect to Conservative Temperature from Absolute Salinity and Conservative Temperature.

Parameters:

SA : array_like

Absolute salinity [g kg -1]

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

p : array_like

pressure [dbar]

Returns:

alpha_CT_exact : array_like

thermal expansion coefficient [K -1] with respect to Conservative Temperature

Notes

This function uses the full Gibbs function. There is an alternative to calling this function, namely alpha_wrt_CT(SA, CT, p) which uses the computationally efficient 48-term expression for density in terms of SA, CT and p (McDougall et al., (2011)).

References

[R72]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).
[R73]McDougall T.J., P.M. Barker, R. Feistel and D.R. Jackett, 2011: A computationally efficient 48-term expression for the density of seawater in terms of Conservative Temperature, and related properties of seawater.

Examples

TODO

gsw.gibbs.density_enthalpy_ct_exact.beta_CT_exact(SA, CT, p)[source]

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

Parameters:

SA : array_like

Absolute salinity [g kg -1]

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

p : array_like

pressure [dbar]

Returns:

beta_CT_exact : array_like

thermal expansion coefficient [K -1]

Notes

This function uses the full Gibbs function. There is an alternative to calling this function, namely beta_const_CT(SA, CT, p) which uses the computationally efficient 48-term expression for density in terms of SA, CT and p (McDougall et al., (2011)).

References

[R74]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).
[R75]McDougall T.J., P.M. Barker, R. Feistel and D.R. Jackett, 2011: A computationally efficient 48-term expression for the density of seawater in terms of Conservative Temperature, and related properties of seawater.

Examples

TODO

gsw.gibbs.density_enthalpy_ct_exact.CT_from_rho_exact(rho, SA, p)[source]

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

Parameters:

rho : array_like

density of a seawater sample [kg/m**3] This input has not had 1000 kg/m^3 subtracted from it (e.g. 1026 kg m**-3), that is, it is density, NOT density anomaly.

SA : array_like

Absolute Salinity [g/kg]

p : array_like

sea pressure [dbar]

Returns:

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

CT_multiple : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

Notes

At low salinities, in brackish water, there are two possible temperatures for a single density. This program will output both valid solutions (t, t_multiple), if there is only one possible solution the second variable will be set to NaN.

References

[R76]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.
gsw.gibbs.density_enthalpy_ct_exact.CT_maxdensity_exact(SA, p)[source]

Calculates the Conservative Temperature of maximum density of seawater. This function returns the Conservative temperature at which the density of seawater is a maximum, at given Absolute Salinity, SA, and sea pressure, p (in dbar).

Parameters:

SA : array_like

Absolute Salinity [g/kg]

p : array_like

sea pressure [dbar]

Returns:

CT_maxdensity_exact : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)] at which the density of seawater is a maximum for given SA and p.

References

[R77]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.42.

Examples

TODO

gsw.gibbs.density_enthalpy_ct_exact.dynamic_enthalpy_CT_exact(SA, CT, p)[source]

Calculates the dynamic enthalpy of seawater from Absolute Salinity and Conservative Temperature and pressure. Dynamic enthalpy is defined as enthalpy minus potential enthalpy (Young, 2010).

Parameters:

SA : array_like

Absolute Salinity [g/kg]

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

p : array_like

sea pressure [dbar]

Returns:

dynamic_enthalpy_CT_exact : array_like

dynamic enthalpy [J/kg]

See also

TODO

Notes

This function uses the full Gibbs function. There is an alternative to calling this function, namely dynamic_enthalpy(SA, CT, p), which uses the computationally efficient 48-term expression for density in terms of SA, CT and p (McDougall et al., 2011).

References

[R78]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 apendix A.30.
[R79]McDougall T.J., P.M. Barker, R. Feistel and D.R. Jackett, 2011: A computationally efficient 48-term expression for the density of seawater in terms of Conservative Temperature, and related properties of seawater.
[R80]Young, W.R., 2010: Dynamic enthalpy, Conservative Temperature, and the seawater Boussinesq approximation. Journal of Physical Oceanography, 40, 394-400.

Examples

TODO

gsw.gibbs.density_enthalpy_ct_exact.enthalpy_CT_exact(SA, CT, p)[source]

Calculates specific enthalpy of seawater from Absolute Salinity and Conservative Temperature and pressure.

Parameters:

SA : array_like

Absolute Salinity [g/kg]

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

p : array_like

sea pressure [dbar]

Returns:

enthalpy_CT_exact : array_like

specific enthalpy [J/kg]

Notes

This function uses the full Gibbs function. There is an alternative to calling this function, namely enthalpy_CT(SA, CT, p), which uses the computationally-efficient 48-term expression for density in terms of SA, CT and p (McDougall et al., 2011).

References

[R81]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.
[R82]McDougall T.J., P.M. Barker, R. Feistel and D.R. Jackett, 2011: A computationally efficient 48-term expression for the density of seawater in terms of Conservative Temperature, and related properties of seawater.

Examples

TODO

gsw.gibbs.density_enthalpy_ct_exact.enthalpy_diff_CT_exact(SA, CT, p_shallow, p_deep)[source]

Calculates the difference of the specific enthalpy of seawater between two different pressures, p_deep (the deeper pressure) and p_shallow (the shallower pressure), at the same values of SA and CT. The output (enthalpy_diff_CT_exact) is the specific enthalpy evaluated at (SA, CT, p_deep) minus the specific enthalpy at (SA,CT,p_shallow).

Parameters:

SA : array_like

Absolute Salinity [g/kg]

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

p_shallow : array_like

lower sea pressure [dbar]

p_deep : array_like

upper sea pressure [dbar]

Returns:

enthalpy_diff_CT_exact : array_like

difference of specific enthalpy [J/kg] (deep minus shallow)

Notes

This function uses the full Gibbs function. There is an alternative to calling this function, namely enthalpy_diff_CT(SA, CT, p), which uses the computationally efficient 48-term expression for density in terms of SA, CT and p (McDougall et al., 2011).

References

[R83]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 (3.32.2).
[R84]McDougall T.J., P.M. Barker, R. Feistel and D.R. Jackett, 2011: A computationally efficient 48-term expression for the density of seawater in terms of Conservative Temperature, and related properties of seawater.

Examples

TODO

gsw.gibbs.density_enthalpy_ct_exact.internal_energy_CT_exact(SA, CT, p)[source]

Calculates the specific internal energy of seawater from Absolute Salinity, Conservative Temperature and pressure.

Parameters:

SA : array_like

Absolute Salinity [g/kg]

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

p : array_like

sea pressure [dbar]

Returns:

internal_energy_CT_exact: array_like

specific internal energy (u) [J/kg]

References

[R85]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).

Examples

TODO

gsw.gibbs.density_enthalpy_ct_exact.rho_CT_exact(SA, CT, p)[source]

Calculates in-situ density from Absolute Salinity and Conservative Temperature.

Parameters:

SA : array_like

Absolute Salinity [g/kg]

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

p : array_like

sea pressure [dbar]

Returns:

rho_CT_exact : array_like

in-situ density [kg/m**3]

Notes

The potential density with respect to reference pressure, p_ref, is obtained by calling this function with the pressure argument being p_ref (i.e. “rho_CT_exact(SA, CT, p_ref)”). This function uses the full Gibbs function. There is an alternative to calling this function, namely rho_CT(SA, CT, p), which uses the computationally efficient 48-term expression for density in terms of SA, CT and p (McDougall et al., 2011).

References

[R86]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.8.2).
[R87]McDougall T.J., P.M. Barker, R. Feistel and D.R. Jackett, 2011: A computationally efficient 48-term expression for the density of seawater in terms of Conservative Temperature, and related properties of seawater.

Examples

TODO

gsw.gibbs.density_enthalpy_ct_exact.rho_alpha_beta_CT_exact(SA, CT, p)[source]

Calculates in-situ density, the appropriate thermal expansion coefficient and the appropriate saline contraction coefficient of seawater from Absolute Salinity and Conservative Temperature. See the individual functions rho_CT_exact, alpha_CT_exact, and beta_CT_exact. Retained for compatibility with the Matlab GSW toolbox.

gsw.gibbs.density_enthalpy_ct_exact.SA_from_rho_CT_exact(rho, CT, p)[source]

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

Parameters:

rho : array_like

density of a seawater sample [kg/m**3] This input has not had 1000 kg/m^3 subtracted from it (e.g. 1026 kg m**-3), that is, it is density, NOT density anomaly.

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

p : array_like

sea pressure [dbar]

Returns:

SA : array_like

Absolute Salinity [g/kg]

Notes

This function uses the full Gibbs function. There is an alternative to calling this function, namely SA_from_rho_CT(rho, CT, p), which uses the computationally efficient 48-term expression for density in terms of SA, CT and p (McDougall et al., 2011).

References

[R88]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.
[R89]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.

Examples

TODO

gsw.gibbs.density_enthalpy_ct_exact.sigma0_CT_exact(SA, CT)[source]

Calculates potential density anomaly with reference pressure of 0 dbar, this being this particular potential density minus 1000 kg/m^3. This function has inputs of Absolute Salinity and Conservative Temperature.

Parameters:

SA : array_like

Absolute Salinity [g/kg]

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

Returns:

sigma0_CT_exact: array_like

Potential density anomaly with [kg/m**3] respect to a reference pressure of 0 dbar that is, this potential density - 1000 kg/m**3.

Notes

Note that this function uses the full Gibbs function. There is an alternative to calling this function, namely gsw_sigma0_CT(SA,CT,p), which uses the computationally efficient 48-term expression for density in terms of SA, CT and p (McDougall et al., 2011).

References

[R90]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. (A.30.1).
[R91]McDougall T.J., P.M. Barker, R. Feistel and D.R. Jackett, 2011: A computationally efficient 48-term expression for the density of seawater in terms of Conservative Temperature, and related properties of seawater.

Examples

TODO

gsw.gibbs.density_enthalpy_ct_exact.sigma1_CT_exact(SA, CT)[source]

Calculates potential density anomaly with reference pressure of 1000 dbar.

gsw.gibbs.density_enthalpy_ct_exact.sigma2_CT_exact(SA, CT)[source]

Calculates potential density anomaly with reference pressure of 2000 dbar.

gsw.gibbs.density_enthalpy_ct_exact.sigma3_CT_exact(SA, CT)[source]

Calculates potential density anomaly with reference pressure of 3000 dbar.

gsw.gibbs.density_enthalpy_ct_exact.sigma4_CT_exact(SA, CT)[source]

Calculates potential density anomaly with reference pressure of 4000 dbar.

gsw.gibbs.density_enthalpy_ct_exact.sound_speed_CT_exact(SA, CT, p)[source]

Calculates the speed of sound in seawater from Absolute Salinity and Conservative Temperature and pressure.

Parameters:

SA : array_like

Absolute salinity [g kg -1]

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

p : array_like

pressure [dbar]

Returns:

sound_speed_CT_exact : array_like

Speed of sound in seawater [m s -1]

References

[R92]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).

Examples

TODO

gsw.gibbs.density_enthalpy_ct_exact.specvol_CT_exact(SA, CT, p)[source]

Calculates specific volume from Absolute Salinity, Conservative Temperature and pressure.

Parameters:

SA : array_like

Absolute Salinity [g/kg]

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

p : array_like

sea pressure [dbar]

Returns:

specvol_CT_exact : array_like

specific volume [m**3/kg]

Notes

This function uses the full Gibbs function. There is an alternative to calling this function, namely specvol_CT(SA, CT, p), which uses the computationally efficient 48-term expression for density in terms of SA, CT and p (McDougall et al., 2011).

References

[R93]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.7.2).
[R94]McDougall T.J., P.M. Barker, R. Feistel and D.R. Jackett, 2011: A computationally efficient 48-term expression for the density of seawater in terms of Conservative Temperature, and related properties of seawater.

Examples

TODO

gsw.gibbs.density_enthalpy_ct_exact.specvol_anom_CT_exact(SA, CT, p)[source]

Calculates specific volume anomaly from Absolute Salinity, Conservative Temperature and pressure. The reference value of Absolute Salinity is SSO and the reference value of Conservative Temperature is equal to 0 deg C.

Parameters:

SA : array_like

Absolute Salinity [g/kg]

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

p : array_like

sea pressure [dbar]

Returns:

specvol_anom_CT_exact : array_like

specific volume anomaly [m**3/kg]

Notes

This function uses the full Gibbs function. There is an alternative to calling this function, namely specvol_anom_CT(SA, CT, p), which uses the computationally efficient 48-term expression for density in terms of SA, CT and p (McDougall et al., 2011).

References

[R95]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).
[R96]McDougall T.J., P.M. Barker, R. Feistel and D.R. Jackett, 2011: A computationally efficient 48-term expression for the density of seawater in terms of Conservative Temperature, and related properties of seawater.

Examples

TODO

gsw.gibbs.derivatives module

gsw.gibbs.derivatives.CT_first_derivatives(SA, pt)[source]

Calculates the following two derivatives of Conservative Temperature (1) CT_SA, the derivative with respect to Absolute Salinity at constant

potential temperature (with pr = 0 dbar), and
  1. CT_pt, the derivative with respect to potential temperature (the regular potential temperature which is referenced to 0 dbar) at constant Absolute Salinity.
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_SA : array_like

The derivative of CT with respect to SA at constant potential temperature reference sea pressure of 0 dbar. [K (g kg -1) -1]

CT_pt : array_like

The derivative of CT with respect to pt at constant SA. [ unitless ]

References

[R97]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. (A.12.3) and (A.12.9a,b).
[R98]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.

Examples

>>> import gsw
>>> SA = 34.7118
>>> pt = 28.7832
>>> gsw.CT_first_derivatives(SA, pt)
(-0.041981092877805957, 1.0028149372966355)
gsw.gibbs.derivatives.CT_second_derivatives(SA, pt)[source]

Calculates the following three, second-order derivatives of Conservative Temperature (1) CT_SA_SA, the second derivative with respect to Absolute Salinity at

constant potential temperature (with p_ref = 0 dbar),
  1. CT_SA_pt, the derivative with respect to potential temperature (the regular potential temperature which is referenced to 0 dbar) and Absolute Salinity, and
  2. CT_pt_pt, the second derivative with respect to potential temperature (the regular potential temperature which is referenced to 0 dbar) at constant Absolute Salinity.
Parameters:

SA : array_like

Absolute salinity [g kg -1]

pt : array_like

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

Returns:

CT_SA_SA : array_like

The second derivative of Conservative Temperature with respect to Absolute Salinity at constant potential temperature (the regular potential temperature which has reference sea pressure of 0 dbar). [K/((g/kg)^2)]

CT_SA_pt : array_like

The derivative of Conservative Temperature with respect to potential temperature (the regular one with p_ref = 0 dbar) and Absolute Salinity. [1/(g/kg)]

CT_pt_pt : array_like

The second derivative of Conservative Temperature with respect to potential temperature (the regular one with p_ref = 0 dbar) at constant SA. [1/K]

References

[R99]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.12.
[R100]McDougall T.J., P.M. Barker, R. Feistel and D.R. Jackett, 2011: A computationally efficient 48-term expression for the density of seawater in terms of Conservative Temperature, and related properties of seawater. To be submitted to Ocean Science Discussions.

Examples

TODO

gsw.gibbs.derivatives.enthalpy_first_derivatives(SA, CT, p)[source]

Calculates the following three derivatives of specific enthalpy (h) (1) h_SA, the derivative with respect to Absolute Salinity at

constant CT and p, and
  1. h_CT, derivative with respect to CT at constant SA and p.
  2. h_P, derivative with respect to pressure (in Pa) at constant SA and CT.
Parameters:

SA : array_like

Absolute salinity [g kg -1]

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

p : array_like

pressure [dbar]

Returns:

h_SA : array_like

The first derivative of specific enthalpy with respect to Absolute Salinity at constant CT and p. [J/(kg (g/kg))] i.e. [J/g]

h_CT : array_like

The first derivative of specific enthalpy with respect to CT at constant SA and p. [J/(kg K)]

h_P : array_like

The first partial derivative of specific enthalpy with respect to pressure (in Pa) at fixed SA and CT. Note that h_P is specific volume (1/rho.)

References

[R101]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. (A.11.18), (A.11.15) and (A.11.12.)

Examples

TODO

gsw.gibbs.derivatives.enthalpy_second_derivatives(SA, CT, p)[source]

Calculates the following three second-order derivatives of specific enthalpy (h) (1) h_SA_SA, second-order derivative with respect to Absolute Salinity

at constant CT & p.
  1. h_SA_CT, second-order derivative with respect to SA & CT at constant p.
  2. h_CT_CT, second-order derivative with respect to CT at constant SA and p.
Parameters:

SA : array_like

Absolute salinity [g kg -1]

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

p : array_like

pressure [dbar]

Returns:

h_SA_SA : array_like

The second derivative of specific enthalpy with respect to Absolute Salinity at constant CT & p. [J/(kg (g/kg)^2)]

h_SA_CT : array_like

The second derivative of specific enthalpy with respect to SA and CT at constant p. [J/(kg K(g/kg))]

h_CT_CT : array_like

The second derivative of specific enthalpy with respect to CT at constant SA and p. [J/(kg K^2)]

References

[R102]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. (A.11.18), (A.11.15) and (A.11.12.)
[R103]McDougall T.J., P.M. Barker, R. Feistel and D.R. Jackett, 2011: A computationally efficient 48-term expression for the density of seawater in terms of Conservative Temperature, and related properties of seawater. To be submitted to Ocean Science Discussions.

Examples

TODO

gsw.gibbs.derivatives.entropy_first_derivatives(SA, CT)[source]

Calculates the following two partial derivatives of specific entropy (eta) (1) eta_SA, the derivative with respect to Absolute Salinity at constant

Conservative Temperature, and
  1. eta_CT, the derivative with respect to Conservative Temperature at constant Absolute Salinity.
Parameters:

SA : array_like

Absolute salinity [g kg -1]

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

Returns:

eta_SA : array_like

The derivative of specific entropy with respect to SA at constant CT [J g -1 K -1]

eta_CT : array_like

The derivative of specific entropy with respect to CT at constant SA [ J (kg K -2) -1 ]

References

[R104]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. (A.12.8) and (P.14a,c).

Examples

>>> import gsw
>>> SA = 34.7118
>>> CT = 28.8099
>>> gsw.entropy_first_derivatives(SA, CT)
(-0.26328680071165517, 13.221031210083824)
gsw.gibbs.derivatives.entropy_second_derivatives(SA, CT)[source]

Calculates the following three second-order partial derivatives of specific entropy (eta) (1) eta_SA_SA, the second derivative with respect to Absolute Salinity at

constant Conservative Temperature, and
  1. eta_SA_CT, the derivative with respect to Absolute Salinity and Conservative Temperature.
  2. eta_CT_CT, the second derivative with respect to Conservative Temperature at constant Absolute Salinity.
Parameters:

SA : array_like

Absolute salinity [g kg -1]

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

Returns:

eta_SA_SA : array_like

The second derivative of specific entropy with respect to SA at constant CT [J (kg K (g kg -1 ) 2) -1]

eta_SA_CT : array_like

The second derivative of specific entropy with respect to SA and CT [J (kg (g kg -1 ) K 2) -1 ]

eta_CT_CT : array_like

The second derivative of specific entropy with respect to CT at constant SA [J (kg K 3) -1 ]

References

[R105]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. (P.14b) and (P.15a,b).

Examples

>>> import gsw
>>> SA = 34.7118
>>> CT = 28.8099
>>> gsw.entropy_second_derivatives(SA, CT)
(-0.0076277189296690626, -0.0018331042167510211, -0.043665023731108685)
gsw.gibbs.derivatives.pt_first_derivatives(SA, CT)[source]

Calculates the following two partial derivatives of potential temperature (the regular potential temperature whose reference sea pressure is 0 dbar) (1) pt_SA, the derivative with respect to Absolute Salinity at

constant Conservative Temperature, and
  1. pt_CT, the derivative with respect to Conservative Temperature at constant Absolute Salinity.
Parameters:

SA : array_like

Absolute salinity [g kg -1]

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

Returns:

pt_SA : array_like

The derivative of potential temperature with respect to Absolute Salinity at constant Conservative Temperature. [K/(g/kg)]

pt_CT : array_like

The derivative of potential temperature with respect to Conservative Temperature at constant Absolute Salinity. [unitless]

References

[R106]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. (A.12.6), (A.12.3), (P.6) and (P.8).
[R107]McDougall T.J., P.M. Barker, R. Feistel and D.R. Jackett, 2011: A computationally efficient 48-term expression for the density of seawater in terms of Conservative Temperature, and related properties of seawater. To be submitted to Ocean Science Discussions.

Examples

TODO

gsw.gibbs.derivatives.pt_second_derivatives(SA, CT)[source]

Calculates the following three second-order derivatives of potential temperature (the regular potential temperature which has a reference sea pressure of 0 dbar), (1) pt_SA_SA, the second derivative with respect to Absolute Salinity at

constant Conservative Temperature,
  1. pt_SA_CT, the derivative with respect to Conservative Temperature and Absolute Salinity, and
  2. pt_CT_CT, the second derivative with respect to Conservative Temperature at constant Absolute Salinity.
Parameters:

SA : array_like

Absolute salinity [g kg -1]

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

Returns:

pt_SA_SA : array_like

The second derivative of potential temperature (the regular potential temperature which has reference sea pressure of 0 dbar) with respect to Absolute Salinity at constant Conservative Temperature. [K/((g/kg)^2)]

pt_SA_CT : array_like

The derivative of potential temperature with respect to Absolute Salinity and Conservative Temperature. [1/(g/kg)]

pt_CT_CT : array_like

The second derivative of potential temperature (the regular one with p_ref = 0 dbar) with respect to Conservative Temperature at constant SA. [1/K]

References

[R108]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. (A.12.9) and (A.12.10).
[R109]McDougall T.J., P.M. Barker, R. Feistel and D.R. Jackett, 2011: A computationally efficient 48-term expression for the density of seawater in terms of Conservative Temperature, and related properties of seawater. To be submitted to Ocean Science Discussions.

Examples

TODO

gsw.gibbs.earth module

gsw.gibbs.earth.distance(lon, lat, p=0)[source]

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)
\[ \begin{align}\begin{aligned}a = \sin^2(\delta ext{lat}/2) + \cos( ext{lat}_1) \cos( ext{lat}_2) \sin^2(\delta ext{lon}/2)\\c = 2 imes ext{atan2}(\sqrt{a}, \sqrt{(1-a)})\\d = R imes c\end{aligned}\end{align} \]
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]

Notes

z is height and is negative in the oceanographic.

Distances are probably good to better than 1% of the “true” distance on the ellipsoidal earth.

References

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

Examples

>>> import gsw
>>> lon = [159, 220]
>>> lat = [-35, 35]
>>> gsw.distance(lon, lat)
array([[ 10030974.652916]])
>>> p = [200, 1000]
>>> gsw.distance(lon, lat, p)
array([[ 10030661.63927133]])
>>> p = [[200], [1000]]
>>> gsw.distance(lon, lat, p)
array([[ 10030661.63927133],
       [ 10029412.58932955]])
gsw.gibbs.earth.f(lat)[source]
Calculates the Coriolis parameter (f) defined by:
f = 2*omega*sin(lat)
where,
omega = 7.292115e-5 (Groten, 2004) [radians s -1]
Parameters:

lat : array_like

latitude [degrees north]

Returns:

f : array_like

Coriolis paramter [s -1]

References

[R111]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.
[R112]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.
gsw.gibbs.earth.grav(lat, p=0)[source]

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]

Notes

In the ocean z is negative.

References

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

Examples

>>> import 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

gsw.gibbs.freezing module

gsw.gibbs.freezing.brineSA_CT(CT, p, saturation_fraction=1)[source]

Calculates the Absolute Salinity of seawater at the freezing temperature. That is, the output is the Absolute Salinity of seawater, with the fraction saturation_fraction of dissolved air, that is in equilibrium with ice at Conservative Temperature CT and pressure p. If the input values are such that there is no positive value of Absolute Salinity for which seawater is frozen, the output, brineSA_CT, is put equal to -99.

Parameters:

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

p : array_like

sea pressure [dbar]

saturation_fraction : fraction between 0, 1. The saturation fraction of

dissolved air in seawater. Default is 0 or completely saturated.

Returns:

brine_SA_CT : array_like

Absolute Salinity of seawater when it freezes [ g/kg ]

References

[R116]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.33.

Examples

TODO

gsw.gibbs.freezing.brineSA_t(t, p, saturation_fraction=1)[source]

Calculates the Absolute Salinity of seawater at the freezing temperature. That is, the output is the Absolute Salinity of seawater, with the fraction saturation_fraction of dissolved air, that is in equilibrium with ice at in-situ temperature t and pressure p. If the input values are such that there is no positive value of Absolute Salinity for which seawater is frozen, the output, brineSA_t, is put equal to -99.

Parameters:

t : array_like

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

p : array_like

sea pressure [dbar]

saturation_fraction : fraction between 0, 1. The saturation fraction of

dissolved air in seawater. Default is 0 or completely saturated.

Returns:

brine_SA_t : array_like

Absolute Salinity of seawater when it freezes [ g/kg ]

References

[R117]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.33.

Examples

TODO

gsw.gibbs.freezing.CT_freezing(SA, p, saturation_fraction=1)[source]

Calculates the Conservative Temperature at which seawater freezes.

Parameters:

SA : array_like

Absolute Salinity [g/kg]

p : array_like

sea pressure [dbar]

saturation_fraction : fraction between 0, 1. The saturation fraction of

dissolved air in seawater. Default is 0 or completely saturated.

Returns:

CT_freezing : array_like

Conservative Temperature at freezing of seawater [\(^\circ\) C (ITS-90)]

References

[R118]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.33 and 3.34.

Examples

TODO

gsw.gibbs.freezing.t_freezing(SA, p, saturation_fraction=1)[source]

Calculates the in-situ temperature at which seawater freezes.

Parameters:

SA : array_like

Absolute Salinity [g/kg]

p : array_like

sea pressure [dbar]

saturation_fraction : fraction between 0, 1. The saturation fraction of

dissolved air in seawater. Default is 0 or completely saturated.

Returns:

t_freezing : array_like

in-situ temperature at which seawater freezes [\(^\circ\) C (ITS-90)]

References

[R119]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.33 and 3.34.

Examples

TODO

gsw.gibbs.geostrophic module

gsw.gibbs.geostrophic_48 module

gsw.gibbs.isobaric module

gsw.gibbs.isobaric.latentheat_evap_t(SA, t)[source]

Calculates latent heat, or enthalpy, of evaporation at p = 0 (the surface). It is defined as a function of Absolute Salinity, SA, and in-situ temperature, t, and is valid in the ranges 0 < SA < 40 g/kg and 0 < CT < 42 deg C. The errors range between -0.4 and 0.6 J/kg.

Parameters:

SA : array_like

Absolute salinity [g kg -1]

t : array_like

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

Returns:

latentheat_evap_t : array_like

latent heat of evaporation [J kg -1]

References

[R120]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.39.

gsw.gibbs.library module

gsw.gibbs.library.deltaSA_atlas(p, lon, lat)[source]

Absolute Salinity anomaly atlas value (excluding the Baltic Sea). Calculates the Absolute Salinity anomaly atlas value, SA - SR, in the open ocean by spatially interpolating the global reference data set of deltaSA_atlas to the location of the seawater sample. This function uses version 3.0 of the deltaSA_ref look up table.

Parameters:

p : array_like

sea pressure (absolute pressure - 10.1325 dbar) [dbar]

lon : array_like

decimal degrees east (will be treated modulo 360)

lat : array_like

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

Returns:

deltaSA_atlas : masked array; masked where no nearby ocean is found in data

Absolute Salinity anomaly atlas value [g/kg]

Notes

The Absolute Salinity anomaly atlas value in the Baltic Sea is evaluated separately, since it is a function of Practical Salinity, not of space. The present function returns a deltaSA_atlas of zero for data in the Baltic Sea. The correct way of calculating Absolute Salinity in the Baltic Sea is by calling SA_from_SP. The mask 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.

The algorithm is taken from the Matlab implementation of the references, but the numpy implementation here differs substantially from the Matlab implementation.

References

[R121]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.
[R122]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
gsw.gibbs.library.enthalpy_SSO_0(p)[source]

Enthalpy at SSO and CT(T=0) (75-term equation)

This function calculates enthalpy at the Standard Ocean Salinty, SSO, and at a Conservative Temperature of zero degrees C, as a function of pressure, p, in dbar, using a streamlined version of the 75-term computationally-efficient expression for specific volume, that is, a streamlined version of the code “enthalpy(SA,CT,p)”.

Parameters:

p : array_like

pressure [dbar]

Returns:

enthalpy_SSO_0 : array_like

enthalpy at (SSO, CT=0, p) [J kg -1]

References

[R123]Roquet, F., G. Madec, T.J. McDougall, P.M. Barker, 2015: Accurate polynomial expressions for the density and specifc volume of seawater using the TEOS-10 standard. Ocean Modelling.

Examples

>>> import gsw
>>> p = np.array([10, 50, 125, 250, 600, 1000])
>>> gsw.library.enthalpy_SSO_0(p)
array([   97.26388583,   486.27439853,  1215.47517122,  2430.24907325,
        5827.90879421,  9704.32030926])
gsw.gibbs.library.enthalpy_SSO_0_p(p)[source]

This function calculates enthalpy at the Standard Ocean Salinty, SSO, and at a Conservative Temperature of zero degrees C, as a function of pressure, p, in dbar, using a streamlined version of the 48-term CT version of the Gibbs function, that is, a streamlined version of the code “enthalpy(SA,CT,p).

References

[R124]McDougall T.J., P.M. Barker, R. Feistel and D.R. Jackett, 2013: A computationally efficient 48-term expression for the density of seawater in terms of Conservative Temperature, and related properties of seawater. To be submitted to J. Atm. Ocean. Technol., xx, yyy-zzz.

Examples

>>> import gsw
>>> p = np.array([10, 50, 125, 250, 600, 1000])
>>> gsw.library.enthalpy_SSO_0_p(p)
array([   97.26388276,   486.27439004,  1215.47518168,  2430.24919716,
        5827.90973888,  9704.32296903])
gsw.gibbs.library.entropy_part(SA, t, p)[source]

Calculates entropy, except that it does not evaluate any terms that are functions of Absolute Salinity alone.

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_part : array_like

entropy minus the terms that due to SA alone [J kg -1 K -1]

Notes

By not calculating these terms, which are a function only of Absolute Salinity, several unnecessary computations are avoided (including saving the computation of a natural logarithm). These terms are a necessary part of entropy, but are not needed when calculating potential temperature from in situ temperature.

gsw.gibbs.library.entropy_part_zerop(SA, pt0)[source]

Calculates entropy at a sea surface (p = 0 dbar), except that it does not evaluate any terms that are functions of Absolute Salinity alone.

Parameters:

SA : array_like

Absolute salinity [g kg -1]

pt0 : array_like

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

Returns:

entropy_part_zerop : array_like

[J kg -1 K -1]

Notes

By not calculating these terms, which are a function only of Absolute Salinity, several unnecessary computations are avoided (including saving the computation of a natural logarithm). These terms are a necessary part of entropy, but are not needed when calculating potential temperature from in situ temperature.

gsw.gibbs.library.Fdelta(p, lon, lat)[source]

Fdelta from the Absolute Salinity Anomaly Ratio (SAAR): Fdelta = (1 + r1)SAAR/(1 - r1*SAAR)

= (SA/Sstar) - 1

with r1 being the constant 0.35 based on the work of Pawlowicz et al. (2011). Note that since SAAR is everywhere less than 0.001 in the global ocean, Fdelta is only slightly different to 1.35*SAAR.

Parameters:

p : array_like

pressure [dbar]

lon : array_like

decimal degrees east (will be treated modulo 360)

lat : array_like

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

Returns:

Fdelta : masked array; masked where no nearby ocean is found in data

Ratio of SA to Sstar, minus 1 [unitless]

Notes

The mask 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

[R125]IOC, SCOR and IAPSO, 2010: The international thermodynamic equation of seawater - 2010: Calculation and use of thermodynamic properties. Intergovernmental Oceanographic Commission, Manuals and Guides No. 56, UNESCO (English), 196 pp. Available from http://www.TEOS-10.org See section 2.5 and appendices A.4 and A.5 of this TEOS-10 Manual.
[R126]McDougall, T.J., D.R. Jackett, F.J. Millero, R. Pawlowicz and P.M. Barker, 2012: A global algorithm for estimating Absolute Salinity. Ocean Science, 8, 1123-1134. http://www.ocean-sci.net/8/1123/2012/os-8-1123-2012.pdf
[R127]Pawlowicz, R., D.G. Wright and F.J. Millero, 2011; The effects of biogeochemical processes on oceanic conductivity/salinty/density relationships and the characterization of real seawater. Ocean Science, 7, 363-387. http://www.ocean-sci.net/7/363/2011/os-7-363-2011.pdf
gsw.gibbs.library.gibbs(ns, nt, npr, SA, t, p)[source]

Calculates specific Gibbs energy and its derivatives up to order 2 for seawater. The Gibbs function approach allows the calculation of internal energy, entropy, enthalpy, potential enthalpy and the chemical potentials of seawater as well as the freezing temperature, and the latent heats of freezing and of evaporation. These quantities were not available from EOS-80 but are essential for the accurate accounting of heat in the ocean and for the consistent and accurate treatment of air-sea and ice-sea heat fluxes.

Parameters:

ns : int

order of SA derivative [0, 1 or 2 ]

nt : int

order of t derivative [0, 1 or 2 ]

npr : int

order of p derivative [0, 1 or 2 ]

SA : array_like

Absolute salinity [g kg -1]

t : array_like

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

p : array_like

pressure [dbar]

Returns:

gibbs : array_like

Specific Gibbs energy or its derivatives. Gibbs energy (ns=nt=npr=0) has units of: [J kg -1] Absolute Salinity derivatives are output in units of: [(J kg -1) (g kg -1) -ns] Temperature derivatives are output in units of: [(J kg -1) K -nt] Pressure derivatives are output in units of: [(J kg -1) Pa -npr] The mixed derivatives are output in units of: [(J kg -1) (g kg -1) -ns K -nt Pa -npr]

Notes

The Gibbs function for seawater is that of TEOS-10 (IOC et al., 2010), being the sum of IAPWS-08 for the saline part and IAPWS-09 for the pure water part. These IAPWS releases are the officially blessed IAPWS descriptions of Feistel (2008) and the pure water part of Feistel (2003). Absolute Salinity, SA, in all of the GSW routines is expressed on the Reference-Composition Salinity Scale of 2008 (RCSS-08) of Millero et al. (2008). The derivatives are taken with respect to pressure in Pa, not withstanding that the pressure input into this routine is in dbar.

References

[R128]Feistel, R., 2003: A new extended Gibbs thermodynamic potential of seawater Progr. Oceanogr., 58, 43-114.
[R129]Feistel, R., 2008: A Gibbs function for seawater thermodynamics for -6 to 80 \(^\circ\) C and salinity up to 120 g kg -1, Deep-Sea Res. I, 55, 1639-1671.
[R130]IAPWS, 2008: Release on the IAPWS Formulation 2008 for the Thermodynamic Properties of Seawater. The International Association for the Properties of Water and Steam. Berlin, Germany, September 2008, available from http://www.iapws.org. This Release is referred to as IAPWS-08.
[R131]IAPWS, 2009: Supplementary Release on a Computationally Efficient Thermodynamic Formulation for Liquid Water for Oceanographic Use. The International Association for the Properties of Water and Steam. Doorwerth, The Netherlands, September 2009, available from http://www.iapws.org. This Release is referred to as IAPWS-09.
[R132]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.6 and appendices A.6, G and H.
[R133]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.
gsw.gibbs.library.gibbs_pt0_pt0(SA, pt0)[source]

Calculates the second derivative of the specific Gibbs function with respect to temperature at zero sea pressure or _gibbs(0,2,0,SA,t,0).

Parameters:

SA : array_like

Absolute salinity [g kg -1]

pt0 : array_like

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

Returns:

gibbs_pt0_pt0 : array_like

TODO: write the eq. for the second derivative of the specific Gibbs function. FIXME: [units]

Notes

This library function is called by both “pt_from_CT(SA,CT)” and “pt0_from_t(SA,t,p)”.

gsw.gibbs.library.Hill_ratio_at_SP2(t)[source]
USAGE:
Hill_ratio = Hill_ratio_at_SP2(t)
DESCRIPTION:
Calculates the Hill ratio, which is the adjustment needed to apply for Practical Salinities smaller than 2. This ratio is defined at a Practical Salinity = 2 and in-situ temperature, t using PSS-78. The Hill ratio is the ratio of 2 to the output of the Hill et al. (1986) formula for Practical Salinity at the conductivity ratio, Rt, at which Practical Salinity on the PSS-78 scale is exactly 2.
INPUT:
t = in-situ temperature (ITS-90) [ deg C ]
OUTPUT:
Hill_ratio = Hill ratio at SP of 2 [ unitless ]
AUTHOR:
Trevor McDougall and Paul Barker

VERSION NUMBER: 3.0 (26th March, 2011)

gsw.gibbs.library.infunnel(SA, CT, p)[source]

Oceanographic funnel check for the 25-term equation

Parameters:

SA : array_like

Absolute Salinity [g/kg]

CT : array_like

Conservative Temperature [°C]

p : array_like

sea pressure [dbar]

(ie. absolute pressure - 10.1325 dbar)

Returns:

in_funnel : boolean ndarray or scalar

True, if SA, CT and p are inside the “funnel” False, if SA, CT and p are outside the “funnel”,

or one of the values are NaN or masked

Note. The term “funnel” describes the range of SA, CT and p over which

the error in the fit of the computationally-efficient 25-term

expression for density in terms of SA, CT and p was calculated

(McDougall et al., 2010).

gsw.gibbs.library.interp_ref_cast(spycnl, A='gn')[source]

Translation of: [SA_iref_cast, CT_iref_cast, p_iref_cast] = interp_ref_cast(spycnl, A) interp_ref_cast linear interpolation of the reference cast ========================================================================== This function interpolates the reference cast with respect to the interpolating variable “spycnl”. This reference cast is at the location 188E,4N from the reference data set which underlies the Jackett & McDougall (1997) Neutral Density computer code. This function finds the values of SA, CT and p on this reference cast which correspond to the value of isopycnal which is passed to this function from the function “geo_strf_isopycnal_CT”. The isopycnal could be either gamma_n or sigma_2. If A is set to any of the following ‘s2’,’S2’,’sigma2’,’sigma_2’ the interpolation will take place in sigma 2 space, any other input will result in the programme working in gamma_n space. VERSION NUMBER: 3.0 (14th April, 2011) REFERENCE: Jackett, D. R. and T. J. McDougall, 1997: A neutral density variable for the world<92>s oceans. Journal of Physical Oceanography, 27, 237-263. FIXME? Do we need argument checking here to handle masked arrays, etc.? I suspect not, since I don’t think this is intended to be user-callable, but is instead used internally by user-callable functions. Note: The v3.03 matlab code is incorrectly using approximate numbers for the gamma_n case, even when the sigma_2 case is in effect. That bug is fixed here.

gsw.gibbs.library.interp_SA_CT(SA, CT, p, p_i)[source]

function [SA_i, CT_i] = interp_SA_CT(SA,CT,p,p_i) interp_SA_CT linear interpolation to p_i on a cast ========================================================================== This function interpolates the cast with respect to the interpolating variable p. This function finds the values of SA, CT at p_i on this cast.

gsw.gibbs.library.SAAR(p, lon, lat)[source]

Absolute Salinity Anomaly Ratio (excluding the Baltic Sea). Calculates the Absolute Salinity Anomaly Ratio, SAAR, in the open ocean by spatially interpolating the global reference data set of SAAR to the location of the seawater sample. This function uses version 3.0 of the SAAR look up table.

Parameters:

p : array_like

pressure [dbar]

lon : array_like

decimal degrees east (will be treated modulo 360)

lat : array_like

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

Returns:

SAAR : array

Absolute Salinity Anomaly Ratio [unitless]

in_ocean : boolean array

Notes

The Absolute Salinity Anomaly Ratio in the Baltic Sea is evaluated separately, since it is a function of Practical Salinity, not of space. The present function returns a SAAR of zero for data in the Baltic Sea. The correct way of calculating Absolute Salinity in the Baltic Sea is by calling SA_from_SP. 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.

The algorithm is taken from the matlab implementation of the references, but the numpy implementation here differs substantially from the matlab implementation.

References

[R134]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.
[R135]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
gsw.gibbs.library.SA_from_SP_Baltic(SP, lon, lat)[source]

Computes absolute salinity from practical in the Baltic Sea.

Parameters:

SP : array_like or masked array

Practical salinity (PSS-78)

lon, lat : array_like or masked arrays

geographical position

Returns:

SA : masked array, at least 1D

Absolute salinity [g/kg] masked where inputs are masked or position outside the Baltic

gsw.gibbs.library.specvol_SSO_0_p(p)[source]

This function calculates specific volume at the Standard Ocean Salinity, SSO, and at a Conservative Temperature of zero degrees C, as a function of pressure, p, in dbar, using a streamlined version of the 48-term CT version of specific volume, that is, a streamlined version of the code “specvol(SA, CT, p)”.

gsw.gibbs.library.SP_from_SA_Baltic(SA, lon, lat)[source]

Calculates Practical Salinity (SP) for the Baltic Sea, from a value computed analytically from Absolute Salinity.

Parameters:

SA : array_like

Absolute salinity [g kg -1]

lon : array_like

decimal degrees east [0..+360]

lat : array_like

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

Returns:

SP_baltic : array_like

salinity [psu (PSS-78)], unitless

See also

SP_from_SA, SP_from_Sstar

Notes

This program will only produce Practical Salinity values for the Baltic Sea.

References

[R136]Feistel, R., S. Weinreben, H. Wolf, S. Seitz, P. Spitzer, B. Adel, G. Nausch, B. Schneider and D. G. Wright, 2010c: Density and Absolute Salinity of the Baltic Sea 2006-2009. Ocean Science, 6, 3-24. http://www.ocean-sci.net/6/3/2010/os-6-3-2010.pdf
[R137]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.
[R138]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

Examples

TODO

gsw.gibbs.neutral_nonlinear_48 module

gsw.gibbs.practical_salinity module

gsw.gibbs.practical_salinity.C_from_SP(SP, t, p)[source]

Calculates conductivity, C, from (SP, t, p) using PSS-78 in the range 2 < SP < 42. If the input Practical Salinity is less than 2 then a modified form of the Hill et al. (1986) fomula is used for Practical Salinity. The modification of the Hill et al. (1986) expression is to ensure that it is exactly consistent with PSS-78 at SP = 2.

The conductivity ratio returned by this function is consistent with the input value of Practical Salinity, SP, to 2x10^-14 psu over the full range of input parameters (from pure fresh water up to SP = 42 psu). This error of 2x10^-14 psu is machine precision at typical seawater salinities. This accuracy is achieved by having four different polynomials for the starting value of Rtx (the square root of Rt) in four different ranges of SP, and by using one and a half iterations of a computationally efficient modified Newton-Raphson technique to find the root of the equation.

Parameters:

SP : array

Practical Salinity [psu (PSS-78), unitless]

t : array

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

p : array

sea pressure [dbar] (i.e. absolute pressure - 10.1325 dbar)

Returns:

C : array

conductivity [mS cm -1]

Notes

Note that strictly speaking PSS-78 (Unesco, 1983) defines Practical Salinity in terms of the conductivity ratio, R, without actually specifying the value of C(35,15,0) (which we currently take to be 42.9140 mS/cm).

References

[R139]Hill, K.D., T.M. Dauphinee and D.J. Woods, 1986: The extension of the Practical Salinity Scale 1978 to low salinities. IEEE J. Oceanic Eng., OE-11, 1, 109 - 112.
[R140]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 E.
[R141]Unesco, 1983: Algorithms for computation of fundamental properties of seawater. Unesco Technical Papers in Marine Science, 44, 53 pp.

Examples

TODO

gsw.gibbs.practical_salinity.R_from_SP(SP, t, p)[source]

Calculates conductivity ratio from (SP, t, p) using PSS-78 in the range 2 < SP < 42. If the input Practical Salinity is less than 2 then a modified form of the Hill et al. (1986) formula is used for Practical Salinity. The modification of the Hill et al. (1986) expression is to ensure that it is exactly consistent with PSS-78 at SP = 2.

The conductivity ratio returned by this function is consistent with the input value of Practical Salinity, SP, to 2x10^-14 psu over the full range of input parameters (from pure fresh water up to SP = 42 psu). This error of 2x10^-14 psu is machine precision at typical seawater salinities. This accuracy is achieved by having four different polynomials for the starting value of Rtx (the square root of Rt) in four different ranges of SP, and by using one and a half iterations of a computationally efficient modified Newton-Raphson technique to find the root of the equation.

Parameters:

SP : array

Practical Salinity [psu (PSS-78), unitless]

t : array_like

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

p : array

sea pressure [dbar] (i.e. absolute pressure - 10.1325 dbar)

Returns:

R : array_like

conductivity ratio [unitless]

Notes

Strictly speaking PSS-78 (Unesco, 1983) defines Practical Salinity in terms of the conductivity ratio, R, without actually specifying the value of C(35, 15, 0) (which we currently take to be 42.9140 mS cm^-1. Culkin and Smith, 1980).

References

[R142]Culkin and Smith, 1980: Determination of the Concentration of Potassium Chloride Solution Having the Same Electrical Conductivity, at 15C and Infinite Frequency, as Standard Seawater of Salinity 35.0000 (Chlorinity 19.37394), IEEE J. Oceanic Eng, 5, 22-23.
[R143]Hill, K.D., T.M. Dauphinee & D.J. Woods, 1986: The extension of the Practical Salinity Scale 1978 to low salinities. IEEE J. Oceanic Eng., 11, 109 - 112.
[R144]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. Appendix E.
[R145]Unesco, 1983: Algorithms for computation of fundamental properties of seawater. Unesco Technical Papers in Marine Science, 44, 53 pp.

Examples

TODO

gsw.gibbs.practical_salinity.SK_from_SP(SP)[source]

Calculates Knudsen Salinity from Practical Salinity.

Parameters:

SP : array

Practical Salinity [psu (PSS-78), unitless]

Returns:

SK : array_like

Knudsen Salinity [parts per thousand, ppt]

References

[R146]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.3.

Examples

TODO

gsw.gibbs.practical_salinity.SP_from_C(C, t, p)[source]

Calculates Practical Salinity, SP, from conductivity, C, primarily using the PSS-78 algorithm. Note that the PSS-78 algorithm for Practical Salinity is only valid in the range 2 < SP < 42. If the PSS-78 algorithm produces a Practical Salinity that is less than 2 then the Practical Salinity is recalculated with a modified form of the Hill et al. (1986) formula. The modification of the Hill et al. (1986) expression is to ensure that it is exactly consistent with PSS-78 at SP = 2. Note that the input values of conductivity need to be in units of mS/cm (not S/m).

Parameters:

C : array

conductivity [mS cm -1]

t : array

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

p : array

sea pressure [dbar] (i.e. absolute pressure - 10.1325 dbar)

Returns:

SP : array

Practical Salinity [psu (PSS-78), unitless]

References

[R147]Culkin and Smith, 1980: Determination of the Concentration of Potassium Chloride Solution Having the Same Electrical Conductivity, at 15C and Infinite Frequency, as Standard Seawater of Salinity 35.0000 (Chlorinity 19.37394), IEEE J. Oceanic Eng, 5, 22-23.
[R148]Hill, K.D., T.M. Dauphinee & D.J. Woods, 1986: The extension of the Practical Salinity Scale 1978 to low salinities. IEEE J. Oceanic Eng., 11, 109 - 112.
[R149]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. Appendix E.
[R150]Unesco, 1983: Algorithms for computation of fundamental properties of seawater. Unesco Technical Papers in Marine Science, 44, 53 pp.

Examples

TODO

gsw.gibbs.practical_salinity.SP_from_R(R, t, p)[source]

Calculates Practical Salinity, SP, from the conductivity ratio, R, primarily using the PSS-78 algorithm. Note that the PSS-78 algorithm for Practical Salinity is only valid in the range 2 < SP < 42. If the PSS-78 algorithm produces a Practical Salinity that is less than 2 then the Practical Salinity is recalculated with a modified form of the Hill et al. (1986) formula. The modification of the Hill et al. (1986) expression are to ensure that it is exactly consistent with PSS-78 at SP = 2.

Parameters:

R : array_like

conductivity ratio [unitless]

t : array_like

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

p : array

sea pressure [dbar] (i.e. absolute pressure - 10.1325 dbar)

Returns:

SP : array

Practical Salinity [psu (PSS-78), unitless]

References

[R151]Culkin and Smith, 1980: Determination of the Concentration of Potassium Chloride Solution Having the Same Electrical Conductivity, at 15C and Infinite Frequency, as Standard Seawater of Salinity 35.0000 (Chlorinity 19.37394), IEEE J. Oceanic Eng, 5, 22-23.
[R152]Hill, K.D., T.M. Dauphinee & D.J. Woods, 1986: The extension of the Practical Salinity Scale 1978 to low salinities. IEEE J. Oceanic Eng., 11, 109 - 112.
[R153]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. Appendix E.
[R154]Unesco, 1983: Algorithms for computation of fundamental properties of seawater. Unesco Technical Papers in Marine Science, 44, 53 pp.

Examples

TODO

gsw.gibbs.practical_salinity.SP_from_SK(SK)[source]

Calculates Practical Salinity from Knudsen Salinity.

Parameters:

SK : array_like

Knudsen Salinity [parts per thousand, ppt]

Returns:

SP : array

Practical Salinity [psu (PSS-78), unitless]

References

[R155]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.3.

Examples

TODO

gsw.gibbs.practical_salinity.SP_salinometer(Rt, t)[source]

Calculates Practical Salinity SP from a salinometer, primarily using the PSS-78 algorithm. Note that the PSS-78 algorithm for Practical Salinity is only valid in the range 2 < SP < 42. If the PSS-78 algorithm produces a Practical Salinity that is less than 2 then the Practical Salinity is recalculated with a modified form of the Hill et al. (1986) formula. The modification of the Hill et al. (1986) expression is to ensure that it is exactly consistent with PSS-78 at SP = 2.

A laboratory salinometer has the ratio of conductivities, Rt, as an output, and the present function uses this conductivity ratio and the temperature t of the salinometer bath as the two input variables.

Parameters:

Rt : array

C(SP,t_68,0)/C(SP=35,t_68,0) [unitless] conductivity ratio :math:`R =

rac{C(S, t_68, 0)}{C(35, 15(IPTS-68),0)} [unitless]

t : array

Temperature of the bath of the salinometer [\(^\circ\) C (ITS-90)]

Returns:

SP : array

Practical Salinity [psu (PSS-78), unitless]

Examples

TODO

gsw.gibbs.steric module

gsw.gibbs.thermodynamics_from_t module

gsw.gibbs.thermodynamics_from_t.adiabatic_lapse_rate_from_t(SA, t, p)[source]

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

sea pressure [dbar]

SA & t need to have the same dimensions. p may have dimensions 1x1 or Mx1 or 1xN or MxN, where SA & t are MxN.

Returns:

adiabatic_lapse_rate : array_like

adiabatic lapse rate [K Pa:sup:-1] Note. The output is in unit of degress Celsius per Pa, (or equivilently K/Pa) not in units of K/dbar.

References

[R159]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.
gsw.gibbs.thermodynamics_from_t.adiabatic_lapse_rate_from_CT(SA, CT, p)[source]

Calculates the adiabatic lapse rate of sea water from Conservative Temperature.

Parameters:

SA : array_like

Absolute salinity [g kg -1]

CT : array_like

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

p : array_like

sea pressure [dbar]

Returns:

adiabatic_lapse_rate : array_like

adiabatic lapse rate [K Pa:sup:-1] Note. The output is in unit of degress Celsius per Pa, (or equivilently K/Pa) not in units of K/dbar.

References

[R160]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.
gsw.gibbs.thermodynamics_from_t.alpha_wrt_CT_t_exact(SA, t, p)[source]

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]

References

[R161]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).

Examples

>>> import 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_t_exact(SA, t, p)
array([ 0.00032471,  0.00032272,  0.00028118,  0.00017314,  0.00014627,
        0.00012943])
gsw.gibbs.thermodynamics_from_t.alpha_wrt_pt_t_exact(SA, t, p)[source]

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]

References

[R162]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).

Examples

>>> import 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_t_exact(SA, t, p)
array([ 0.00032562,  0.00032355,  0.00028164,  0.00017314,  0.00014623,
        0.00012936])
gsw.gibbs.thermodynamics_from_t.alpha_wrt_t_exact(SA, t, p)[source]

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]

References

[R163]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)
[R164]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.

Examples

>>> import 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_exact(SA, t, p)
array([ 0.0003256 ,  0.00032345,  0.00028141,  0.00017283,  0.00014557,
        0.00012836])
gsw.gibbs.thermodynamics_from_t.beta_const_CT_t_exact(SA, t, p)[source]

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

References

[R165]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)

Examples

>>> import 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_t_exact(SA, t, p)
array([ 0.00071749,  0.00071765,  0.00072622,  0.00075051,  0.00075506,
        0.00075707])
gsw.gibbs.thermodynamics_from_t.beta_const_pt_t_exact(SA, t, p)[source]

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]

References

[R166]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)

Examples

>>> import 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_t_exact(SA, t, p)
array([ 0.00073112,  0.00073106,  0.00073599,  0.00075375,  0.00075712,
        0.00075843])
gsw.gibbs.thermodynamics_from_t.beta_const_t_exact(SA, t, p)[source]

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]

References

[R167]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)

Examples

>>> import 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_exact(SA, t, p)
array([ 0.00073112,  0.00073107,  0.00073602,  0.00075381,  0.00075726,
        0.00075865])
gsw.gibbs.thermodynamics_from_t.chem_potential_relative_t_exact(SA, t, p)[source]

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_relative : array_like

relative chemical potential [J g:sup:-1]

References

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

Examples

>>> import 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_t_exact(SA, t, p)
array([ 79.4254481 ,  79.25989214,  74.69154859,  65.64063719,
        61.22685656,  57.21298557])
gsw.gibbs.thermodynamics_from_t.chem_potential_salt_t_exact(SA, t, p)[source]

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]

References

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

Examples

>>> import 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_t_exact(SA, t, p)
array([ 70.87988696,  71.25180659,  69.58756845,  65.00656941,
        64.56242337,  64.76842002])
gsw.gibbs.thermodynamics_from_t.chem_potential_water_t_exact(SA, t, p)[source]

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/g]

References

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

Examples

>>> import 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_t_exact(SA, t, p)
array([-8.54556115, -8.00808555, -5.10398014, -0.63406778,  3.3355668 ,
        7.55543445])
gsw.gibbs.thermodynamics_from_t.cp_t_exact(SA, t, p)[source]

Calculates the isobaric 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:

cp_t_exact : array_like

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

References

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

Examples

>>> import 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_t_exact(SA, t, p)
array([ 4002.88800396,  4000.98028393,  3995.54646889,  3985.07676902,
        3973.59384348,  3960.18408479])
gsw.gibbs.thermodynamics_from_t.dynamic_enthalpy_t_exact(SA, t, p)[source]

Calculates the dynamic enthalpy of seawater from Absolute Salinity, in situ temperature and pressure. Dynamic enthalpy was defined by Young (2010) as the difference between enthalpy and potential enthalpy. Note that this function uses the full TEOS-10 Gibbs function (i.e. the sum of the IAPWS-09 and IAPWS-08 Gibbs functions, see the TEOS-10 Manual, IOC et al. (2010)).

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:

dynamic_enthalpy_t_exact : array_like

dynamic enthalpy [J -1]

References

[R172]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.
[R173]Young, W.R., 2010: Dynamic enthalpy, Conservative Temperature, and the seawater. Boussinesq approximation. Journal of Physical Oceanography, 40, 394-400.

Examples

TODO

gsw.gibbs.thermodynamics_from_t.enthalpy_t_exact(SA, t, p)[source]

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) \]

rac{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]

Examples

>>> import 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([ 115006.17814391,  113989.54266878,   92276.22539967,
         43390.30953352,   33323.12790127,   27290.12626901])
gsw.gibbs.thermodynamics_from_t.entropy_t_exact(SA, t, p)[source]

Calculates specific entropy of seawater. The specific entropy of seawater \(\eta\) is given by:

\[\eta(SA, t, p) = -g_T = \]

rac{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]

Examples

>>> import 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_t_exact(SA, t, p)
array([ 400.38942528,  395.43817843,  319.8664982 ,  146.79088159,
         98.64734087,   62.79150873])
gsw.gibbs.thermodynamics_from_t.Helmholtz_energy_t_exact(SA, t, p)[source]

Calculates the Helmholtz energy of seawater. The specific Helmholtz energy of seawater \(f\) is given by:

\[f(SA, t, p) = g - (p + P_0) \]
u =
g - (p + P_0)

rac{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]

Examples

>>> import 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_t_exact(SA, t, p)
array([-5985.58288209, -5830.81845224, -3806.96617841,  -877.66369421,
        -462.17033905,  -245.50407205])
gsw.gibbs.thermodynamics_from_t.internal_energy_t_exact(SA, t, p)[source]

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)\]
u =
g - (T_0 + t)
rac{partial g}{partial T}Big|_{SA,p} -
(p + P_0)

rac{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]

Examples

>>> import 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_t_exact(SA, t, p)
array([ 114906.23847309,  113426.57417062,   90860.81858842,
         40724.34005719,   27162.66600185,   17182.50522667])
gsw.gibbs.thermodynamics_from_t.isochoric_heat_cap_t_exact(SA, t, p)[source]

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]

References

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

Examples

>>> import 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_t_exact(SA, t, p)
array([ 3928.13708702,  3927.27381633,  3941.36418525,  3966.26126146,
        3960.50903222,  3950.13901342])
gsw.gibbs.thermodynamics_from_t.kappa_const_t_exact(SA, t, p)[source]

Calculates isothermal compressibility of seawater at constant in situ temperature.

\[\kappa^t(SA, t, p) =\]

ho^{-1} rac{partial ho}{partial P}Big|_{SA,T} =

u^{-1} rac{partial u}{partial P}Big|_{SA,T} =

rac{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]

Notes

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

Examples

>>> import 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_exact(SA, t, p)
array([  4.19071646e-10,   4.18743202e-10,   4.22265764e-10,
         4.37735100e-10,   4.40373818e-10,   4.41156577e-10])
gsw.gibbs.thermodynamics_from_t.kappa_t_exact(SA, t, p)[source]

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) =\]

ho^{-1} rac{partial ho}{partial P}Big|_{SA,eta} =

u^{-1} rac{partial u}{partial P}Big|_{SA,eta} =

ho^{-1} rac{partial ho}{partial P}Big|_{SA, heta} =

u^{-1} rac{partial u}{partial P}Big|_{SA, heta} =

rac{ (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]

Notes

The output is Pascal and not dbar.

Examples

>>> import 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_t_exact(SA, t, p)
array([  4.11245799e-10,   4.11029072e-10,   4.16539558e-10,
         4.35668338e-10,   4.38923693e-10,   4.40037576e-10])
gsw.gibbs.thermodynamics_from_t.osmotic_coefficient_t_exact(SA, t, p)[source]

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]

References

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

Examples

>>> import 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_t_exact(SA,t , p)
array([ 0.90284718,  0.90298624,  0.90238866,  0.89880927,  0.89801054,
        0.89767912])
gsw.gibbs.thermodynamics_from_t.osmotic_pressure_t_exact(SA, t, pw)[source]

Calculates the osmotic pressure of seawater.

Parameters:

SA : array_like

Absolute salinity [g kg -1]

t : array_like

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

pw : array_like

sea pressure of the pure water side [dbar]

Returns:

osmotic_pressure_t_exact : array_like

dynamic osmotic pressure of seawater [dbar]

References

[R182]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.41.

Examples

TODO

gsw.gibbs.thermodynamics_from_t.pot_rho_t_exact(SA, t, p, p_ref=0)[source]

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]

p_ref : int, float, optional

reference pressure, default = 0

Returns:

pot_rho : array_like

potential density [kg m -3]

References

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

Examples

>>> import 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_t_exact(SA, t, p)
array([ 1021.79814581,  1022.05248442,  1023.89358365,  1026.66762112,
        1027.10723087,  1027.40963126])
>>> gsw.pot_rho_t_exact(SA, t, p, p_ref=1000)
array([ 1025.95554512,  1026.21306986,  1028.12563226,  1031.1204547 ,
        1031.63768355,  1032.00240412])
gsw.gibbs.thermodynamics_from_t.rho_t_exact(SA, t, p)[source]

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_t_exact : array_like

in situ density [kg m -3]

References

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

Examples

>>> import 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_t_exact(SA, t, p)
array([ 1021.84017319,  1022.26268993,  1024.42771594,  1027.79020181,
        1029.83771473,  1032.00240412])
gsw.gibbs.thermodynamics_from_t.SA_from_rho_t(rho, t, p)[source]
gsw.gibbs.thermodynamics_from_t.SA_from_rho_t_exact(rho, t, p)[source]
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 :math:`
ho` 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]

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 \(^ imes\) 10 -13 [g kg -1]

Examples

>>> import 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_t_exact(rho, t, p)
array([ 34.71022966,  34.89057683,  35.02332421,  34.84952096,
        34.73824809,  34.73188384])
gsw.gibbs.thermodynamics_from_t.sigma0_pt0_exact(SA, pt0)[source]

Calculates potential density anomaly with reference sea pressure of zero (0) dbar. The temperature input to this function is potential temperature referenced to zero dbar.

Parameters:

SA : array_like

Absolute salinity [g kg -1]

pt0 : array_like

potential temperature [\(^\circ\) C (ITS-90)] with respect to a reference sea pressure of 0 dbar

Returns:

sigma0_pt0_exact : array_like

potential density anomaly [kg m -3] respect to a reference pressure of 0 dbar

References

[R187]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.6.1).

Examples

>>> import gsw
>>> SA = [34.7118, 34.8915, 35.0256, 34.8472, 34.7366, 34.7324]
>>> pt0 = [28.7832, 28.4209, 22.7850, 10.2305, 6.8292, 4.3245]
>>> gsw.sigma0_pt0_exact(SA, pt0)
array([ 21.79814475,  22.05251193,  23.89356369,  26.66762521,
        27.10723499,  27.4096324 ])
gsw.gibbs.thermodynamics_from_t.sound_speed_t_exact(SA, t, p)[source]

Calculates the speed of sound in seawater. The speed of sound in seawater \(c\) is given by: .. math:

c(SA, t, p) = \sqrt{ \partial P  / \partial 
ho |_{SA,eta}} =
sqrt{(
hokappa)^{-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]

Examples

>>> import 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_t_exact(SA, t, p)
array([ 1542.61580359,  1542.70353407,  1530.84497914,  1494.40999692,
        1487.37710252,  1483.93460908])
gsw.gibbs.thermodynamics_from_t.specvol_anom_t_exact(SA, t, p)[source]

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_t_exact : array_like

specific volume anomaly [m 3 kg -1]

References

[R189]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)

Examples

>>> import 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_t_exact(SA, t, p)
array([  6.01044463e-06,   5.78602432e-06,   4.05564999e-06,
         1.42198662e-06,   1.04351837e-06,   7.63964850e-07])
gsw.gibbs.thermodynamics_from_t.specvol_t_exact(SA, t, p)[source]

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]

References

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

Examples

>>> import 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.00097862,  0.00097822,  0.00097616,  0.00097297,  0.00097104,
        0.000969  ])
gsw.gibbs.thermodynamics_from_t.t_from_rho_exact(rho, SA, p)[source]

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

Parameters:

rho : array_like

in situ density [kg m -3]

SA : array_like

Absolute salinity [g kg -1]

p : array_like

pressure [dbar]

Returns:

t : array_like

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

t_multiple : array_like

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

Notes

At low salinities, in brackish water, there are two possible temperatures for a single density. This program will output both valid solutions (t, t_multiple), if there is only one possible solution the second variable will be set to NaN.

References

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

Examples

TODO

gsw.gibbs.thermodynamics_from_t.t_maxdensity_exact(SA, p)[source]

Calculates the in-situ temperature of maximum density of seawater. This function returns the in-situ temperature at which the density of seawater is a maximum, at given Absolute Salinity, SA, and sea pressure, p (in dbar).

Parameters:

SA : array_like

Absolute salinity [g kg -1]

p : array_like

pressure [dbar]

Returns:

t_maxdensity_exact : array_like

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

References

[R192]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.42.

Examples

TODO

gsw.gibbs.water_column_48 module

gsw.gibbs.water_column_48.IPV_vs_fNsquared_ratio(SA, CT, p, p_ref=0)[source]

Calculates the ratio of the vertical gradient of potential density to the vertical gradient of locally-referenced potential density. This ratio is also the ratio of the planetary Isopycnal Potential Vorticity (IPV) to f times N^2, hence the name for this variable, IPV_vs_fNsquared_ratio (see Eqn. (3.20.5) of IOC et al. (2010)). The reference sea pressure, p_ref, of the potential density surface must have a constant value.

IPV_vs_fNsquared_ratio is evaluated at the mid pressure between the individual data points in the vertical. This function uses the computationally-efficient 48-term expression for density in terms of SA, CT and p (McDougall et al., 2011).

Parameters:

SA : array_like

Absolute Salinity [g/kg]

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

p : array_like

sea pressure [dbar]

p_ref : int, float, optional

reference pressure, default = 0

Returns:

IPV_vs_fNsquared_ratio : array_like

The ratio of the vertical gradient of potential density, on the same (M-1)xN grid as p_mid. [unitless] referenced to p_ref, to the vertical gradient of locally- referenced potential density.

p_mid : array_like

Mid pressure between p grid [dbar]

Notes

The 48-term equation has been fitted in a restricted range of parameter space, and is most accurate inside the “oceanographic funnel” described in McDougall et al. (2011). The GSW library function “infunnel(SA, CT, p)” is available to be used if one wants to test if some of one’s data lies outside this “funnel”.

References

[R193]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.20.5).
[R194]McDougall T.J., P.M. Barker, R. Feistel and D.R. Jackett, 2011: A computationally efficient 48-term expression for the density of seawater in terms of Conservative Temperature, and related properties of seawater.

Examples

TODO

gsw.gibbs.water_column_48.Nsquared(SA, CT, p, lat=None)[source]

Calculates the buoyancy frequency squared (N^2) (i.e. the Brunt-Väisälä frequency squared) at the mid pressure from the equation:

  2      2             beta.d(SA) - alpha.d(CT)
N   =  g   rho_local   ------------------------
                                 dP

dP in the above formula is in Pascals

Parameters:

SA : array_like

Absolute Salinity [g/kg]

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

p : array_like

sea pressure [dbar]

lat : array_like, optional

latitude in decimal degrees north [-90..+90] If lat is not supplied, a default gravitational acceleration of 9.7963 m/s^2 (Griffies, 2004) will be used.

Returns:

N2 : array_like

Brunt-Väisälä frequency squared [1 s \(-2\)]

p_mid : array_like

Mid pressure between p grid [dbar]

Notes

This routine uses rho from the computationally efficient 48-term expression for density in terms of SA, CT and p.

The 48-term equation has been fitted in a restricted range of parameter space, and is most accurate inside the “oceanographic funnel” described in IOC et al. (2010). The GSW library function “infunnel(SA, CT, p)” is available to be used if one wants to test if some of one’s data lies outside this “funnel”.

References

[R195]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.10 and Eqn. (3.10.2).
[R196]Griffies, S. M., 2004: Fundamentals of Ocean Climate Models. Princeton, NJ: Princeton University Press, 518 pp + xxxiv.

Examples

TODO

gsw.gibbs.water_column_48.Turner_Rsubrho(SA, CT, p)[source]

Calculates the Turner angle and the Rsubrho as a function of pressure down a vertical water column. These quantities express the relative contributions of the vertical gradients of Conservative Temperature and Absolute Salinity to the vertical stability (the square of the Brunt-Väisälä Frequency squared, N^2). Tu and Rsubrho are evaluated at the mid pressure between the individual data points in the vertical. This function uses computationally-efficient 48-term expression for density in terms of SA, CT and p (McDougall et al., 2011). Note that in the double-diffusive literature, papers concerned with the “diffusive” form of double-diffusive convection often define the stability ratio as the reciprocal of what is defined here as the stability ratio.

Parameters:

SA : array_like

Absolute Salinity [g/kg]

CT : array_like

Conservative Temperature [\(^\circ\) C (ITS-90)]

p : array_like

sea pressure [dbar]

Returns:

Tu : array_like

Turner angle, on the same (M-1)xN grid as p_mid. [degrees of rotation]

Rsubrho : array_like

Stability Ratio, on the same (M-1)xN grid as p_mid. [unitless]

p_mid : array_like

Mid pressure between p grid [dbar]

Notes

The 48-term equation has been fitted in a restricted range of parameter space, and is most accurate inside the “oceanographic funnel” described in McDougall et al. (2011). The GSW library function “infunnel(SA, CT, p)” is available to be used if one wants to test if some of one’s data lies outside this “funnel”.

References

[R197]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. (3.15.1) and (3.16.1).
[R198]McDougall T.J., P.M. Barker, R. Feistel and D.R. Jackett, 2011: A computationally efficient 48-term expression for the density of seawater in terms of Conservative Temperature, and related properties of seawater.

Examples

TODO

Module contents