clouddrift.kinematics#

Functions for kinematic computations.

Functions

inertial_oscillation_from_position(...[, ...])

Extract inertial oscillations from consecutive geographical positions.

kinetic_energy(u[, v])

Compute kinetic energy from zonal and meridional velocities.

position_from_velocity(u, v, time, x_origin, ...)

Compute positions from arrays of velocities and time and a pair of origin coordinates.

residual_position_from_displacement(...)

Return residual longitudes and latitudes along a trajectory on the spherical Earth after correcting for zonal and meridional displacements x and y in meters.

spin(u, v, time[, difference_scheme, time_axis])

Compute spin continuously from velocities and times.

velocity_from_position(x, y, time[, ...])

Compute velocity from arrays of positions and time.

clouddrift.kinematics.inertial_oscillation_from_position(longitude: ndarray, latitude: ndarray, relative_bandwidth: float | None = None, wavelet_duration: float | None = None, time_step: float | None = 3600.0, relative_vorticity: float | ndarray | None = 0.0) ndarray[source]#

Extract inertial oscillations from consecutive geographical positions.

This function acts by performing a time-frequency analysis of horizontal displacements with analytic Morse wavelets. It extracts the portion of the wavelet transform signal that follows the inertial frequency (opposite of Coriolis frequency) as a function of time, potentially shifted in frequency by a measure of relative vorticity. The result is a pair of zonal and meridional relative displacements in meters.

This function is equivalent to a bandpass filtering of the horizontal displacements. The characteristics of the filter are defined by the relative bandwidth of the wavelet transform or by the duration of the wavelet, see the parameters below.

Parameters#

longitudearray-like

Longitude sequence. Unidimensional array input.

latitudearray-like

Latitude sequence. Unidimensional array input.

relative_bandwidthfloat, optional

Bandwidth of the frequency-domain equivalent filter for the extraction of the inertial oscillations; a number less or equal to one which is a fraction of the inertial frequency. A value of 0.1 leads to a bandpass filter equivalent of +/- 10 percent of the inertial frequency.

wavelet_durationfloat, optional

Duration of the wavelet, or inverse of the relative bandwidth, which can be passed instead of the relative bandwidth.

time_stepfloat, optional

The constant time interval between data points in seconds. Default is 3600.

relative_vorticity: Optional, float or array-like

Relative vorticity adding to the local Coriolis frequency. If “f” is the Coriolis frequency then “f” + relative_vorticity will be the effective Coriolis frequency as defined by Kunze (1985). Positive values correspond to cyclonic vorticity, irrespectively of the latitudes of the data points.

Returns#

xhatarray-like

Zonal relative displacement in meters from inertial oscillations.

yhatarray-like

Meridional relative displacement in meters from inertial oscillations.

Examples#

To extract displacements from inertial oscillations from sequences of longitude and latitude values, equivalent to bandpass around 20 percent of the local inertial frequency:

>>> xhat, yhat = inertial_oscillation_from_position(longitude, latitude, relative_bandwidth=0.2)

The same result can be obtained by specifying the wavelet duration instead of the relative bandwidth:

>>> xhat, yhat = inertial_oscillation_from_position(longitude, latitude, wavelet_duration=5)

Next, the residual positions from the inertial displacements can be obtained with another function:

>>> residual_longitudes, residual_latitudes = residual_position_from_displacement(longitude, latitude, xhat, yhat)

Raises#

ValueError

If longitude and latitude arrays do not have the same shape. If relative_vorticity is an array and does not have the same shape as longitude and latitude. If time_step is not a float. If both relative_bandwidth and wavelet_duration are specified. If neither relative_bandwidth nor wavelet_duration are specified. If the absolute value of relative_bandwidth is not in the range (0,1]. If the wavelet duration is not greater than or equal to 1.

See Also#

residual_position_from_displacement(), wavelet_transform, morse_wavelet

clouddrift.kinematics.kinetic_energy(u: float | list | ndarray | DataArray | Series, v: float | list | ndarray | DataArray | Series | None = None) float | ndarray | DataArray[source]#

Compute kinetic energy from zonal and meridional velocities.

Parameters#

ufloat or array-like

Zonal velocity.

vfloat or array-like, optional.

Meridional velocity. If not provided, the flow is assumed one-dimensional in time and defined by u.

Returns#

kefloat or array-like

Kinetic energy.

Examples#

