Rotate Functions
Contains functions for rotating data through frames of reference (FoR):
‘beam’: Follows the acoustic beam FoR, where velocity data is organized by beam number 1-3 or 1-4.
‘inst’: The instrument’s XYZ Cartesian directions. For ADVs, this orientation is from the mark on the ADV body/battery canister, not the sensor head. For TRDI 4-beam instruments, the fourth velocity term is the error velocity (aka XYZE). For Nortek 4-beam instruments, this is XYZ1 Z2, where E=Z2-Z1.
‘earth’: East North UP (ENU) FoR. Based on either magnetic or true North, depending on whether or not DOLfYN has a magnetic declination associated with the dataset. Instruments do not internally record magnetic declination, unless it has been supplied via external software like TRDI’s VMDAS.
‘principal’: Rotates velocity data into a streamwise, cross-stream, and vertical FoR based on the principal flow direction. One must calculate principal heading first.
Rotate a dataset to a new coordinate system. |
|
Set the magnetic declination |
|
Compute the principal angle of the horizontal velocity. |
|
Set the instrument to head rotation matrix for the Nortek ADV if it hasn't already been set through a '.userdata.json' file. |
|
Calculate the orientation matrix from DOLfYN-defined euler angles. |
|
Calculate DOLfYN-defined euler angles from the orientation matrix. |
|
Calculate orientation from Nortek AHRS quaternions, where q = [W, X, Y, Z] instead of the standard q = [X, Y, Z, W] = [q1, q2, q3, q4] |
|
Calculate "tilt", the vertical inclination, from pitch and roll. |
These functions pertain to both ADCPs and ADVs:
>> import dolfyn
>> dat = dolfyn.read_example('burst_mode01.VEC')
>> dolfyn.set_declination(dat, 12)
>> dolfyn.rotate2(dat, 'earth')
>> dat.attrs['principal_heading'] = dolfyn.calc_principal_heading(dat['vel'])
>> dolfyn.rotate2(dat, 'principal')
- dolfyn.rotate.api.rotate2(ds, out_frame='earth', inplace=True)[source]
Rotate a dataset to a new coordinate system.
- Parameters:
ds (xr.Dataset) – The dolfyn dataset (ADV or ADCP) to rotate.
out_frame (string {'beam', 'inst', 'earth', 'principal'}) – The coordinate system to rotate the data into.
inplace (bool (default: True)) – When True
ds
is modified. When False a copy is returned.
- Returns:
ds (xarray.Dataset or None) – Returns a new rotated dataset when ``inplace=False``, otherwise returns None.
Notes
This function rotates all variables in
ds.attrs['rotate_vars']
.In order to rotate to the ‘principal’ frame, a value should exist for
ds.attrs['principal_heading']
. The functioncalc_principal_heading
is recommended for this purpose, e.g.:ds.attrs[‘principal_heading’] = dolfyn.calc_principal_heading(ds[‘vel’].mean(range))
where here we are using the depth-averaged velocity to calculate the principal direction.
- dolfyn.rotate.api.calc_principal_heading(vel, tidal_mode=True)[source]
Compute the principal angle of the horizontal velocity.
- Parameters:
vel (np.ndarray (2,...,Nt), or (3,...,Nt)) – The 2D or 3D velocity array (3rd-dim is ignored in this calculation)
tidal_mode (bool (default: True)) – If true, range is set from 0 to +/-180 degrees. If false, range is 0 to 360 degrees
- Returns:
p_heading (float or ndarray) – The principal heading in degrees clockwise from North.
Notes
When tidal_mode=True, this tool calculates the heading that is aligned with the bidirectional flow. It does so following these steps:
rotates vectors with negative velocity by 180 degrees
then doubles those angles to make a complete circle again
computes a mean direction from this, and halves that angle (to undo the doubled-angles in step 2)
The returned angle is forced to be between 0 and 180. So, you may need to add 180 to this if you want your positive direction to be in the western-half of the plane.
Otherwise, this function simply computes the average direction using a vector method.
- dolfyn.rotate.api.set_declination(ds, declin, inplace=True)[source]
Set the magnetic declination
- Parameters:
ds (xarray.Dataset or
dolfyn.velocity.Velocity
) – The input dataset or velocity classdeclination (float) – The value of the magnetic declination in degrees (positive values specify that Magnetic North is clockwise from True North)
inplace (bool (default: True)) – When True
ds
is modified. When False a copy is returned.
- Returns:
ds (xarray.Dataset or None) – Returns a new dataset with declination set when ``inplace=False``, otherwise returns None.
Notes
This function modifies the data object in the following ways:
If the dataset is in the earth reference frame at the time of setting declination, it will be rotated into the “True-East, True-North, Up” (hereafter, ETU) coordinate system
dat['orientmat']
is modified to be an ETU to instrument (XYZ) rotation matrix (rather than the magnetic-ENU to XYZ rotation matrix). Therefore, all rotations to/from the ‘earth’ frame will now be to/from this ETU coordinate system.The value of the specified declination will be stored in
dat.attrs['declination']
dat['heading']
is adjusted for declination (i.e., it is relative to True North).If
dat.attrs['principal_heading']
is set, it is adjusted to account for the orientation of the new ‘True’ earth coordinate system (i.e., calling set_declination on a data object in the principal coordinate system, then calling dat.rotate2(‘earth’) will yield a data object in the new ‘True’ earth coordinate system)
- dolfyn.rotate.api.set_inst2head_rotmat(ds, rotmat, inplace=True)[source]
Set the instrument to head rotation matrix for the Nortek ADV if it hasn’t already been set through a ‘.userdata.json’ file.
- Parameters:
ds (xarray.Dataset) – The data set to assign inst2head_rotmat
rotmat (float) – 3x3 rotation matrix
inplace (bool (default: True)) – When True
ds
is modified. When False a copy is returned.
- Returns:
ds (xarray.Dataset or None) – Returns a new dataset with inst2head_rotmat set when ``inplace=False``, otherwise returns None.
Notes
If the data object is in earth or principal coords, it is first rotated to ‘inst’ before assigning inst2head_rotmat, it is then rotated back to the coordinate system in which it was input. This way the inst2head_rotmat gets applied correctly (in inst coordinate system).
- dolfyn.rotate.base.euler2orient(heading, pitch, roll, units='degrees')[source]
Calculate the orientation matrix from DOLfYN-defined euler angles.
This function is not likely to be called during data processing since it requires DOLfYN-defined euler angles. It is intended for testing DOLfYN.
The matrices H, P, R are the transpose of the matrices for rotation about z, y, x as shown here https://en.wikipedia.org/wiki/Rotation_matrix. The transpose is used because in DOLfYN the orientation matrix is organized for rotation from EARTH –> INST, while the wiki’s matrices are organized for rotation from INST –> EARTH.
- Parameters:
heading (
numpy.ndarray
(Nt)) – The heading angle.pitch (
numpy.ndarray
(Nt)) – The pitch angle.roll (
numpy.ndarray
(Nt)) – The roll angle.units (str {'degrees' (default), 'radians'}) –
- Returns:
omat (|np.ndarray| (3x3xNt)) – The orientation matrix of the data. The returned orientation matrix obeys the following conventions:
a “ZYX” rotation order. That is, these variables are computed assuming that rotation from the earth -> instrument frame happens by rotating around the z-axis first (heading), then rotating around the y-axis (pitch), then rotating around the x-axis (roll). Note this requires matrix multiplication in the reverse order.
heading is defined as the direction the x-axis points, positive clockwise from North (this is opposite the right-hand-rule around the Z-axis), range 0-360 degrees.
pitch is positive when the x-axis pitches up (this is opposite the right-hand-rule around the Y-axis)
roll is positive according to the right-hand-rule around the instrument’s x-axis
- dolfyn.rotate.base.orient2euler(omat)[source]
Calculate DOLfYN-defined euler angles from the orientation matrix.
- Parameters:
omat (
numpy.ndarray
) – The orientation matrix- Returns:
heading (|np.ndarray|) – The heading angle. Heading is defined as the direction the x-axis points, positive clockwise from North (this is opposite the right-hand-rule around the Z-axis), range 0-360 degrees.
pitch (np.ndarray) – The pitch angle (degrees). Pitch is positive when the x-axis pitches up (this is opposite the right-hand-rule around the Y-axis).
roll (np.ndarray) – The roll angle (degrees). Roll is positive according to the right-hand-rule around the instrument’s x-axis.
- dolfyn.rotate.base.quaternion2orient(quaternions)[source]
Calculate orientation from Nortek AHRS quaternions, where q = [W, X, Y, Z] instead of the standard q = [X, Y, Z, W] = [q1, q2, q3, q4]
- Parameters:
quaternions (xarray.DataArray) – Quaternion dataArray from the raw dataset
- Returns:
orientmat (|np.ndarray|) – The earth2inst rotation maxtrix as calculated from the quaternions
See also
scipy.spatial.transform.Rotation
- dolfyn.rotate.base.calc_tilt(pitch, roll)[source]
Calculate “tilt”, the vertical inclination, from pitch and roll.
- Parameters:
roll (numpy.ndarray or xarray.DataArray) – Instrument roll
pitch (numpy.ndarray or xarray.DataArray) – Instrument pitch
- Returns:
tilt (numpy.ndarray) – Vertical inclination of the instrument