Skip to content

Capacitance Matrix Method

Extends fast spectral solvers to irregular / masked domains using the Sherman-Morrison capacitance correction (Buzbee, Golub & Nielson 1970).

Solver

finitevolx.CapacitanceSolver

Bases: Module

Spectral Poisson/Helmholtz solver for masked irregular domains.

Uses the capacitance matrix method (Buzbee, Golub & Nielson 1970) to extend a fast rectangular spectral solver to a domain defined by a binary mask.

The algorithm (Buzbee, Golub & Nielson 1970):

  1. Solve the PDE on the full rectangle using a fast spectral solver (DST/DCT/FFT), ignoring the mask. Call this u.
  2. u generally violates ψ = 0 at inner-boundary points. Correct it: ψ = u − Σ_k α_k g_k, where g_k are precomputed Green's functions (rectangular-domain response to δ-sources at each boundary point b_k).
  3. The coefficients α are found by requiring ψ(b_k) = 0 at all N_b boundary points, giving the linear system C α = u[B] where C[k,l] = g_l(b_k) is the capacitance matrix.

Construct with :func:build_capacitance_solver.

Attributes:

Name Type Description
_C_inv Float[Array, 'Nb Nb']

Pre-inverted capacitance matrix.

_green_flat Float[Array, 'Nb NyNx']

Green's functions (one row per boundary point), stored flat.

_mask Float[Array, 'Ny Nx']

Domain mask (1.0 = interior, 0.0 = exterior). Applied to the output so that values outside the physical domain are exactly zero.

_j_b Array

Row indices of inner-boundary points.

_i_b Array

Column indices of inner-boundary points.

dx float

Grid spacing in x.

dy float

Grid spacing in y.

lambda_ float

Helmholtz parameter.

base_bc str

Spectral solver used as the rectangular base.

Source code in .venv/lib/python3.12/site-packages/spectraldiffx/_src/fourier/capacitance.py
class CapacitanceSolver(eqx.Module):
    """Spectral Poisson/Helmholtz solver for masked irregular domains.

    Uses the **capacitance matrix method** (Buzbee, Golub & Nielson 1970) to
    extend a fast rectangular spectral solver to a domain defined by a binary
    mask.

    The algorithm (Buzbee, Golub & Nielson 1970):

    1. Solve the PDE on the **full rectangle** using a fast spectral solver
       (DST/DCT/FFT), ignoring the mask.  Call this ``u``.
    2. ``u`` generally violates ψ = 0 at inner-boundary points.  Correct it:
       ``ψ = u − Σ_k α_k g_k``, where ``g_k`` are precomputed Green's
       functions (rectangular-domain response to δ-sources at each
       boundary point b_k).
    3. The coefficients α are found by requiring ψ(b_k) = 0 at all
       N_b boundary points, giving the linear system ``C α = u[B]``
       where ``C[k,l] = g_l(b_k)`` is the **capacitance matrix**.

    Construct with :func:`build_capacitance_solver`.

    Attributes
    ----------
    _C_inv : Float[Array, "Nb Nb"]
        Pre-inverted capacitance matrix.
    _green_flat : Float[Array, "Nb NyNx"]
        Green's functions (one row per boundary point), stored flat.
    _mask : Float[Array, "Ny Nx"]
        Domain mask (1.0 = interior, 0.0 = exterior).  Applied to the output
        so that values outside the physical domain are exactly zero.
    _j_b : Array
        Row indices of inner-boundary points.
    _i_b : Array
        Column indices of inner-boundary points.
    dx : float
        Grid spacing in x.
    dy : float
        Grid spacing in y.
    lambda_ : float
        Helmholtz parameter.
    base_bc : str
        Spectral solver used as the rectangular base.
    """

    _C_inv: Float[Array, "Nb Nb"]
    _green_flat: Float[Array, "Nb NyNx"]
    _mask: Float[Array, "Ny Nx"]
    _j_b: Array
    _i_b: Array
    dx: float
    dy: float
    lambda_: float = eqx.field(static=True)
    base_bc: str = eqx.field(static=True)

    def __call__(
        self,
        rhs: Float[Array, "Ny Nx"],
    ) -> Float[Array, "Ny Nx"]:
        """Solve (∇² − λ)ψ = rhs on the masked domain.

        Given the precomputed capacitance matrix C⁻¹ and Green's functions
        G, the online solve proceeds in four steps:

        1. Rectangular solve:   u = L_rect⁻¹ rhs              [Ny, Nx]
        2. Sample boundary:     u_B = u[j_b, i_b]             [N_b]
        3. Correction weights:  α = C⁻¹ · u_B                 [N_b]
        4. Subtract correction: ψ = u − Σ_k α_k g_k          [Ny, Nx]
        5. Mask exterior:       ψ = ψ * mask                  [Ny, Nx]

        The result satisfies ψ(b_k) ≈ 0 at all N_b inner-boundary points,
        (∇² − λ)ψ = rhs at interior points, and ψ = 0 outside the mask.

        Parameters
        ----------
        rhs : Float[Array, "Ny Nx"]
            Right-hand side on the full rectangular grid.
            Values outside the physical domain (mask = False) are ignored.

        Returns
        -------
        Float[Array, "Ny Nx"]
            Solution ψ on the full rectangular grid.  Exactly zero outside
            the mask, and ψ ≈ 0 at inner-boundary points.
        """
        Ny, Nx = rhs.shape
        # Step 1: rectangular spectral solve
        u = _spectral_solve(rhs, self.dx, self.dy, self.lambda_, self.base_bc)
        # Step 2: values of u at inner-boundary points
        u_b = u[self._j_b, self._i_b]  # [Nb]
        # Step 3: correction coefficients  alpha = C^{-1} u_b
        alpha = self._C_inv @ u_b  # [Nb]
        # Step 4: correction field  sum_k alpha_k g_k
        correction = (self._green_flat.T @ alpha).reshape(Ny, Nx)  # [Ny, Nx]
        # Step 5: zero out exterior cells
        return (u - correction) * self._mask