>>> import numpy as np
>>> from clouddrift.kinematics import kinetic_energy
>>> u = np.array([1., 2., 3., 4.])
>>> v = np.array([1., 1., 1., 1.])
>>> kinetic_energy(u, v)
array([1. , 2.5, 5. , 8.5])
>>> u = np.reshape(np.tile([1., 2., 3., 4.], 2), (2, 4))
>>> v = np.reshape(np.tile([1., 1., 1., 1.], 2), (2, 4))
>>> kinetic_energy(u, v)
array([[1. , 2.5, 5. , 8.5],
       [1. , 2.5, 5. , 8.5]])
clouddrift.kinematics.position_from_velocity(u: ndarray, v: ndarray, time: ndarray, x_origin: float, y_origin: float, coord_system: str | None = 'spherical', integration_scheme: str | None = 'forward', time_axis: int | None = -1) tuple[ndarray, ndarray][source]#

Compute positions from arrays of velocities and time and a pair of origin coordinates.

The units of the result are degrees if coord_system == "spherical" (default). If coord_system == "cartesian", the units of the result are equal to the units of the input velocities multiplied by the units of the input time. For example, if the input velocities are in meters per second and the input time is in seconds, the units of the result will be meters.

Integration scheme can take one of three values:

  1. “forward” (default): integration from x[i] to x[i+1] is performed

    using the velocity at x[i].

  2. “backward”: integration from x[i] to x[i+1] is performed using the

    velocity at x[i+1].

  3. “centered”: integration from x[i] to x[i+1] is performed using the

    arithmetic average of the velocities at x[i] and x[i+1]. Note that this method introduces some error due to the averaging.

u, v, and time can be multi-dimensional arrays. If the time axis, along which the finite differencing is performed, is not the last one (i.e. x.shape[-1]), use the time_axis optional argument to specify along which axis should the differencing be done. x, y, and time must have the same shape.

This function will not do any special handling of longitude ranges. If the integrated trajectory crosses the antimeridian (dateline) in either direction, the longitude values will not be adjusted to stay in any specific range such as [-180, 180] or [0, 360]. If you need your longitudes to be in a specific range, recast the resulting longitude from this function using the function clouddrift.sphere.recast_lon().

Parameters#

unp.ndarray

An array of eastward velocities.

vnp.ndarray

An array of northward velocities.

timenp.ndarray

An array of time values.

x_originfloat

Origin x-coordinate or origin longitude.

y_originfloat

Origin y-coordinate or origin latitude.

coord_systemstr, optional

The coordinate system of the input. Can be “spherical” or “cartesian”. Default is “spherical”.

integration_schemestr, optional

The difference scheme to use for computing the position. Can be “forward” or “backward”. Default is “forward”.

time_axisint, optional

The axis of the time array. Default is -1, which corresponds to the last axis.

Returns#

xnp.ndarray

An array of zonal displacements or longitudes.

ynp.ndarray

An array of meridional displacements or latitudes.

Examples#

Simple integration on a plane, using the forward scheme by default:

>>> import numpy as np
>>> from clouddrift.analysis import position_from_velocity
>>> u = np.array([1., 2., 3., 4.])
>>> v = np.array([1., 1., 1., 1.])
>>> time = np.array([0., 1., 2., 3.])
>>> x, y = position_from_velocity(u, v, time, 0, 0, coord_system="cartesian")
>>> x
array([0., 1., 3., 6.])
>>> y
array([0., 1., 2., 3.])

As above, but using centered scheme:

>>> x, y = position_from_velocity(u, v, time, 0, 0, coord_system="cartesian", integration_scheme="centered")
>>> x
array([0., 1.5, 4., 7.5])
>>> y
array([0., 1., 2., 3.])

Simple integration on a sphere (default):

>>> u = np.array([1., 2., 3., 4.])
>>> v = np.array([1., 1., 1., 1.])
>>> time = np.array([0., 1., 2., 3.]) * 1e5
>>> x, y = position_from_velocity(u, v, time, 0, 0)
>>> x
array([0.        , 0.89839411, 2.69584476, 5.39367518])
>>> y
array([0.        , 0.89828369, 1.79601515, 2.69201609])

Integrating across the antimeridian (dateline) by default does not recast the resulting longitude:

>>> u = np.array([1., 1.])
>>> v = np.array([0., 0.])
>>> time = np.array([0, 1e5])
>>> x, y = position_from_velocity(u, v, time, 179.5, 0)
>>> x
array([179.5      , 180.3983205])
>>> y
array([0., 0.])

