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):
- Solve the PDE on the full rectangle using a fast spectral solver
(DST/DCT/FFT), ignoring the mask. Call this
u. ugenerally violates ψ = 0 at inner-boundary points. Correct it:ψ = u − Σ_k α_k g_k, whereg_kare precomputed Green's functions (rectangular-domain response to δ-sources at each boundary point b_k).- The coefficients α are found by requiring ψ(b_k) = 0 at all
N_b boundary points, giving the linear system
C α = u[B]whereC[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
61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 | |
__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:
- Rectangular solve: u = L_rect⁻¹ rhs [Ny, Nx]
- Sample boundary: u_B = u[j_b, i_b] [N_b]
- Correction weights: α = C⁻¹ · u_B [N_b]
- Subtract correction: ψ = u − Σ_k α_k g_k [Ny, Nx]
- 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
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. When an :class: |
required |
dx
|
float
|
Grid spacing in x. |
required |
dy
|
float
|
Grid spacing in y. |
required |
lambda_
|
float
|
Helmholtz parameter λ. Use |
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
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: |
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. |