__call__(rhs)

Solve (∇² − λ)ψ = rhs on the masked domain.

Given the precomputed capacitance matrix C⁻¹ and Green's functions G, the online solve proceeds in four steps:

  1. Rectangular solve: u = L_rect⁻¹ rhs [Ny, Nx]
  2. Sample boundary: u_B = u[j_b, i_b] [N_b]
  3. Correction weights: α = C⁻¹ · u_B [N_b]
  4. Subtract correction: ψ = u − Σ_k α_k g_k [Ny, Nx]
  5. Mask exterior: ψ = ψ * mask [Ny, Nx]

The result satisfies ψ(b_k) ≈ 0 at all N_b inner-boundary points, (∇² − λ)ψ = rhs at interior points, and ψ = 0 outside the mask.

Parameters:

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

Right-hand side on the full rectangular grid. Values outside the physical domain (mask = False) are ignored.

required

Returns:

Type Description
Float[Array, 'Ny Nx']

Solution ψ on the full rectangular grid. Exactly zero outside the mask, and ψ ≈ 0 at inner-boundary points.

Source code in .venv/lib/python3.12/site-packages/spectraldiffx/_src/fourier/capacitance.py
def __call__(
    self,
    rhs: Float[Array, "Ny Nx"],
) -> Float[Array, "Ny Nx"]:
    """Solve (∇² − λ)ψ = rhs on the masked domain.

    Given the precomputed capacitance matrix C⁻¹ and Green's functions
    G, the online solve proceeds in four steps:

    1. Rectangular solve:   u = L_rect⁻¹ rhs              [Ny, Nx]
    2. Sample boundary:     u_B = u[j_b, i_b]             [N_b]
    3. Correction weights:  α = C⁻¹ · u_B                 [N_b]
    4. Subtract correction: ψ = u − Σ_k α_k g_k          [Ny, Nx]
    5. Mask exterior:       ψ = ψ * mask                  [Ny, Nx]

    The result satisfies ψ(b_k) ≈ 0 at all N_b inner-boundary points,
    (∇² − λ)ψ = rhs at interior points, and ψ = 0 outside the mask.

    Parameters
    ----------
    rhs : Float[Array, "Ny Nx"]
        Right-hand side on the full rectangular grid.
        Values outside the physical domain (mask = False) are ignored.

    Returns
    -------
    Float[Array, "Ny Nx"]
        Solution ψ on the full rectangular grid.  Exactly zero outside
        the mask, and ψ ≈ 0 at inner-boundary points.
    """
    Ny, Nx = rhs.shape
    # Step 1: rectangular spectral solve
    u = _spectral_solve(rhs, self.dx, self.dy, self.lambda_, self.base_bc)
    # Step 2: values of u at inner-boundary points
    u_b = u[self._j_b, self._i_b]  # [Nb]
    # Step 3: correction coefficients  alpha = C^{-1} u_b
    alpha = self._C_inv @ u_b  # [Nb]
    # Step 4: correction field  sum_k alpha_k g_k
    correction = (self._green_flat.T @ alpha).reshape(Ny, Nx)  # [Ny, Nx]
    # Step 5: zero out exterior cells
    return (u - correction) * self._mask

Factory

finitevolx.build_capacitance_solver(mask, dx, dy, lambda_=0.0, base_bc='fft')

Pre-compute the capacitance matrix and return a ready-to-use solver.

This is an offline function that performs N_b rectangular spectral solves (N_b = number of inner-boundary points). The result is a :class:CapacitanceSolver whose __call__ method is JIT-compilable.

Parameters:

Name Type Description Default
mask np.ndarray of bool shape (Ny, Nx), or ArakawaCGridMask

Physical domain mask. True = interior (ocean/fluid), False = exterior (land/walls).

When an :class:ArakawaCGridMask is passed, the psi staggering mask is extracted automatically.

required
dx float

Grid spacing in x.

required
dy float

Grid spacing in y.

required
lambda_ float

Helmholtz parameter λ. Use 0.0 for pure Poisson.

