Skip to content

Raw Stencil Primitives

Pure index-arithmetic stencils for Arakawa C-grids. No metric scaling, no ghost-ring padding. See the theory page for the mathematical background.

1-D Difference Stencils

finitevolx.diff_x_fwd_1d(h)

Forward difference in x (centre → east face), 1-D.

Δx h[i+½] = h[i+1] − h[i]

Maps T-points → U-points.

Parameters:

Name Type Description Default
h Float[Array, Nx]

Field on T-points, including ghost cells.

required

Returns:

Type Description
Float[Array, Nx - 2]

Raw difference at interior U-points.

Notes

No metric scaling is applied. The caller divides by dx, R·cos(φ)·dλ, or whatever metric is appropriate.

Source code in finitevolx/_src/operators/stencils.py
def diff_x_fwd_1d(h: Float[Array, "Nx"]) -> Float[Array, "Nx-2"]:
    """Forward difference in x (centre → east face), 1-D.

    Δx h[i+½] = h[i+1] − h[i]

    Maps T-points → U-points.

    Parameters
    ----------
    h : Float[Array, "Nx"]
        Field on T-points, including ghost cells.

    Returns
    -------
    Float[Array, "Nx-2"]
        Raw difference at interior U-points.

    Notes
    -----
    No metric scaling is applied.  The caller divides by ``dx``,
    ``R·cos(φ)·dλ``, or whatever metric is appropriate.
    """
    return h[2:] - h[1:-1]

finitevolx.diff_x_bwd_1d(h)

Backward difference in x (east face → centre), 1-D.

Δx h[i] = h[i+½] − h[i−½]

Maps U-points → T-points.

Parameters:

Name Type Description Default
h Float[Array, Nx]

Field on U-points, including ghost cells.

required

Returns:

Type Description
Float[Array, Nx - 2]

Raw difference at interior T-points.

Notes

No metric scaling is applied. The caller divides by dx, R·cos(φ)·dλ, or whatever metric is appropriate.

Source code in finitevolx/_src/operators/stencils.py
def diff_x_bwd_1d(h: Float[Array, "Nx"]) -> Float[Array, "Nx-2"]:
    """Backward difference in x (east face → centre), 1-D.

    Δx h[i] = h[i+½] − h[i−½]

    Maps U-points → T-points.

    Parameters
    ----------
    h : Float[Array, "Nx"]
        Field on U-points, including ghost cells.

    Returns
    -------
    Float[Array, "Nx-2"]
        Raw difference at interior T-points.

    Notes
    -----
    No metric scaling is applied.  The caller divides by ``dx``,
    ``R·cos(φ)·dλ``, or whatever metric is appropriate.
    """
    return h[1:-1] - h[:-2]

2-D Difference Stencils

finitevolx.diff_x_fwd(h)

Forward difference in x (centre → east face).

Δx h[j, i+½] = h[j, i+1] − h[j, i]

Maps T-points → U-points (or V-points → X-points).

Parameters:

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

Field on source points (T or V), including ghost ring.

required

Returns:

Type Description
Float[Array, 'Ny-2 Nx-2']

Raw difference at interior destination points (U or X).

Notes

No metric scaling is applied. The caller divides by dx, R·cos(φ)·dλ, or whatever metric is appropriate.

Source code in finitevolx/_src/operators/stencils.py
def diff_x_fwd(h: Float[Array, "Ny Nx"]) -> Float[Array, "Ny-2 Nx-2"]:
    """Forward difference in x (centre → east face).

    Δx h[j, i+½] = h[j, i+1] − h[j, i]

    Maps T-points → U-points (or V-points → X-points).

    Parameters
    ----------
    h : Float[Array, "Ny Nx"]
        Field on source points (T or V), including ghost ring.

    Returns
    -------
    Float[Array, "Ny-2 Nx-2"]
        Raw difference at interior destination points (U or X).

    Notes
    -----
    No metric scaling is applied.  The caller divides by ``dx``,
    ``R·cos(φ)·dλ``, or whatever metric is appropriate.
    """
    return h[1:-1, 2:] - h[1:-1, 1:-1]

finitevolx.diff_y_fwd(h)

Forward difference in y (centre → north face).

Δy h[j+½, i] = h[j+1, i] − h[j, i]

Maps T-points → V-points (or U-points → X-points).

Parameters:

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

Field on source points (T or U), including ghost ring.

required

Returns:

