API - CertifScripts

Pseudospectra certification module.

BallArithmetic.CertifScripts._evaluate_sample_ogita_bigfloatMethod
_evaluate_sample_ogita_bigfloat(T_matrix::BallMatrix, z::Number, idx::Int;
                                 max_iterations::Int=4,
                                 target_precision::Int=256,
                                 distance_threshold::Real=1e-4,
                                 use_cache::Bool=true)

Evaluate the smallest singular value at z using Ogita refinement with caching.

This function uses a worker-local cache to store the last refined BigFloat SVD. For nearby points (distance < distance_threshold), it uses the cached SVD as a starting point for Ogita refinement, which requires fewer iterations than starting from Float64.

Cache Strategy

  • If cache hit: use cached BigFloat SVD → 1-2 Ogita iterations
  • If cache miss: use Float64 SVD → 3-4 Ogita iterations

Arguments

  • T_matrix: The Schur matrix T (BallMatrix)
  • z: Sample point on the circle
  • idx: Sample index for logging
  • max_iterations: Number of Ogita iterations for cache miss (3-4 for 256-bit)
  • target_precision: Target precision in bits
  • distance_threshold: Maximum distance for cache reuse
  • use_cache: Whether to use caching (default: true)
source
BallArithmetic.CertifScripts._evaluate_sample_with_ogita_cacheMethod
_evaluate_sample_with_ogita_cache(T, z, idx; ogita_distance_threshold, ogita_quality_threshold)

Evaluate sample with Ogita optimization: try to refine from cached SVD if available.

If the cached SVD is from a nearby point (distance < ogitadistancethreshold), attempt Ogita refinement. If the refined SVD has acceptable quality (residual < ogitaqualitythreshold), use it; otherwise fall back to full SVD.

source
BallArithmetic.CertifScripts.adaptive_arcs!Method
adaptive_arcs!(arcs, cache, pending, η; kwargs...)

Drive the adaptive refinement routine. When job channels are provided the refinement uses asynchronous workers; otherwise the evaluation is carried out serially using the supplied evaluator.

source
BallArithmetic.CertifScripts.compute_schur_and_errorMethod
compute_schur_and_error(A; polynomial = nothing)

Compute the Schur decomposition of A and certified bounds for the orthogonality defect, the reconstruction error, and the norms of Z and Z⁻¹. When polynomial is provided (as coefficients in ascending order), additional bounds are computed for p(A) and p(T).

Supports both Float64 and BigFloat precision based on the element type of A. For BigFloat, the Schur decomposition is computed in Float64 and then refined to higher precision using iterative refinement.

source
BallArithmetic.CertifScripts.dowork_ogitaMethod
dowork_ogita(jobs, results; ogita_distance_threshold=1e-4, ogita_quality_threshold=1e-10)

Process tasks with Ogita optimization enabled. Similar to dowork but tries to use cached SVD from previous evaluations to speed up nearby points.

The worker maintains a cache of the last computed SVD. When a new point z arrives, if it's within ogita_distance_threshold of the cached point, Ogita refinement is attempted. If the refined SVD has acceptable quality (relative residual < ogita_quality_threshold), it's used; otherwise a full SVD is computed.

This is beneficial when the adaptive bisection sends consecutive jobs for nearby points, which happens naturally as arcs get smaller during refinement.

source
BallArithmetic.CertifScripts.dowork_ogita_bigfloatMethod
dowork_ogita_bigfloat(jobs, results; target_precision=256, max_ogita_iterations=3,
                      distance_threshold=1e-4)

Process tasks with BigFloat precision using Ogita SVD refinement with caching.

For each job (id, z), this worker:

  1. Checks if a cached BigFloat SVD exists from a nearby point
  2. If cache hit: refines from cached SVD (1-2 iterations)
  3. If cache miss: computes Float64 SVD and refines to BigFloat (3-4 iterations)
  4. Certifies with Miyajima bounds
  5. Caches the result for future reuse

This is the distributed equivalent of run_certification_ogita for parallel execution. The Schur matrix must be a BigFloat BallMatrix registered with set_schur_matrix!.

Performance

  • Cache miss: Ogita from Float64 requires ~4 iterations (10^-16 → 10^-64)
  • Cache hit: Ogita from BigFloat requires ~1-2 iterations (already at high precision)
  • For adaptive bisection with nearby points, expect 50-90% cache hit rate
source
BallArithmetic.CertifScripts.poly_from_rootsMethod
poly_from_roots(roots::AbstractVector)

Given a list of roots r₁, r₂, …, rₙ, returns the coefficients [a₀, a₁, …, aₙ] of the monic polynomial p(x) = (x - r₁)(x - r₂)…(x - rₙ) so that p(x) = a₀ + a₁x + a₂x² + … + aₙ*xⁿ.

