Skip to content

Interpolation Operators

Operators for interpolating fields between staggered grid points (T, U, V, X).

finitevolx.Interpolation1D

Bases: Module

1-D averaging operators on an Arakawa C-grid.

Parameters:

Name Type Description Default
grid ArakawaCGrid1D

The underlying 1-D grid.

required
Source code in finitevolx/_src/operators/interpolation.py
class Interpolation1D(eqx.Module):
    """1-D averaging operators on an Arakawa C-grid.

    Parameters
    ----------
    grid : ArakawaCGrid1D
        The underlying 1-D grid.
    """

    grid: ArakawaCGrid1D

    def T_to_U(self, h: Float[Array, "Nx"]) -> Float[Array, "Nx"]:
        """Interpolate T-point -> U-point (east face).

        h_on_u[i+1/2] = 1/2 * (h[i] + h[i+1])

        Parameters
        ----------
        h : Float[Array, "Nx"]
            Scalar at T-points.

        Returns
        -------
        Float[Array, "Nx"]
            Scalar interpolated to U-points.
        """
        out = interior(avg_x_fwd_1d(h), h)
        return out

    def U_to_T(self, u: Float[Array, "Nx"]) -> Float[Array, "Nx"]:
        """Interpolate U-point -> T-point (cell centre).

        u_on_h[i] = 1/2 * (u[i+1/2] + u[i-1/2])

        Parameters
        ----------
        u : Float[Array, "Nx"]
            Velocity at U-points.

        Returns
        -------
        Float[Array, "Nx"]
            Velocity interpolated to T-points.
        """
        out = interior(avg_x_bwd_1d(u), u)
        return out

T_to_U(h)

Interpolate T-point -> U-point (east face).

h_on_u[i+1/2] = 1/2 * (h[i] + h[i+1])

Parameters:

Name Type Description Default
h Float[Array, Nx]

Scalar at T-points.

required

Returns:

Type Description
Float[Array, Nx]

Scalar interpolated to U-points.

Source code in finitevolx/_src/operators/interpolation.py
def T_to_U(self, h: Float[Array, "Nx"]) -> Float[Array, "Nx"]:
    """Interpolate T-point -> U-point (east face).

    h_on_u[i+1/2] = 1/2 * (h[i] + h[i+1])

    Parameters
    ----------
    h : Float[Array, "Nx"]
        Scalar at T-points.

    Returns
    -------
    Float[Array, "Nx"]
        Scalar interpolated to U-points.
    """
    out = interior(avg_x_fwd_1d(h), h)
    return out

U_to_T(u)

Interpolate U-point -> T-point (cell centre).

u_on_h[i] = 1/2 * (u[i+1/2] + u[i-1/2])

Parameters:

Name Type Description Default
u Float[Array, Nx]

Velocity at U-points.

required

Returns:

Type Description
Float[Array, Nx]

Velocity interpolated to T-points.

Source code in finitevolx/_src/operators/interpolation.py
def U_to_T(self, u: Float[Array, "Nx"]) -> Float[Array, "Nx"]:
    """Interpolate U-point -> T-point (cell centre).

    u_on_h[i] = 1/2 * (u[i+1/2] + u[i-1/2])

    Parameters
    ----------
    u : Float[Array, "Nx"]
        Velocity at U-points.

    Returns
    -------
    Float[Array, "Nx"]
        Velocity interpolated to T-points.
    """
    out = interior(avg_x_bwd_1d(u), u)
    return out

finitevolx.Interpolation2D

Bases: Module

2-D averaging operators on an Arakawa C-grid.

Parameters:

Name Type Description Default
grid ArakawaCGrid2D

The underlying 2-D grid.

