CTD Functions

Background

The Ocean Observatories Initiative deploys the Sea-Bird Electronics conductivity-temperature-depth (CTD) family of instruments across its moored, profiling, and mobile platforms. The table below lists the OOI instrument classes covered by this module.

Class Hardware Platform type Designator meaning
CTDBP SBE 16Plus V2 Moored (fixed depth) CTD, Bottom Pumped
CTDMO SBE 37IM Moored (fixed depth) CTD, Modem (Inductive)
CTDPF SBE 16Plus V2 (A/B) or SBE 52MP (C/K/L) Profiling CTD, Profiler
CTDGV SBE GPCTD (Seabird Payload CTD) Gliders CTD, Glider Vehicle

ctd_functions.py converts raw Sea-Bird Electronics CTD data into calibrated L1 engineering products (TEMPWAT, PRESWAT, CONDWAT) and computes the L2 derived products Practical Salinity (PRACSAL) and in-situ Density (DENSITY) using the TEOS-10 GSW library. All calibration coefficients are from factory calibration values supplied with individual instruments.

The L1 engineering products (TEMPWAT, PRESWAT, CONDWAT) may be reported directly by the instrument (computed onboard the sensor using vendor firmware) or they may be reported as L0 values in raw units (counts, frequency, or linearly scaled integers). Depending on the instrument, telemetered data may be reported directly in L1 units while the recovered instrument data is reported in L0 units requiring the conversions below.


TEMPWAT_L1 — Seawater Temperature

The raw seawater temperature (TEMPWAT_L0, \(t_0\)), is reported in either counts or in L1 units (\(^\circ\)C), but scaled to an integer to compress the data for internal storage or transmission. The converted seawater temperature, TEMPWAT_L1, is reported in \(^\circ\)C. Conversions from L0 to L1 depend on the instrument class and/or data delivery method.

SBE 16Plus (CTDBP, CTDPF-A/B): Raw counts \(t_0\) are converted to seawater temperature (\(^\circ\)C) via an intermediate resistance quantity \(R\):

\[MV = \frac{t_0 - 524288}{1.6 \times 10^7}\]
\[R = \frac{MV \times 2.9 \times 10^9 + 1.024 \times 10^8}{2.048 \times 10^4 - MV \times 2.0 \times 10^5}\]
\[T_{L1} = \frac{1}{a_0 + a_1 \ln R + a_2 \ln^2 R + a_3 \ln^3 R} - 273.15\]

where \(a_0\), \(a_1\), \(a_2\), \(a_3\) are factory calibration coefficients.

SBE 37IM telemetered and recovered_host (CTDMO): The instrument outputs engineering units scaled to an integer. The L1 conversion is:

\[T_{L1} = t_0 / 10000 - 10\]

SBE 37IM instrument-recovered (CTDMO): Raw decimal counts \(t_0\) are retained in the instrument-stored file and require the calibration-coefficient equation:

\[T_{L1} = \frac{1}{a_0 + a_1 \ln t_0 + a_2 \ln^2 t_0 + a_3 \ln^3 t_0} - 273.15\]

SBE 52MP (CTDPF-C/K/L): The instrument outputs engineering units scaled to an integer. The L1 conversion is:

\[T_{L1} = t_0 / 10000 - 5\]

CTDGV (glider): The SBE GPCTD computes temperature onboard the vehicle using vendor software and transmits the result already in \(^\circ\)C; ion-functions is not invoked for those deployments.

Output accuracy: SBE 16Plus V2 \(\pm 0.005\ ^\circ\)C; SBE 37IM \(\pm 0.002\ ^\circ\)C.


PRESWAT_L1 — Seawater Pressure

The raw seawater pressure (PRESWAT_L0, \(p_0\)), is reported in either counts, frequency or in scaled L1 units (dbar). The converted seawater pressure, PRESWAT_L1, is reported in dbar relative to one standard atmosphere (10.1325 dbar). As with temperature, conversions from L0 to L1 depend on the instrument class and/or data delivery method.

SBE 16Plus — strain-gauge pressure sensor (CTDBP except N/O, CTDPF-A/B):

\[t_v = $pt_0$ \div 13107\]
\[t = \text{PTEMPA0} + \text{PTEMPA1} \times t_v + \text{PTEMPA2} \times t_v^2\]
\[x = p_0 - \text{PTCA0} - \text{PTCA1 }\times t - \text{PTCA2} \times t^2\]
\[n = \frac{x \times \text{PTCB0}}{\text{PTCB0} + \text{PTCB1} \times t + \text{PTCB2} \times t^2}\]
\[p_{\text{psi}} = \text{PA0} + \text{PA1} \times n + \text{PA2} \times n^2\]
\[P_{L1}\ = p_{\text{psi}} \times 0.689475729 - 10.1325 + \delta\]

where \(pt_0\) is the pressure sensor thermistor in counts, \(\delta\) is an optional Druck sensor offset correction (default 0 dbar) and all calibration coefficients are from factory calibration sheets.

SBE 16Plus — digiquartz pressure sensor (CTDBP-N and CTDBP-O only):

\[p_f = p_0 \div 256\]
\[t_v = t_0 \div 13107\]
\[U = 23.7(t_v + 9.7917) - 273.15\]
\[C = C_1 + C_2 U + C_3 U^2\]
\[D = D_1 + D_2 U\]
\[T_0 = T_1 + T_2 U + T_3 U^2 + T_4 U^3 + T_5 U^4\]
\[\tau = (1\div p_f) \times 10^6\]
\[p_{\text{psi}} = C\left(1 - \frac{T_0^2}{\tau^2}\right)\left(1 - D\left(1 - \frac{T_0^2}{\tau^2}\right)\right)\]
\[P_{L1} = p_{\text{psi}} \times 0.689475729 - 10.1325\]

All calibration coefficients are from factory calibration sheets.

SBE 37IM telemetered and recovered_host (CTDMO): Raw pressure is a scaled integer relative to a factory-set full-scale pressure range \(P_{rng}\) (in psi):

\[P_{rng,dbar} = (P_{rng,psi} - 14.7) \times 0.6894757\]
\[P_{L1} = \frac{p_0 \times P_{rng,dbar}}{0.85 \times 65536} - 0.05 \times P_{rng,dbar}\]

SBE 37IM instrument-recovered (CTDMO): Uses the same strain-gauge polynomial as the SBE 16Plus strain-gauge path, but with the raw thermistor count \(pt_0\) used directly in the polynomial rather than being first converted to voltage.

SBE 52MP (CTDPF-C/K/L):

\[P_{L1} = p_0 \div 100 - 10\]

CTDGV (glider): The SBE GPCTD reports pressure in bar; ion-functions converts to dbar:

\[P_{L1}\ [\text{dbar}] = p_{\text{bar}} \times 10\]

Output accuracy: SBE 16Plus V2 and SBE 37IM 0.1 % of full-scale range.


CONDWAT_L1 — Seawater Conductivity

The raw seawater conductivity (CONDWAT_L0, \(c_0\)), is reported in either counts or scaled L1 units. The converted seawater conductivity, CONDWAT_L1, is reported in S m\(^{-1}\). As with temperature, conversions from L0 to L1 depend on the instrument class and/or data delivery method.

SBE 16Plus (CTDBP, CTDPF-A/B): The raw count is converted to a frequency \(f\) in kHz, then evaluated with a polynomial corrected for temperature and pressure:

\[f = \frac{c_0 \div 256}{1000}\]
\[C_{L1} = \frac{g + h f^2 + i f^3 + j f^4}{1 + \text{CTcor}\times T + \text{CPcor}\times P}\]

where \(T\) is TEMPWAT_L1 (\(^\circ\text{C}\)), \(P\) is PRESWAT_L1 (dbar), and \(g\), \(h\), \(i\), \(j\), CTcor, CPcor are factory calibration coefficients.

SBE 37IM instrument-recovered (CTDMO): Same polynomial as the SBE 16Plus path, but includes an additional wbotc correction to the frequency:

\[f = \frac{c_0 \div 256}{1000} \times \sqrt{1 + \text{wbotc}\times T}\]

SBE 37IM telemetered and recovered_host (CTDMO):

\[C_{L1} = c_0 \div 100000 - 0.5\]

SBE 52MP (CTDPF-C/K/L): Linear scaling with a unit conversion from mS cm\(^{-1}\) to S m\(^{-1}\):

\[C_{L1} = (c_0 \div 10000 - 0.5) \times 0.1\]

CTDGV (glider): The SBE GPCTD computes conductivity onboard the vehicle using vendor software and transmits the result already in S m\(^{-1}\); ion-functions is not invoked for those deployments.

Output accuracy: SBE 16Plus V2 \(\pm 0.0005\) S m\(^{-1}\); SBE 37IM \(\pm 0.0003\) S m\(^{-1}\).


PRACSAL_L2 — Seawater Practical Salinity

Seawater practical salinity is computed from L1 conductivity, temperature, and pressure using the TEOS-10 GSW library function gsw.SP_from_C, which implements the Practical Salinity Scale 1978 (PSS-78) algorithm. Conductivity must be converted from S m\(^{-1}\) to mS cm\(^{-1}\) (multiply by 10) before calling the GSW function.

Seawater practical salinity is dimensionless and reported without units on the PSS-78 scale.

\[PS_{L2} =\text{gsw.SP_from_C}(C_{L1} * 10, T_{L1}, P_{L1})\]

DENSITY_L2 — In-situ Density

In-situ seawater density is computed via a three-step chain using the TEOS-10 GSW library:

Step 1 — Absolute Salinity:

\[S_A = \text{gsw.SA_from_SP}(PS_{L2}, P_{L1}, \text{lon}, \text{lat})\]

Absolute Salinity \(S_A\) (g kg\(^{-1}\)) is derived from Practical Salinity using a lookup table of the Absolute Salinity Anomaly (SAAR) as a function of location and pressure. For moored instruments, latitude and longitude are the mooring position metadata; for gliders, they are the vehicle position at each sample.

Step 2 — Conservative Temperature:

\[\Theta = \text{gsw.CT_from_t}(S_A, T_{L1}, P_{L1})\]

Step 3 — In-situ density:

\[\rho = \text{gsw.rho}(S_A, \Theta, p)\]

Density is computed using the computationally-efficient 48-term expression described in McDougall et al. (2011).


Core functions

ctd_sbe16plus_tempwat(t0, a0, a1, a2, a3)

Compute TEMPWAT_L1 from SBE 16Plus raw temperature counts.

Applies the Sea-Bird ITS-90 thermistor conversion for the SBE 16Plus V2 (OutputFormat 0). Covers all CTDBP instrument series and CTDPF series A and B.

Parameters:
  • t0 (array_like) –

    Raw temperature counts (TEMPWAT_L0) [counts].

  • a0 (float) –

    Temperature calibration coefficient a0.

  • a1 (float) –

    Temperature calibration coefficient a1.

  • a2 (float) –

    Temperature calibration coefficient a2.

  • a3 (float) –

    Temperature calibration coefficient a3.

Returns:
  • t( ndarray ) –

    Seawater temperature (TEMPWAT_L1) [degC, ITS-90].

Notes

Algorithm converts the 6-character hex integer t0 (already decoded to decimal counts by the CTD driver) to ITS-90 temperature via:

MV = (t0 - 524288) / 1.6e7
R  = (MV * 2.9e9 + 1.024e8) / (2.048e4 - MV * 2.0e5)
T  = 1 / (a0 + a1*ln(R) + a2*ln(R)^2 + a3*ln(R)^3) - 273.15

Calibration coefficients a0-a3 are from factory calibration sheets.

Source code in ion_functions/data/ctd_functions.py
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
def ctd_sbe16plus_tempwat(t0, a0, a1, a2, a3):
    """
    Compute TEMPWAT_L1 from SBE 16Plus raw temperature counts.

    Applies the Sea-Bird ITS-90 thermistor conversion for the SBE 16Plus
    V2 (OutputFormat 0). Covers all CTDBP instrument series and CTDPF
    series A and B.

    Parameters
    ----------
    t0 : array_like
        Raw temperature counts (TEMPWAT_L0) [counts].
    a0 : float
        Temperature calibration coefficient a0.
    a1 : float
        Temperature calibration coefficient a1.
    a2 : float
        Temperature calibration coefficient a2.
    a3 : float
        Temperature calibration coefficient a3.

    Returns
    -------
    t : ndarray
        Seawater temperature (TEMPWAT_L1) [degC, ITS-90].

    Notes
    -----
    Algorithm converts the 6-character hex integer t0 (already decoded
    to decimal counts by the CTD driver) to ITS-90 temperature via:

        MV = (t0 - 524288) / 1.6e7
        R  = (MV * 2.9e9 + 1.024e8) / (2.048e4 - MV * 2.0e5)
        T  = 1 / (a0 + a1*ln(R) + a2*ln(R)^2 + a3*ln(R)^3) - 273.15

    Calibration coefficients a0-a3 are from factory calibration sheets.
    """
    mv = (t0 - 524288) / 1.6e7
    r = (mv * 2.9e9 + 1.024e8) / (2.048e4 - mv * 2.0e5)
    t = 1 / (a0 + a1 * np.log(r) + a2 * np.log(r)**2 + a3 * np.log(r)**3) - 273.15
    return t

History

Date Author Change
2013-04-12 Luke Campbell Initial code
2013-04-12 Christopher Wingard Minor edits
2013-05-10 Christopher Wingard Minor edits to comments
2014-01-31 Russell Desiderio Standardized comment format
2023-08-15 Samuel Dahlberg Removed use of numexpr
2025-04-17 Christopher Wingard Converted to NumPy docstring format; updated documentation

ctd_sbe37im_tempwat_instrument_recovered(t0, a0, a1, a2, a3)

Compute TEMPWAT_L1 from SBE 37IM instrument-recovered counts.

Applies the Sea-Bird ITS-90 thermistor conversion for data recovered directly from an SBE 37IM instrument (all series), where the raw value is stored as a decimal count rather than scaled engineering units.

Parameters:
  • t0 (array_like) –

    Raw temperature counts (TEMPWAT_L0) recovered directly from the CTD instrument [counts].

  • a0 (float) –

    Temperature calibration coefficient a0.

  • a1 (float) –

    Temperature calibration coefficient a1.

  • a2 (float) –

    Temperature calibration coefficient a2.

  • a3 (float) –

    Temperature calibration coefficient a3.

Returns:
  • t( ndarray ) –

    Seawater temperature (TEMPWAT_L1) [degC, ITS-90].

Notes

Unlike telemetered and recovered_host data (see ctd_sbe37im_tempwat), instrument-recovered data retains the raw count t0, which is used directly in the temperature equation:

T = 1 / (a0 + a1*ln(t0) + a2*ln(t0)^2 + a3*ln(t0)^3) - 273.15

As of June 2016 this processing path was not included in the TEMPWAT DPS (DCN 1341-00010). Calibration coefficients a0--a3 are from factory calibration sheets.

