Skip to content

Exported Functions

RadialBasisFunctions.PHS Method
julia
function PHS(n::T=3; poly_deg::T=2) where {T<:Int}

Convienience contructor for polyharmonic splines.

source
RadialBasisFunctions.classify_stencil Method
julia
classify_stencil(is_boundary, boundary_conditions, eval_idx,
                neighbors, global_to_boundary)

Classify stencil type for dispatch in kernel execution.

source
RadialBasisFunctions.directional Method
julia
function directional(data, eval_points, v, basis, is_boundary, boundary_conditions, normals; k=autoselect_k(data, basis))

Builds a Hermite-compatible RadialBasisOperator where the operator is the directional derivative, Directional. The additional boundary information enables Hermite interpolation with proper boundary condition handling.

source
RadialBasisFunctions.directional Method
julia
function directional(data, eval_points, v, basis; k=autoselect_k(data, basis))

Builds a RadialBasisOperator where the operator is the directional derivative, Directional.

source
RadialBasisFunctions.directional Method
julia
function directional(data, v, basis; k=autoselect_k(data, basis))

Builds a RadialBasisOperator where the operator is the directional derivative, Directional.

source
RadialBasisFunctions.gradient Method
julia
function gradient(data, basis; k=autoselect_k(data, basis))

Builds a RadialBasisOperator where the operator is the gradient, Gradient.

source
RadialBasisFunctions.gradient Method
julia
function gradient(data, eval_points, basis, is_boundary, boundary_conditions, normals; k=autoselect_k(data, basis))

Builds a Hermite-compatible RadialBasisOperator where the operator is the gradient, Gradient. The additional boundary information enables Hermite interpolation with proper boundary condition handling.

source
RadialBasisFunctions.gradient Method
julia
function gradient(data, eval_points, basis; k=autoselect_k(data, basis))

Builds a RadialBasisOperator where the operator is the gradient, Gradient. The resulting operator will only evaluate at eval_points.

source
RadialBasisFunctions.laplacian Method
julia
function laplacian(data, eval_points, basis, is_boundary, boundary_conditions, normals; k=autoselect_k(data, basis))

Builds a Hermite-compatible RadialBasisOperator where the operator is the Laplacian, Laplacian. The additional boundary information enables Hermite interpolation with proper boundary condition handling.

source
RadialBasisFunctions.partial Method
julia
function partial(data, eval_points, order, dim, basis, is_boundary, boundary_conditions, normals; k=autoselect_k(data, basis))

Builds a Hermite-compatible RadialBasisOperator where the operator is the partial derivative, Partial. The additional boundary information enables Hermite interpolation with proper boundary condition handling.

source
RadialBasisFunctions.partial Method
julia
function partial(data, eval_points, order, dim, basis; k=autoselect_k(data, basis))

Builds a RadialBasisOperator where the operator is the partial derivative, Partial. The resulting operator will only evaluate at eval_points.

source
RadialBasisFunctions.partial Method
julia
function partial(data, order, dim, basis; k=autoselect_k(data, basis))

Builds a RadialBasisOperator where the operator is the partial derivative, Partial, of order with respect to dim.

source
RadialBasisFunctions.regrid Method
julia
function regrid(data, eval_points, basis=PHS(3; poly_deg=2); k=autoselect_k(data, basis))

Builds a RadialBasisOperator where the operator is an regrid from one set of points to another, data -> eval_points.

source
RadialBasisFunctions.update_hermite_stencil_data! Method
julia
update_hermite_stencil_data!(hermite_data, global_data, neighbors,
                             is_boundary, boundary_conditions, normals,
                             global_to_boundary)

Populate local Hermite stencil data from global arrays. Used within kernels to extract boundary info for specific neighborhoods.

source
RadialBasisFunctions.∂virtual Method
julia
function ∂virtual(data, eval_points, dim, Δ, basis; k=autoselect_k(data, basis))

Builds a virtual RadialBasisOperator whichi will be evaluated at eval_points where the operator is the partial derivative with respect to dim. Virtual operators interpolate the data to structured points at a distance Δ for which standard finite difference formulas can be applied.

source
RadialBasisFunctions.∂virtual Method
julia
function ∂virtual(data, dim, Δ, basis; k=autoselect_k(data, basis))