required
Source code in finitevolx/_src/operators/interpolation.py
class Interpolation2D(eqx.Module):
    """2-D averaging operators on an Arakawa C-grid.

    Parameters
    ----------
    grid : ArakawaCGrid2D
        The underlying 2-D grid.
    """

    grid: ArakawaCGrid2D

    # ------------------------------------------------------------------
    # T-point -> faces / corners
    # ------------------------------------------------------------------

    def T_to_U(self, h: Float[Array, "Ny Nx"]) -> Float[Array, "Ny Nx"]:
        """T-point -> U-point (east face), x-average.

        h_on_u[j, i+1/2] = 1/2 * (h[j, i] + h[j, i+1])

        Parameters
        ----------
        h : Float[Array, "Ny Nx"]
            Scalar at T-points.

        Returns
        -------
        Float[Array, "Ny Nx"]
            Scalar interpolated to U-points.
        """
        out = interior(avg_x_fwd(h), h)
        return out

    def T_to_V(self, h: Float[Array, "Ny Nx"]) -> Float[Array, "Ny Nx"]:
        """T-point -> V-point (north face), y-average.

        h_on_v[j+1/2, i] = 1/2 * (h[j, i] + h[j+1, i])

        Parameters
        ----------
        h : Float[Array, "Ny Nx"]
            Scalar at T-points.

        Returns
        -------
        Float[Array, "Ny Nx"]
            Scalar interpolated to V-points.
        """
        out = interior(avg_y_fwd(h), h)
        return out

    def T_to_X(self, h: Float[Array, "Ny Nx"]) -> Float[Array, "Ny Nx"]:
        """T-point -> X-point (NE corner), bilinear average.

        h_on_q[j+1/2, i+1/2] = 1/4 * (h[j,i] + h[j,i+1] + h[j+1,i] + h[j+1,i+1])

        Parameters
        ----------
        h : Float[Array, "Ny Nx"]
            Scalar at T-points.

        Returns
        -------
        Float[Array, "Ny Nx"]
            Scalar interpolated to X-points (corners).
        """
        out = interior(avg_xy_fwd(h), h)
        return out

    # ------------------------------------------------------------------
    # X-point -> face points
    # ------------------------------------------------------------------

    def X_to_U(self, q: Float[Array, "Ny Nx"]) -> Float[Array, "Ny Nx"]:
        """X-point (corner) -> U-point (east face), y-average.

        q_on_u[j, i+1/2] = 1/2 * (q[j+1/2, i+1/2] + q[j-1/2, i+1/2])

        Parameters
        ----------
        q : Float[Array, "Ny Nx"]
            Scalar at X-points.

        Returns
        -------
        Float[Array, "Ny Nx"]
            Scalar interpolated to U-points.
        """
        out = interior(avg_y_bwd(q), q)
        return out

    def X_to_V(self, q: Float[Array, "Ny Nx"]) -> Float[Array, "Ny Nx"]:
        """X-point (corner) -> V-point (north face), x-average.

        q_on_v[j+1/2, i] = 1/2 * (q[j+1/2, i+1/2] + q[j+1/2, i-1/2])

        Parameters
        ----------
        q : Float[Array, "Ny Nx"]
            Scalar at X-points.

        Returns
        -------
        Float[Array, "Ny Nx"]
            Scalar interpolated to V-points.
        """
        out = interior(avg_x_bwd(q), q)
        return out

    # ------------------------------------------------------------------
    # Face points -> T-point
    # ------------------------------------------------------------------

    def U_to_T(self, u: Float[Array, "Ny Nx"]) -> Float[Array, "Ny Nx"]:
        """U-point -> T-point, x-average.

        u_on_h[j, i] = 1/2 * (u[j, i+1/2] + u[j, i-1/2])

        Parameters
        ----------
        u : Float[Array, "Ny Nx"]
            Velocity at U-points.

        Returns
        -------
        Float[Array, "Ny Nx"]
            Velocity interpolated to T-points.
        """
        out = interior(avg_x_bwd(u), u)
        return out

    def V_to_T(self, v: Float[Array, "Ny Nx"]) -> Float[Array, "Ny Nx"]:
        """V-point -> T-point, y-average.

        v_on_h[j, i] = 1/2 * (v[j+1/2, i] + v[j-1/2, i])

        Parameters
        ----------
        v : Float[Array, "Ny Nx"]
            Velocity at V-points.

        Returns
        -------
        Float[Array, "Ny Nx"]
            Velocity interpolated to T-points.
        """
        out = interior(avg_y_bwd(v), v)
        return out

    def X_to_T(self, q: Float[Array, "Ny Nx"]) -> Float[Array, "Ny Nx"]:
        """X-point (corner) -> T-point, bilinear average.

        q_on_h[j, i] = 1/4 * (q[j+1/2,i+1/2] + q[j-1/2,i+1/2]
                             + q[j+1/2,i-1/2] + q[j-1/2,i-1/2])

        Parameters
        ----------
        q : Float[Array, "Ny Nx"]
            Scalar at X-points.

        Returns
        -------
        Float[Array, "Ny Nx"]
            Scalar interpolated to T-points.
        """
        out = interior(avg_xy_bwd(q), q)
        return out

    # ------------------------------------------------------------------
    # Face/center -> X-point (corner)
    # ------------------------------------------------------------------

    def U_to_X(self, u: Float[Array, "Ny Nx"]) -> Float[Array, "Ny Nx"]:
        """U-point -> X-point (corner), y-average.

        u_on_q[j+1/2, i+1/2] = 1/2 * (u[j, i+1/2] + u[j+1, i+1/2])

        Parameters
        ----------
        u : Float[Array, "Ny Nx"]
            Velocity at U-points.

        Returns
        -------
        Float[Array, "Ny Nx"]
            Velocity interpolated to X-points.
        """
        out = interior(avg_y_fwd(u), u)
        return out

    def V_to_X(self, v: Float[Array, "Ny Nx"]) -> Float[Array, "Ny Nx"]:
        """V-point -> X-point (corner), x-average.

        v_on_q[j+1/2, i+1/2] = 1/2 * (v[j+1/2, i] + v[j+1/2, i+1])

        Parameters
        ----------
        v : Float[Array, "Ny Nx"]
            Velocity at V-points.

        Returns
        -------
        Float[Array, "Ny Nx"]
            Velocity interpolated to X-points.
        """
        out = interior(avg_x_fwd(v), v)
        return out

    # ------------------------------------------------------------------
    # Cross-face (bilinear 4-point)
    # ------------------------------------------------------------------

    def U_to_V(self, u: Float[Array, "Ny Nx"]) -> Float[Array, "Ny Nx"]:
        """U-point -> V-point (cross-face bilinear, 4-point).

        u_on_v[j+1/2, i] = 1/4 * (u[j,   i+1/2] + u[j+1, i+1/2]
                                 + u[j,   i-1/2] + u[j+1, i-1/2])

        Parameters
        ----------
        u : Float[Array, "Ny Nx"]
            Velocity at U-points.

        Returns
        -------
        Float[Array, "Ny Nx"]
            Velocity interpolated to V-points.
        """
        out = interior(avg_xbwd_yfwd(u), u)
        return out

    def V_to_U(self, v: Float[Array, "Ny Nx"]) -> Float[Array, "Ny Nx"]:
        """V-point -> U-point (cross-face bilinear, 4-point).

        v_on_u[j, i+1/2] = 1/4 * (v[j+1/2, i] + v[j-1/2, i]
                                 + v[j+1/2, i+1] + v[j-1/2, i+1])

        Parameters
        ----------
        v : Float[Array, "Ny Nx"]
            Velocity at V-points.

        Returns
        -------
        Float[Array, "Ny Nx"]
            Velocity interpolated to U-points.
        """
        out = interior(avg_xfwd_ybwd(v), v)
        return out