0.0
base_bc ('fft', 'dst', 'dct')

Rectangular spectral solver used as the base.

"fft"

Returns:

Type Description
CapacitanceSolver

A callable equinox Module with all precomputed arrays baked in.

Source code in finitevolx/_src/solvers/elliptic.py
def build_capacitance_solver(
    mask: np.ndarray | ArakawaCGridMask,
    dx: float,
    dy: float,
    lambda_: float = 0.0,
    base_bc: str = "fft",
) -> CapacitanceSolver:
    """Pre-compute the capacitance matrix and return a ready-to-use solver.

    This is an **offline** function that performs *N_b* rectangular spectral
    solves (``N_b`` = number of inner-boundary points).  The result is a
    :class:`CapacitanceSolver` whose ``__call__`` method is JIT-compilable.

    Parameters
    ----------
    mask : np.ndarray of bool shape (Ny, Nx), or ArakawaCGridMask
        Physical domain mask.  ``True`` = interior (ocean/fluid),
        ``False`` = exterior (land/walls).

        When an :class:`ArakawaCGridMask` is passed, the ``psi``
        staggering mask is extracted automatically.
    dx : float
        Grid spacing in x.
    dy : float
        Grid spacing in y.
    lambda_ : float
        Helmholtz parameter λ.  Use ``0.0`` for pure Poisson.
    base_bc : {"fft", "dst", "dct"}
        Rectangular spectral solver used as the base.

    Returns
    -------
    CapacitanceSolver
        A callable equinox Module with all precomputed arrays baked in.
    """
    if isinstance(mask, ArakawaCGridMask):
        mask = np.asarray(mask.psi, dtype=bool)
    return _build_capacitance_solver_base(mask, dx, dy, lambda_, base_bc)

Masked Laplacian

finitevolx.masked_laplacian(psi, mask, dx, dy, lambda_=0.0)

Apply the masked discrete Helmholtz operator (∇² − λ)·ψ.

Enforces homogeneous Dirichlet conditions at the mask boundary by zeroing psi outside the mask before applying the 5-point stencil. The output is also zeroed outside the mask.

Neighbors at the rectangle edges wrap around (periodic roll), which is consistent with using the FFT as a preconditioner.

Parameters:

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

Field to which the operator is applied.

required
mask Float[Array, 'Ny Nx'] or ArakawaCGridMask

Binary mask: 1 inside the physical domain, 0 outside (land/exterior). When an :class:ArakawaCGridMask is passed, the psi staggering mask is used (all four surrounding h-cells must be wet).

required
dx float

Grid spacing in x.

required
dy float

Grid spacing in y.

required
lambda_ float

Helmholtz parameter. Default: 0.0 (pure Laplacian).

0.0

Returns:

Type Description
Float[Array, 'Ny Nx']

Result of (∇² − λ)·(ψ·mask), zeroed outside the mask.

Source code in finitevolx/_src/solvers/iterative.py
def masked_laplacian(
    psi: Float[Array, "Ny Nx"],
    mask: Float[Array, "Ny Nx"] | ArakawaCGridMask,
    dx: float,
    dy: float,
    lambda_: float = 0.0,
) -> Float[Array, "Ny Nx"]:
    """Apply the masked discrete Helmholtz operator (∇² − λ)·ψ.

    Enforces homogeneous Dirichlet conditions at the mask boundary by zeroing
    *psi* outside the mask before applying the 5-point stencil.  The output
    is also zeroed outside the mask.

    Neighbors at the rectangle edges wrap around (periodic roll), which is
    consistent with using the FFT as a preconditioner.

    Parameters
    ----------
    psi : Float[Array, "Ny Nx"]
        Field to which the operator is applied.
    mask : Float[Array, "Ny Nx"] or ArakawaCGridMask
        Binary mask: 1 inside the physical domain, 0 outside (land/exterior).
        When an :class:`ArakawaCGridMask` is passed, the ``psi`` staggering
        mask is used (all four surrounding h-cells must be wet).
    dx : float
        Grid spacing in x.
    dy : float
        Grid spacing in y.
    lambda_ : float
        Helmholtz parameter.  Default: 0.0 (pure Laplacian).

    Returns
    -------
    Float[Array, "Ny Nx"]
        Result of (∇² − λ)·(ψ·mask), zeroed outside the mask.
    """
    if isinstance(mask, ArakawaCGridMask):
        mask_arr = mask.psi.astype(psi.dtype)
    else:
        mask_arr = mask
    psi_m = psi * mask_arr  # enforce zero outside domain
    # 5-point finite-difference stencil with periodic roll at edges
    lap = (
        jnp.roll(psi_m, 1, axis=1) + jnp.roll(psi_m, -1, axis=1) - 2.0 * psi_m
    ) / dx**2 + (
        jnp.roll(psi_m, 1, axis=0) + jnp.roll(psi_m, -1, axis=0) - 2.0 * psi_m
    ) / dy**2
    return (lap - lambda_ * psi_m) * mask_arr