API Reference

High-Level API

AdaptiveCrossApproximation.H.assembleFunction
H.assemble(op, testspace, trialspace; tree=..., kwargs...)

Assemble a hierarchical matrix with automatic tree-based blocking.

High-level convenience function that automatically constructs a hierarchical clustering tree and assembles the complete HMatrix with near and far-field blocks. This is the recommended entry point for most applications.

Arguments

  • op: Operator/kernel for matrix entry evaluation

  • testspace: Test space (row basis/points)

  • trialspace: Trial space (column basis/points)

  • tree: Hierarchical tree structure (auto-generated if not provided)

  • kwargs...: Additional parameters passed to HMatrix:

    • tol: Convergence tolerance (default 1e-4)
    • maxrank: Maximum rank for compression (default 40)
    • compressor: ACA or iACA instance (default ACA(tol=tol))
    • isnear: Admissibility predicate (default isnear())
    • spaceordering: Space ordering strategy (default PermuteSpaceInPlace())

Returns

  • HMatrix: Assembled hierarchical matrix ready for matrix-vector products

Notes

Requires H2Trees.jl to be loaded for automatic tree generation. If providing a custom tree, ensure it implements the tree interface expected by HMatrix assembly.

Examples

# With BEAST boundary integral operators
hmat = H.assemble(op, testspace, trialspace; tol=1e-5, maxrank=100)
y = hmat * x  # Efficient matrix-vector product

See also

HMatrix

source
AdaptiveCrossApproximation.AdaptiveCrossApproximationModule
AdaptiveCrossApproximation

Adaptive Cross Approximation (ACA) algorithms for hierarchical low-rank matrix compression.

Provides full-rank and incomplete adaptive cross approximation algorithms optimized for boundary integral operators and other kernel-based matrices. Key features:

  • ACA: Standard adaptive cross approximation selecting pivots by alternating rows/columns
  • iACA: Incomplete variant for geometric pivoting in hierarchical matrix construction
  • HMatrix: Hierarchical matrix representation combining dense near-field and low-rank far-field blocks
  • Pivoting strategies: Maximum value, geometric (Leja, Fill Distance, Mimicry), random sampling
  • Convergence criteria: Frobenius norm estimation, random sampling, extrapolation, combined criteria
  • Kernel matrices: Support for BEAST boundary integral operators and point-based kernels

Main API

Compressors:

  • ACA: Standard row-first variant
  • iACA: Incomplete variant for geometric pivoting and hierarchical matrices
  • aca: Convenience function for matrix compression

Hierarchical matrices:

Pivoting strategies:

Convergence criteria:

Kernel matrices:

Example

using AdaptiveCrossApproximation

# Compress a matrix with default settings
U, V = aca(M; tol=1e-4, maxrank=50)
approx = U * V

# Or build an HMatrix with automatic tree-based blocking
hmat = H.assemble(operator, testspace, trialspace)
y = hmat * x  # Matrix-vector product

See also

  • BEAST.jl for boundary integral operators
  • H2Trees.jl for hierarchical clustering
  • BlockSparseMatrices.jl for sparse block storage
source
AdaptiveCrossApproximation.ACAType
ACA{RowPivType,ColPivType,ConvCritType}

Adaptive Cross Approximation (ACA) compressor for low-rank matrix approximation.

Computes M ≈ U * V by iteratively selecting rows and columns (pivots) until a convergence criterion is met. The algorithm starts with row, samples it to select a column pivot, then alternates between row and column selection.

Fields

  • rowpivoting::RowPivType: Strategy for selecting row pivots
  • columnpivoting::ColPivType: Strategy for selecting column pivots
  • convergence::ConvCritType: Convergence criterion to stop iterations
source
AdaptiveCrossApproximation.ACAMethod
(aca::ACA)(A, colbuffer, rowbuffer, rows, cols, rowidcs, colidcs, maxrank)

Compute ACA approximation with preallocated buffers (main computational routine).

Fills colbuffer and rowbuffer with low-rank factors U and V such that A[rowidcs, colidcs] ≈ U * V. Uses deflation to ensure orthogonality of pivots.

Arguments

  • A: Matrix or matrix-like object (must support nextrc! interface)
  • colbuffer::AbstractArray{K}: Buffer for U factors, size (length(rowidcs), maxrank)
  • rowbuffer::AbstractArray{K}: Buffer for V factors, size (maxrank, length(colidcs))
  • rows::Vector{Int}: Storage for selected row indices
  • cols::Vector{Int}: Storage for selected column indices
  • rowidcs::Vector{Int}: Global row indices of the block to compress
  • colidcs::Vector{Int}: Global column indices of the block to compress
  • maxrank::Int: Maximum number of pivots (hard limit on rank)

