ADV Module
This module contains routines for reading and working with ADV data. It contains:
Read a binary Nortek (e.g., .VEC, .wpr, .ad2cp, etc.) or RDI (.000, .PD0, .ENX, etc.) data file. |
|
Load xarray dataset from netCDF (.nc) |
|
Rotate a dataset to a new coordinate system. |
|
Set the instrument to head rotation matrix for the Nortek ADV if it hasn't already been set through a '.userdata.json' file. |
|
Compute the principal angle of the horizontal velocity. |
|
Module containing functions to clean data |
|
This function performs motion correction on an IMU-ADV data object. |
|
This is the base binning (averaging) tool. |
|
A class that builds upon VelBinner for calculating turbulence statistics and velocity spectra from ADV data |
|
Functional version of ADVBinner that computes a suite of turbulence statistics for the input dataset, and returns a binned data object. |
Quick Example
# Start by importing DOLfYN:
import dolfyn
import dolfyn.adv.api as api
# Then read a file containing adv data:
dat = dolfyn.read_example('vector_data01.VEC')
# Clean the file using the Goring+Nikora method:
mask = api.clean.GN2002(dat.vel)
dat['vel'] = api.clean.clean_fill(dat.vel, mask, npt=12, method='cubic')
# Rotate that data from the instrument to earth frame:
# First set the magnetic declination
dolfyn.set_declination(dat, 10) # 10 degrees East
dolfyn.rotate2(dat, 'earth')
# Rotate it into a 'principal axes frame':
# First calculate the principal heading
dat.attrs['principal_heading'] = dolfyn.calc_principal_heading(dat.vel)
dolfyn.rotate2(dat, 'principal')
# Define an averaging object, and create an 'ensembled' data set:
binner = api.ADVBinner(n_bin=9600, fs=dat.fs, n_fft=4096)
dat_binned = binner(dat)
# At any point you can save the data:
#dolfyn.save(dat_binned, 'adv_data.nc')
# And reload the data:
#dat_bin_copy = dolfyn.load('adv_data.nc')
Cleaning Data
Interpolate over mask values in timeseries data using the specified method |
|
Fill missing values with the ensemble mean. |
|
Returns a logical vector where a spike in u of magnitude greater than thresh occurs. |
|
Returns a logical vector that is True where the values of u are outside of range. |
|
The Goring & Nikora 2002 'despiking' method, with Wahl2003 correction. |
Module containing functions to clean data
- dolfyn.adv.clean.clean_fill(u, mask, npt=12, method='cubic', maxgap=6)[source]
Interpolate over mask values in timeseries data using the specified method
- Parameters:
u (xarray.DataArray) – The dataArray to clean.
mask (bool) – Logical tensor of elements to “nan” out (from spikeThresh, rangeLimit, or GN2002) and replace
npt (int) – The number of points on either side of the bad values that interpolation occurs over
method (string (default: 'cubic')) – Interpolation scheme to use (linear, cubic, pchip, etc)
maxgap (int (default: 6)) – Max number of consective nan’s to interpolate across
- Returns:
da (xarray.DataArray) – The dataArray with nan’s filled in
See also
xarray.DataArray.interpolate_na
- dolfyn.adv.clean.fill_nan_ensemble_mean(u, mask, fs, window)[source]
Fill missing values with the ensemble mean.
- Parameters:
u (xarray.DataArray (..., time)) – The dataArray to clean. Can be 1D or 2D.
mask (bool) – Logical tensor of elements to “nan” out (from spikeThresh, rangeLimit, or GN2002) and replace
fs (int) – Instrument sampling frequency
window (int) – Size of window in seconds used to calculate ensemble means
- Returns:
da (xarray.DataArray) – The dataArray with nan’s filled in
Notes
Gaps larger than the ensemble size will not get filled in.
- dolfyn.adv.clean.spike_thresh(u, thresh=10)[source]
Returns a logical vector where a spike in u of magnitude greater than thresh occurs. Both ‘Negative’ and ‘positive’ spikes are found.
- Parameters:
u (xarray.DataArray) – The timeseries data to clean.
thresh (int (default: 10)) – Magnitude of velocity spike, must be positive.
- Returns:
mask (|np.ndarray|) – Logical vector with spikes labeled as ‘True’
- dolfyn.adv.clean.range_limit(u, range=[-5, 5])[source]
Returns a logical vector that is True where the values of u are outside of range.
- Parameters:
u (xarray.DataArray) – The timeseries data to clean.
range (list (default: [-5, 5])) – Min and max magnitudes beyond which are masked
- Returns:
mask (|np.ndarray|) – Logical vector with spikes labeled as ‘True’
- dolfyn.adv.clean.GN2002(u, npt=5000)[source]
The Goring & Nikora 2002 ‘despiking’ method, with Wahl2003 correction. Returns a logical vector that is true where spikes are identified.
- Parameters:
u (xarray.DataArray) – The velocity array (1D or 3D) to clean.
npt (int (default: 5000)) – The number of points over which to perform the method.
- Returns:
mask (|np.ndarray|) – Logical vector with spikes labeled as ‘True’
- class dolfyn.adv.motion.CalcMotion(ds, accel_filtfreq=None, vel_filtfreq=None, to_earth=True)[source]
Bases:
object
A ‘calculator’ for computing the velocity of points that are rigidly connected to an ADV-body with an IMU.
- Parameters:
ds (xarray.Dataset) – The IMU-adv data that will be used to compute motion.
accel_filtfreq (float) – the frequency at which to high-pass filter the acceleration sigal to remove low-frequency drift. (default: 0.03 Hz)
vel_filtfreq (float (optional)) – a second frequency to high-pass filter the integrated acceleration. (default: 1/3 of accel_filtfreq)
- calc_velacc()[source]
Calculates the translational velocity from the high-pass filtered acceleration signal.
- Returns:
velacc (numpy.ndarray (3 x n_time)) – The acceleration-induced velocity array (3, n_time).
- calc_velrot(vec, to_earth=None)[source]
Calculate the induced velocity due to rotations of the instrument about the IMU center.
- Parameters:
vec (numpy.ndarray (len(3) or 3 x M)) – The vector in meters (or vectors) from the body-origin (center of head end-cap) to the point of interest (in the body coord-sys).
- Returns:
velrot (numpy.ndarray (3 x M x N_time)) – The rotation-induced velocity array (3, n_time).
- dolfyn.adv.motion.correct_motion(ds, accel_filtfreq=None, vel_filtfreq=None, to_earth=True, separate_probes=False)[source]
This function performs motion correction on an IMU-ADV data object. The IMU and ADV data should be tightly synchronized and contained in a single data object.
- Parameters:
ds (xarray.Dataset) – Cleaned ADV dataset in “inst” coordinates
accel_filtfreq (float) – the frequency at which to high-pass filter the acceleration sigal to remove low-frequency drift.
vel_filtfreq (float (optional)) – a second frequency to high-pass filter the integrated acceleration. (default: 1/3 of accel_filtfreq)
to_earth (bool (optional, default: True)) – All variables in the ds.props[‘rotate_vars’] list will be rotated into either the earth frame (to_earth=True) or the instrument frame (to_earth=False).
separate_probes (bool (optional, default: False)) – a flag to perform motion-correction at the probe tips, and perform motion correction in beam-coordinates, then transform back into XYZ/earth coordinates. This correction seems to be lower than the noise levels of the ADV, so the default is to not use it (False).
- Returns:
This function returns None, it operates on the input data object,
ds
. The following attributes are added to ds –velraw
is the uncorrected velocityvelrot
is the rotational component of the head motion (fromangrt)
velacc
is the translational component of the head motion (fromaccel, the high-pass filtered accel sigal)
acclow
is the low-pass filtered accel sigal (i.e.,The primary velocity vector attribute,
vel
, is motion correctedsuch that – vel = velraw + velrot + velacc
The sigs are correct in this equation. The measured velocity
induced by head-motion is *in the opposite direction of the head*
motion. i.e. when the head moves one way in stationary flow, it
measures a velocity in the opposite direction. Therefore, to
remove the motion from the raw sigal we *add the head motion.*
Notes
Acceleration signals from inertial sensors are notorious for having a small bias that can drift slowly in time. When integrating these sigals to estimate velocity the bias is amplified and leads to large errors in the estimated velocity. There are two methods for removing these errors,
high-pass filter the acceleration sigal prior and/or after integrating. This implicitly assumes that the low-frequency translational velocity is zero.
provide a slowly-varying reference position (often from a GPS) to an IMU that can use the sigal (usually using Kalman filters) to debias the acceleration sigal.
Because method (1) removes real low-frequency acceleration, method (2) is more accurate. However, providing reference position estimates to undersea instruments is practically challenging and expensive. Therefore, lacking the ability to use method (2), this function utilizes method (1).
For deployments in which the ADV is mounted on a mooring, or other semi-fixed structure, the assumption of zero low-frequency translational velocity is a reasonable one. However, for deployments on ships, gliders, or other moving objects it is not. The measured velocity, after motion-correction, will still hold some of this contamination and will be a sum of the ADV motion and the measured velocity on long time scales. If low-frequency motion is known separate from the ADV (e.g. from a bottom-tracking ADP, or from a ship’s GPS), it may be possible to remove that sigal from the ADV sigal in post-processing.