Source code in ion_functions/data/ctd_functions.py
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
def ctd_sbe37im_tempwat_instrument_recovered(t0, a0, a1, a2, a3):
    """
    Compute TEMPWAT_L1 from SBE 37IM instrument-recovered counts.

    Applies the Sea-Bird ITS-90 thermistor conversion for data recovered
    directly from an SBE 37IM instrument (all series), where the raw
    value is stored as a decimal count rather than scaled engineering
    units.

    Parameters
    ----------
    t0 : array_like
        Raw temperature counts (TEMPWAT_L0) recovered directly from
        the CTD instrument [counts].
    a0 : float
        Temperature calibration coefficient a0.
    a1 : float
        Temperature calibration coefficient a1.
    a2 : float
        Temperature calibration coefficient a2.
    a3 : float
        Temperature calibration coefficient a3.

    Returns
    -------
    t : ndarray
        Seawater temperature (TEMPWAT_L1) [degC, ITS-90].

    Notes
    -----
    Unlike telemetered and recovered_host data (see
    ctd_sbe37im_tempwat), instrument-recovered data retains the raw
    count t0, which is used directly in the temperature equation:

        T = 1 / (a0 + a1*ln(t0) + a2*ln(t0)^2 + a3*ln(t0)^3) - 273.15

    As of June 2016 this processing path was not included in the
    TEMPWAT DPS (DCN 1341-00010). Calibration coefficients a0--a3 are
    from factory calibration sheets.
    """
    t = 1 / (a0 + a1 * np.log(t0) + a2 * np.log(t0)**2 + a3 * np.log(t0)**3) - 273.15
    return t

History

Date Author Change
2016-06-16 Russell Desiderio Initial code
2023-08-15 Samuel Dahlberg Removed use of numexpr
2025-04-17 Christopher Wingard Converted to NumPy docstring format; updated documentation

ctd_sbe37im_tempwat(t0)

Compute TEMPWAT_L1 from SBE 37IM telemetered/recovered_host counts.

Applies the linear scaling for SBE 37IM (all series) telemetered and recovered_host data streams (CTDMO instrument class). The instrument pre-scales raw counts to engineering units before transmission; no calibration coefficients are required.

Parameters:
  • t0 (array_like) –

    Raw temperature counts (TEMPWAT_L0) [counts].

Returns:
  • t( ndarray ) –

    Seawater temperature (TEMPWAT_L1) [degC, ITS-90].

Notes

The SBE 37IM encodes temperature as engineering units in hexadecimal (OutputFormat 0). The conversion is:

T = t0 / 10000 - 10

This function does not apply to instrument-recovered data; use ctd_sbe37im_tempwat_instrument_recovered for that stream.

Source code in ion_functions/data/ctd_functions.py
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
def ctd_sbe37im_tempwat(t0):
    """
    Compute TEMPWAT_L1 from SBE 37IM telemetered/recovered_host counts.

    Applies the linear scaling for SBE 37IM (all series) telemetered
    and recovered_host data streams (CTDMO instrument class). The
    instrument pre-scales raw counts to engineering units before
    transmission; no calibration coefficients are required.

    Parameters
    ----------
    t0 : array_like
        Raw temperature counts (TEMPWAT_L0) [counts].

    Returns
    -------
    t : ndarray
        Seawater temperature (TEMPWAT_L1) [degC, ITS-90].

    Notes
    -----
    The SBE 37IM encodes temperature as engineering units in
    hexadecimal (OutputFormat 0). The conversion is:

        T = t0 / 10000 - 10

    This function does not apply to instrument-recovered data; use
    ctd_sbe37im_tempwat_instrument_recovered for that stream.
    """
    t = t0 / 10000.0 - 10.0
    return t

History

Date Author Change
2014-02-05 Russell Desiderio Initial code
2025-04-17 Christopher Wingard Converted to NumPy docstring format; updated documentation

ctd_sbe52mp_tempwat(t0)

Compute TEMPWAT_L1 from SBE 52MP raw temperature counts.

Applies the linear scaling for SBE 52MP instruments (CTDPF series C, K, and L).

Parameters:
  • t0 (array_like) –

    Raw temperature counts (TEMPWAT_L0) [counts].

Returns:
  • t( ndarray ) –

    Seawater temperature (TEMPWAT_L1) [degC, ITS-90].

Notes

The conversion is:

T = t0 / 10000 - 5
Source code in ion_functions/data/ctd_functions.py
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
def ctd_sbe52mp_tempwat(t0):
    """
    Compute TEMPWAT_L1 from SBE 52MP raw temperature counts.

    Applies the linear scaling for SBE 52MP instruments (CTDPF series
    C, K, and L).

    Parameters
    ----------
    t0 : array_like
        Raw temperature counts (TEMPWAT_L0) [counts].

    Returns
    -------
    t : ndarray
        Seawater temperature (TEMPWAT_L1) [degC, ITS-90].

    Notes
    -----
    The conversion is:

        T = t0 / 10000 - 5
    """
    t = t0 / 10000.0 - 5.0
    return t

History

Date Author Change
2014-02-17 Russell Desiderio Initial code
2025-04-17 Christopher Wingard Converted to NumPy docstring format; updated documentation

ctd_sbe16plus_preswat(p0, t0, ptempa0, ptempa1, ptempa2, ptca0, ptca1, ptca2, ptcb0, ptcb1, ptcb2, pa0, pa1, pa2, offset=0)

Compute PRESWAT_L1 from SBE 16Plus strain-gauge pressure counts.

Applies the Sea-Bird strain-gauge pressure conversion for SBE 16Plus instruments equipped with an internal strain-gauge pressure sensor (PType=1). Covers most CTDBP instrument series (exceptions: N and O) and CTDPF series A and B.

Parameters:
  • p0 (array_like) –

    Raw pressure counts (PRESWAT_L0) [counts].

  • t0 (array_like) –

    Raw temperature counts from the pressure-sensor thermistor [counts].

  • ptempa0 (float) –

    Strain-gauge pressure calibration coefficient PTEMPA0.

  • ptempa1 (float) –

    Strain-gauge pressure calibration coefficient PTEMPA1.

  • ptempa2 (float) –

    Strain-gauge pressure calibration coefficient PTEMPA2.

  • ptca0 (float) –

    Strain-gauge pressure calibration coefficient PTCA0.

  • ptca1 (float) –

    Strain-gauge pressure calibration coefficient PTCA1.

  • ptca2 (float) –

    Strain-gauge pressure calibration coefficient PTCA2.

  • ptcb0 (float) –

    Strain-gauge pressure calibration coefficient PTCB0.

  • ptcb1 (float) –

    Strain-gauge pressure calibration coefficient PTCB1.

  • ptcb2 (float) –

    Strain-gauge pressure calibration coefficient PTCB2.

  • pa0 (float) –

    Strain-gauge pressure calibration coefficient PA0.

  • pa1 (float) –

    Strain-gauge pressure calibration coefficient PA1.

  • pa2 (float) –

    Strain-gauge pressure calibration coefficient PA2.

  • offset (float, default: 0 ) –

    Druck sensor offset correction [dbar]. Default is 0.

Returns:
  • p( ndarray ) –

    Seawater pressure (PRESWAT_L1) [dbar].

Notes

Algorithm (from OOI DPS DCN 1341-00020, PType=1):

t_v = t0 / 13107
t   = PTEMPA0 + PTEMPA1*t_v + PTEMPA2*t_v^2
x   = p0 - PTCA0 - PTCA1*t - PTCA2*t^2
n   = x * PTCB0 / (PTCB0 + PTCB1*t + PTCB2*t^2)
p_psi = PA0 + PA1*n + PA2*n^2
P_L1  = p_psi * 0.689475729 - 10.1325 + offset

All calibration coefficients are from factory calibration sheets. The optional offset corrects a known Druck sensor bias.