Returns

  • npivot::Int: Number of pivots computed (≤ maxrank). The approximation is A[rowidcs, colidcs] ≈ colbuffer[:, 1:npivot] * rowbuffer[1:npivot, :]
source
AdaptiveCrossApproximation.AbstractKernelMatrixType
AbstractKernelMatrix{T}

Abstract matrix-like interface for kernel-based entry evaluation used by ACA-style compressors.

Arguments

  • T: scalar element type returned by kernel evaluations

Returns

A subtype that supports lazy matrix entry access through the kernel matrix interface.

Notes

Implement this type when matrix entries are computed on demand from geometric/operator data.

See also

AbstractKernelMatrix(operator, testspace, trialspace; args...)

source
AdaptiveCrossApproximation.AbstractKernelMatrixMethod
AbstractKernelMatrix(operator, testspace::AbstractVector, trialspace::AbstractVector; args...)

Create a PointMatrix kernel wrapper from vector point data and operator.

Constructs a PointMatrix when given vector collections of test and trial points. Matches the exported AbstractKernelMatrix factory interface to provide point-based kernel evaluation.

Arguments

  • operator: kernel function or operator
  • testspace::AbstractVector: test point locations (indexed collection)
  • trialspace::AbstractVector: trial point locations (indexed collection)
  • args...: additional keyword arguments (unused for point-based evaluation)

Returns

  • PointMatrix with element type inferred from operator
source
AdaptiveCrossApproximation.AbstractKernelMatrixMethod
AbstractKernelMatrix(operator, testspace, trialspace; args...)

Construct a concrete kernel matrix wrapper from operator and space data.

Arguments

  • operator: operator or kernel definition
  • testspace: space for row evaluation points or basis data
  • trialspace: space for column evaluation points or basis data
  • args...: backend-specific keyword arguments

Returns

A concrete subtype of AbstractKernelMatrix provided by method dispatch.

Notes

This declaration defines the interface entry point. Concrete backends provide specialized methods for specific operator/space types.

See also

AbstractKernelMatrix, nextrc!

source
AdaptiveCrossApproximation.AbstractKernelMatrixMethod
AbstractKernelMatrix(operator::Function, testspace::AbstractVector, trialspace::AbstractVector; args...)

Create a PointMatrix from a plain Julia function kernel (with warning).

Provides fallback support for plain function kernels. Issues a warning as operators with structure (kernels, boundary integral operators) are typically preferred for better performance and correctness.

Arguments

  • operator::Function: plain Julia function with signature operator(testpoint, trialpoint)
  • testspace::AbstractVector: test point locations (indexed collection)
  • trialspace::AbstractVector: trial point locations (indexed collection)
  • args...: additional keyword arguments (unused)

Returns

  • PointMatrix with element type inferred from function evaluation
source
AdaptiveCrossApproximation.BEASTKernelMatrixType
BEASTKernelMatrix{T,NearBlockAssemblerType} <: AbstractKernelMatrix{T}

Kernel matrix wrapper for BEAST operator assembly.

Provides lazy matrix entry evaluation through a BEAST near-field block assembler, which computes matrix entries on demand from operator and basis function data.

Fields

  • nearassembler::NearBlockAssemblerType: BEAST assembler providing matrix entries

Type parameters

  • T: scalar element type returned by kernel evaluations
  • NearBlockAssemblerType: type of the underlying BEAST assembler
source
AdaptiveCrossApproximation.CombinedConvCritType
CombinedConvCrit

Composite convergence criterion combining multiple criteria. Converges when any constituent criterion is satisfied.

Fields

  • crits::Vector{ConvCrit}: Vector of convergence criteria to combine
source
AdaptiveCrossApproximation.CombinedPivStratType
CombinedPivStrat

Composite pivoting strategy that switches between multiple strategies based on convergence.

Combines multiple pivoting strategies with a combined convergence criterion, allowing the pivot selection method to change as different convergence criteria are satisfied. For example, can start with geometric pivoting and switch to value-based pivoting once a certain accuracy is reached.

Fields

  • strats::Vector{PivStrat}: Ordered list of pivoting strategies to use