source
BallArithmetic.CertifScripts.run_certificationMethod
run_certification(A, circle; schur_data = nothing, polynomial = nothing, kwargs...)

Run the adaptive certification routine on circle using a serial evaluator.

Arguments

  • A: matrix to certify. Converted to BallMatrix when required.
  • circle: CertificationCircle describing the contour used for the adaptive refinement.

Keyword Arguments

  • schur_data = nothing: pre-computed Schur data as the 5-tuple (S, errF, errT, norm_Z, norm_Z_inv) returned by compute_schur_and_error. When provided, the expensive Schur computation is skipped entirely. This is useful for reusing the same Schur decomposition across multiple circles or when the default Float64-seeded refinement fails (e.g., matrices with eigenvalues below ~10⁻¹⁶).
  • polynomial = nothing: optional coefficients (ascending order) describing a polynomial p. When provided the certification is carried out on p(T) and the returned error corresponds to the reconstruction error of p(A).
  • η = 0.5: admissible threshold for the adaptive refinement. Must lie in the open unit interval.
  • check_interval = 100: number of processed arcs between progress reports and consistency checks.
  • log_io = stdout: destination IO for log messages.
  • Cbound = 1.0: constant used by bound_res_original when lifting resolvent bounds back to the original matrix.

The return value is a named tuple containing the computed Schur form, the accumulated certification log, and the resolvent bounds for both the Schur factor and the original matrix.

source
BallArithmetic.CertifScripts.run_certification_ogitaMethod
run_certification_ogita(A, circle; target_precision=256, kwargs...)

Optimized BigFloat certification using Ogita SVD refinement.

This function is specifically designed for BigFloat precision certification where computing fresh BigFloat SVDs at each sample point is expensive. Instead, it:

  1. Computes Schur decomposition with BigFloat refinement
  2. At each sample point z, computes Float64 SVD of (T - zI)
  3. Refines the Float64 SVD to BigFloat using Ogita's algorithm (3-4 iterations)

Due to quadratic convergence, 4 Ogita iterations from Float64 (~10^-16 error) achieve ~10^-64 error, saturating 256-bit precision.

Performance

Typically 10-100x faster than computing fresh BigFloat SVDs at each point.

Arguments

  • A: matrix to certify (will be converted to BigFloat BallMatrix)
  • circle: CertificationCircle describing the contour
  • target_precision::Int=256: precision in bits for BigFloat
  • max_ogita_iterations::Int=3: Ogita iterations (3 is optimal for 256-bit precision)
  • Other kwargs passed to standard certification
source
BallArithmetic.CertifScripts.run_certification_parametricMethod
run_certification_parametric(A, circle; k=nothing, config=config_v2(), kwargs...)

Run the adaptive certification routine using the parametric Sylvester-based certifier.

This method uses block-diagonalization via Sylvester equation to certify the resolvent norm, which can be more efficient than full SVD for large matrices.

Arguments

  • A: matrix to certify. Converted to BallMatrix when required.
  • circle: CertificationCircle describing the contour.

Keyword Arguments

  • k::Union{Nothing, Int}=nothing: Split index for the Schur block decomposition. If nothing, automatically selects k ≈ n/4.
  • config::ResolventBoundConfig=config_v2(): Configuration specifying estimators. Options: config_v1(), config_v2(), config_v2p5(), config_v3().
  • polynomial = nothing: optional coefficients for polynomial certification.
  • η = 0.5: admissible threshold for adaptive refinement.
  • check_interval = 100: number of arcs between progress reports.
  • log_io = stdout: destination for log messages.
  • Cbound = 1.0: constant for resolvent bound lifting.

Example

using BallArithmetic
A = randn(ComplexF64, 50, 50)
circle = CertificationCircle(0.0, 1.5; samples=64)

# Use V2 configuration (default)
result = run_certification_parametric(A, circle)

# Use V3 configuration with Neumann bounds
result = run_certification_parametric(A, circle; config=config_v3())

# Specify custom split
result = run_certification_parametric(A, circle; k=10)
source
BallArithmetic.CertifScripts.schur_to_original_resolventMethod
schur_to_original_resolvent(resolvent_schur, ϵ; Cbound=1.0)

Transform a Schur resolvent bound to the original matrix using the Schur perturbation theorem. Returns Inf when the bound is not computable (ε·resolvent too large for the perturbation estimate).

source
BallArithmetic.CertifScripts.set_parametric_config!Method
set_parametric_config!(precomp, R, config; k=0)

Set the parametric certifier configuration.

Arguments

  • precomp: Precomputed Sylvester quantities from sylvester_resolvent_precompute
  • R: Sylvester residual matrix
  • config: ResolventBoundConfig specifying estimators
  • k: Split index (for reference)
source