Builds a virtual RadialBasisOperator whichi will be evaluated at the input points (data) where the operator is the partial derivative with respect to dim. Virtual operators interpolate the data to structured points at a distance Δ for which standard finite difference formulas can be applied.

source
RadialBasisFunctions.AbstractPHS Type

abstract type AbstractPHS <: AbstractRadialBasis

Supertype of all Polyharmonic Splines.

source
RadialBasisFunctions.AbstractRadialBasis Type
julia
abstract type AbstractRadialBasis <: AbstractBasis end
source
RadialBasisFunctions.BoundaryCondition Type
julia
BoundaryCondition{T}

Unified boundary condition representation: Bu = α_u + β_∂ₙu

Special cases:

  • Dirichlet: α=1, β=0

  • Neumann: α=0, β=1

  • Robin: α≠0, β≠0

  • Internal: α=0, β=0 (sentinel for interior points)

source
RadialBasisFunctions.Custom Type
julia
Custom <: ScalarValuedOperator

Builds an operator for a first order partial derivative.

source
RadialBasisFunctions.Directional Type
julia
Directional{Dim,T} <: ScalarValuedOperator

Operator for the directional derivative, or the inner product of the gradient and a direction vector.

source
RadialBasisFunctions.Gaussian Type
julia
struct Gaussian{T,D<:Int} <: AbstractRadialBasis

Gaussian radial basis function:ϕ(r)=e(εr)2

source
RadialBasisFunctions.Gradient Type
julia
Gradient{Dim} <: VectorValuedOperator

Builds an operator for the gradient of a function.

source
RadialBasisFunctions.HermiteStencilData Type
julia
HermiteStencilData{T}

Local stencil data for Hermite interpolation with boundary conditions.

Fields:

  • data: Coordinates of k stencil points

  • is_boundary: Boolean flags for each point

  • boundary_conditions: BC for each point (use Internal() for interior)

  • normals: Normal vectors (zero for interior points)

Note: For interior points (is_boundary[i] == false), boundary_conditions[i] and normals[i] contain sentinel values and should not be accessed.

source
RadialBasisFunctions.HermiteStencilData Method

Pre-allocation constructor for HermiteStencilData

source
RadialBasisFunctions.Interpolator Type
julia
struct Interpolator

Construct a radial basis interpolation.

source
RadialBasisFunctions.Interpolator Method
julia
function Interpolator(x, y, basis::B=PHS())

Construct a radial basis interpolator.

source
RadialBasisFunctions.Laplacian Type
julia
Laplacian <: ScalarValuedOperator

Operator for the sum of the second derivatives w.r.t. each independent variable.

source
RadialBasisFunctions.MonomialBasis Type
julia
struct MonomialBasis{Dim,Deg} <: AbstractBasis

Dim dimensional monomial basis of order Deg.

source
RadialBasisFunctions.PHS1 Type
julia
struct PHS1{T<:Int} <: AbstractPHS

Polyharmonic spline radial basis function:ϕ(r)=r

source
RadialBasisFunctions.PHS3 Type
julia
struct PHS3{T<:Int} <: AbstractPHS

Polyharmonic spline radial basis function:ϕ(r)=r3

source
RadialBasisFunctions.PHS5 Type
julia
struct PHS5{T<:Int} <: AbstractPHS

Polyharmonic spline radial basis function:ϕ(r)=r5

source
RadialBasisFunctions.PHS7 Type
julia
struct PHS7{T<:Int} <: AbstractPHS

Polyharmonic spline radial basis function:ϕ(r)=r7

source
RadialBasisFunctions.Partial Type
julia
Partial <: ScalarValuedOperator

Builds an operator for a first order partial derivative.

source
RadialBasisFunctions.RadialBasisOperator Type
julia
struct RadialBasisOperator

Operator of data using a radial basis with potential monomial augmentation.

source
RadialBasisFunctions.RadialBasisOperator Method
julia
function (op::RadialBasisOperator)(y, x)

Evaluate the operator at x in-place and store the result in y.

source
RadialBasisFunctions.RadialBasisOperator Method
julia
function (op::RadialBasisOperator)(x)

Evaluate the operator at x.

