API Reference
High-Level API
AdaptiveCrossApproximation.H.assemble — Function
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 evaluationtestspace: Test space (row basis/points)trialspace: Trial space (column basis/points)tree: Hierarchical tree structure (auto-generated if not provided)kwargs...: Additional parameters passed toHMatrix:tol: Convergence tolerance (default1e-4)maxrank: Maximum rank for compression (default40)compressor: ACA or iACA instance (defaultACA(tol=tol))isnear: Admissibility predicate (defaultisnear())spaceordering: Space ordering strategy (defaultPermuteSpaceInPlace())
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 productSee also
AdaptiveCrossApproximation.AdaptiveCrossApproximation — Module
AdaptiveCrossApproximationAdaptive 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 variantiACA: Incomplete variant for geometric pivoting and hierarchical matricesaca: Convenience function for matrix compression
Hierarchical matrices:
HMatrix: Hierarchical matrix typeAdaptiveCrossApproximation.H.assemble: Automatic assembly with tree-based blockingnearmatrix,farmatrix: Extract system components
Pivoting strategies:
MaximumValue: Partial pivoting (largest absolute value)Leja2: Modified Leja point selectionFillDistance: Geometric fill distance minimizationMimicryPivoting: Distributional mimicry without full matrix accessTreeMimicryPivoting: Tree-aware variant for H²-matrices
Convergence criteria:
FNormEstimator,iFNormEstimator: Frobenius norm-based stoppingFNormExtrapolator: Extrapolation-enhanced convergence detectionRandomSampling: Random sampling-based convergence (convergence module)
Kernel matrices:
AbstractKernelMatrix: Interface for lazy matrix entry evaluationPointMatrix: Point cloud kernel matricesBEASTKernelMatrix: BEAST boundary integral operator 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 productSee also
- BEAST.jl for boundary integral operators
- H2Trees.jl for hierarchical clustering
- BlockSparseMatrices.jl for sparse block storage
AdaptiveCrossApproximation.ACA — Type
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 pivotscolumnpivoting::ColPivType: Strategy for selecting column pivotsconvergence::ConvCritType: Convergence criterion to stop iterations
AdaptiveCrossApproximation.ACA — Method
(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 supportnextrc!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 indicescols::Vector{Int}: Storage for selected column indicesrowidcs::Vector{Int}: Global row indices of the block to compresscolidcs::Vector{Int}: Global column indices of the block to compressmaxrank::Int: Maximum number of pivots (hard limit on rank)
Returns
npivot::Int: Number of pivots computed (≤ maxrank). The approximation isA[rowidcs, colidcs] ≈ colbuffer[:, 1:npivot] * rowbuffer[1:npivot, :]
AdaptiveCrossApproximation.AbstractKernelMatrix — Type
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...)
AdaptiveCrossApproximation.AbstractKernelMatrix — Method
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 operatortestspace::AbstractVector: test point locations (indexed collection)trialspace::AbstractVector: trial point locations (indexed collection)args...: additional keyword arguments (unused for point-based evaluation)
Returns
PointMatrixwith element type inferred fromoperator
AdaptiveCrossApproximation.AbstractKernelMatrix — Method
AbstractKernelMatrix(operator, testspace, trialspace; args...)Construct a concrete kernel matrix wrapper from operator and space data.
Arguments
operator: operator or kernel definitiontestspace: space for row evaluation points or basis datatrialspace: space for column evaluation points or basis dataargs...: 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!
AdaptiveCrossApproximation.AbstractKernelMatrix — Method
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 signatureoperator(testpoint, trialpoint)testspace::AbstractVector: test point locations (indexed collection)trialspace::AbstractVector: trial point locations (indexed collection)args...: additional keyword arguments (unused)
Returns
PointMatrixwith element type inferred from function evaluation
AdaptiveCrossApproximation.BEASTKernelMatrix — Type
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 evaluationsNearBlockAssemblerType: type of the underlying BEAST assembler
AdaptiveCrossApproximation.CombinedConvCrit — Type
CombinedConvCritComposite convergence criterion combining multiple criteria. Converges when any constituent criterion is satisfied.
Fields
crits::Vector{ConvCrit}: Vector of convergence criteria to combine
AdaptiveCrossApproximation.CombinedPivStrat — Type
CombinedPivStratComposite 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
AdaptiveCrossApproximation.ConvCrit — Type
ConvCritAbstract 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
AdaptiveCrossApproximation.ConvCritFunctor — Type
ConvCritFunctorAbstract 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!
AdaptiveCrossApproximation.ConvPivStrat — Type
ConvPivStrat <: PivStratAbstract type for convergence-driven pivoting strategies.
These strategies adapt their behavior based on convergence information or use randomization to improve robustness.
Concrete Types
RandomSampling: Random pivot selection for convergence estimation
AdaptiveCrossApproximation.ConvPivStratFunctor — Type
ConvPivStratFunctor <: PivStratFunctorAbstract type for stateful convergence-driven pivoting functors.
AdaptiveCrossApproximation.FNormEstimator — Type
FNormEstimator{F} <: ConvCritFrobenius norm-based convergence criterion for standard ACA.
Fields
tol::F: Relative tolerance threshold
AdaptiveCrossApproximation.FNormExtrapolator — Type
FNormExtrapolator{F} <: ConvCritConvergence 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
AdaptiveCrossApproximation.FNormExtrapolator — Method
FNormExtrapolator(tol::F)Construct extrapolator with Frobenius norm estimator.
Arguments
tol::F: Convergence tolerance
AdaptiveCrossApproximation.FillDistance — Type
FillDistance{D,F<:Real} <: GeoPivStratGeometric 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 dimensionF: Floating point type for coordinates
AdaptiveCrossApproximation.GeoPivStrat — Type
GeoPivStrat <: PivStratAbstract 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 pointsLeja2: Maximizes product of distances to selected points
AdaptiveCrossApproximation.GeoPivStratFunctor — Type
GeoPivStratFunctor <: PivStratFunctorAbstract type for stateful geometric pivoting functors.
AdaptiveCrossApproximation.HMatrix — Type
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 contributionfarinteractions: collection of compressed far-field interaction blocksdim::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
AdaptiveCrossApproximation.HMatrix — Method
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 evaluationtestspace: test space used for row indexingtrialspace: trial space used for column indexingtree: hierarchical clustering/tree structure controlling block partitioningtol: compression tolerance (default1e-4)compressor: ACA-style compressor, e.g.ACA(; tol=tol)isnear: near-field predicate controlling admissibilitymaxrank: maximum rank for far-field block compressionspaceordering: strategy for applying tree permutations to spacesnearmatrixdata: optional data passed to near-field assemblyfarmatrixdata: optional data passed to far-field assemblyscheduler: 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
AdaptiveCrossApproximation.Leja2 — Type
Leja2{D,F<:Real} <: GeoPivStratGeometric 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 dimensionF: Floating point type for coordinates
AdaptiveCrossApproximation.LowRankMatrix — Type
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)whereris the rankV::Matrix{K}: Right factor, size(r, n)z::Vector{K}: Temporary buffer for intermediate productsV*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)).
AdaptiveCrossApproximation.MaximumValue — Type
MaximumValue <: ValuePivStratPivoting 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.
AdaptiveCrossApproximation.MimicryPivoting — Type
MimicryPivoting{D,F<:Real} <: GeoPivStratGeometric 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 domainpos::Vector{SVector{D,F}}: Positions from which to select pivots
Type Parameters
D: Spatial dimensionF: Floating point type for coordinates
AdaptiveCrossApproximation.PivStrat — Type
PivStratAbstract 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 forlenindices(strategy::MyPivStrat)(idcs::AbstractArray{Int}): Create functor for index array
AdaptiveCrossApproximation.PivStratFunctor — Type
PivStratFunctorAbstract 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
AdaptiveCrossApproximation.PointMatrix — Type
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 evaluatetestpoints::PointCollectionType: Test space point locationstrialpoints::PointCollectionType: Trial space point locations
Type parameters
T: scalar element type returned by kernel evaluationsFunctorType: type of the kernel functorPointCollectionType: type of point collections (e.g., Vector, AbstractArray)
AdaptiveCrossApproximation.RandomSampling — Type
RandomSampling{F<:Real} <: ConvCritConvergence criterion based on random matrix entry sampling. Monitors approximation error at randomly sampled positions.
Fields
nsamples::Int: Number of random samples to takefactor::F: Factor for automatic sample count (nsamples = factor * (nrows + ncols))tol::F: Convergence tolerance
AdaptiveCrossApproximation.RandomSampling — Method
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)
AdaptiveCrossApproximation.RandomSamplingPivoting — Type
RandomSamplingPivoting <: ConvPivStratPivoting 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
AdaptiveCrossApproximation.TreeMimicryPivoting — Type
TreeMimicryPivoting{D,T,TreeType} <: GeoPivStratTree-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 positionstree::TreeType: Tree structure providing cluster centers, children and values
Type parameters
D: spatial dimensionT: numeric type for coordinatesTreeType: type of the tree adapter (must implementcenter,values,children,firstchild)
AdaptiveCrossApproximation.ValuePivStrat — Type
ValuePivStrat <: PivStratAbstract 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)
AdaptiveCrossApproximation.ValuePivStratFunctor — Type
ValuePivStratFunctor <: PivStratFunctorAbstract type for stateful value-based pivoting functors.
AdaptiveCrossApproximation.iACA — Type
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 pivotsconvergence::ConvCritType: Convergence criterion
AdaptiveCrossApproximation.iACA — Method
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 indicesspos::Vector{SVector{D,F}}: geometric positions for trial indices
Returns
An iACA instance using MaximumValue/MimicryPivoting with FNormExtrapolator.
See also
iACA, MimicryPivoting, FNormExtrapolator
AdaptiveCrossApproximation.iACA — Method
(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 compresscolbuffer::AbstractArray{K}: Buffer for column datarowbuffer::AbstractArray{K}: Buffer for row datamaxrank::Int: Maximum rankrows::Vector{Int}: Row indices storagecols::Vector{Int}: Column indices storagerowidcs::Vector{Int}: Row index range
Returns
npivot::Int: Number of pivots computedrows::Vector{Int}: Selected row indices (global)cols::Vector{Int}: Selected column indices
AdaptiveCrossApproximation.iACA — Method
(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 compresscolbuffer::AbstractMatrix{K}: Buffer for column datarowbuffer::AbstractMatrix{K}: Buffer for row datamaxrank::Int: Maximum rankrows::Vector{Int}: Row indices storagecols::Vector{Int}: Column indices storagecolidcs::Vector{Int}: Column index range
Returns
npivot::Int: Number of pivots computedrows::Vector{Int}: Selected row indicescols::Vector{Int}: Selected column indices (global)
AdaptiveCrossApproximation.iFNormEstimator — Type
iFNormEstimator{F} <: ConvCritFrobenius norm-based convergence criterion for incomplete ACA (iACA). Uses moving average norm estimate for geometric pivoting scenarios.
Fields
tol::F: Relative tolerance threshold
AdaptiveCrossApproximation.aca — Method
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 tolerancerowpivoting = MaximumValue(): Row pivot selection strategycolumnpivoting = MaximumValue(): Column pivot selection strategyconvergence = FNormEstimator(tol): Convergence criterionmaxrank::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)wherer ≤ maxrankV::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.
AdaptiveCrossApproximation.assemblefars — Method
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 evaluationtestspace: Test space (row basis) — will be reorderedtrialspace: Trial space (column basis) — will be reorderedtree: Hierarchical tree structure::PermuteSpaceInPlace: Space ordering strategy markercompressor = ACA(): Low-rank compressor algorithmisnear: Admissibility predicate (default:isnear())matrixdata: Assembly data passed to kernel matrix (optional)maxrank = 50: Maximum rank for compressed blocksscheduler: 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.
AdaptiveCrossApproximation.assemblefars — Method
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 evaluationtestspace: Test space (row basis/evaluation points)trialspace: Trial space (column basis/evaluation points)tree: Hierarchical tree structure controlling block layout::PreserveSpaceOrder: Space ordering strategy markercompressor = ACA(): Low-rank compressor algorithmisnear: Admissibility predicate (default:isnear())matrixdata: Assembly data passed to kernel matrix (optional)maxrank = 50: Maximum rank for compressed blocksscheduler: 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.
AdaptiveCrossApproximation.assemblenears — Method
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 evaluationtestspace: Test space (row basis) — will be reorderedtrialspace: Trial space (column basis) — will be reorderedtree: Hierarchical tree structure::PermuteSpaceInPlace: Space ordering strategy markerisnear: 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.
AdaptiveCrossApproximation.assemblenears — Method
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 evaluationtestspace: Test space (row basis/evaluation points)trialspace: Trial space (column basis/evaluation points)tree: Hierarchical tree structure::PreserveSpaceOrder: Space ordering strategy markerisnear: 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.
AdaptiveCrossApproximation.farinteractions — Method
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 interfaceargs...: 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.
AdaptiveCrossApproximation.farinteractions_consecutive — Method
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 objectargs...: 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.
AdaptiveCrossApproximation.farmatrix — Method
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
AdaptiveCrossApproximation.isnear — Function
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.
AdaptiveCrossApproximation.nearinteractions — Method
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 interfaceargs...: 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.
AdaptiveCrossApproximation.nearinteractions_consecutive — Method
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 objectargs...: 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.
AdaptiveCrossApproximation.nearmatrix — Method
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
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
AdaptiveCrossApproximation.storage — Method
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 GBsummary size: Total memory footprint including Julia object overhead (GB)compression ratio: Ratio of hierarchical storage to dense storage
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.