Type Description
Float[Array, 'Ny-2 Nx-2']

Raw difference at interior destination points (V or X).

Notes

No metric scaling is applied. The caller divides by dy, R·dφ, or whatever metric is appropriate.

Source code in finitevolx/_src/operators/stencils.py
def diff_y_fwd(h: Float[Array, "Ny Nx"]) -> Float[Array, "Ny-2 Nx-2"]:
    """Forward difference in y (centre → north face).

    Δy h[j+½, i] = h[j+1, i] − h[j, i]

    Maps T-points → V-points (or U-points → X-points).

    Parameters
    ----------
    h : Float[Array, "Ny Nx"]
        Field on source points (T or U), including ghost ring.

    Returns
    -------
    Float[Array, "Ny-2 Nx-2"]
        Raw difference at interior destination points (V or X).

    Notes
    -----
    No metric scaling is applied.  The caller divides by ``dy``,
    ``R·dφ``, or whatever metric is appropriate.
    """
    return h[2:, 1:-1] - h[1:-1, 1:-1]

finitevolx.diff_x_bwd(h)

Backward difference in x (east face → centre).

Δx h[j, i] = h[j, i+½] − h[j, i−½]

Maps U-points → T-points (or X-points → V-points).

Parameters:

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

Field on source points (U or X), including ghost ring.

required

Returns:

Type Description
Float[Array, 'Ny-2 Nx-2']

Raw difference at interior destination points (T or V).

Notes

No metric scaling is applied. The caller divides by dx, R·cos(φ)·dλ, or whatever metric is appropriate.

Source code in finitevolx/_src/operators/stencils.py
def diff_x_bwd(h: Float[Array, "Ny Nx"]) -> Float[Array, "Ny-2 Nx-2"]:
    """Backward difference in x (east face → centre).

    Δx h[j, i] = h[j, i+½] − h[j, i−½]

    Maps U-points → T-points (or X-points → V-points).

    Parameters
    ----------
    h : Float[Array, "Ny Nx"]
        Field on source points (U or X), including ghost ring.

    Returns
    -------
    Float[Array, "Ny-2 Nx-2"]
        Raw difference at interior destination points (T or V).

    Notes
    -----
    No metric scaling is applied.  The caller divides by ``dx``,
    ``R·cos(φ)·dλ``, or whatever metric is appropriate.
    """
    return h[1:-1, 1:-1] - h[1:-1, :-2]

finitevolx.diff_y_bwd(h)

Backward difference in y (north face → centre).

Δy h[j, i] = h[j+½, i] − h[j−½, i]

Maps V-points → T-points (or X-points → U-points).

Parameters:

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

Field on source points (V or X), including ghost ring.

required

Returns:

Type Description
Float[Array, 'Ny-2 Nx-2']

Raw difference at interior destination points (T or U).

Notes

No metric scaling is applied. The caller divides by dy, R·dφ, or whatever metric is appropriate.

Source code in finitevolx/_src/operators/stencils.py
def diff_y_bwd(h: Float[Array, "Ny Nx"]) -> Float[Array, "Ny-2 Nx-2"]:
    """Backward difference in y (north face → centre).

    Δy h[j, i] = h[j+½, i] − h[j−½, i]

    Maps V-points → T-points (or X-points → U-points).

    Parameters
    ----------
    h : Float[Array, "Ny Nx"]
        Field on source points (V or X), including ghost ring.

    Returns
    -------
    Float[Array, "Ny-2 Nx-2"]
        Raw difference at interior destination points (T or U).

    Notes
    -----
    No metric scaling is applied.  The caller divides by ``dy``,
    ``R·dφ``, or whatever metric is appropriate.
    """
    return h[1:-1, 1:-1] - h[:-2, 1:-1]

3-D Difference Stencils

finitevolx.diff_x_fwd_3d(h)

Forward difference in x over all z-levels (centre → east face).

Δx h[k, j, i+½] = h[k, j, i+1] − h[k, j, i]

Maps T-points → U-points (or V-points → X-points).

Parameters:

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

Field on source points (T or V), including ghost ring.

required

Returns:

Type Description
Float[Array, 'Nz-2 Ny-2 Nx-2']

Raw difference at interior destination points (U or X).

Notes

No metric scaling is applied. The caller divides by dx, R·cos(φ)·dλ, or whatever metric is appropriate.

