API reference

Assembly

StencilAssembly.assembleFunction
assemble(st, row, col) -> SparseMatrixCSC

Build the sparsity pattern (colptr, rowval) of the operator induced by st between row and col and allocate nzval uninitialised. Call update! to populate nzval, or use build to do both in one shot.

row and col are interpreted on a shared integer mesh: row[d], col[d] :: AbstractUnitRange{Int} cover the d-th coordinate of the operator's row / column index space.

Each stencil type adds its own methods — see LinearStencil and StarStencil.

source
StencilAssembly.assembleMethod
assemble(st::LinearStencil{1, O, L, T, A, ColumnAccess},
         row::NTuple{1, AbstractUnitRange{Int}},
         col::NTuple{1, AbstractUnitRange{Int}}) -> SparseMatrixCSC{T, Int}

Build the sparsity pattern (colptr, rowval) of the operator induced by st between row and col and allocate nzval uninitialised. Call update! to populate nzval, or use build to do both in one shot.

row and col are interpreted on a shared integer mesh; the entry at column c ∈ col[1] and offset δ lands on row c − δ iff it lies in row[1].

Dispatch pins D = 1, N = 1, and S = ColumnAccess (misuse → MethodError); the L − 1 ≤ length(row[1]) guard is the three-phase kernel's exact correctness boundary.

source
StencilAssembly.assembleMethod
assemble(st::StarStencil{L, 1, M, T, A, ColumnAccess}, row, col) -> SparseMatrixCSC{T, Int}

1-D StarStencil assembly delegates to the equivalent LinearStencil{1}.

source
StencilAssembly.assembleMethod
assemble(st::StarStencil{L, N, M, T, A, ColumnAccess}, row, col) -> SparseMatrixCSC{T, Int}

N-D entry (N ≥ 2, S = ColumnAccess). Enforces the per-axis guard 2L ≤ length(row[d]). Builds colptr / rowval and allocates uninitialised nzval; call update! to populate, or use build.

source
StencilAssembly.buildMethod
build(st, row, col) -> SparseMatrixCSC

Equivalent to update!(assemble(st, row, col), st, row, col). Use when you don't need the pattern / values split — the returned matrix is fully populated and ready to use.

source
StencilAssembly.update!Function
update!(mat, st, row, col) -> mat

Write mat.nzval in place by re-walking row / col with st. mat must have been produced by a matching assemble (same st, row, col) so its colptr / rowval align with the kernel's sweep.

Each stencil type adds its own methods — see LinearStencil and StarStencil.

source
StencilAssembly.update!Method
update!(mat::SparseMatrixCSC{T, Int},
        st::LinearStencil{1, O, L, T, A, ColumnAccess},
        row::NTuple{1, AbstractUnitRange{Int}},
        col::NTuple{1, AbstractUnitRange{Int}}) -> mat

Write mat.nzval in place by re-walking row/col with st. mat must have been produced by a matching assemble (same st, row, col) so its colptr/rowval align with the kernel's sweep. Carries the same L − 1 ≤ length(row[1]) guard.

source
StencilAssembly.update!Method
update!(mat, st::StarStencil{L, N, M, T, A, ColumnAccess}, row, col) -> mat

N-D in-place value update (N ≥ 2, S = ColumnAccess); same guard as assemble. mat must come from a matching assemble.

source
StencilAssembly.update!Method
update!(mat, st::StarStencil{L, 1, M, T, A, ColumnAccess}, row, col) -> mat

1-D StarStencil update delegates to the equivalent LinearStencil{1}.

source

Stencil types (re-exported from StencilCore)

Full documentation lives in StencilCore; they are re-exported here for convenience.

StencilCore.LinearStencilType
LinearStencil{D, O, L, T, A<:ArrayOrTermLike{SVector{L, T}}, S<:AccessStyle}
    <: AbstractStencil{S, T}

Variable-coefficient stencil with contiguous offsets, aligned with mesh dimension D. Offsets are diagonal indices: for a column j and a row i the diagonal number is k = j − i. For column c at mesh position p_c and offset δ, the matrix entry lands on row p_c − δ with coefficient term[p_c][δ − O + 1] (column-anchored under S = ColumnAccess); each element term[p_c] is the SVector{L, T} of all L per-offset coefficients on that column in ascending offset order (term[p_c][1] ↦ offset O).