Source code in ion_functions/data/ctd_functions.py
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
def ctd_sbe16plus_preswat(p0, t0, ptempa0, ptempa1, ptempa2,
                          ptca0, ptca1, ptca2, ptcb0, ptcb1, ptcb2,
                          pa0, pa1, pa2, offset=0):
    """
    Compute PRESWAT_L1 from SBE 16Plus strain-gauge pressure counts.

    Applies the Sea-Bird strain-gauge pressure conversion for SBE 16Plus
    instruments equipped with an internal strain-gauge pressure sensor
    (PType=1). Covers most CTDBP instrument series (exceptions: N and O)
    and CTDPF series A and B.

    Parameters
    ----------
    p0 : array_like
        Raw pressure counts (PRESWAT_L0) [counts].
    t0 : array_like
        Raw temperature counts from the pressure-sensor thermistor
        [counts].
    ptempa0 : float
        Strain-gauge pressure calibration coefficient PTEMPA0.
    ptempa1 : float
        Strain-gauge pressure calibration coefficient PTEMPA1.
    ptempa2 : float
        Strain-gauge pressure calibration coefficient PTEMPA2.
    ptca0 : float
        Strain-gauge pressure calibration coefficient PTCA0.
    ptca1 : float
        Strain-gauge pressure calibration coefficient PTCA1.
    ptca2 : float
        Strain-gauge pressure calibration coefficient PTCA2.
    ptcb0 : float
        Strain-gauge pressure calibration coefficient PTCB0.
    ptcb1 : float
        Strain-gauge pressure calibration coefficient PTCB1.
    ptcb2 : float
        Strain-gauge pressure calibration coefficient PTCB2.
    pa0 : float
        Strain-gauge pressure calibration coefficient PA0.
    pa1 : float
        Strain-gauge pressure calibration coefficient PA1.
    pa2 : float
        Strain-gauge pressure calibration coefficient PA2.
    offset : float, optional
        Druck sensor offset correction [dbar]. Default is 0.

    Returns
    -------
    p : ndarray
        Seawater pressure (PRESWAT_L1) [dbar].

    Notes
    -----
    Algorithm (from OOI DPS DCN 1341-00020, PType=1):

        t_v = t0 / 13107
        t   = PTEMPA0 + PTEMPA1*t_v + PTEMPA2*t_v^2
        x   = p0 - PTCA0 - PTCA1*t - PTCA2*t^2
        n   = x * PTCB0 / (PTCB0 + PTCB1*t + PTCB2*t^2)
        p_psi = PA0 + PA1*n + PA2*n^2
        P_L1  = p_psi * 0.689475729 - 10.1325 + offset

    All calibration coefficients are from factory calibration sheets.
    The optional offset corrects a known Druck sensor bias.
    """
    # compute calibration parameters
    tv = t0 / 13107.0
    t = ptempa0 + ptempa1 * tv + ptempa2 * tv**2
    x = p0 - ptca0 - ptca1 * t - ptca2 * t**2
    n = x * ptcb0 / (ptcb0 + ptcb1 * t + ptcb2 * t**2)

    # compute pressure in psi, rescale and compute in dbar and return
    p_psi = pa0 + pa1 * n + pa2 * n**2
    p_dbar = (p_psi * 0.689475729) - 10.1325
    return p_dbar + offset

History

Date Author Change
2013-04-12 Christopher Wingard Initial code
2013-05-10 Christopher Wingard Minor edits to comments
2014-01-31 Russell Desiderio Standardized comment format
2017-03-31 Dan Mergens Added Druck offset correction
2025-04-17 Christopher Wingard Converted to NumPy docstring format; updated documentation

ctd_sbe16digi_preswat(p0, t0, C1, C2, C3, D1, D2, T1, T2, T3, T4, T5)

Compute PRESWAT_L1 from SBE 16Plus digiquartz pressure counts.

Applies the Sea-Bird digiquartz pressure conversion for SBE 16Plus instruments equipped with a digiquartz pressure sensor (PType=3). Applies exclusively to CTDBP-N and CTDBP-O instruments.

Parameters:
  • p0 (array_like) –

    Raw pressure counts (PRESWAT_L0) [counts].

  • t0 (array_like) –

    Raw temperature counts from the pressure-sensor thermistor [counts].

  • C1 (float) –

    Digiquartz pressure calibration coefficient C1.

  • C2 (float) –

    Digiquartz pressure calibration coefficient C2.

  • C3 (float) –

    Digiquartz pressure calibration coefficient C3.

  • D1 (float) –

    Digiquartz pressure calibration coefficient D1.

  • D2 (float) –

    Digiquartz pressure calibration coefficient D2.

  • T1 (float) –

    Digiquartz pressure calibration coefficient T1.

  • T2 (float) –

    Digiquartz pressure calibration coefficient T2.

  • T3 (float) –

    Digiquartz pressure calibration coefficient T3.

  • T4 (float) –

    Digiquartz pressure calibration coefficient T4.

  • T5 (float) –

    Digiquartz pressure calibration coefficient T5.

Returns:
  • p( ndarray ) –

    Seawater pressure (PRESWAT_L1) [dbar].

Notes

Algorithm (from OOI DPS DCN 1341-00020, PType=3):

pf  = p0 / 256          (pressure frequency in Hz)
t_v = t0 / 13107        (thermistor voltage)
U   = 23.7*(t_v + 9.7917) - 273.15
C   = C1 + C2*U + C3*U^2
D   = D1 + D2*U
T0  = T1 + T2*U + T3*U^2 + T4*U^3 + T5*U^4
T   = (1/pf) * 1e6      (pressure period in microseconds)
p_psi = C*(1 - T0^2/T^2)*(1 - D*(1 - T0^2/T^2))
P_L1  = p_psi * 0.689475729 - 10.1325

All calibration coefficients are from factory calibration sheets.

Source code in ion_functions/data/ctd_functions.py
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
def ctd_sbe16digi_preswat(p0, t0, C1, C2, C3, D1, D2, T1, T2, T3, T4, T5):
    """
    Compute PRESWAT_L1 from SBE 16Plus digiquartz pressure counts.

    Applies the Sea-Bird digiquartz pressure conversion for SBE 16Plus
    instruments equipped with a digiquartz pressure sensor (PType=3).
    Applies exclusively to CTDBP-N and CTDBP-O instruments.

    Parameters
    ----------
    p0 : array_like
        Raw pressure counts (PRESWAT_L0) [counts].
    t0 : array_like
        Raw temperature counts from the pressure-sensor thermistor
        [counts].
    C1 : float
        Digiquartz pressure calibration coefficient C1.
    C2 : float
        Digiquartz pressure calibration coefficient C2.
    C3 : float
        Digiquartz pressure calibration coefficient C3.
    D1 : float
        Digiquartz pressure calibration coefficient D1.
    D2 : float
        Digiquartz pressure calibration coefficient D2.
    T1 : float
        Digiquartz pressure calibration coefficient T1.
    T2 : float
        Digiquartz pressure calibration coefficient T2.
    T3 : float
        Digiquartz pressure calibration coefficient T3.
    T4 : float
        Digiquartz pressure calibration coefficient T4.
    T5 : float
        Digiquartz pressure calibration coefficient T5.

    Returns
    -------
    p : ndarray
        Seawater pressure (PRESWAT_L1) [dbar].

    Notes
    -----
    Algorithm (from OOI DPS DCN 1341-00020, PType=3):

        pf  = p0 / 256          (pressure frequency in Hz)
        t_v = t0 / 13107        (thermistor voltage)
        U   = 23.7*(t_v + 9.7917) - 273.15
        C   = C1 + C2*U + C3*U^2
        D   = D1 + D2*U
        T0  = T1 + T2*U + T3*U^2 + T4*U^3 + T5*U^4
        T   = (1/pf) * 1e6      (pressure period in microseconds)
        p_psi = C*(1 - T0^2/T^2)*(1 - D*(1 - T0^2/T^2))
        P_L1  = p_psi * 0.689475729 - 10.1325

    All calibration coefficients are from factory calibration sheets.
    """
    # Convert raw pressure input to frequency [Hz]
    pf = p0 / 256.0

    # Convert raw temperature input to voltage
    tv = t0 / 13107.0

    # Calculate U (thermistor temp):
    U = (23.7 * (tv + 9.7917)) - 273.15

    # Calculate calibration parameters
    C = C1 + C2 * U + C3 * U**2
    D = D1 + D2 * U
    T0 = T1 + T2 * U + T3 * U**2 + T4 * U**3 + T5 * U**4

    # Calculate T (pressure period, in microseconds):
    T = (1.0 / pf) * 1.0e6

    # compute pressure in psi, rescale and compute in dbar and return
    p_psi = C * (1.0 - T0**2 / T**2) * (1.0 - D * (1.0 - T0**2 / T**2))
    p_dbar = (p_psi * 0.689475729) - 10.1325
    return p_dbar