T_to_U(h)

T-point -> U-point (east face), x-average.

h_on_u[j, i+1/2] = 1/2 * (h[j, i] + h[j, i+1])

Parameters:

Name Type Description Default
h Float[Array, 'Ny Nx']

Scalar at T-points.

required

Returns:

Type Description
Float[Array, 'Ny Nx']

Scalar interpolated to U-points.

Source code in finitevolx/_src/operators/interpolation.py
def T_to_U(self, h: Float[Array, "Ny Nx"]) -> Float[Array, "Ny Nx"]:
    """T-point -> U-point (east face), x-average.

    h_on_u[j, i+1/2] = 1/2 * (h[j, i] + h[j, i+1])

    Parameters
    ----------
    h : Float[Array, "Ny Nx"]
        Scalar at T-points.

    Returns
    -------
    Float[Array, "Ny Nx"]
        Scalar interpolated to U-points.
    """
    out = interior(avg_x_fwd(h), h)
    return out

T_to_V(h)

T-point -> V-point (north face), y-average.

h_on_v[j+1/2, i] = 1/2 * (h[j, i] + h[j+1, i])

Parameters:

Name Type Description Default
h Float[Array, 'Ny Nx']

Scalar at T-points.

required

Returns:

Type Description
Float[Array, 'Ny Nx']