source
AdaptiveCrossApproximation.ConvCritType
ConvCrit

Abstract base type for convergence criteria used by ACA and iACA compressors.

Notes

Concrete subtypes define how stopping decisions are made and are converted into stateful ConvCritFunctor objects during block compression.

See also

ConvCritFunctor, FNormEstimator, iFNormEstimator, FNormExtrapolator, RandomSampling

source
AdaptiveCrossApproximation.ConvCritFunctorType
ConvCritFunctor

Abstract base type for stateful convergence criterion functors.

Notes

Instances are called during ACA iterations and return (npivot, continue::Bool). Subtypes should implement reset! to reinitialize internal state for a new block.

See also

ConvCrit, reset!, normF!

source
AdaptiveCrossApproximation.ConvPivStratType
ConvPivStrat <: PivStrat

Abstract type for convergence-driven pivoting strategies.

These strategies adapt their behavior based on convergence information or use randomization to improve robustness.

Concrete Types

source
AdaptiveCrossApproximation.FNormExtrapolatorType
FNormExtrapolator{F} <: ConvCrit

Convergence criterion using polynomial extrapolation of pivot norms. Combines norm estimation with quadratic extrapolation to predict convergence.

Fields

  • estimator::Union{FNormEstimator{F},iFNormEstimator{F}}: Underlying norm estimator
source
AdaptiveCrossApproximation.FillDistanceType
FillDistance{D,F<:Real} <: GeoPivStrat

Geometric pivoting strategy based on fill distance minimization.

Selects pivots to minimize the fill distance, promoting well-distributed sampling in geometric space.

Fields

  • pos::Vector{SVector{D,F}}: Geometric positions of all points (D-dimensional)

Type Parameters

  • D: Spatial dimension
  • F: Floating point type for coordinates
source
AdaptiveCrossApproximation.GeoPivStratType
GeoPivStrat <: PivStrat

Abstract type for geometric/spatial pivoting strategies.

These strategies select pivots based on spatial/geometric properties rather than matrix values. Useful when geometric information about rows/columns is available.

Concrete Types

  • FillDistance: Maximizes minimum distance to already selected points
  • Leja2: Maximizes product of distances to selected points
source
AdaptiveCrossApproximation.HMatrixType
HMatrix{K,NearInteractionType,FarInteractionType}

Hierarchical matrix that stores near-field interactions explicitly and far-field interactions as low-rank block data.

Arguments

  • nearinteractions: block-sparse near-field contribution
  • farinteractions: collection of compressed far-field interaction blocks
  • dim::Tuple{Int,Int}: matrix dimensions (m, n)

Returns

An HMatrix linear map that supports matrix-vector products and conversion to a dense matrix via Matrix.

Notes

HMatrix is typically created through the high-level constructor HMatrix(operator, testspace, trialspace, tree; kwargs...) or H.assemble(...).

See also

HMatrix, H.assemble, farmatrix, nearmatrix

source
AdaptiveCrossApproximation.HMatrixMethod
HMatrix(operator, testspace, trialspace, tree; kwargs...)

Assemble a hierarchical matrix approximation of an operator on test and trial spaces.

Arguments

  • operator: bilinear form or kernel operator used for matrix entry evaluation
  • testspace: test space used for row indexing
  • trialspace: trial space used for column indexing
  • tree: hierarchical clustering/tree structure controlling block partitioning
  • tol: compression tolerance (default 1e-4)
  • compressor: ACA-style compressor, e.g. ACA(; tol=tol)
  • isnear: near-field predicate controlling admissibility
  • maxrank: maximum rank for far-field block compression
  • spaceordering: strategy for applying tree permutations to spaces
  • nearmatrixdata: optional data passed to near-field assembly
  • farmatrixdata: optional data passed to far-field assembly
  • scheduler: thread scheduler used for assembly

Returns

An HMatrix containing assembled near and far interactions.

Notes

Use this constructor as the main entry point when you already have a tree. For a convenience entry point, use H.assemble.

See also

H.assemble, ACA, iACA, farmatrix, nearmatrix

source
AdaptiveCrossApproximation.Leja2Type
Leja2{D,F<:Real} <: GeoPivStrat

Geometric pivoting strategy based on Leja points (product of distances).

A modified more efficient version of the fill distance approach. This leads to well-separated point sequences. These points have been introduced as modified leja points and will, therefore, be referred to as Leja2 points within this package.

Fields

  • pos::Vector{SVector{D,F}}: Geometric positions of all points (D-dimensional)