Additional Notes

Per the SBE 16Plus V2 User Manual, the raw pressure input p0 is in units of Hz (counts divided by 256), not raw A/D counts; the code implements this correctly.

History

Date Author Change
2013-05-10 Christopher Wingard Initial code
2014-01-31 Russell Desiderio Standardized comment format; corrected pressure period calculation to use Hz input
2025-04-17 Christopher Wingard Converted to NumPy docstring format; updated documentation

ctd_sbe37im_preswat_instrument_recovered(p0, pt0, ptempa0, ptempa1, ptempa2, ptca0, ptca1, ptca2, ptcb0, ptcb1, ptcb2, pa0, pa1, pa2)

Compute PRESWAT_L1 from SBE 37IM instrument-recovered pressure counts.

Applies the Sea-Bird strain-gauge pressure conversion for data recovered directly from an SBE 37IM instrument (all series).

Parameters:
  • p0 (array_like) –

    Raw pressure counts (PRESWAT_L0) recovered directly from the CTD instrument [counts].

  • pt0 (array_like) –

    Raw temperature counts from the pressure-sensor thermistor [counts].

  • ptempa0 (float) –

    Strain-gauge pressure calibration coefficient PTEMPA0.

  • ptempa1 (float) –

    Strain-gauge pressure calibration coefficient PTEMPA1.

  • ptempa2 (float) –

    Strain-gauge pressure calibration coefficient PTEMPA2.

  • ptca0 (float) –

    Strain-gauge pressure calibration coefficient PTCA0.

  • ptca1 (float) –

    Strain-gauge pressure calibration coefficient PTCA1.

  • ptca2 (float) –

    Strain-gauge pressure calibration coefficient PTCA2.

  • ptcb0 (float) –

    Strain-gauge pressure calibration coefficient PTCB0.

  • ptcb1 (float) –

    Strain-gauge pressure calibration coefficient PTCB1.

  • ptcb2 (float) –

    Strain-gauge pressure calibration coefficient PTCB2.

  • pa0 (float) –

    Strain-gauge pressure calibration coefficient PA0.

  • pa1 (float) –

    Strain-gauge pressure calibration coefficient PA1.

  • pa2 (float) –

    Strain-gauge pressure calibration coefficient PA2.

Returns:
  • p( ndarray ) –

    Seawater pressure (PRESWAT_L1) [dbar].

Notes

Unlike telemetered and recovered_host data (see ctd_sbe37im_preswat), instrument-recovered data retains raw strain-gauge counts. The algorithm matches the SBE 16Plus strain-gauge path but uses pt0 directly (not scaled to voltage):

t   = PTEMPA0 + PTEMPA1*pt0 + PTEMPA2*pt0^2
x   = p0 - PTCA0 - PTCA1*t - PTCA2*t^2
n   = x * PTCB0 / (PTCB0 + PTCB1*t + PTCB2*t^2)
p_psi = PA0 + PA1*n + PA2*n^2
P_L1  = p_psi * 0.689475729 - 10.1325

As of June 2016 this processing path was not included in the PRESWAT DPS (DCN 1341-00020). Calibration coefficients are from factory calibration sheets.

Source code in ion_functions/data/ctd_functions.py
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
def ctd_sbe37im_preswat_instrument_recovered(p0, pt0, ptempa0, ptempa1, ptempa2,
                                             ptca0, ptca1, ptca2, ptcb0, ptcb1, ptcb2,
                                             pa0, pa1, pa2):
    """
    Compute PRESWAT_L1 from SBE 37IM instrument-recovered pressure counts.

    Applies the Sea-Bird strain-gauge pressure conversion for data
    recovered directly from an SBE 37IM instrument (all series).

    Parameters
    ----------
    p0 : array_like
        Raw pressure counts (PRESWAT_L0) recovered directly from the
        CTD instrument [counts].
    pt0 : array_like
        Raw temperature counts from the pressure-sensor thermistor
        [counts].
    ptempa0 : float
        Strain-gauge pressure calibration coefficient PTEMPA0.
    ptempa1 : float
        Strain-gauge pressure calibration coefficient PTEMPA1.
    ptempa2 : float
        Strain-gauge pressure calibration coefficient PTEMPA2.
    ptca0 : float
        Strain-gauge pressure calibration coefficient PTCA0.
    ptca1 : float
        Strain-gauge pressure calibration coefficient PTCA1.
    ptca2 : float
        Strain-gauge pressure calibration coefficient PTCA2.
    ptcb0 : float
        Strain-gauge pressure calibration coefficient PTCB0.
    ptcb1 : float
        Strain-gauge pressure calibration coefficient PTCB1.
    ptcb2 : float
        Strain-gauge pressure calibration coefficient PTCB2.
    pa0 : float
        Strain-gauge pressure calibration coefficient PA0.
    pa1 : float
        Strain-gauge pressure calibration coefficient PA1.
    pa2 : float
        Strain-gauge pressure calibration coefficient PA2.

    Returns
    -------
    p : ndarray
        Seawater pressure (PRESWAT_L1) [dbar].

    Notes
    -----
    Unlike telemetered and recovered_host data (see ctd_sbe37im_preswat),
    instrument-recovered data retains raw strain-gauge counts. The
    algorithm matches the SBE 16Plus strain-gauge path but uses pt0
    directly (not scaled to voltage):

        t   = PTEMPA0 + PTEMPA1*pt0 + PTEMPA2*pt0^2
        x   = p0 - PTCA0 - PTCA1*t - PTCA2*t^2
        n   = x * PTCB0 / (PTCB0 + PTCB1*t + PTCB2*t^2)
        p_psi = PA0 + PA1*n + PA2*n^2
        P_L1  = p_psi * 0.689475729 - 10.1325

    As of June 2016 this processing path was not included in the
    PRESWAT DPS (DCN 1341-00020). Calibration coefficients are from
    factory calibration sheets.
    """
    # compute calibration parameters
    t = ptempa0 + ptempa1 * pt0 + ptempa2 * pt0**2
    x = p0 - ptca0 - ptca1 * t - ptca2 * t**2
    n = x * ptcb0 / (ptcb0 + ptcb1 * t + ptcb2 * t**2)

    # compute pressure in psi, rescale and compute in dbar and return
    p_psi = pa0 + pa1 * n + pa2 * n**2
    p_dbar = (p_psi * 0.689475729) - 10.1325
    return p_dbar

History

Date Author Change
2016-06-16 Russell Desiderio Initial code
2025-04-17 Christopher Wingard Converted to NumPy docstring format; updated documentation

ctd_sbe37im_preswat(p0, p_range_psia)

Compute PRESWAT_L1 from SBE 37IM telemetered/recovered_host counts.

Applies the linear pressure scaling for SBE 37IM (all series) telemetered and recovered_host data streams (CTDMO instrument class).

Parameters:
  • p0 (array_like) –

    Raw pressure counts (PRESWAT_L0) [counts].

  • p_range_psia (float) –

    Pressure range calibration coefficient [psia], a factory-set value stored in the instrument metadata.

Returns:
  • p( ndarray ) –

    Seawater pressure (PRESWAT_L1) [dbar].

Notes

Algorithm (from OOI DPS DCN 1341-00020, SBE 37IM):

P_range_dbar = (p_range_psia - 14.7) * 0.6894757
P_L1 = p0 * P_range_dbar / (0.85 * 65536) - 0.05 * P_range_dbar

The pressure range is a factory-set calibration coefficient. This function does not apply to instrument-recovered data; use ctd_sbe37im_preswat_instrument_recovered for that stream.