Scalar interpolated to V-points.

Source code in finitevolx/_src/operators/interpolation.py
def T_to_V(self, h: Float[Array, "Ny Nx"]) -> Float[Array, "Ny Nx"]:
    """T-point -> V-point (north face), y-average.

    h_on_v[j+1/2, i] = 1/2 * (h[j, i] + h[j+1, i])

    Parameters
    ----------
    h : Float[Array, "Ny Nx"]
        Scalar at T-points.

    Returns
    -------
    Float[Array, "Ny Nx"]
        Scalar interpolated to V-points.
    """
    out = interior(avg_y_fwd(h), h)
    return out

T_to_X(h)

T-point -> X-point (NE corner), bilinear average.

h_on_q[j+1/2, i+1/2] = 1/4 * (h[j,i] + h[j,i+1] + h[j+1,i] + h[j+1,i+1])

Parameters:

Name Type Description Default
h Float[Array, 'Ny Nx']

Scalar at T-points.

required

Returns:

Type Description
Float[Array, 'Ny Nx']

Scalar interpolated to X-points (corners).

Source code in finitevolx/_src/operators/interpolation.py
def T_to_X(self, h: Float[Array, "Ny Nx"]) -> Float[Array, "Ny Nx"]:
    """T-point -> X-point (NE corner), bilinear average.

    h_on_q[j+1/2, i+1/2] = 1/4 * (h[j,i] + h[j,i+1] + h[j+1,i] + h[j+1,i+1])

    Parameters
    ----------
    h : Float[Array, "Ny Nx"]
        Scalar at T-points.

    Returns
    -------
    Float[Array, "Ny Nx"]
        Scalar interpolated to X-points (corners).
    """
    out = interior(avg_xy_fwd(h), h)
    return out

U_to_T(u)

U-point -> T-point, x-average.

u_on_h[j, i] = 1/2 * (u[j, i+1/2] + u[j, i-1/2])

Parameters:

Name Type Description Default
u Float[Array, 'Ny Nx']

Velocity at U-points.

required

Returns:

Type Description
Float[Array, 'Ny Nx']

Velocity interpolated to T-points.

Source code in finitevolx/_src/operators/interpolation.py
def U_to_T(self, u: Float[Array, "Ny Nx"]) -> Float[Array, "Ny Nx"]:
    """U-point -> T-point, x-average.

    u_on_h[j, i] = 1/2 * (u[j, i+1/2] + u[j, i-1/2])

    Parameters
    ----------
    u : Float[Array, "Ny Nx"]
        Velocity at U-points.

    Returns
    -------
    Float[Array, "Ny Nx"]
        Velocity interpolated to T-points.
    """
    out = interior(avg_x_bwd(u), u)
    return out

U_to_V(u)

U-point -> V-point (cross-face bilinear, 4-point).

u_on_v[j+1/2, i] = 1/4 * (u[j, i+1/2] + u[j+1, i+1/2] + u[j, i-1/2] + u[j+1, i-1/2])

Parameters:

Name Type Description Default
u Float[Array, 'Ny Nx']

Velocity at U-points.

required

Returns:

Type Description
Float[Array, 'Ny Nx']

Velocity interpolated to V-points.

Source code in finitevolx/_src/operators/interpolation.py
def U_to_V(self, u: Float[Array, "Ny Nx"]) -> Float[Array, "Ny Nx"]:
    """U-point -> V-point (cross-face bilinear, 4-point).

    u_on_v[j+1/2, i] = 1/4 * (u[j,   i+1/2] + u[j+1, i+1/2]
                             + u[j,   i-1/2] + u[j+1, i-1/2])

    Parameters
    ----------
    u : Float[Array, "Ny Nx"]
        Velocity at U-points.

    Returns
    -------
    Float[Array, "Ny Nx"]
        Velocity interpolated to V-points.
    """
    out = interior(avg_xbwd_yfwd(u), u)
    return out