Type Parameters

  • D: Spatial dimension
  • F: Floating point type for coordinates
source
AdaptiveCrossApproximation.LowRankMatrixType
LowRankMatrix{K} <: LinearMaps.LinearMap{K}

Efficient storage and evaluation of low-rank matrix blocks as U*V.

Stores left and right factors separately, performing matrix-vector products as (U*V)*x = U*(V*x) to avoid forming the dense U*V product. Critical for efficient hierarchical matrix computations where blocks are stored in low-rank form.

Fields

  • U::Matrix{K}: Left factor, size (m, r) where r is the rank
  • V::Matrix{K}: Right factor, size (r, n)
  • z::Vector{K}: Temporary buffer for intermediate products V*x

Type parameters

  • K: scalar element type

Notes

Create via LowRankMatrix(U, V) which allocates the temporary buffer automatically. Block dimensions are (m, n) = (size(U, 1), size(V, 2)).

source
AdaptiveCrossApproximation.MaximumValueType
MaximumValue <: ValuePivStrat

Pivoting strategy that selects the index with maximum absolute value.

This is the standard pivoting strategy used in classical ACA algorithms also referred to as partial pivoting. At each iteration, it chooses the row or column with the largest absolute value among the unused indices, ensuring numerical stability and good approximation quality.

source
AdaptiveCrossApproximation.MimicryPivotingType
MimicryPivoting{D,F<:Real} <: GeoPivStrat

Geometric pivoting strategy that mimics point distribution of a fully pivoted ACA geometrically.

Selects pivots to reproduce the spatial distribution of a fully pivoted ACA. The strategy balances three objectives: geometric separation (Leja-like behavior), proximity to the reference distribution, and fill distance maximization. Particularly useful for H²–matrix compression where incomplete factorizations are sufficient.

Fields

  • refpos::Vector{SVector{D,F}}: Positions of test or expansion domain
  • pos::Vector{SVector{D,F}}: Positions from which to select pivots

Type Parameters

  • D: Spatial dimension
  • F: Floating point type for coordinates
source
AdaptiveCrossApproximation.PivStratType
PivStrat

Abstract base type for pivoting strategies in cross approximation algorithms.

Pivoting strategies determine which row or column to select at each iteration of the ACA algorithm. Concrete subtypes are typically stateless and callable, creating stateful PivStratFunctor instances when invoked with index information.

Subtypes by Category

  • GeoPivStrat: Geometric strategies (e.g., fill distance, Leja points)
  • ValuePivStrat: Value-based strategies (e.g., maximum absolute value)
  • ConvPivStrat: Convergence-driven strategies (e.g., random sampling)

Interface

Concrete subtypes should implement:

  • (strategy::MyPivStrat)(len::Int): Create functor for len indices
  • (strategy::MyPivStrat)(idcs::AbstractArray{Int}): Create functor for index array
source
AdaptiveCrossApproximation.PivStratFunctorType
PivStratFunctor

Abstract base type for stateful pivoting functors.

Functors maintain state during the pivot selection process (e.g., tracking which indices have been used). Created by calling a PivStrat instance with index information.

Interface

Concrete subtypes should implement:

  • (functor::MyPivStratFunctor)(): Select initial pivot (no data available)
  • (functor::MyPivStratFunctor)(rc::AbstractArray): Select next pivot based on row/column data
source
AdaptiveCrossApproximation.PointMatrixType
PointMatrix{T,FunctorType,PointCollectionType} <: AbstractKernelMatrix{T}

Kernel matrix evaluating entries from a kernel functor and point collections.

Stores point sample locations and a kernel function, computing matrix entries on demand by evaluating the kernel at test and trial point pairs. Useful for point cloud interactions or problems with explicit coordinate data.

Fields

  • functor::FunctorType: Kernel function or operator to evaluate
  • testpoints::PointCollectionType: Test space point locations
  • trialpoints::PointCollectionType: Trial space point locations

Type parameters

  • T: scalar element type returned by kernel evaluations
  • FunctorType: type of the kernel functor
  • PointCollectionType: type of point collections (e.g., Vector, AbstractArray)
source
AdaptiveCrossApproximation.RandomSamplingType
RandomSampling{F<:Real} <: ConvCrit

Convergence criterion based on random matrix entry sampling. Monitors approximation error at randomly sampled positions.