Source code in ion_functions/data/ctd_functions.py
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
def ctd_sbe37im_preswat(p0, p_range_psia):
    """
    Compute PRESWAT_L1 from SBE 37IM telemetered/recovered_host counts.

    Applies the linear pressure scaling for SBE 37IM (all series)
    telemetered and recovered_host data streams (CTDMO instrument class).

    Parameters
    ----------
    p0 : array_like
        Raw pressure counts (PRESWAT_L0) [counts].
    p_range_psia : float
        Pressure range calibration coefficient [psia], a factory-set
        value stored in the instrument metadata.

    Returns
    -------
    p : ndarray
        Seawater pressure (PRESWAT_L1) [dbar].

    Notes
    -----
    Algorithm (from OOI DPS DCN 1341-00020, SBE 37IM):

        P_range_dbar = (p_range_psia - 14.7) * 0.6894757
        P_L1 = p0 * P_range_dbar / (0.85 * 65536) - 0.05 * P_range_dbar

    The pressure range is a factory-set calibration coefficient.
    This function does not apply to instrument-recovered data; use
    ctd_sbe37im_preswat_instrument_recovered for that stream.
    """
    # compute pressure range in units of dbar
    p_range_dbar = (p_range_psia - 14.7) * 0.6894757

    # compute pressure in dbar and return
    p_dbar = p0 * p_range_dbar / (0.85 * 65536.0) - 0.05 * p_range_dbar
    return p_dbar

History

Date Author Change
2014-02-05 Russell Desiderio Initial code
2025-04-17 Christopher Wingard Converted to NumPy docstring format; updated documentation

ctd_glider_preswat(pr_bar)

Compute PRESWAT_L1 from glider CTD pressure in bar.

Converts pressure reported by a Seabird CTD installed on a glider (CTDGV instrument class) from bar to dbar.

Parameters:
  • pr_bar (array_like) –

    Seawater pressure from glider (PRESWAT_L0) [bar].

Returns:
  • pr_dbar( ndarray ) –

    Seawater pressure (PRESWAT_L1) [dbar].

Notes

Conversion:

P_L1 [dbar] = pr_bar * 10
Source code in ion_functions/data/ctd_functions.py
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
def ctd_glider_preswat(pr_bar):
    """
    Compute PRESWAT_L1 from glider CTD pressure in bar.

    Converts pressure reported by a Seabird CTD installed on a glider
    (CTDGV instrument class) from bar to dbar.

    Parameters
    ----------
    pr_bar : array_like
        Seawater pressure from glider (PRESWAT_L0) [bar].

    Returns
    -------
    pr_dbar : ndarray
        Seawater pressure (PRESWAT_L1) [dbar].

    Notes
    -----
    Conversion:

        P_L1 [dbar] = pr_bar * 10
    """
    pr_dbar = pr_bar * 10.0
    return pr_dbar

History

Date Author Change
2015-10-28 Russell Desiderio Initial code
2025-04-17 Christopher Wingard Converted to NumPy docstring format; updated documentation

ctd_sbe52mp_preswat(p0)

Compute PRESWAT_L1 from SBE 52MP raw pressure counts.

Applies the linear pressure scaling for SBE 52MP instruments (CTDPF series C, K, and L).

Parameters:
  • p0 (array_like) –

    Raw pressure counts (PRESWAT_L0) [counts].

Returns:
  • p( ndarray ) –

    Seawater pressure (PRESWAT_L1) [dbar].

Notes

Conversion:

P_L1 = p0 / 100 - 10
Source code in ion_functions/data/ctd_functions.py
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
def ctd_sbe52mp_preswat(p0):
    """
    Compute PRESWAT_L1 from SBE 52MP raw pressure counts.

    Applies the linear pressure scaling for SBE 52MP instruments
    (CTDPF series C, K, and L).

    Parameters
    ----------
    p0 : array_like
        Raw pressure counts (PRESWAT_L0) [counts].

    Returns
    -------
    p : ndarray
        Seawater pressure (PRESWAT_L1) [dbar].

    Notes
    -----
    Conversion:

        P_L1 = p0 / 100 - 10
    """
    p_dbar = p0 / 100.0 - 10.0
    return p_dbar

History

Date Author Change
2014-02-17 Russell Desiderio Initial code
2025-04-17 Christopher Wingard Converted to NumPy docstring format; updated documentation

ctd_sbe16plus_condwat(c0, t1, p1, g, h, i, j, cpcor, ctcor)

Compute CONDWAT_L1 from SBE 16Plus raw conductivity counts.

Applies the Sea-Bird conductivity frequency conversion for SBE 16Plus instruments. Covers all CTDBP instrument series and CTDPF series A and B.

Parameters:
  • c0 (array_like) –

    Raw conductivity counts (CONDWAT_L0) [counts].

  • t1 (array_like) –

    Seawater temperature (TEMPWAT_L1) [degC].

  • p1 (array_like) –

    Seawater pressure (PRESWAT_L1) [dbar].

  • g (float) –

    Conductivity calibration coefficient g.

  • h (float) –

    Conductivity calibration coefficient h.

  • i (float) –

    Conductivity calibration coefficient i.

  • j (float) –

    Conductivity calibration coefficient j.

  • cpcor (float) –

    Conductivity calibration coefficient CPcor.

  • ctcor (float) –

    Conductivity calibration coefficient CTcor.

Returns:
  • c( ndarray ) –

    Seawater conductivity (CONDWAT_L1) [S m^-1].

Notes

Algorithm (from OOI DPS DCN 1341-00030, SBE 16Plus):

f    = (c0 / 256) / 1000    (conductivity frequency in kHz)
C_L1 = (g + h*f^2 + i*f^3 + j*f^4) / (1 + CTcor*T + CPcor*P)

where T is TEMPWAT_L1 (degC) and P is PRESWAT_L1 (dbar). All calibration coefficients are from factory calibration sheets.

Source code in ion_functions/data/ctd_functions.py
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
def ctd_sbe16plus_condwat(c0, t1, p1, g, h, i, j, cpcor, ctcor):
    """
    Compute CONDWAT_L1 from SBE 16Plus raw conductivity counts.

    Applies the Sea-Bird conductivity frequency conversion for SBE
    16Plus instruments. Covers all CTDBP instrument series and CTDPF
    series A and B.

    Parameters
    ----------
    c0 : array_like
        Raw conductivity counts (CONDWAT_L0) [counts].
    t1 : array_like
        Seawater temperature (TEMPWAT_L1) [degC].
    p1 : array_like
        Seawater pressure (PRESWAT_L1) [dbar].
    g : float
        Conductivity calibration coefficient g.
    h : float
        Conductivity calibration coefficient h.
    i : float
        Conductivity calibration coefficient i.
    j : float
        Conductivity calibration coefficient j.
    cpcor : float
        Conductivity calibration coefficient CPcor.
    ctcor : float
        Conductivity calibration coefficient CTcor.

    Returns
    -------
    c : ndarray
        Seawater conductivity (CONDWAT_L1) [S m^-1].

    Notes
    -----
    Algorithm (from OOI DPS DCN 1341-00030, SBE 16Plus):

        f    = (c0 / 256) / 1000    (conductivity frequency in kHz)
        C_L1 = (g + h*f^2 + i*f^3 + j*f^4) / (1 + CTcor*T + CPcor*P)

    where T is TEMPWAT_L1 (degC) and P is PRESWAT_L1 (dbar).
    All calibration coefficients are from factory calibration sheets.
    """
    # convert raw conductivity measurement to frequency
    f = (c0 / 256.0) / 1000.0

    # calculate conductivity [S m-1]
    c = (g + h * f**2 + i * f**3 + j * f**4) / (1 + ctcor * t1 + cpcor * p1)
    return c

History

Date Author Change
2013-04-12 Christopher Wingard Initial code
2013-05-10 Christopher Wingard Minor edits to comments
2014-01-31 Russell Desiderio Standardized comment format
2025-04-17 Christopher Wingard Converted to NumPy docstring format; updated documentation

ctd_sbe37im_condwat_instrument_recovered(c0, t1, p1, g, h, i, j, cpcor, ctcor, wbotc)

Compute CONDWAT_L1 from SBE 37IM instrument-recovered counts.