U_to_X(u)

U-point -> X-point (corner), y-average.

u_on_q[j+1/2, i+1/2] = 1/2 * (u[j, i+1/2] + u[j+1, i+1/2])

Parameters:

Name Type Description Default
u Float[Array, 'Ny Nx']

Velocity at U-points.

required

Returns:

Type Description
Float[Array, 'Ny Nx']

Velocity interpolated to X-points.

Source code in finitevolx/_src/operators/interpolation.py
def U_to_X(self, u: Float[Array, "Ny Nx"]) -> Float[Array, "Ny Nx"]:
    """U-point -> X-point (corner), y-average.

    u_on_q[j+1/2, i+1/2] = 1/2 * (u[j, i+1/2] + u[j+1, i+1/2])

    Parameters
    ----------
    u : Float[Array, "Ny Nx"]
        Velocity at U-points.

    Returns
    -------
    Float[Array, "Ny Nx"]
        Velocity interpolated to X-points.
    """
    out = interior(avg_y_fwd(u), u)
    return out

V_to_T(v)

V-point -> T-point, y-average.

v_on_h[j, i] = 1/2 * (v[j+1/2, i] + v[j-1/2, i])

Parameters:

Name Type Description Default
v Float[Array, 'Ny Nx']

Velocity at V-points.

required

Returns:

Type Description
Float[Array, 'Ny Nx']

Velocity interpolated to T-points.

Source code in finitevolx/_src/operators/interpolation.py
def V_to_T(self, v: Float[Array, "Ny Nx"]) -> Float[Array, "Ny Nx"]:
    """V-point -> T-point, y-average.

    v_on_h[j, i] = 1/2 * (v[j+1/2, i] + v[j-1/2, i])

    Parameters
    ----------
    v : Float[Array, "Ny Nx"]
        Velocity at V-points.

    Returns
    -------
    Float[Array, "Ny Nx"]
        Velocity interpolated to T-points.
    """
    out = interior(avg_y_bwd(v), v)
    return out

V_to_U(v)

V-point -> U-point (cross-face bilinear, 4-point).

v_on_u[j, i+1/2] = 1/4 * (v[j+1/2, i] + v[j-1/2, i] + v[j+1/2, i+1] + v[j-1/2, i+1])

Parameters:

Name Type Description Default
v Float[Array, 'Ny Nx']

Velocity at V-points.

required

Returns:

Type Description
Float[Array, 'Ny Nx']

Velocity interpolated to U-points.

Source code in finitevolx/_src/operators/interpolation.py
def V_to_U(self, v: Float[Array, "Ny Nx"]) -> Float[Array, "Ny Nx"]:
    """V-point -> U-point (cross-face bilinear, 4-point).

    v_on_u[j, i+1/2] = 1/4 * (v[j+1/2, i] + v[j-1/2, i]
                             + v[j+1/2, i+1] + v[j-1/2, i+1])

    Parameters
    ----------
    v : Float[Array, "Ny Nx"]
        Velocity at V-points.

    Returns
    -------
    Float[Array, "Ny Nx"]
        Velocity interpolated to U-points.
    """
    out = interior(avg_xfwd_ybwd(v), v)
    return out

V_to_X(v)

V-point -> X-point (corner), x-average.

v_on_q[j+1/2, i+1/2] = 1/2 * (v[j+1/2, i] + v[j+1/2, i+1])

Parameters:

Name Type Description Default
v Float[Array, 'Ny Nx']

Velocity at V-points.

required

Returns:

Type Description
Float[Array, 'Ny Nx']

Velocity interpolated to X-points.

Source code in finitevolx/_src/operators/interpolation.py
def V_to_X(self, v: Float[Array, "Ny Nx"]) -> Float[Array, "Ny Nx"]:
    """V-point -> X-point (corner), x-average.

    v_on_q[j+1/2, i+1/2] = 1/2 * (v[j+1/2, i] + v[j+1/2, i+1])

    Parameters
    ----------
    v : Float[Array, "Ny Nx"]
        Velocity at V-points.

    Returns
    -------
    Float[Array, "Ny Nx"]
        Velocity interpolated to X-points.
    """
    out = interior(avg_x_fwd(v), v)
    return out