Fields

  • nsamples::Int: Number of random samples to take
  • factor::F: Factor for automatic sample count (nsamples = factor * (nrows + ncols))
  • tol::F: Convergence tolerance
source
AdaptiveCrossApproximation.RandomSamplingMethod
RandomSampling(; factor=1.0, nsamples=0, tol=1e-4)

Construct random sampling convergence criterion.

Arguments

  • factor::F: Multiplier for automatic sample count (default: 1.0)
  • nsamples::Int: Fixed sample count (default: 0, use factor instead)
  • tol::F: Convergence tolerance (default: 1e-4)
source
AdaptiveCrossApproximation.RandomSamplingPivotingType
RandomSamplingPivoting <: ConvPivStrat

Pivoting strategy that uses the error of random sampling from the convergence estimation.

Instead of selecting pivots based on maximum values or geometric properties, this strategy chooses pivots from randomly sampled indices used by a random sampling convergence criterion. Works in conjunction with the random-sampling convergence functor to provide statistically based pivot selection.

Fields

  • rc::Int: Index indicating which coordinate (row=1 or column=2) to select from
source
AdaptiveCrossApproximation.TreeMimicryPivotingType
TreeMimicryPivoting{D,T,TreeType} <: GeoPivStrat

Tree-aware mimicry pivoting strategy.

This strategy adapts the MimicryPivoting idea to a hierarchical tree of clusters. Instead of selecting individual points directly, it navigates the tree to pick clusters and then nodes within those clusters so that the selected pivots mimic a reference distribution at multiple scales.

Fields

  • refpos::Vector{SVector{D,T}}: Reference positions to mimic (e.g., parent pivots)
  • pos::Vector{SVector{D,T}}: Candidate point positions
  • tree::TreeType: Tree structure providing cluster centers, children and values

Type parameters

  • D: spatial dimension
  • T: numeric type for coordinates
  • TreeType: type of the tree adapter (must implement center, values, children, firstchild)
source
AdaptiveCrossApproximation.ValuePivStratType
ValuePivStrat <: PivStrat

Abstract type for value-based pivoting strategies.

These strategies select pivots based on matrix element values sampled during the ACA algorithm. Most common approach for general matrices.

Concrete Types

  • MaximumValue: Selects index with maximum absolute value (standard ACA)
  • RandomSampling: Random selection (for statistical approaches)
source
AdaptiveCrossApproximation.iACAType
iACA{RowPivType,ColPivType,ConvCritType}

Incomplete Adaptive Cross Approximation (iACA) compressor.

Unlike standard ACA, iACA computes only one side per iteration and relies on geometric pivoting strategies (for example mimicry or tree mimicry) to select pivots from spatial information. This reduces matrix entry evaluations in hierarchical matrix construction, where only selected row or column samples are required.

Fields

  • rowpivoting::RowPivType: Strategy for selecting row pivots (geometric)
  • columnpivoting::ColPivType: Strategy for selecting column pivots
  • convergence::ConvCritType: Convergence criterion
source
AdaptiveCrossApproximation.iACAMethod
iACA(tpos, spos)

Create a default incomplete ACA compressor for geometrically indexed row/column sets.

Arguments

  • tpos::Vector{SVector{D,F}}: geometric positions for test indices
  • spos::Vector{SVector{D,F}}: geometric positions for trial indices

Returns

An iACA instance using MaximumValue/MimicryPivoting with FNormExtrapolator.

See also

iACA, MimicryPivoting, FNormExtrapolator

source
AdaptiveCrossApproximation.iACAMethod
(iaca::iACA{ValuePivStratFunctor,GeoPivStratFunctor,ConvCritFunctor})(A, colbuffer, rowbuffer, maxrank, rows, cols, rowidcs)

Main computational routine for column matrix iACA (value-based row pivoting, geometric column pivoting). Performs incomplete ACA compression where columns are selected geometrically and rows by maximum value.

Arguments

  • A: Matrix to compress
  • colbuffer::AbstractArray{K}: Buffer for column data
  • rowbuffer::AbstractArray{K}: Buffer for row data
  • maxrank::Int: Maximum rank
  • rows::Vector{Int}: Row indices storage
  • cols::Vector{Int}: Column indices storage
  • rowidcs::Vector{Int}: Row index range

Returns

  • npivot::Int: Number of pivots computed
  • rows::Vector{Int}: Selected row indices (global)
  • cols::Vector{Int}: Selected column indices