Applies the Sea-Bird conductivity frequency conversion for data recovered directly from an SBE 37IM instrument (all series).

Parameters:
  • c0 (array_like) –

    Raw conductivity counts (CONDWAT_L0) recovered directly from the CTD instrument [counts].

  • t1 (array_like) –

    Seawater temperature (TEMPWAT_L1) [degC].

  • p1 (array_like) –

    Seawater pressure (PRESWAT_L1) [dbar].

  • g (float) –

    Conductivity calibration coefficient g.

  • h (float) –

    Conductivity calibration coefficient h.

  • i (float) –

    Conductivity calibration coefficient i.

  • j (float) –

    Conductivity calibration coefficient j.

  • cpcor (float) –

    Conductivity calibration coefficient CPcor.

  • ctcor (float) –

    Conductivity calibration coefficient CTcor.

  • wbotc (float) –

    Conductivity calibration coefficient wbotc.

Returns:
  • c( ndarray ) –

    Seawater conductivity (CONDWAT_L1) [S m^-1].

Notes

Algorithm applies a wbotc correction to the frequency before the standard polynomial evaluation:

f    = (c0/256)/1000 * sqrt(1 + wbotc*T)
C_L1 = (g + h*f^2 + i*f^3 + j*f^4) / (1 + CTcor*T + CPcor*P)

where T is TEMPWAT_L1 (degC) and P is PRESWAT_L1 (dbar). As of June 2016 this processing path was not included in the CONDWAT DPS (DCN 1341-00030). Calibration coefficients are from factory calibration sheets.

Source code in ion_functions/data/ctd_functions.py
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
def ctd_sbe37im_condwat_instrument_recovered(c0, t1, p1, g, h, i, j, cpcor, ctcor, wbotc):
    """
    Compute CONDWAT_L1 from SBE 37IM instrument-recovered counts.

    Applies the Sea-Bird conductivity frequency conversion for data
    recovered directly from an SBE 37IM instrument (all series).

    Parameters
    ----------
    c0 : array_like
        Raw conductivity counts (CONDWAT_L0) recovered directly from
        the CTD instrument [counts].
    t1 : array_like
        Seawater temperature (TEMPWAT_L1) [degC].
    p1 : array_like
        Seawater pressure (PRESWAT_L1) [dbar].
    g : float
        Conductivity calibration coefficient g.
    h : float
        Conductivity calibration coefficient h.
    i : float
        Conductivity calibration coefficient i.
    j : float
        Conductivity calibration coefficient j.
    cpcor : float
        Conductivity calibration coefficient CPcor.
    ctcor : float
        Conductivity calibration coefficient CTcor.
    wbotc : float
        Conductivity calibration coefficient wbotc.

    Returns
    -------
    c : ndarray
        Seawater conductivity (CONDWAT_L1) [S m^-1].

    Notes
    -----
    Algorithm applies a wbotc correction to the frequency before the
    standard polynomial evaluation:

        f    = (c0/256)/1000 * sqrt(1 + wbotc*T)
        C_L1 = (g + h*f^2 + i*f^3 + j*f^4) / (1 + CTcor*T + CPcor*P)

    where T is TEMPWAT_L1 (degC) and P is PRESWAT_L1 (dbar).
    As of June 2016 this processing path was not included in the
    CONDWAT DPS (DCN 1341-00030). Calibration coefficients are from
    factory calibration sheets.
    """
    # convert raw conductivity measurement to frequency
    f = (c0 / 256.0) / 1000.0 * np.sqrt(1.0 + wbotc * t1)

    # calculate conductivity [S m-1]
    c = (g + h * f**2 + i * f**3 + j * f**4) / (1 + ctcor * t1 + cpcor * p1)
    return c

History

Date Author Change
2016-06-16 Russell Desiderio Initial code
2023-08-15 Samuel Dahlberg Removed use of numexpr
2025-04-17 Christopher Wingard Converted to NumPy docstring format; updated documentation

ctd_sbe37im_condwat(c0)

Compute CONDWAT_L1 from SBE 37IM telemetered/recovered_host counts.

Applies the linear conductivity scaling for SBE 37IM (all series) telemetered and recovered_host data streams (CTDMO instrument class).

Parameters:
  • c0 (array_like) –

    Raw conductivity counts (CONDWAT_L0) [counts].

Returns:
  • c( ndarray ) –

    Seawater conductivity (CONDWAT_L1) [S m^-1].

Notes

The SBE 37IM encodes conductivity as engineering units in hexadecimal (OutputFormat 0). The conversion is:

C_L1 = c0 / 100000 - 0.5

This function does not apply to instrument-recovered data; use ctd_sbe37im_condwat_instrument_recovered for that stream.

Source code in ion_functions/data/ctd_functions.py
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
def ctd_sbe37im_condwat(c0):
    """
    Compute CONDWAT_L1 from SBE 37IM telemetered/recovered_host counts.

    Applies the linear conductivity scaling for SBE 37IM (all series)
    telemetered and recovered_host data streams (CTDMO instrument class).

    Parameters
    ----------
    c0 : array_like
        Raw conductivity counts (CONDWAT_L0) [counts].

    Returns
    -------
    c : ndarray
        Seawater conductivity (CONDWAT_L1) [S m^-1].

    Notes
    -----
    The SBE 37IM encodes conductivity as engineering units in
    hexadecimal (OutputFormat 0). The conversion is:

        C_L1 = c0 / 100000 - 0.5

    This function does not apply to instrument-recovered data; use
    ctd_sbe37im_condwat_instrument_recovered for that stream.
    """
    c = c0 / 100000.0 - 0.5
    return c

History

Date Author Change
2014-02-05 Russell Desiderio Initial code
2025-04-17 Christopher Wingard Converted to NumPy docstring format; updated documentation

ctd_sbe52mp_condwat(c0)

Compute CONDWAT_L1 from SBE 52MP raw conductivity counts.

Applies the linear conductivity scaling for SBE 52MP instruments (CTDPF series C, K, and L).

Parameters:
  • c0 (array_like) –

    Raw conductivity counts (CONDWAT_L0) [counts].

Returns:
  • c( ndarray ) –

    Seawater conductivity (CONDWAT_L1) [S m^-1].

Notes

Two-step conversion:

C [mmho/cm] = c0 / 10000 - 0.5
C [S m^-1]  = C [mmho/cm] * 0.1
Source code in ion_functions/data/ctd_functions.py
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
def ctd_sbe52mp_condwat(c0):
    """
    Compute CONDWAT_L1 from SBE 52MP raw conductivity counts.

    Applies the linear conductivity scaling for SBE 52MP instruments
    (CTDPF series C, K, and L).

    Parameters
    ----------
    c0 : array_like
        Raw conductivity counts (CONDWAT_L0) [counts].

    Returns
    -------
    c : ndarray
        Seawater conductivity (CONDWAT_L1) [S m^-1].

    Notes
    -----
    Two-step conversion:

        C [mmho/cm] = c0 / 10000 - 0.5
        C [S m^-1]  = C [mmho/cm] * 0.1
    """
    c_mmho_cm = c0 / 10000.0 - 0.5
    c_S_m = 0.1 * c_mmho_cm
    return c_S_m

History

Date Author Change
2014-02-17 Russell Desiderio Initial code
2025-04-17 Christopher Wingard Converted to NumPy docstring format; updated documentation

ctd_pracsal(c, t, p)

Compute PRACSAL_L2 from L1 conductivity, temperature, and pressure.

Calculates Practical Salinity on the PSS-78 scale using the TEOS-10 GSW library function gsw.SP_from_C, which implements the UNESCO 1983 PSS-78 algorithm with the Hill et al. (1986) extension for SP < 2.

Parameters:
  • c (array_like) –

    Seawater conductivity (CONDWAT_L1) [S m^-1].

  • t (array_like) –

    Seawater temperature (TEMPWAT_L1) [degC, ITS-90].

  • p (array_like) –

    Seawater pressure (PRESWAT_L1) [dbar].