Source code in finitevolx/_src/operators/stencils.py
def diff_x_fwd_3d(
    h: Float[Array, "Nz Ny Nx"],
) -> Float[Array, "Nz-2 Ny-2 Nx-2"]:
    """Forward difference in x over all z-levels (centre → east face).

    Δx h[k, j, i+½] = h[k, j, i+1] − h[k, j, i]

    Maps T-points → U-points (or V-points → X-points).

    Parameters
    ----------
    h : Float[Array, "Nz Ny Nx"]
        Field on source points (T or V), including ghost ring.

    Returns
    -------
    Float[Array, "Nz-2 Ny-2 Nx-2"]
        Raw difference at interior destination points (U or X).

    Notes
    -----
    No metric scaling is applied.  The caller divides by ``dx``,
    ``R·cos(φ)·dλ``, or whatever metric is appropriate.
    """
    return h[1:-1, 1:-1, 2:] - h[1:-1, 1:-1, 1:-1]

finitevolx.diff_y_fwd_3d(h)

Forward difference in y over all z-levels (centre → north face).

Δy h[k, j+½, i] = h[k, j+1, i] − h[k, j, i]

Maps T-points → V-points (or U-points → X-points).

Parameters:

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

Field on source points (T or U), including ghost ring.

required

Returns:

Type Description
Float[Array, 'Nz-2 Ny-2 Nx-2']

Raw difference at interior destination points (V or X).

Notes

No metric scaling is applied. The caller divides by dy, R·dφ, or whatever metric is appropriate.

Source code in finitevolx/_src/operators/stencils.py
def diff_y_fwd_3d(
    h: Float[Array, "Nz Ny Nx"],
) -> Float[Array, "Nz-2 Ny-2 Nx-2"]:
    """Forward difference in y over all z-levels (centre → north face).

    Δy h[k, j+½, i] = h[k, j+1, i] − h[k, j, i]

    Maps T-points → V-points (or U-points → X-points).

    Parameters
    ----------
    h : Float[Array, "Nz Ny Nx"]
        Field on source points (T or U), including ghost ring.

    Returns
    -------
    Float[Array, "Nz-2 Ny-2 Nx-2"]
        Raw difference at interior destination points (V or X).

    Notes
    -----
    No metric scaling is applied.  The caller divides by ``dy``,
    ``R·dφ``, or whatever metric is appropriate.
    """
    return h[1:-1, 2:, 1:-1] - h[1:-1, 1:-1, 1:-1]

finitevolx.diff_x_bwd_3d(h)

Backward difference in x over all z-levels (east face → centre).

Δx h[k, j, i] = h[k, j, i+½] − h[k, j, i−½]

Maps U-points → T-points (or X-points → V-points).

Parameters:

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

Field on source points (U or X), including ghost ring.

required

Returns:

Type Description
Float[Array, 'Nz-2 Ny-2 Nx-2']

Raw difference at interior destination points (T or V).

Notes

No metric scaling is applied. The caller divides by dx, R·cos(φ)·dλ, or whatever metric is appropriate.

Source code in finitevolx/_src/operators/stencils.py
def diff_x_bwd_3d(
    h: Float[Array, "Nz Ny Nx"],
) -> Float[Array, "Nz-2 Ny-2 Nx-2"]:
    """Backward difference in x over all z-levels (east face → centre).

    Δx h[k, j, i] = h[k, j, i+½] − h[k, j, i−½]

    Maps U-points → T-points (or X-points → V-points).

    Parameters
    ----------
    h : Float[Array, "Nz Ny Nx"]
        Field on source points (U or X), including ghost ring.

    Returns
    -------
    Float[Array, "Nz-2 Ny-2 Nx-2"]
        Raw difference at interior destination points (T or V).

    Notes
    -----
    No metric scaling is applied.  The caller divides by ``dx``,
    ``R·cos(φ)·dλ``, or whatever metric is appropriate.
    """
    return h[1:-1, 1:-1, 1:-1] - h[1:-1, 1:-1, :-2]

finitevolx.diff_y_bwd_3d(h)

Backward difference in y over all z-levels (north face → centre).

Δy h[k, j, i] = h[k, j+½, i] − h[k, j−½, i]

Maps V-points → T-points (or X-points → U-points).

Parameters:

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

Field on source points (V or X), including ghost ring.

required

Returns:

Type Description
Float[Array, 'Nz-2 Ny-2 Nx-2']

Raw difference at interior destination points (T or U).

Notes