source
RadialBasisFunctions.Regrid Type
julia
Regrid

Builds an operator for interpolating from one set of points to another.

source

Private

RadialBasisFunctions._build_collocation_matrix! Method
julia
_build_collocation_matrix!(A, data, basis, mon, k)

Build RBF collocation matrix for interior stencil (no boundary conditions).

Matrix structure:

julia
┌─────────────────┬─────────┐
Φ(xᵢ, xⱼ)      │ P(xᵢ)   │  k×k RBF + k×nmon polynomial
├─────────────────┼─────────┤
P(xⱼ)ᵀ         │   0     │  nmon×k poly + nmon×nmon zero
└─────────────────┴─────────┘
source
RadialBasisFunctions._build_collocation_matrix! Method
julia
_build_collocation_matrix!(A, data::HermiteStencilData, basis, mon, k)

Build RBF collocation matrix for Hermite stencil (with boundary conditions).

For boundary points with Neumann/Robin conditions, basis functions are modified:

  • Instead of Φ(·,xⱼ), we use B₂Φ(·,xⱼ) where B is the boundary operator

  • This maintains matrix symmetry

source
RadialBasisFunctions._build_rhs! Method
julia
_build_rhs!(b, ℒrbf::Tuple, ℒmon::Tuple, data, eval_point, basis, k)

Build RHS vector for interior stencil, multiple operators.

source
RadialBasisFunctions._build_rhs! Method
julia
_build_rhs!(b, ℒrbf::Tuple, ℒmon::Tuple, data::HermiteStencilData, eval_point, basis, mon, k)

Build RHS vector for Hermite stencil, multiple operators.

source
RadialBasisFunctions._build_rhs! Method
julia
_build_rhs!(b, ℒrbf, ℒmon, data, eval_point, basis, k)

Build RHS vector for interior stencil, single operator.

source
RadialBasisFunctions._build_rhs! Method
julia
_build_rhs!(b, ℒrbf, ℒmon, data::HermiteStencilData, eval_point, basis, mon, k)

Build RHS vector for Hermite stencil, single operator. Applies boundary conditions to both stencil points and evaluation point.

source
RadialBasisFunctions._build_stencil! Method
julia
_build_stencil!(A, b, ℒrbf, ℒmon, data, eval_point, basis, mon, k)

Assemble complete stencil: build collocation matrix, build RHS, solve for weights. Works for both interior and Hermite stencils via multiple dispatch on data type.

Returns: weights (first k rows of solution, size k×num_ops)

source
RadialBasisFunctions._build_weights Method
julia
_build_weights(ℒ, data, eval_points, adjl, basis)

Apply operator to basis functions and route to weight computation.

source
RadialBasisFunctions._build_weights Method
julia
_build_weights(data, eval_points, adjl, basis, ℒrbf, ℒmon, mon;
              batch_size=10, device=CPU())

Build weights for interior-only problems (no boundary conditions). Creates empty BoundaryData to indicate all interior points.

source
RadialBasisFunctions._build_weights Method
julia
_build_weights(ℒ, data, eval_points, adjl, basis,
              is_boundary, boundary_conditions, normals)

Generic Hermite dispatcher for operators. Applies operator to basis and routes to Hermite weight computation.

This eliminates repetitive _build_weights methods across operator files. Note: Type constraint removed to avoid circular dependency with operators.jl

source
RadialBasisFunctions._build_weights Method
julia
_build_weights(data, eval_points, adjl, basis, ℒrbf, ℒmon, mon,
              is_boundary, boundary_conditions, normals;
              batch_size=10, device=CPU())

Build weights for problems with boundary conditions using Hermite interpolation. Exact allocation: Dirichlet points get single entry, others get full stencil.

source
RadialBasisFunctions._build_weights Method
julia
_build_weights(ℒ, op)

Entry point from operator construction. Extracts configuration from operator and routes to appropriate implementation.

source
RadialBasisFunctions._build_weights Method
julia
_build_weights(ℒ::Directional, data, eval_points, adjl, basis, is_boundary, boundary_conditions, normals)

Hermite-compatible method for building directional derivative weights with boundary condition support.

source
RadialBasisFunctions._num_ops Method

Get number of operators (type-stable)

source
RadialBasisFunctions._prepare_buffer Method