Use the clouddrift.sphere.recast_lon function to recast the longitudes to the desired range:

>>> from clouddrift.sphere import recast_lon
>>> recast_lon(x, -180)
array([ 179.5      , -179.6016795])

Raises#

ValueError

If u and v do not have the same shape. If the time axis is outside of the valid range ([-1, N-1]). If lengths of x, y, and time along time_axis are not equal. If the input coordinate system is not “spherical” or “cartesian”. If the input integration scheme is not “forward”, “backward”, or “centered”

See Also#

velocity_from_position()

clouddrift.kinematics.residual_position_from_displacement(longitude: float | ndarray | DataArray, latitude: float | ndarray | DataArray, x: float | ndarray, y: float | ndarray) tuple[float] | tuple[ndarray][source]#

Return residual longitudes and latitudes along a trajectory on the spherical Earth after correcting for zonal and meridional displacements x and y in meters.

This is applicable as an example when one seeks to correct a trajectory for horizontal oscillations due to inertial motions, tides, etc.

Parameters#

longitudefloat or array-like

Longitude in degrees.

latitudefloat or array-like

Latitude in degrees.

xfloat or np.ndarray

Zonal displacement in meters.

yfloat or np.ndarray

Meridional displacement in meters.

Returns#

residual_longitudefloat or np.ndarray

Residual longitude after correcting for zonal displacement, in degrees.

residual_latitudefloat or np.ndarray

Residual latitude after correcting for meridional displacement, in degrees.

Examples#

Obtain the new geographical position for a displacement of 1/360-th of the circumference of the Earth from original position (longitude,latitude) = (1,0):

>>> from clouddrift.sphere import EARTH_RADIUS_METERS
>>> residual_position_from_displacement(1,0,2 * np.pi * EARTH_RADIUS_METERS / 360,0)
(0.0, 0.0)
clouddrift.kinematics.spin(u: ndarray, v: ndarray, time: ndarray, difference_scheme: str | None = 'forward', time_axis: int | None = -1) float | ndarray[source]#

Compute spin continuously from velocities and times.

Spin is traditionally (Sawford, 1999; Veneziani et al., 2005) defined as (<u’dv’ - v’du’>) / (2 dt EKE) where u’ and v’ are eddy-perturbations of the velocity field, EKE is eddy kinetic energy, dt is the time step, and du’ and dv’ are velocity component increments during dt, and < > denotes ensemble average.

To allow computing spin based on full velocity fields, this function does not do any demeaning of the velocity fields. If you need the spin based on velocity anomalies, ensure to demean the velocity fields before passing them to this function. This function also returns instantaneous spin values, so the rank of the result is not reduced relative to the input.

u, v, and time can be multi-dimensional arrays. If the time axis, along which the finite differencing is performed, is not the last one (i.e. u.shape[-1]), use the time_axis optional argument to specify along which the spin should be calculated. u, v, and time must either have the same shape, or time must be a 1-d array with the same length as u.shape[time_axis].

Difference scheme can be one of three values:

  1. “forward” (default): finite difference is evaluated as dx[i] = dx[i+1] - dx[i];

  2. “backward”: finite difference is evaluated as dx[i] = dx[i] - dx[i-1];

  3. “centered”: finite difference is evaluated as dx[i] = (dx[i+1] - dx[i-1]) / 2.

Forward and backward schemes are effectively the same except that the position at which the velocity is evaluated is shifted one element down in the backward scheme relative to the forward scheme. In the case of a forward or backward difference scheme, the last or first element of the velocity, respectively, is extrapolated from its neighboring point. In the case of a centered difference scheme, the start and end boundary points are evaluated using the forward and backward difference scheme, respectively.

Parameters#

unp.ndarray

Zonal velocity

vnp.ndarray

Meridional velocity

timearray-like

Time

difference_schemestr, optional

Difference scheme to use; possible values are “forward”, “backward”, and “centered”.

time_axisint, optional

Axis along which the time varies (default is -1)

Returns#

sfloat or np.ndarray

Spin

Raises#

ValueError

If u and v do not have the same shape. If the time axis is outside of the valid range ([-1, N-1]). If lengths of u, v, and time along time_axis are not equal. If difference_scheme is not “forward”, “backward”, or “centered”.

Examples#