No metric scaling is applied. The caller divides by dy, R·dφ, or whatever metric is appropriate.

Source code in finitevolx/_src/operators/stencils.py
def diff_y_bwd_3d(
    h: Float[Array, "Nz Ny Nx"],
) -> Float[Array, "Nz-2 Ny-2 Nx-2"]:
    """Backward difference in y over all z-levels (north face → centre).

    Δy h[k, j, i] = h[k, j+½, i] − h[k, j−½, i]

    Maps V-points → T-points (or X-points → U-points).

    Parameters
    ----------
    h : Float[Array, "Nz Ny Nx"]
        Field on source points (V or X), including ghost ring.

    Returns
    -------
    Float[Array, "Nz-2 Ny-2 Nx-2"]
        Raw difference at interior destination points (T or U).

    Notes
    -----
    No metric scaling is applied.  The caller divides by ``dy``,
    ``R·dφ``, or whatever metric is appropriate.
    """
    return h[1:-1, 1:-1, 1:-1] - h[1:-1, :-2, 1:-1]

1-D Averaging Stencils

finitevolx.avg_x_fwd_1d(h)

Forward average in x (centre → east face), 1-D.

h̄[i+½] = ½ (h[i] + h[i+1])

Maps T-points → U-points.

Parameters:

Name Type Description Default
h Float[Array, Nx]

Field on T-points, including ghost cells.

required

Returns:

Type Description
Float[Array, Nx - 2]

Averaged values at interior U-points.

Source code in finitevolx/_src/operators/stencils.py
def avg_x_fwd_1d(h: Float[Array, "Nx"]) -> Float[Array, "Nx-2"]:
    """Forward average in x (centre → east face), 1-D.

    h̄[i+½] = ½ (h[i] + h[i+1])

    Maps T-points → U-points.

    Parameters
    ----------
    h : Float[Array, "Nx"]
        Field on T-points, including ghost cells.

    Returns
    -------
    Float[Array, "Nx-2"]
        Averaged values at interior U-points.
    """
    return 0.5 * (h[1:-1] + h[2:])

finitevolx.avg_x_bwd_1d(h)

Backward average in x (east face → centre), 1-D.

h̄[i] = ½ (h[i−½] + h[i+½])

Maps U-points → T-points.

Parameters:

Name Type Description Default
h Float[Array, Nx]

Field on U-points, including ghost cells.

required

Returns:

Type Description
Float[Array, Nx - 2]

Averaged values at interior T-points.

Source code in finitevolx/_src/operators/stencils.py
def avg_x_bwd_1d(h: Float[Array, "Nx"]) -> Float[Array, "Nx-2"]:
    """Backward average in x (east face → centre), 1-D.

    h̄[i] = ½ (h[i−½] + h[i+½])

    Maps U-points → T-points.

    Parameters
    ----------
    h : Float[Array, "Nx"]
        Field on U-points, including ghost cells.

    Returns
    -------
    Float[Array, "Nx-2"]
        Averaged values at interior T-points.
    """
    return 0.5 * (h[1:-1] + h[:-2])

2-D Averaging Stencils (2-point)

finitevolx.avg_x_fwd(h)

Forward average in x (centre → east face).

h̄[j, i+½] = ½ (h[j, i] + h[j, i+1])

Maps T-points → U-points (or V-points → X-points).

Parameters:

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

Field on source points (T or V), including ghost ring.

required

Returns:

Type Description
Float[Array, 'Ny-2 Nx-2']

Averaged values at interior destination points (U or X).

Source code in finitevolx/_src/operators/stencils.py
def avg_x_fwd(h: Float[Array, "Ny Nx"]) -> Float[Array, "Ny-2 Nx-2"]:
    """Forward average in x (centre → east face).

    h̄[j, i+½] = ½ (h[j, i] + h[j, i+1])

    Maps T-points → U-points (or V-points → X-points).

    Parameters
    ----------
    h : Float[Array, "Ny Nx"]
        Field on source points (T or V), including ghost ring.

    Returns
    -------
    Float[Array, "Ny-2 Nx-2"]
        Averaged values at interior destination points (U or X).
    """
    return 0.5 * (h[1:-1, 1:-1] + h[1:-1, 2:])

finitevolx.avg_y_fwd(h)

Forward average in y (centre → north face).

h̄[j+½, i] = ½ (h[j, i] + h[j+1, i])

Maps T-points → V-points (or U-points → X-points).

Parameters:

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