source
AdaptiveCrossApproximation.iACAMethod
(iaca::iACA{GeoPivStratFunctor,ValuePivStratFunctor,ConvCritFunctor})(A, colbuffer, rowbuffer, maxrank, rows, cols, colidcs)

Main computational routine for row matrix iACA (geometric row pivoting, value-based column pivoting). Performs incomplete ACA compression where rows are selected geometrically and columns by maximum value.

Arguments

  • A: Matrix to compress
  • colbuffer::AbstractMatrix{K}: Buffer for column data
  • rowbuffer::AbstractMatrix{K}: Buffer for row data
  • maxrank::Int: Maximum rank
  • rows::Vector{Int}: Row indices storage
  • cols::Vector{Int}: Column indices storage
  • colidcs::Vector{Int}: Column index range

Returns

  • npivot::Int: Number of pivots computed
  • rows::Vector{Int}: Selected row indices
  • cols::Vector{Int}: Selected column indices (global)
source
AdaptiveCrossApproximation.iFNormEstimatorType
iFNormEstimator{F} <: ConvCrit

Frobenius norm-based convergence criterion for incomplete ACA (iACA). Uses moving average norm estimate for geometric pivoting scenarios.

Fields

  • tol::F: Relative tolerance threshold
source
AdaptiveCrossApproximation.acaMethod
aca(M; tol=1e-4, rowpivoting=MaximumValue(), columnpivoting=MaximumValue(),
    convergence=FNormEstimator(tol), maxrank=40, svdrecompress=false)

Compute adaptive cross approximation of matrix M returning low-rank factors.

High-level convenience function that automatically allocates buffers and returns U, V such that M ≈ U * V.

Arguments

  • M::AbstractMatrix{K}: Matrix to approximate

Keyword Arguments

  • tol::Real = 1e-4: Approximation tolerance
  • rowpivoting = MaximumValue(): Row pivot selection strategy
  • columnpivoting = MaximumValue(): Column pivot selection strategy
  • convergence = FNormEstimator(tol): Convergence criterion
  • maxrank::Int = 40: Maximum rank (hard limit)
  • svdrecompress::Bool = false: Apply SVD-based recompression to reduce rank further

Returns

  • U::Matrix{K}: Left factor, size (size(M,1), r) where r ≤ maxrank
  • V::Matrix{K}: Right factor, size (r, size(M,2))

Satisfies M ≈ U * V with norm(M - U*V) / norm(M) ≲ tol (if maxrank sufficient).

SVD Recompression

When svdrecompress=true, performs QR-SVD recompression: computes M ≈ U*V, then U = Q*R, R*V = Û*Σ*V̂ᵀ, truncates small singular values, and returns optimal rank factors at the cost of additional computation.

source
AdaptiveCrossApproximation.assemblefarsMethod
assemblefars(operator, testspace, trialspace, tree, ::PermuteSpaceInPlace; kwargs...)

Assemble far-field blocks with tree-aligned space reordering.

Compresses inadmissible cluster pairs using ACA-style algorithms, with both spaces permuted to align with the hierarchical tree structure. Produces a specialized storage format optimized for the reordered layout.

Arguments

  • operator: Operator/kernel for matrix entry evaluation
  • testspace: Test space (row basis) — will be reordered
  • trialspace: Trial space (column basis) — will be reordered
  • tree: Hierarchical tree structure
  • ::PermuteSpaceInPlace: Space ordering strategy marker
  • compressor = ACA(): Low-rank compressor algorithm
  • isnear: Admissibility predicate (default: isnear())
  • matrixdata: Assembly data passed to kernel matrix (optional)
  • maxrank = 50: Maximum rank for compressed blocks
  • scheduler: Thread scheduler (default: SerialScheduler())

Returns

  • Vector{VariableBlockCompressedRowStorage}: Level-wise compressed blocks in permuted layout

Notes

Called by HMatrix constructor when spaceordering=PermuteSpaceInPlace(). The spaces are reordered to match tree leaf ordering before compression, improving cache locality and reducing fragmentation. Returns specialized block storage optimized for permuted layout.

source
AdaptiveCrossApproximation.assemblefarsMethod
assemblefars(operator, testspace, trialspace, tree, ::PreserveSpaceOrder; kwargs...)

Assemble far-field blocks without reordering test and trial spaces.

Compresses inadmissible cluster pairs using ACA-style algorithms, maintaining the original space ordering. Produces a collection of level-by-level block-sparse matrices.

