ctdproc package

Subpackages

Submodules

ctdproc.calcs module

ctdproc.calcs.calc_allsal(c, t, p, lon, lat)[source]

Calculate absolute and practical salinity.

Wrapper for gsw functions. Converts conductivity from S/m to mS/cm if output salinity is less than 5.

Parameters:
  • c (array-like) – Conductivity. See notes on units above.

  • t (array-like) – In-situ temperature (ITS-90), degrees C

  • p (array-like) – Sea pressure, dbar

  • lon (array-like) – Longitude, -360 to 360 degrees

  • lat (array-like) – Latitude, -90 to 90 degrees

Returns:

  • SA (array-like, g/kg) – Absolute Salinity

  • SP (array-like) – Practical Salinity

ctdproc.calcs.calc_depth(data)[source]
ctdproc.calcs.calc_sal(data)[source]
ctdproc.calcs.calc_sigma(data)[source]
ctdproc.calcs.calc_temp(data)[source]
ctdproc.calcs.swcalcs(data)[source]

Calculate SA, s, CT, theta, sigma, oxygen, depth, z

Using gsw package. Variable names for practical salinity and potential temperature are kept consistent with the Matlab processing package.

Salinity is calculated both as absolute salinity (SA) and practical salinity (s). Note: Per TEOS-10 instructions, the variable to be archived is practical salinity.

Temperature is converted to conservative temperature (CT) and potential temperature (th). Potential temperature is referenced to 0 dbar.

Potential density anomaly is referenced to 0 dbar pressure. The anomaly is calculated as potential density minus 1000 kg/m^3.

Parameters:

data (array-like) – CTD time series data structure.

Returns:

data – CTD time series data structure.

Return type:

array-like

ctdproc.calcs.wsink(p, Ts, Fs)[source]

Compute sinking velocity from pressure record.

Computes the sinking (or rising) velocity from the pressure signal p by first differencing. The pressure signal is smoothed with a low-pass filter for differentiation. If the input signal is shorter than the smoothing time scale, w is taken as the slope of the linear regression of p.

Adapted from wsink.m - Fabian Wolk, Rockland Oceanographic Services Inc.

Parameters:
  • p (array-like) – Pressure [dbar]

  • Ts (float) – Smoothing time scale [s]

  • Fs (float) – Sampling frequency [Hz]

Returns:

w – Sinking velocity [dbar/s]

Return type:

array-like

ctdproc.helpers module

ctdproc.helpers.atanfit(x, f, Phi, W)[source]
ctdproc.helpers.datetime2mtlb(dt)[source]
ctdproc.helpers.findsegments(ibad)[source]

Find contiguous segments in an array of indices.

Parameters:

ibad (np.array) – Array with indices

Returns:

  • istart (np.array) – Segment start indices

  • istop (np.array) – Segment stop indices

  • seglength (np.array) – Segment length

ctdproc.helpers.glitchcorrect(x, diffx, prodx, ibefore=0, iafter=0)[source]

Remove glitches/spikes in array.

Adapted from tms_tc_glitchcorrect.m

Parameters:
  • x (np.array) – Input array

  • diffx (int) – Threshold for differences

  • prodx (int) – Threshold for products

  • ibefore (int, optional) – Number of elements to interpolate on left side of glitch

  • after (int, optional) – Number of elements to interpolate on right side of glitch

Returns:

y – Interpolated array

Return type:

np.array

ctdproc.helpers.inearby(ibad, inearm, inearp, n)[source]

Find indices left and right of given indices and add them to the index array.

Differs from inearby.m as it breaks up ibad into segments and then finds nearby indices for each segment.

Parameters:
  • ibad (np.array) – Array with indices

  • inearm (int) – Number of elements to include on left side

  • inearp (int) – Number of elements to include on right side

  • n (int) – Maximum index in output array and size of array being indexed

Returns:

k – New array of indices

Return type:

np.array

ctdproc.helpers.interpbadsegments(x, ibad)[source]

Interpolate over segments of bad data.

Parameters:
  • x (np.array) – Input data

  • ibad (np.array) – Indices of bad data

Returns:

y – Interpolated array

Return type:

np.array

ctdproc.helpers.mtlb2datetime(matlab_datenum, strip_microseconds=False, strip_seconds=False)[source]

Convert Matlab datenum format to python datetime. This version also works for vector input and strips milliseconds if desired.

Parameters:
  • matlab_datenum (float or np.array) – Matlab time vector.

  • strip_microseconds (bool) – Get rid of microseconds (optional)

  • strip_seconds (bool) – Get rid of seconds (optional)

Returns:

t – Time in numpy’s datetime64 format.

Return type:

np.datetime64