>>> from clouddrift.kinematics import spin
>>> import numpy as np
>>> u = np.array([1., 2., -1., 4.])
>>> v = np.array([1., 3., -2., 1.])
>>> time = np.array([0., 1., 2., 3.])
>>> spin(u, v, time)
array([ 0.5       , -0.07692308,  1.4       ,  0.41176471])

Use difference_scheme to specify an alternative finite difference scheme for the velocity differences:

>>> spin(u, v, time, difference_scheme="centered")
array([0.5       , 0.        , 0.6       , 0.41176471])
>>> spin(u, v, time, difference_scheme="backward")
array([ 0.5       ,  0.07692308, -0.2       ,  0.41176471])

References#

clouddrift.kinematics.velocity_from_position(x: ndarray, y: ndarray, time: ndarray, coord_system: str | None = 'spherical', difference_scheme: str | None = 'forward', time_axis: int | None = -1) tuple[DataArray, DataArray][source]#

Compute velocity from arrays of positions and time.

x and y can be provided as longitude and latitude in degrees if coord_system == “spherical” (default), or as easting and northing if coord_system == “cartesian”.

The units of the result are meters per unit of time if coord_system == “spherical”. For example, if the time is provided in the units of seconds, the resulting velocity is in the units of meters per second. Otherwise, if coord_system == “cartesian”, the units of the resulting velocity correspond to the units of the input. For example, if zonal and meridional displacements are in the units of kilometers and time is in the units of hours, the resulting velocity is in the units of kilometers per hour.

x, y, and time can be multi-dimensional arrays. If the time axis, along which the finite differencing is performed, is not the last one (i.e. x.shape[-1]), use the time_axis optional argument to specify along which axis should the differencing be done. x, y, and time must have the same shape.

Difference scheme can take one of three values:

  1. “forward” (default): finite difference is evaluated as dx[i] = dx[i+1] - dx[i];

  2. “backward”: finite difference is evaluated as dx[i] = dx[i] - dx[i-1];

  3. “centered”: finite difference is evaluated as dx[i] = (dx[i+1] - dx[i-1]) / 2.

Forward and backward schemes are effectively the same except that the position at which the velocity is evaluated is shifted one element down in the backward scheme relative to the forward scheme. In the case of a forward or backward difference scheme, the last or first element of the velocity, respectively, is extrapolated from its neighboring point. In the case of a centered difference scheme, the start and end boundary points are evaluated using the forward and backward difference scheme, respectively.

Parameters#

xarray_like

An N-d array of x-positions (longitude in degrees or zonal displacement in any unit)

yarray_like

An N-d array of y-positions (latitude in degrees or meridional displacement in any unit)

timearray_like

An N-d array of times as floating point values (in any unit)

coord_systemstr, optional

Coordinate system that x and y arrays are in; possible values are “spherical” (default) or “cartesian”.

difference_schemestr, optional

Difference scheme to use; possible values are “forward”, “backward”, and “centered”.

time_axisint, optional

Axis along which to differentiate (default is -1)

Returns#

unp.ndarray

Zonal velocity

vnp.ndarray

Meridional velocity

Raises#

ValueError

If x and y do not have the same shape. If time_axis is outside of the valid range. If lengths of x, y, and time along time_axis are not equal. If coord_system is not “spherical” or “cartesian”. If difference_scheme is not “forward”, “backward”, or “centered”.

Examples#

Simple integration on a sphere, using the forward scheme by default:

>>> import numpy as np
>>> from clouddrift.kinematics import velocity_from_position
>>> lon = np.array([0., 1., 3., 6.])
>>> lat = np.array([0., 1., 2., 3.])
>>> time = np.array([0., 1., 2., 3.]) * 1e5
>>> u, v = velocity_from_position(lon, lat, time)
>>> u
array([1.11307541, 2.22513331, 3.33515501, 3.33515501])
>>> v
array([1.11324496, 1.11409224, 1.1167442 , 1.1167442 ])

Integration on a Cartesian plane, using the forward scheme by default:

>>> x = np.array([0., 1., 3., 6.])
>>> y = np.array([0., 1., 2., 3.])
>>> time = np.array([0., 1., 2., 3.])
>>> u, v = velocity_from_position(x, y, time, coord_system="cartesian")
>>> u
array([1., 2., 3., 3.])
>>> v
array([1., 1., 1., 1.])

See Also#

position_from_velocity()