Arguments

  • operator: Operator/kernel for matrix entry evaluation
  • testspace: Test space (row basis/evaluation points)
  • trialspace: Trial space (column basis/evaluation points)
  • tree: Hierarchical tree structure controlling block layout
  • ::PreserveSpaceOrder: Space ordering strategy marker
  • compressor = ACA(): Low-rank compressor algorithm
  • isnear: Admissibility predicate (default: isnear())
  • matrixdata: Assembly data passed to kernel matrix (optional)
  • maxrank = 50: Maximum rank for compressed blocks
  • scheduler: Thread scheduler (default: SerialScheduler())

Returns

  • Vector{BlockSparseMatrix}: Collection of level-wise compressed blocks

Notes

Called by HMatrix constructor when spaceordering=PreserveSpaceOrder(). Iterates over tree levels, compressing each level's far-field interactions independently. Returns one sparse matrix per level, indexed by tree depth.

source
AdaptiveCrossApproximation.assemblenearsMethod
assemblenears(operator, testspace, trialspace, tree, ::PermuteSpaceInPlace; kwargs...)

Assemble near-field blocks with tree-aligned space reordering.

Computes dense matrix blocks for admissible cluster pairs, with both spaces permuted to align with the hierarchical tree structure. Produces a specialized block storage format optimized for the reordered layout.

Arguments

  • operator: Operator/kernel for matrix entry evaluation
  • testspace: Test space (row basis) — will be reordered
  • trialspace: Trial space (column basis) — will be reordered
  • tree: Hierarchical tree structure
  • ::PermuteSpaceInPlace: Space ordering strategy marker
  • isnear: Admissibility predicate (default: isnear())
  • scheduler: Thread scheduler (default: SerialScheduler())
  • matrixdata: Assembly data passed to kernel matrix (optional)

Returns

  • VariableBlockCompressedRowStorage: Permuted block-sparse near-field storage

Notes

Called by HMatrix constructor when spaceordering=PermuteSpaceInPlace(). The spaces are reordered to match tree leaf ordering before block assembly, improving cache locality and reducing block fragmentation.

source
AdaptiveCrossApproximation.assemblenearsMethod
assemblenears(operator, testspace, trialspace, tree, ::PreserveSpaceOrder; kwargs...)

Assemble near-field blocks without reordering test and trial spaces.

Computes dense (non-low-rank) matrix blocks for all admissible cluster pairs, retaining the original ordering of test and trial spaces. The resulting block-sparse matrix stores these near-field interactions.

Arguments

  • operator: Operator/kernel for matrix entry evaluation
  • testspace: Test space (row basis/evaluation points)
  • trialspace: Trial space (column basis/evaluation points)
  • tree: Hierarchical tree structure
  • ::PreserveSpaceOrder: Space ordering strategy marker
  • isnear: Admissibility predicate (default: isnear())
  • scheduler: Thread scheduler (default: SerialScheduler())
  • matrixdata: Assembly data passed to kernel matrix (optional)

Returns

  • BlockSparseMatrix: Block-sparse storage of near-field blocks

Notes

Called by HMatrix constructor when spaceordering=PreserveSpaceOrder(). Each near-field block is evaluated directly via the kernel matrix, avoiding compression.

source
AdaptiveCrossApproximation.farinteractionsMethod
farinteractions(tree; args...)

Extract far-field (inadmissible) block pairs requiring low-rank compression.

Should be implemented by tree backends (e.g., H2Trees) to return index ranges for pairs of clusters that do NOT satisfy the admissibility criterion. These blocks are compressed via ACA or similar algorithms.

Arguments

  • tree: Hierarchical tree object implementing the tree interface
  • args...: Additional keyword arguments passed by the assembly routine (e.g., isnear)

Returns

  • (values, farptr, farvalues): Triple of row indices, far-field pointers, and column ranges

Notes

Implemented by tree backend extensions. The pointer array farptr enables efficient lookup of far-field couples for each node, reducing access overhead during level-by-level assembly.

source
AdaptiveCrossApproximation.farinteractions_consecutiveMethod
farinteractions_consecutive(tree; args...)

Extract far-field blocks with consecutive storage layout.

Variant used when spaces are permuted to align with tree ordering. Returns far-field block structure compatible with permuted space layout.

Arguments

  • tree: Hierarchical tree object
  • args...: Additional keyword arguments (e.g., isnear)

Returns

  • (values, farptr, farvalues): Index structure for consecutive far-field layout

Notes