X_to_T(q)

X-point (corner) -> T-point, bilinear average.

q_on_h[j, i] = 1/4 * (q[j+1/2,i+1/2] + q[j-1/2,i+1/2] + q[j+1/2,i-1/2] + q[j-1/2,i-1/2])

Parameters:

Name Type Description Default
q Float[Array, 'Ny Nx']

Scalar at X-points.

required

Returns:

Type Description
Float[Array, 'Ny Nx']

Scalar interpolated to T-points.

Source code in finitevolx/_src/operators/interpolation.py
def X_to_T(self, q: Float[Array, "Ny Nx"]) -> Float[Array, "Ny Nx"]:
    """X-point (corner) -> T-point, bilinear average.

    q_on_h[j, i] = 1/4 * (q[j+1/2,i+1/2] + q[j-1/2,i+1/2]
                         + q[j+1/2,i-1/2] + q[j-1/2,i-1/2])

    Parameters
    ----------
    q : Float[Array, "Ny Nx"]
        Scalar at X-points.

    Returns
    -------
    Float[Array, "Ny Nx"]
        Scalar interpolated to T-points.
    """
    out = interior(avg_xy_bwd(q), q)
    return out

X_to_U(q)

X-point (corner) -> U-point (east face), y-average.

q_on_u[j, i+1/2] = 1/2 * (q[j+1/2, i+1/2] + q[j-1/2, i+1/2])

Parameters:

Name Type Description Default
q Float[Array, 'Ny Nx']

Scalar at X-points.

required

Returns:

Type Description
Float[Array, 'Ny Nx']

Scalar interpolated to U-points.

Source code in finitevolx/_src/operators/interpolation.py
def X_to_U(self, q: Float[Array, "Ny Nx"]) -> Float[Array, "Ny Nx"]:
    """X-point (corner) -> U-point (east face), y-average.

    q_on_u[j, i+1/2] = 1/2 * (q[j+1/2, i+1/2] + q[j-1/2, i+1/2])

    Parameters
    ----------
    q : Float[Array, "Ny Nx"]
        Scalar at X-points.

    Returns
    -------
    Float[Array, "Ny Nx"]
        Scalar interpolated to U-points.
    """
    out = interior(avg_y_bwd(q), q)
    return out

X_to_V(q)

X-point (corner) -> V-point (north face), x-average.

q_on_v[j+1/2, i] = 1/2 * (q[j+1/2, i+1/2] + q[j+1/2, i-1/2])

Parameters:

Name Type Description Default
q Float[Array, 'Ny Nx']

Scalar at X-points.

required

Returns:

Type Description
Float[Array, 'Ny Nx']

Scalar interpolated to V-points.

Source code in finitevolx/_src/operators/interpolation.py
def X_to_V(self, q: Float[Array, "Ny Nx"]) -> Float[Array, "Ny Nx"]:
    """X-point (corner) -> V-point (north face), x-average.

    q_on_v[j+1/2, i] = 1/2 * (q[j+1/2, i+1/2] + q[j+1/2, i-1/2])

    Parameters
    ----------
    q : Float[Array, "Ny Nx"]
        Scalar at X-points.

    Returns
    -------
    Float[Array, "Ny Nx"]
        Scalar interpolated to V-points.
    """
    out = interior(avg_x_bwd(q), q)
    return out

finitevolx.Interpolation3D

Bases: Module

3-D averaging operators on an Arakawa C-grid.

Operates on the horizontal (y, x) plane for each z-level. Array shape is [Nz, Ny, Nx].

Parameters:

Name Type Description Default
grid ArakawaCGrid3D

The underlying 3-D grid.