ctdproc.helpers.pad_lr(p, nPad)[source]

Pad array left and right

ctdproc.helpers.preen(x, xmin, xmax)[source]

Eliminate values outside given range and interpolate.

Parameters:
  • x (np.array) – Input array

  • xmin (float) – Lower limit

  • xmax (float) – Upper limit

Returns:

xp – Cleaned array

Return type:

np.array

ctdproc.helpers.unique_arrays(*arrays)[source]

Find unique elements in more than one numpy array.

Parameters:

arrays (A tuple of np.array.ndarray)

Returns:

unique – numpy array with unique elements from the input arrays

Return type:

np.array

ctdproc.io module

class ctdproc.io.CTDHex(filename)[source]

Bases: object

Converter for Seabird CTD data in hex format. Initialize with full path to hex file. xml config file needs to be located in the same directory.

TODO:
  • Add oxygen hysteresis

  • Convert fluorometer voltage

mattime_to_datetime64(dtnum)[source]

Convert Matlab time format to numpy datetime64 time format.

mattime_to_sbetime(dt)[source]

Convert matlab time format to SBE time format.

parse_hex()[source]
physicalunits()[source]
read_xml_config()[source]

Read xml config file.

sbetime_to_mattime(dt)[source]

Convert SBE time format to matlab time format.

to_mat(matname)[source]

Save data in Matlab format.

to_xarray()[source]

Convert data into xarray.Dataset.

Returns:

ds – CTD time series in Dataset format.

Return type:

xarray.Dataset

ctdproc.io.CTDx(filename)[source]
ctdproc.io.add_default_proc_params(ds)[source]
ctdproc.io.prof_to_mat(matname, datad, datau)[source]

ctdproc.proc module

ctdproc.proc.add_tcfit_default(ds)[source]

Get default values for tc fit range depending on depth of cast.

Range for tc fit is 200dbar to maximum pressure if the cast is shallower than 1000dbar, 500dbar to max pressure otherwise.

Parameters:

ds (xarray.Dataset) – CTD time series data structure

Returns:

tcfit – Upper and lower limit for tc fit in phase_correct.

Return type:

tuple

ctdproc.proc.bincast(ds, dz, zmin, zmax)[source]

Depth-bin CTD time series.

Parameters:
  • ds (xr.Dataset) – CTD time series

  • dz (float) – Bin size [m]

  • zmin (float) – Minimum bin depth center [m]

  • zmax (float) – Maximum bin depth center [m]

Returns:

ds – Depth-binned CTD profile

Return type:

xr.Dataset

ctdproc.proc.cleanup(ds)[source]

Clean up CTD raw time series.

  • despike pressure

  • eliminate data near surface

  • remove spikes in other data

  • remove smaller T, C, glitches

ctdproc.proc.cleanup_ud(ds)[source]

More cleaning and calculation of derived variables.

ctdproc.proc.despike(da, spike_threshold)[source]

Set spikes to NaN.

ctdproc.proc.phase_correct(ds)[source]

Bring temperature and conductivity in phase.

Parameters:

ds (dtype) – description

Returns:

ds – description

Return type:

dtype

ctdproc.proc.plot_TS(ds)[source]
ctdproc.proc.plot_profile(ds, var_list)[source]
ctdproc.proc.preen_ctd(ds)[source]

Remove spikes in p, t1, t2, c1, c2.

Parameters:

data (xarray.Dataset) – Dataset with CTD time series.

Returns:

data – Cleaned dataset.

Return type:

xarray.Dataset

ctdproc.proc.remove_out_of_bounds(da, bmin, bmax)[source]

Remove out of bounds data.

ctdproc.proc.rmloops(ds)[source]

Eliminate depth loops in CTD data based on sinking velocity.

All data variables are set to NaN within depth loops.

Parameters:

ds (xarray.Dataset) – CTD time series dataset

Returns:

ds – CTD time series dataset

Return type:

xarray.Dataset

ctdproc.proc.run_all(ds)[source]

Run all standard processing steps on raw CTD time series.

Parameters:

data (xarray.Dataset) – CTD time series data structure

Returns:

  • datad (xarray.Dataset) – Downcast time series

  • datau (xarray.Dataset) – Upcast time series

ctdproc.proc.split_updn(ds)[source]

Separate into down/up-casts and apply corrections.

## TODO: (should this really be here?) - tc lag correction - thermal mass correction - lowpass-filter T, C, oxygen

Parameters:

data (xarray.Dataset) – CTD time series

Returns:

  • datad (xarray.Dataset) – CTD time series for downcast

  • datau (xarray.Dataset) – CTD time series for upcast

Module contents

Library for CTD data processing.