Prepare RHS buffer with correct type (type-stable)

source
RadialBasisFunctions.allocate_sparse_arrays Method
julia
allocate_sparse_arrays(TD, k, N_eval, num_ops, adjl, boundary_data)

Allocate sparse matrix arrays for COO format sparse matrix construction. Exactly counts non-zeros: interior points get k entries, Dirichlet points get 1 entry.

source
RadialBasisFunctions.autoselect_k Method
julia
autoselect_k(data::Vector, basis<:AbstractRadialBasis)

See Bayona, 2017 - https://doi.org/10.1016/j.jcp.2016.12.008

source
RadialBasisFunctions.build_weights_kernel Method
julia
build_weights_kernel(data, eval_points, adjl, basis, ℒrbf, ℒmon, mon,
                    boundary_data; batch_size, device)

Main orchestrator: allocate memory, launch kernel, construct sparse matrix.

source
RadialBasisFunctions.calculate_batch_range Method

Calculate batch index range for kernel execution

source
RadialBasisFunctions.compute_hermite_poly_entry! Method
julia
compute_hermite_poly_entry!(a, i, data, mon)

Compute polynomial entries for Hermite interpolation. Applies boundary operators to polynomial basis at boundary points.

source
RadialBasisFunctions.compute_hermite_rbf_entry Method
julia
compute_hermite_rbf_entry(i, j, data, basis)

Compute single RBF matrix entry with Hermite boundary modifications. Dispatches based on point types (Interior/Dirichlet/NeumannRobin).

source
RadialBasisFunctions.construct_global_to_boundary Method
julia
construct_global_to_boundary(is_boundary)

Construct mapping from global indices to boundary-only indices. For boundary points: global_to_boundary[i] = boundary array index For interior points: global_to_boundary[i] = 0 (sentinel)

source
RadialBasisFunctions.count_nonzeros Method
julia
count_nonzeros(adjl, is_boundary, boundary_conditions)

Count exact number of non-zero entries for optimized allocation. Returns (total_nnz, row_offsets) where row_offsets[i] is the starting position for row i.

source
RadialBasisFunctions.fill_dirichlet_entry! Method

Fill Dirichlet identity row for optimized allocation

source
RadialBasisFunctions.fill_entries! Method

Fill sparse matrix entries using indexed storage (row_offsets)

source
RadialBasisFunctions.hermite_mono_rhs! Method
julia
hermite_mono_rhs!(bmono, ℒmon, mon, eval_point, is_bound, bc, normal, T)

Apply boundary conditions to monomial operator evaluation.

  • Interior/Dirichlet: apply operator ℒ to monomials

  • Neumann/Robin: α_ℒ(P) + β_ℒ(∂ₙP)

source
RadialBasisFunctions.hermite_rbf_rhs Method
julia
hermite_rbf_rhs(ℒrbf, eval_point, data_point, is_bound, bc, normal)

Apply boundary conditions to RBF operator evaluation.

  • Interior/Dirichlet: standard evaluation ℒΦ(x_eval, x_data)

  • Neumann/Robin: α_ℒΦ + β_ℒ(∂ₙΦ)

source
RadialBasisFunctions.launch_kernel! Method
julia
launch_kernel!(I, J, V, data, eval_points, adjl, basis, ℒrbf, ℒmon, mon,
               boundary_data, row_offsets, batch_size, N_eval, n_batches,
               k, nmon, num_ops, device)

Launch parallel kernel for weight computation. Handles Dirichlet/Interior/Hermite stencil classification via dispatch.

source
RadialBasisFunctions.operator_arity Method

Extract operator arity at compile time

source
RadialBasisFunctions.point_type Method

Determine boundary type of a single point

source
RadialBasisFunctions.AbstractBasis Type
julia
abstract type AbstractBasis end
source
RadialBasisFunctions.BoundaryData Type
julia
BoundaryData{T}

Wrapper for global boundary information (replaces fragile tuples).

source
RadialBasisFunctions.BoundaryPointType Type

Trait types for individual point boundary classification

source
RadialBasisFunctions.OperatorArity Type

Operator arity traits for compile-time dispatch

source
RadialBasisFunctions.StencilType Type

Trait types for stencil classification

source