Capacitance Solver¶
Capacitance matrix method for elliptic PDEs on masked/irregular domains. See the theory page for the algorithm.
build_capacitance_solver(mask, dx, dy, lambda_=0.0, base_bc='fft')
¶
Pre-compute the capacitance matrix and return a ready-to-use solver.
Offline algorithm (Buzbee, Golub & Nielson 1970):
- Detect inner boundary — find the N_b mask-interior cells that are
4-connected to at least one exterior cell (using
scipy.ndimage.binary_dilationwith a cross-shaped structuring element). - Compute Green's functions — for each boundary point b_k, solve
L_rect g_k = e_{b_k}on the full rectangle using the base spectral solver. Stores G as a [N_b, Ny*Nx] matrix. - Build capacitance matrix —
C[k, l] = g_l(b_k), i.e. the response at boundary point k due to a unit source at boundary point l. Shape: [N_b, N_b]. - Invert —
C⁻¹ = inv(C)via dense NumPy (offline, not JIT-traced).
Complexity¶
- Offline (this function): O(N_b · Ny · Nx · log(Ny·Nx)) time, O(N_b · Ny · Nx) memory for the Green's function matrix.
- Online (
CapacitanceSolver.__call__): O(N_b² + Ny · Nx · log(Ny·Nx)) time per solve.
Parameters¶
mask : np.ndarray of bool, shape (Ny, Nx)
Physical domain mask. True = interior (ocean/fluid),
False = exterior (land/walls).
Inner-boundary points are computed as wet (True) cells that are
4-connected to at least one dry (False) cell.
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.
Raises¶
ValueError If the mask has no inner-boundary points (e.g. all-ones mask).
Source code in spectraldiffx/_src/fourier/capacitance.py
158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 | |
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¶
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 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 | |
Functions¶
__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¶
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.