API reference

StencilCore.ArrayOrTermLikeType
ArrayOrTermLike{T} = Union{AbstractArray{T}, AbstractPointwise{T}}

A coefficient that is either a concrete array (assemblable) or a symbolic term (must be materialized first), both with element type T. Stencil types parameterise their coefficient field over this union.

source
StencilCore.ôConstant
ô, ê₁ … ê₉

Predefined StaticShift constants: ô is the zero shift (identity), and êᵢ is the unit offset along axis i. Combine them with the +/-/*Int algebra to write lattice offsets, e.g. -2ê₁, 3ê₁ + ê₂.

source
StencilCore.AbstractPointwiseType
AbstractPointwise{T} <: AbstractStencil{T}

Supertype of every symbolic grid expression — and, by virtue of being a subtype of AbstractStencil, of every diagonal stencil. An AbstractPointwise{T} behaves like a dimension- and size-less array whose eltype is T: its grid rank N is unknown until it is materialized against concrete arrays, but its element type T (the value each cell will hold once materialized) is fixed at construction.

A pointwise term acts as a diagonal stencil under *(stencil, pointwise): applying it to another pointwise term multiplies element-wise. Concrete subtypes (Slot, Pointwise, Shifted, Fill, IdentityStencil, DiagonalStencil) live in the StencilCalculus package.

source
StencilCore.AbstractStencilType
AbstractStencil{T}

Abstract supertype for every stencil — both the diagonal stencils that live in pointwise-land (AbstractPointwise and its subtypes) and the neighborhood stencils that carry an AccessStyle and a tuple of offsets (NeighborhoodStencil). The single type parameter T is the linear-map space — the element type of the per-cell coefficient, mirroring how Unity{T} carries the multiplicative-identity space in scalar-land.

For a LinearStencil / StarStencil whose coefficient stores SVector{L, F} per cell, T === F (not SVector{L, F}). The match rule for applying a stencil to an AbstractPointwise{U} is T === _unity_space(U): scalar-on- scalar for U <: Number, SMatrix{N, N, F}-on-SVector{N, F} for vector-valued fields.

Provides a Base.eltype accessor — subtypes inherit it without redefining. The AccessStyle trait is defined on NeighborhoodStencil, not on AbstractStencil itself, since diagonal stencils have no offsets.

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.ConstantType
Constant{T}(val) / Constant(val)

Literal scalar leaf carrying a value val::T as-is. Materializes to val. T is any concrete type — Number, SArray, etc. — making Constant the right carrier for boundary literals such as x + 1, v + SVector(1).

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

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.NeighborhoodStencilType
NeighborhoodStencil{T, S<:AccessStyle} <: AbstractStencil{T}

Abstract supertype for stencils that carry off-diagonal offsets and therefore require an AccessStyle anchor. Concrete subtypes — Stencil, LinearStencil, StarStencil — declare <: NeighborhoodStencil{T, S}, carrying the coefficient element type T as their second-to-last type parameter and the access style S as their last type parameter.

Provides the AccessStyle trait accessor; subtypes inherit it without redefining.

source
StencilCore.ScalarType
Scalar(fn, args::Tuple{Vararg{AbstractScalar}})

Interior node of a scalar-tree: applies fn to scalar args component-wise. The element type T = Base.promote_op(fn, eltype.(args)...) is computed at construction; a Union{} result throws (the node is unconstructable). Pointwise-side analogue: AbstractPointwise / Pointwise (lives in StencilCalculus).

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

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.StaticShiftType
StaticShift{P<:Tuple{Vararg{StaticPair}}}

A lattice offset: a normalized collection of StaticPairs. Invariants (enforced by the constructor): pairs are sorted ascending by dimension, no two pairs share a dimension (same-dimension pairs are summed), and no pair has zero offset (dropped). The empty shift StaticShift{Tuple{}} is the identity.

Construct from pairs, or via the basis symbols ê₁ … ê₉ and the +/-/*Int algebra:

3ê₁ + ê₂            # StaticShift{Tuple{StaticPair{1,3}, StaticPair{2,1}}}

Alias: SShift.

source
StencilCore.StencilType
Stencil{M, C<:NTuple{M,StaticShift}, A<:NTuple{M,ArrayOrTermLike}, T, S<:AccessStyle}
    <: NeighborhoodStencil{T, S}

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.UnityType
Unity{T}()

Type-level multiplicative identity / structural one for AbstractScalar. Materializes to one(T); requires T to be a square scalar — a type with one(T) defined (Number, square SMatrix{N,N,F}, …). Construction rejects T lacking one(T) (e.g. SVector, non-square SMatrix).

Lets the scalar simplify rules collapse multiplicative identities by dispatch — structurally, with no .val inspection — mirroring how Null collapses additive identities.

source
StencilCore.VarType
Var{S, T}()

Named, runtime-substituted scalar parameter S (a Symbol) of concrete type T (default Float64). Materializes to the value supplied at the keyword S of the pairs NamedTuple — the scalar-side analogue of a slot variable, but without per-cell indexing.

source
StencilCore._interlaceFunction
_interlace(terms::NTuple{M, ArrayOrTermLike}) -> ArrayOrTermLike{<:SVector{M}}

Combine M per-offset (structure-of-arrays) coefficients into the single array-of-structs SVector{M}-valued coefficient that LinearStencil / StarStencil store — the representation switch performed by narrowing.

The concrete-array method (stack element-wise into an array of SVector{M}) lives here; the symbolic-term method (Term(SVector, terms)) is added by the StencilCalculus package.

source
StencilCore._var_as_fieldMethod
_scalar_body_expr(s::AbstractScalar) -> Expr or literal

Lower a scalar tree to an Expr that evaluates it in the StencilCalculus codegen args::NamedTuple context (no per-cell indices — scalars are position-independent). Used by _body_expr(::Fill{<:AbstractScalar}, …) on the term side. Returns an Expr suitable for embedding in a larger AST.

source
StencilCore.as_linearMethod
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_starMethod
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
StencilCore.derivativeFunction
derivative(f, ::Val{i}, args...) -> AbstractScalar

The symbolic partial derivative ∂f/∂(argᵢ) of a primitive f applied to scalar args (ChainRules frule style). The 0/1 cases use the structural Null / Unity carriers via _unity (which routes through _unity_space so a Number leaf stays Number and an SVector leaf lands in the square SMatrix Jacobian space).

source
StencilCore.differentiateMethod
differentiate(s::AbstractScalar, v::Var{S}) -> AbstractScalar

Differentiate s with respect to the named scalar parameter v (matched on the symbol only). The result is an AbstractScalar whose eltype is the Jacobian element type J = _jacobian_type(eltype(s), eltype(v)): promote_type(Tout, Tin) for Number/Number, SMatrix{N, N, promote_type(F1, F2)} for matching-N SVector{N, F1} / SVector{N, F2}. Mixed shape-classes throw. Returns Null{J}() when s does not depend on v; otherwise the chain-rule expression, simplified. Scalar-side analogue of StencilCalculus.differentiate.

source
StencilCore.materializeMethod
materialize(s::AbstractScalar, pairs::NamedTuple = (;))

Substitute the Var leaves named in pairs into s and reduce to a single value of type eltype(s). The scalar-side analogue of StencilCalculus.materialize.

source
StencilCore.saturate_scalarMethod
saturate_scalar(s::AbstractScalar; theory, timeout, kwargs...) -> AbstractScalar

Simplify a scalar expression via Metatheory e-graph equality saturation.

Requires loading Metatheory:

using StencilCore, Metatheory
saturate_scalar(expr)

simplify is the fast, type-stable alternative for small trees. Use saturate_scalar when you need equational completeness or want to compose with external Metatheory theories via StencilCore.MetatheoryExt.SCALAR_DEFAULT_RULES.

source
StencilCore.scalar_default_rulesFunction
scalar_default_rules() -> Vector{RewriteRule}

Return the default Metatheory rewrite theory for scalar simplification (SCALAR_DEFAULT_RULES in the MetatheoryExt extension).

Requires the MetatheoryExt extension to be loaded (using Metatheory).

source
StencilCore.simplifyMethod
simplify(s::AbstractScalar) -> AbstractScalar

Structurally simplify a scalar expression tree (post-order, single pass):

  • Coefficient fold: Constant(v₁) ⊕ Constant(v₂) → Constant(v₁ ⊕ v₂) for any binary or unary op (type-stable: result eltype is promote_op); also folds mixed combinations where each operand carries a <:Number coefficient (Constant{<:Number}, Unity as 1, Null as 0)
  • Null additive identity: 0+x→x, x+0→x, x−0→x
  • Null subtraction: 0−x→−x (with double-negation folding)
  • Null multiplicative annihilator: 0*x→0, x*0→0, 0/x→0
  • Unity multiplicative identity: 1*x→x, x*1→x, x/1→x (shape-guarded)
  • Unity self-interaction: 1*1→1, 1/1→1
  • Double negation: −(−x)→x

All rules are type-dispatched. The return type is fully inferred at compile time.

source
StencilCore.@addmethodMacro
@addmethod fname scalar (v1, v2, ...)
@addmethod fname scalar (n1 = T1, n2 = T2, ...)

Generate and register a method whose body is derived from the AbstractScalar expression scalar, with each Var{S} replaced by a direct reference to the argument named S. The argument list is a tuple; its shape selects between two forms:

Tuple-of-Var form(v1, v2, ...)fname(v1::T1, v2::T2, ...) = <body>. Each vi must be a Var instance; its type parameter provides the declared argument type and its symbol parameter the argument name.

NamedTuple override form(n1 = T1, ...)fname(n1::T1, ...) = <body>. Argument names and (possibly widened) types come from the NamedTuple, in key order. The scalar tree is still built and validated at the Vars' concrete types — only the emitted signature widens. Each override must be a supertype of the Var's declared type (widening only). Because StaticArrays element parameters are invariant (SVector{2,Float64} <: SVector{2,Number} is false), shape widening must use the UnionAll form SVector{2,<:Number}.

A single-element list needs a trailing comma — (v,) / (n = T,) — like any Julia tuple literal.

If fname already has methods, the new method is added via ordinary Julia multiple-dispatch — no existing method is overwritten, so both forms can register methods for the same fname.

An ArgumentError is raised at call time if:

  • any Var appearing in scalar is not listed among the arguments,
  • scalar calls fname itself (the generated method would recurse infinitely; e.g. getindex cannot be defined by indexing the same argument),
  • (Var form) the argument list contains duplicate names or a listed argument is not a Var instance,
  • (override form) a listed name is not a Var in the tree, a value is not a Type, or an override is not a supertype of the declared type.

Examples

@var x Float64
@var y Int

foo = x ^ y
@addmethod myfun foo (y, x)                       # myfun(y::Int, x::Float64) = x ^ y
@addmethod myfun foo (y = Integer, x = AbstractFloat)  # myfun(y::Integer, x::AbstractFloat) = x ^ y

using StaticArrays
@var A SMatrix{2, 2, Float64, 4}
@var u SVector{2, Float64}

bar = A * u
@addmethod myfun bar (A, u)    # adds  myfun(A::SMatrix{2,2,Float64,4}, u::SVector{2,Float64}) = A * u

@var V SVector{2, Float64}
vv = V + V
@addmethod vfun vv (V = SVector{2, <:Number},)   # vfun(V::SVector{2,<:Number}) = V + V
source
StencilCore.@varMacro
@var name [T = Float64]

Bind name to Var{:name, T}(). @var ττ = Var{:τ, Float64}(); @var τ Float32τ = Var{:τ, Float32}(). Term-side analogue: @slot (lives in StencilCalculus).

source