required
Source code in finitevolx/_src/operators/interpolation.py
class Interpolation3D(eqx.Module):
    """3-D averaging operators on an Arakawa C-grid.

    Operates on the horizontal (y, x) plane for each z-level.
    Array shape is [Nz, Ny, Nx].

    Parameters
    ----------
    grid : ArakawaCGrid3D
        The underlying 3-D grid.
    """

    grid: ArakawaCGrid3D

    def T_to_U(self, h: Float[Array, "Nz Ny Nx"]) -> Float[Array, "Nz Ny Nx"]:
        """T -> U (x-average) over all z-levels.

        h_on_u[k, j, i+1/2] = 1/2 * (h[k, j, i] + h[k, j, i+1])
        """
        out = interior(avg_x_fwd_3d(h), h)
        return out

    def T_to_V(self, h: Float[Array, "Nz Ny Nx"]) -> Float[Array, "Nz Ny Nx"]:
        """T -> V (y-average) over all z-levels.

        h_on_v[k, j+1/2, i] = 1/2 * (h[k, j, i] + h[k, j+1, i])
        """
        out = interior(avg_y_fwd_3d(h), h)
        return out

    def U_to_T(self, u: Float[Array, "Nz Ny Nx"]) -> Float[Array, "Nz Ny Nx"]:
        """U -> T (x-average) over all z-levels.

        u_on_h[k, j, i] = 1/2 * (u[k, j, i+1/2] + u[k, j, i-1/2])
        """
        out = interior(avg_x_bwd_3d(u), u)
        return out

    def V_to_T(self, v: Float[Array, "Nz Ny Nx"]) -> Float[Array, "Nz Ny Nx"]:
        """V -> T (y-average) over all z-levels.

        v_on_h[k, j, i] = 1/2 * (v[k, j+1/2, i] + v[k, j-1/2, i])
        """
        out = interior(avg_y_bwd_3d(v), v)
        return out

T_to_U(h)

T -> U (x-average) over all z-levels.

h_on_u[k, j, i+1/2] = 1/2 * (h[k, j, i] + h[k, j, i+1])

Source code in finitevolx/_src/operators/interpolation.py
def T_to_U(self, h: Float[Array, "Nz Ny Nx"]) -> Float[Array, "Nz Ny Nx"]:
    """T -> U (x-average) over all z-levels.

    h_on_u[k, j, i+1/2] = 1/2 * (h[k, j, i] + h[k, j, i+1])
    """
    out = interior(avg_x_fwd_3d(h), h)
    return out

T_to_V(h)

T -> V (y-average) over all z-levels.

h_on_v[k, j+1/2, i] = 1/2 * (h[k, j, i] + h[k, j+1, i])

Source code in finitevolx/_src/operators/interpolation.py
def T_to_V(self, h: Float[Array, "Nz Ny Nx"]) -> Float[Array, "Nz Ny Nx"]:
    """T -> V (y-average) over all z-levels.

    h_on_v[k, j+1/2, i] = 1/2 * (h[k, j, i] + h[k, j+1, i])
    """
    out = interior(avg_y_fwd_3d(h), h)
    return out

U_to_T(u)

U -> T (x-average) over all z-levels.

u_on_h[k, j, i] = 1/2 * (u[k, j, i+1/2] + u[k, j, i-1/2])

Source code in finitevolx/_src/operators/interpolation.py
def U_to_T(self, u: Float[Array, "Nz Ny Nx"]) -> Float[Array, "Nz Ny Nx"]:
    """U -> T (x-average) over all z-levels.

    u_on_h[k, j, i] = 1/2 * (u[k, j, i+1/2] + u[k, j, i-1/2])
    """
    out = interior(avg_x_bwd_3d(u), u)
    return out

V_to_T(v)

V -> T (y-average) over all z-levels.

v_on_h[k, j, i] = 1/2 * (v[k, j+1/2, i] + v[k, j-1/2, i])

Source code in finitevolx/_src/operators/interpolation.py
def V_to_T(self, v: Float[Array, "Nz Ny Nx"]) -> Float[Array, "Nz Ny Nx"]:
    """V -> T (y-average) over all z-levels.

    v_on_h[k, j, i] = 1/2 * (v[k, j+1/2, i] + v[k, j-1/2, i])
    """
    out = interior(avg_y_bwd_3d(v), v)
    return out