Field on source points (T or U), including ghost ring.

required

Returns:

Type Description
Float[Array, 'Ny-2 Nx-2']

Averaged values at interior destination points (V or X).

Source code in finitevolx/_src/operators/stencils.py
def avg_y_fwd(h: Float[Array, "Ny Nx"]) -> Float[Array, "Ny-2 Nx-2"]:
    """Forward average in y (centre → north face).

    h̄[j+½, i] = ½ (h[j, i] + h[j+1, i])

    Maps T-points → V-points (or U-points → X-points).

    Parameters
    ----------
    h : Float[Array, "Ny Nx"]
        Field on source points (T or U), including ghost ring.

    Returns
    -------
    Float[Array, "Ny-2 Nx-2"]
        Averaged values at interior destination points (V or X).
    """
    return 0.5 * (h[1:-1, 1:-1] + h[2:, 1:-1])

finitevolx.avg_x_bwd(h)

Backward average in x (east face → centre).

h̄[j, i] = ½ (h[j, i−½] + h[j, i+½])

Maps U-points → T-points (or X-points → V-points).

Parameters:

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

Field on source points (U or X), including ghost ring.

required

Returns:

Type Description
Float[Array, 'Ny-2 Nx-2']

Averaged values at interior destination points (T or V).

Source code in finitevolx/_src/operators/stencils.py
def avg_x_bwd(h: Float[Array, "Ny Nx"]) -> Float[Array, "Ny-2 Nx-2"]:
    """Backward average in x (east face → centre).

    h̄[j, i] = ½ (h[j, i−½] + h[j, i+½])

    Maps U-points → T-points (or X-points → V-points).

    Parameters
    ----------
    h : Float[Array, "Ny Nx"]
        Field on source points (U or X), including ghost ring.

    Returns
    -------
    Float[Array, "Ny-2 Nx-2"]
        Averaged values at interior destination points (T or V).
    """
    return 0.5 * (h[1:-1, 1:-1] + h[1:-1, :-2])

finitevolx.avg_y_bwd(h)

Backward average in y (north face → centre).

h̄[j, i] = ½ (h[j−½, i] + h[j+½, i])

Maps V-points → T-points (or X-points → U-points).

Parameters:

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

Field on source points (V or X), including ghost ring.

required

Returns:

Type Description
Float[Array, 'Ny-2 Nx-2']

Averaged values at interior destination points (T or U).

Source code in finitevolx/_src/operators/stencils.py
def avg_y_bwd(h: Float[Array, "Ny Nx"]) -> Float[Array, "Ny-2 Nx-2"]:
    """Backward average in y (north face → centre).

    h̄[j, i] = ½ (h[j−½, i] + h[j+½, i])

    Maps V-points → T-points (or X-points → U-points).

    Parameters
    ----------
    h : Float[Array, "Ny Nx"]
        Field on source points (V or X), including ghost ring.

    Returns
    -------
    Float[Array, "Ny-2 Nx-2"]
        Averaged values at interior destination points (T or U).
    """
    return 0.5 * (h[1:-1, 1:-1] + h[:-2, 1:-1])

2-D Averaging Stencils (4-point bilinear)

finitevolx.avg_xy_fwd(h)

Bilinear average to NE corner (forward in both x and y).

h̄[j+½, i+½] = ¼ (h[j,i] + h[j,i+1] + h[j+1,i] + h[j+1,i+1])

Maps T-points → X-points.

Parameters:

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

Field on T-points, including ghost ring.

required

Returns:

Type Description
Float[Array, 'Ny-2 Nx-2']

Averaged values at interior X-points.

Source code in finitevolx/_src/operators/stencils.py
def avg_xy_fwd(h: Float[Array, "Ny Nx"]) -> Float[Array, "Ny-2 Nx-2"]:
    """Bilinear average to NE corner (forward in both x and y).

    h̄[j+½, i+½] = ¼ (h[j,i] + h[j,i+1] + h[j+1,i] + h[j+1,i+1])

    Maps T-points → X-points.

    Parameters
    ----------
    h : Float[Array, "Ny Nx"]
        Field on T-points, including ghost ring.

    Returns
    -------
    Float[Array, "Ny-2 Nx-2"]
        Averaged values at interior X-points.
    """
    return 0.25 * (h[1:-1, 1:-1] + h[1:-1, 2:] + h[2:, 1:-1] + h[2:, 2:])

finitevolx.avg_xy_bwd(h)

