clouddrift.kinematics.position_from_velocity

clouddrift.kinematics.position_from_velocity#

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()