Implemented by tree backend extensions. Must account for space permutation applied by PermuteSpaceInPlace() ordering strategy.

source
AdaptiveCrossApproximation.farmatrixMethod
farmatrix(hmat::HMatrix)

Extract the far-field (low-rank compressed) contributions from a hierarchical matrix.

Returns a new HMatrix containing only the low-rank far-field blocks, with all near-field blocks removed. Useful for analyzing the rank structure or applying operations that exploit low-rank properties.

Arguments

  • hmat::HMatrix: Source hierarchical matrix

Returns

  • HMatrix: New matrix with identical far-field interactions but empty near-field
source
AdaptiveCrossApproximation.isnearFunction
isnear(η::Real=Float64(1.0))

Create an admissibility predicate for near-field block detection.

Clusters with scaled geometric distance less than η are considered admissible for near-field assembly (direct evaluation). The scaling factor η controls the geometric separation required for far-field low-rank approximation.

Arguments

  • η::Real = 1.0: Admissibility parameter (dimensionless geometric distance threshold)

Returns

  • IsNearFunctor: Admissibility predicate functor

Notes

For well-separated clusters (distance > η × cluster diameter), the interaction is computed via ACA compression; otherwise, via direct near-field assembly. Typical values: η ≈ 1.0 to 3.0 depending on required accuracy.

source
AdaptiveCrossApproximation.nearinteractionsMethod
nearinteractions(tree; args...)

Extract near-field (admissible) block pairs from a hierarchical tree structure.

Should be implemented by tree backends (e.g., H2Trees) to return index ranges for pairs of clusters that satisfy the admissibility criterion. This is the primary method when space ordering is not modified.

Arguments

  • tree: Hierarchical tree object implementing the tree interface
  • args...: Additional keyword arguments passed by the assembly routine (e.g., isnear)

Returns

  • (values, nearvalues): Tuple of index vectors for row/column near-field blocks

Notes

Implemented by tree backend extensions (e.g., ACAH2Trees for H2Trees). Should return rows and corresponding near-field column indices that are admissible for direct (non-low-rank) assembly.

source
AdaptiveCrossApproximation.nearinteractions_consecutiveMethod
nearinteractions_consecutive(tree; args...)

Extract near-field blocks with consecutive storage layout.

Variant used when spaces are permuted to align with tree ordering. Returns near-field block structure compatible with permuted space layout.

Arguments

  • tree: Hierarchical tree object
  • args...: Additional keyword arguments (e.g., isnear)

Returns

  • (values, nearvalues): Index vectors for consecutive near-field layout

Notes

Implemented by tree backend extensions. Must account for space permutation applied by PermuteSpaceInPlace() ordering strategy.

source
AdaptiveCrossApproximation.nearmatrixMethod
nearmatrix(hmat::HMatrix)

Extract the near-field (dense block) contributions from a hierarchical matrix.

Returns only the block-sparse near-field interactions, removing all low-rank far-field blocks. Useful for analyzing near-field structure or applying operations specific to dense blocks.

Arguments

  • hmat::HMatrix: Source hierarchical matrix

Returns

  • BlockSparseMatrix: Near-field interactions as block-sparse matrix
source
AdaptiveCrossApproximation.reset!Method
reset!(convcrit::ConvCritFunctor)

Reset a convergence functor before starting compression of a new block.

Notes

Concrete subtypes should overload this method. The default fallback throws ArgumentError.

See also

ConvCritFunctor

source
AdaptiveCrossApproximation.storageMethod
storage(hmat::HMatrix)

Analyze and report memory storage requirements of a hierarchical matrix.

Prints detailed statistics (to stdout) comparing the actual storage needed for the hierarchical representation versus dense storage, including compression ratio. Also reports total memory footprint via Julia's summarysize.

Arguments

  • hmat::HMatrix: Hierarchical matrix to analyze

Returns

  • Float64: Total storage in GB used by hierarchical blocks

Output

Prints to stdout:

  • storage: Total block storage in GB
  • summary size: Total memory footprint including Julia object overhead (GB)
  • compression ratio: Ratio of hierarchical storage to dense storage
source
AdaptiveCrossApproximation.update_refcentroid!Method
update_refcentroid!(functor::PivStratFunctor, refidcs)

Update the stored reference-domain centroid in pivoting functors that expose refcentroid. For those functors, the reference positions are taken from functor.pivoting.refpos. For functors without refcentroid, this is a no-op.

source