Bilinear average from corners to centre (backward in both x and y).

h̄[j, i] = ¼ (h[j+½,i+½] + h[j−½,i+½] + h[j+½,i−½] + h[j−½,i−½])

Maps X-points → T-points.

Parameters:

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

Field on X-points, including ghost ring.

required

Returns:

Type Description
Float[Array, 'Ny-2 Nx-2']

Averaged values at interior T-points.

Source code in finitevolx/_src/operators/stencils.py
def avg_xy_bwd(h: Float[Array, "Ny Nx"]) -> Float[Array, "Ny-2 Nx-2"]:
    """Bilinear average from corners to centre (backward in both x and y).

    h̄[j, i] = ¼ (h[j+½,i+½] + h[j−½,i+½] + h[j+½,i−½] + h[j−½,i−½])

    Maps X-points → T-points.

    Parameters
    ----------
    h : Float[Array, "Ny Nx"]
        Field on X-points, including ghost ring.

    Returns
    -------
    Float[Array, "Ny-2 Nx-2"]
        Averaged values at interior T-points.
    """
    return 0.25 * (h[1:-1, 1:-1] + h[:-2, 1:-1] + h[1:-1, :-2] + h[:-2, :-2])

finitevolx.avg_xbwd_yfwd(h)

Bilinear average: U-points → V-points (backward-x, forward-y).

h̄[j+½, i] = ¼ (h[j,i+½] + h[j+1,i+½] + h[j,i−½] + h[j+1,i−½])

Parameters:

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

Field on U-points, including ghost ring.

required

Returns:

Type Description
Float[Array, 'Ny-2 Nx-2']

Averaged values at interior V-points.

Source code in finitevolx/_src/operators/stencils.py
def avg_xbwd_yfwd(h: Float[Array, "Ny Nx"]) -> Float[Array, "Ny-2 Nx-2"]:
    """Bilinear average: U-points → V-points (backward-x, forward-y).

    h̄[j+½, i] = ¼ (h[j,i+½] + h[j+1,i+½] + h[j,i−½] + h[j+1,i−½])

    Parameters
    ----------
    h : Float[Array, "Ny Nx"]
        Field on U-points, including ghost ring.

    Returns
    -------
    Float[Array, "Ny-2 Nx-2"]
        Averaged values at interior V-points.
    """
    return 0.25 * (h[1:-1, 1:-1] + h[2:, 1:-1] + h[1:-1, :-2] + h[2:, :-2])

finitevolx.avg_xfwd_ybwd(h)

Bilinear average: V-points → U-points (forward-x, backward-y).

h̄[j, i+½] = ¼ (h[j+½,i] + h[j−½,i] + h[j+½,i+1] + h[j−½,i+1])

Parameters:

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

Field on V-points, including ghost ring.

required

Returns:

Type Description
Float[Array, 'Ny-2 Nx-2']

Averaged values at interior U-points.

Source code in finitevolx/_src/operators/stencils.py
def avg_xfwd_ybwd(h: Float[Array, "Ny Nx"]) -> Float[Array, "Ny-2 Nx-2"]:
    """Bilinear average: V-points → U-points (forward-x, backward-y).

    h̄[j, i+½] = ¼ (h[j+½,i] + h[j−½,i] + h[j+½,i+1] + h[j−½,i+1])

    Parameters
    ----------
    h : Float[Array, "Ny Nx"]
        Field on V-points, including ghost ring.

    Returns
    -------
    Float[Array, "Ny-2 Nx-2"]
        Averaged values at interior U-points.
    """
    return 0.25 * (h[1:-1, 1:-1] + h[:-2, 1:-1] + h[1:-1, 2:] + h[:-2, 2:])

3-D Averaging Stencils

finitevolx.avg_x_fwd_3d(h)

Forward average in x over all z-levels (centre → east face).

h̄[k, j, i+½] = ½ (h[k, j, i] + h[k, j, i+1])

Maps T-points → U-points (or V-points → X-points).

Parameters:

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

Field on source points (T or V), including ghost ring.

required

Returns:

Type Description
Float[Array, 'Nz-2 Ny-2 Nx-2']

Averaged values at interior destination points (U or X).