Returns:
  • SP( ndarray ) –

    Practical salinity (PRACSAL_L2) [PSS-78, unitless].

Notes

Conductivity is converted from S m^-1 to mS cm^-1 (multiply by 10) before passing to gsw.SP_from_C(C, t, p), which expects mS cm^-1.

Source code in ion_functions/data/ctd_functions.py
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
def ctd_pracsal(c, t, p):
    """
    Compute PRACSAL_L2 from L1 conductivity, temperature, and pressure.

    Calculates Practical Salinity on the PSS-78 scale using the TEOS-10
    GSW library function gsw.SP_from_C, which implements the UNESCO 1983
    PSS-78 algorithm with the Hill et al. (1986) extension for SP < 2.

    Parameters
    ----------
    c : array_like
        Seawater conductivity (CONDWAT_L1) [S m^-1].
    t : array_like
        Seawater temperature (TEMPWAT_L1) [degC, ITS-90].
    p : array_like
        Seawater pressure (PRESWAT_L1) [dbar].

    Returns
    -------
    SP : ndarray
        Practical salinity (PRACSAL_L2) [PSS-78, unitless].

    Notes
    -----
    Conductivity is converted from S m^-1 to mS cm^-1 (multiply by 10)
    before passing to gsw.SP_from_C(C, t, p), which expects mS cm^-1.
    """
    # Convert L1 Conductivity from S/m to mS/cm
    C10 = c * 10.0

    # Calculate the Practical Salinity (PSS-78) [unitless]
    SP = gsw.SP_from_C(C10, t, p)
    return SP

History

Date Author Change
2013-03-13 Christopher Wingard Initial code
2013-05-10 Christopher Wingard Minor edits to comments
2014-01-31 Russell Desiderio Standardized comment format
2023-08-15 Samuel Dahlberg Replaced incompatible pygsw with GSW library
2025-04-17 Christopher Wingard Converted to NumPy docstring format; updated documentation

ctd_density(SP, t, p, lat, lon)

Compute DENSITY_L2 from practical salinity, temperature, pressure, and position using TEOS-10.

Calculates in-situ seawater density using the TEOS-10 GSW library via a three-step chain: Practical Salinity to Absolute Salinity, in-situ temperature to Conservative Temperature, then density from the computationally-efficient 48-term expression.

Parameters:
  • SP (array_like) –

    Practical salinity (PRACSAL_L2) [PSS-78, unitless].

  • t (array_like) –

    Seawater temperature (TEMPWAT_L1) [degC, ITS-90].

  • p (array_like) –

    Seawater pressure (PRESWAT_L1) [dbar].

  • lat (array_like) –

    Latitude of observation [decimal degrees north].

  • lon (array_like) –

    Longitude of observation [decimal degrees east].

Returns:
  • rho( ndarray ) –

    In-situ seawater density (DENSITY_L2) [kg m^-3].

Notes

TEOS-10 processing chain:

SA  = gsw.SA_from_SP(SP, p, lon, lat)   (Absolute Salinity)
CT  = gsw.CT_from_t(SA, t, p)           (Conservative Temp)
rho = gsw.rho(SA, CT, p)                (density, 48-term)

For moored instruments, latitude (lat) and longitude (lon) are from the mooring position metadata; for gliders, they are the vehicle's position at the time of each measurement.

Source code in ion_functions/data/ctd_functions.py
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
def ctd_density(SP, t, p, lat, lon):
    """
    Compute DENSITY_L2 from practical salinity, temperature, pressure,
    and position using TEOS-10.

    Calculates in-situ seawater density using the TEOS-10 GSW library
    via a three-step chain: Practical Salinity to Absolute Salinity,
    in-situ temperature to Conservative Temperature, then density from
    the computationally-efficient 48-term expression.

    Parameters
    ----------
    SP : array_like
        Practical salinity (PRACSAL_L2) [PSS-78, unitless].
    t : array_like
        Seawater temperature (TEMPWAT_L1) [degC, ITS-90].
    p : array_like
        Seawater pressure (PRESWAT_L1) [dbar].
    lat : array_like
        Latitude of observation [decimal degrees north].
    lon : array_like
        Longitude of observation [decimal degrees east].

    Returns
    -------
    rho : ndarray
        In-situ seawater density (DENSITY_L2) [kg m^-3].

    Notes
    -----
    TEOS-10 processing chain:

        SA  = gsw.SA_from_SP(SP, p, lon, lat)   (Absolute Salinity)
        CT  = gsw.CT_from_t(SA, t, p)           (Conservative Temp)
        rho = gsw.rho(SA, CT, p)                (density, 48-term)

    For moored instruments, latitude (lat) and longitude (lon)
    are from the mooring position metadata; for gliders, they are
    the vehicle's position at the time of each measurement.
    """
    # Calculate the density [kg m-3]
    sa = gsw.SA_from_SP(SP, p, lon, lat)  # absolute salinity
    ct = gsw.CT_from_t(sa, t, p)  # conservative temperature
    rho = gsw.rho(sa, ct, p)  # density
    return rho

History

Date Author Change
2013-03-11 Christopher Mueller Initial code
2013-03-13 Christopher Wingard Added commenting; moved to ctd_functions
2013-05-10 Christopher Wingard Minor edits to comments
2014-01-31 Russell Desiderio Standardized comment format
2023-08-15 Samuel Dahlberg Replaced incompatible pygsw with GSW library
2025-04-17 Christopher Wingard Converted to NumPy docstring format; updated documentation

References

Feistel, R. (2008). A Gibbs function for seawater thermodynamics for -6 to 80 degrees C and salinity up to 120 g/kg. Deep Sea Research I, 55, 1639-1671.

Fofonoff, N. P. and Millard, R. C. (1983). Algorithms for computation of fundamental properties of seawater. UNESCO Technical Papers in Marine Science, 44, 1-53.

Hill, K. D., Dauphinee, T. M., and Woods, D. J. (1986). The extension of the Practical Salinity Scale 1978 to low salinities. IEEE Journal of Oceanic Engineering, OE-11(1), 109-112.

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, 196 pp.

McDougall, T. J., Barker, P. M., Feistel, R., and Jackett, D. R. (2011). A computationally efficient 48-term expression for the density of seawater in terms of Conservative Temperature, and related properties of seawater. Journal of Atmospheric and Oceanic Technology, 28, 1464-1477.

McDougall, T. J., Jackett, D. R., and Millero, F. J. (2010). An algorithm for estimating Absolute Salinity in the global ocean. Ocean Science Discussions, 6, 215-242.

Millero, F. J., et al. (2008). The composition of Standard Seawater and the definition of the Reference-Composition Salinity Scale. Deep Sea Research I, 55, 50-72.

Pawlowicz, R. (2010). What every oceanographer needs to know about TEOS-10 (The TEOS-10 Primer). Thermodynamic Equation Of Seawater - 2010 (TEOS-10) website: http://www.teos-10.org/.

Sea-Bird Electronics (2008). Application Note 10. Compressibility Compensation of Sea-Bird Conductivity Sensors. Revision March 2008. Sea-Bird Electronics, Inc.

Sea-Bird Electronics (2009). SBE 16plus V2 SEACAT User's Manual. Manual Version 005. Sea-Bird Electronics, Inc.

Sea-Bird Electronics (2010). Application Note 42. ITS-90 Temperature Scale. Revision February 2010. Sea-Bird Electronics, Inc.

Sea-Bird Electronics (2011). SBE 37-IM MicroCAT User's Manual. Manual Version 027. Sea-Bird Electronics, Inc.

OOI (2013). Data Product Specification for Water Temperature. Document Control Number 1341-00010.

OOI (2013). Data Product Specification for Pressure (Depth). Document Control Number 1341-00020.

OOI (2013). Data Product Specification for Conductivity. Document Control Number 1341-00030.

OOI (2013). Data Product Specification for Salinity. Document Control Number 1341-00040.

OOI (2012). Data Product Specification for Density. Document Control Number 1341-00050.