Type parameters:

  • D: mesh dim along which the stencil acts (D ≥ 1; for a concrete coefficient also D ≤ ndims(term)).
  • O = δ_min, L = δ_max − δ_min + 1 — static via SUnitRange.
  • T: scalar coefficient element type (the per-cell coefficient is SVector{L, T}). This is the linear-map space surfaced as the first parameter of AbstractStencil / NeighborhoodStencil.
  • A<:ArrayOrTermLike{SVector{L, T}}: the coefficient container — a concrete AbstractArray{SVector{L, T}} (assemblable) or a symbolic AbstractPointwise{SVector{L, T}}. N (grid rank) is not a stencil parameter: it is ndims(A) for a concrete coefficient, and unknown (resolved at materialize) for a symbolic one.
  • S<:AccessStyle: coefficient anchoring (ColumnAccess for CSC).

The default outer constructor defaults S to ColumnAccess. A catch-all reports ArgumentErrors for non-SUnitRange offsets and ill-typed coefficients.

source
StencilCore.StarStencilType
StarStencil{L, N, M, T, A<:ArrayOrTermLike{SVector{M, T}}, S<:AccessStyle}
    <: AbstractStencil{S, T}

N-D variable-coefficient star-shaped stencil with symmetric reach −L … +L along every mesh dimension, stored interlaced: a single coefficient term::A whose element at each cell is the SVector{M, T} of the whole star, M = 2NL + 1. The entries are in reverse-lexicographic offset order (axis N most significant — the CartesianIndex order) with the diagonal as the explicit middle slot (M+1)/2. Unlike a per-axis decomposition, the diagonal is a single free coefficient (Helmholtz , parabolic ∂ₜ), not a sum of per-axis centers.

For L = 2, N = 3 (M = 13) the slot ↦ offset map is:

1:(d3,-2) 2:(d3,-1) 3:(d2,-2) 4:(d2,-1) 5:(d1,-2) 6:(d1,-1) 7:(diag)
8:(d1,+1) 9:(d1,+2) 10:(d2,+1) 11:(d2,+2) 12:(d3,+1) 13:(d3,+2)

Type parameters:

  • L ≥ 1 per-axis reach; M = 2NL + 1 whole-star offset count.
  • N grid rank, kept explicit; checked to equal (M-1)/(2L) (and the coefficient array's ndims, when concrete).
  • T scalar coefficient element type (the per-cell coefficient is SVector{M, T}). The linear-map space surfaced as the first parameter of AbstractStencil / NeighborhoodStencil.
  • A<:ArrayOrTermLike{SVector{M, T}}: the (single) coefficient container — concrete array (assemblable) or symbolic term.
  • S<:AccessStyle.
source
StencilCore.StencilType
Stencil{M, C<:NTuple{M,StaticShift}, A<:NTuple{M,ArrayOrTermLike}, T, S<:AccessStyle}
    <: AbstractStencil{S, T}

General offset-list stencil in structure-of-arrays form: M offsets shifts (reverse-lexicographically ordered) and a parallel M-tuple terms of per-offset coefficients (terms[k] is the coefficient at offset shifts[k]). Each coefficient is an ArrayOrTermLike — a concrete array or a symbolic term — with a scalar element type, and all coefficients must share the same eltype === T.

The form produced by symbolic differentiation; not assembled directly — narrowed to a LinearStencil / StarStencil via as_linear / as_star, which is where the layout switches to array-of-structs (one SVector{M}-valued coefficient) by _interlaceing terms.

source
StencilCore.AccessStyleType
AccessStyle

Holy trait reporting how a stencil's coefficient arrays are anchored:

  • ColumnAccess: term[c_idx...][k] is the value at column mesh position c_idx. Required for assembly into SparseMatrixCSC.
  • RowAccess: term[r_idx...][k] is the value at row mesh position r_idx. Required for assembly into a row-major format (CSR; not yet implemented).

The trait is a type parameter on the stencil, not a runtime field. Assemblers dispatch on it; mismatching the access style and the target sparse format is a dispatch-time error (MethodError).

source
StencilCore.as_linearFunction
as_linear(st::Stencil{…,S}) -> LinearStencil{D, …, S}

Narrow a Stencil whose offsets are single-axis (same axis D) and contiguous to the equivalent LinearStencil{D}, interlacing terms into the single SVector{D}-valued coefficient. Throws if the offsets are multi-axis, span several axes, or are not contiguous-ascending.

source
StencilCore.as_starFunction
as_star(st::Stencil{…,S}) -> StarStencil{L, …, S}

Narrow a Stencil whose offsets form the canonical reverse-lex star pattern (every axis 1..N with symmetric reach −L..L, plus the diagonal) to the equivalent StarStencil{L}, interlacing terms into the single SVector{M}- valued coefficient. Throws otherwise.

source