Source code in finitevolx/_src/operators/stencils.py
def avg_x_fwd_3d(
    h: Float[Array, "Nz Ny Nx"],
) -> Float[Array, "Nz-2 Ny-2 Nx-2"]:
    """Forward average in x over all z-levels (centre → east face).

    h̄[k, j, i+½] = ½ (h[k, j, i] + h[k, j, i+1])

    Maps T-points → U-points (or V-points → X-points).

    Parameters
    ----------
    h : Float[Array, "Nz Ny Nx"]
        Field on source points (T or V), including ghost ring.

    Returns
    -------
    Float[Array, "Nz-2 Ny-2 Nx-2"]
        Averaged values at interior destination points (U or X).
    """
    return 0.5 * (h[1:-1, 1:-1, 1:-1] + h[1:-1, 1:-1, 2:])

finitevolx.avg_y_fwd_3d(h)

Forward average in y over all z-levels (centre → north face).

h̄[k, j+½, i] = ½ (h[k, j, i] + h[k, j+1, i])

Maps T-points → V-points (or U-points → X-points).

Parameters:

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

Field on source points (T or U), including ghost ring.

required

Returns:

Type Description
Float[Array, 'Nz-2 Ny-2 Nx-2']

Averaged values at interior destination points (V or X).

Source code in finitevolx/_src/operators/stencils.py
def avg_y_fwd_3d(
    h: Float[Array, "Nz Ny Nx"],
) -> Float[Array, "Nz-2 Ny-2 Nx-2"]:
    """Forward average in y over all z-levels (centre → north face).

    h̄[k, j+½, i] = ½ (h[k, j, i] + h[k, j+1, i])

    Maps T-points → V-points (or U-points → X-points).

    Parameters
    ----------
    h : Float[Array, "Nz Ny Nx"]
        Field on source points (T or U), including ghost ring.

    Returns
    -------
    Float[Array, "Nz-2 Ny-2 Nx-2"]
        Averaged values at interior destination points (V or X).
    """
    return 0.5 * (h[1:-1, 1:-1, 1:-1] + h[1:-1, 2:, 1:-1])

finitevolx.avg_x_bwd_3d(h)

Backward average in x over all z-levels (east face → centre).

h̄[k, j, i] = ½ (h[k, j, i−½] + h[k, j, i+½])

Maps U-points → T-points (or X-points → V-points).

Parameters:

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

Field on source points (U or X), including ghost ring.

required

Returns:

Type Description
Float[Array, 'Nz-2 Ny-2 Nx-2']

Averaged values at interior destination points (T or V).

Source code in finitevolx/_src/operators/stencils.py
def avg_x_bwd_3d(
    h: Float[Array, "Nz Ny Nx"],
) -> Float[Array, "Nz-2 Ny-2 Nx-2"]:
    """Backward average in x over all z-levels (east face → centre).

    h̄[k, j, i] = ½ (h[k, j, i−½] + h[k, j, i+½])

    Maps U-points → T-points (or X-points → V-points).

    Parameters
    ----------
    h : Float[Array, "Nz Ny Nx"]
        Field on source points (U or X), including ghost ring.

    Returns
    -------
    Float[Array, "Nz-2 Ny-2 Nx-2"]
        Averaged values at interior destination points (T or V).
    """
    return 0.5 * (h[1:-1, 1:-1, 1:-1] + h[1:-1, 1:-1, :-2])

finitevolx.avg_y_bwd_3d(h)

Backward average in y over all z-levels (north face → centre).

h̄[k, j, i] = ½ (h[k, j−½, i] + h[k, j+½, i])

Maps V-points → T-points (or X-points → U-points).

Parameters:

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

Field on source points (V or X), including ghost ring.

required

Returns:

Type Description
Float[Array, 'Nz-2 Ny-2 Nx-2']

Averaged values at interior destination points (T or U).

Source code in finitevolx/_src/operators/stencils.py
def avg_y_bwd_3d(
    h: Float[Array, "Nz Ny Nx"],
) -> Float[Array, "Nz-2 Ny-2 Nx-2"]:
    """Backward average in y over all z-levels (north face → centre).

    h̄[k, j, i] = ½ (h[k, j−½, i] + h[k, j+½, i])

    Maps V-points → T-points (or X-points → U-points).

    Parameters
    ----------
    h : Float[Array, "Nz Ny Nx"]
        Field on source points (V or X), including ghost ring.

    Returns
    -------
    Float[Array, "Nz-2 Ny-2 Nx-2"]
        Averaged values at interior destination points (T or U).
    """
    return 0.5 * (h[1:-1, 1:-1, 1:-1] + h[1:-1, :-2, 1:-1])