Documentation for Polfed.jl
This page documents POLFED solver entry points and the public Hamiltonian constructors used throughout the tutorials.
Main Solver Entry Point
Use polfed with:
- matrix input:
polfed(mat, x0, howmany, target; ...) - mapping callback:
polfed(f!, x0, howmany, target; ...)
Polfed.PolfedCore.polfed — Functionpolfed(mat, x0, howmany, target; kwargs...)
polfed(f!, x0, howmany, target; kwargs...)Solve an eigenvalue problem using the Polynomial Filtering Eigenvalue Decomposition (POLFED) algorithm with a given matrix input.
This method employs Chebyshev polynomial spectral filtering to efficiently extract eigenvalues and eigenvectors within a specified spectral region of a dense or sparse matrix. When optimize_mapping = true, the operator is internally rescaled and optimized to reduce memory access, particularly effective for Hamiltonians with a small number of unique off-diagonal elements.
Arguments
mat::AbstractMatrix: The Hamiltonian or operator matrix (dense or sparse).f!::Function: In-place linear operator callback with signaturef!(Y, X).x0::AbstractVecOrMat: The initial vector or block of vectors used to start the iteration. If a vector is provided, a standard Lanczos factorization is performed; if a matrix is provided, the Block Lanczos variant is used.howmany::Integer: Number of eigenvalues to compute.target: Spectral region to target (in unscaled units). Accepts a numeric value (treated as unrescaled) or a spec such as:maxdos,:middle,(:offset, frac),(:unrescaled, value), or(:rescaled, value).
Keyword Arguments
produce_report::Bool = PolfedDefaults.produce_reportWhether to return detailed diagnostic information about the run. If enabled, aReportstruct is returned alongside eigenvalues and eigenvectors.mapping::MappingConfig = MappingConfig()Configuration for mapping, rescaling, parallelization, and optional optimized kernels.transform::TransformConfig = TransformConfig()Configuration for polynomial coefficients, normalization, interval, and order.fact::FactorizationConfig = FactorizationConfig()Configuration for the underlying Lanczos/Krylov factorization (seeFactorizationConfig). Includes settings for numerical tolerance, reorthogonalization strategy, and eigenvalue convergence criteria.dos::DoSConfig = DoSConfig()Configuration for the stochastic density-of-states (DoS) estimation (seeDoSConfig). Parameters include the number of random vectors, number of Chebyshev moments, and statistical averaging options.
Returns
(vals, vecs)— Eigenvalues and eigenvectors.(vals, vecs, report)— Ifproduce_report = true, also returns aReportcontaining:- Spectral transform diagnostics (mapping, scaling, polynomial order)
- Factorization statistics (iterations, reorthogonalization)
- Benchmark metrics (wall time, CPU time, and device information)
Notes
- If
x0is not normalized (vector input) or orthonormal (matrix input), POLFED emits a warning and uses a normalized/orthonormalized copy without mutating the caller's array. - The polynomial order is chosen automatically unless explicitly set via
TransformConfig.order. - When
mapping.optimize_mapping = true, a rescaled operator and optimized Clenshaw recurrence are used to reduce memory bandwidth usage, especially beneficial for structured Hamiltonians. - GPU execution is automatically enabled when
is_gpu_array(x0).
See also
MappingConfig, TransformConfig, FactorizationConfig, DoSConfig, display_report
polfed(f!::Function, x0::AbstractVecOrMat, howmany::Integer, target; kwargs...)Operator-callback overload of polfed.
Use this when the matrix is represented implicitly through f!(Y, X).
Model Constructors
Public Hamiltonian builders live under the common Polfed.Models namespace. Their full signatures, keyword arguments, return types, and conventions are documented in Models.
Mapping-Optimization Helpers
These helpers are available when constructing high-performance custom mapping workflows:
Polfed.PolfedCore.optimized_mapping! — Functionoptimized_mapping!(diagonals::AbstractVector, offdiagonals::Offdiagonals, parallel_strategy::Parallelization) -> FunctionBuild an optimized in-place mapping callback f!(Y, X) from packed diagonal and off-diagonal data.
Polfed.PolfedCore.optimized_clenshaw_recurrence_relation! — Functionoptimized_clenshaw_recurrence_relation!(diagonals, offdiagonals, parallel_strategy) -> FunctionBuild an optimized Clenshaw recurrence callback compatible with POLFED's internal recurrence signature.
Polfed.PolfedCore.optimized_clenshaw_final_sum! — Functionoptimized_clenshaw_final_sum!(diagonals, offdiagonals, parallel_strategy) -> FunctionBuild an optimized callback for the terminal Clenshaw combination step.
Low-Level Factorization Utility
Polfed.Lanczos.lanczos_extrema — Functionlanczos_extrema(mat::AbstractMatrix; x0=nothing, maxdim=1000, tol=1e-14, eigentol=1e-8, kwargs...)
lanczos_extrema(f!::Function, x0::AbstractVector; maxdim=1000, tol=1e-14, eigentol=1e-8, kwargs...)Estimate the smallest and largest eigenvalues of a Hermitian/symmetric matrix or linear operator using two short Lanczos probes.
The returned tuple is (Emin, Emax), where Emin is obtained with which=:SR and Emax with which=:LR. If x0 is provided, it is copied and normalized before use. If omitted for matrix input, a random normalized vector with a floating element type compatible with mat is generated.
Extra keyword arguments are forwarded to lanczos, except that which and howmany are fixed internally.
Polfed.Lanczos.lanczos — Functionlanczos(mat::AbstractMatrix, x0::AbstractVecOrMat, howmany::Int;
basistype::Type{<:OrthonormalBasis}=MatrixBasis,
rot::ReOrthTechnique=FullRO(),
maxdim::Int=10howmany,
which::Symbol=:SR,
tol::Real=1e-14,
eigentol::Real=1e-8,
mapvals::Union{Function,Nothing}=nothing)
lanczos(f!::Function, x0::AbstractVecOrMat, howmany::Int;
basistype::Type{<:OrthonormalBasis}=MatrixBasis,
rot::ReOrthTechnique=FullRO(),
maxdim::Int=10howmany,
which::Symbol=:SR,
tol::Real=1e-14,
eigentol::Real=1e-8,
mapvals::Function=f!)Performs the block Lanczos factorization to compute a specified number of extremal eigenvalues and eigenvectors of a matrix or a linear operator. This method automatically detects whether the input data lives on the CPU or GPU, and dispatches to the appropriate optimized backend.
Arguments
mat::AbstractMatrix: The input matrix for which the spectral decomposition is sought.f!::Function: Matrix–vector multiplication routine, defined asf!(Y, X) = A * X. Used when the matrix is represented implicitly or stored on GPU.x0::AbstractVecOrMat: Initial block of orthonormal starting vectors. The number of columns determines the block size.howmany::Int: Number of desired eigenpairs.
Keyword Arguments
basistype::Type{<:OrthonormalBasis}=MatrixBasis: Type of orthonormal basis used (MatrixBasis,HybridMatrixBasis, orVectorBasis).maxdim::Int=500howmany: Maximum dimension of the Krylov subspace.which::Symbol=:SR: Determines which eigenvalues to converge (:SR,:LR, or:amplitude).tol::Real=1e-14: Numerical tolerance for orthogonality loss in the Lanczos basis.eigentol::Real=1e-8: Convergence tolerance for eigenvalue residuals.mapvals::Union{Function,Nothing}=nothing: Optional mapping function applied to the spectral transformation. Defaults to standard matrix–vector multiplication.
Description
The Lanczos algorithm is a Krylov subspace method used to approximate the smallest (or largest) eigenvalues and corresponding eigenvectors of large Hermitian or symmetric operators.
Two calling conventions are supported:
- Matrix interface: Pass a concrete matrix
mat— a multiplication functionf!is generated internally. - Operator interface: Pass a user-defined function
f!(Y, X)for custom or GPU-resident operators.
The function constructs:
- A
LanczosIteratorto manage Krylov basis vectors and reorthogonalization. - A
ConvergenceInfoobject to track eigenvalue convergence and stopping criteria. - A
FactorizationReportfor diagnostic output and benchmarking.
GPU arrays are automatically detected via is_gpu_array, and computations are performed on the GPU via CUDA.
Returns
A tuple containing:
vals: Approximated eigenvalues (sorted according towhich).vecs: Corresponding approximate eigenvectors (as columns).report::FactorizationReport: Detailed factorization report containing timing, convergence, and iteration statistics.
lanczos(f!::Function, x0::AbstractVecOrMat, howmany::Int; kwargs...)Operator-callback overload of lanczos.
Example: Passing Full Configuration into polfed
using Polfed
using Polfed.Models: qsun_hamiltonian
using LinearAlgebra
L_loc = 12
L_grain = 2
g0 = 1.0
α = 0.5
mat = qsun_hamiltonian(L_loc, L_grain, g0, α; use_sparse=true)
x0 = rand(size(mat, 1)); x0 ./= norm(x0)
howmany = 100
target = 0.0
mapping_cfg = MappingConfig(
parallel_strategy = MulColsParallel(),
optimize_mapping = false,
)
transform_cfg = TransformConfig(
normalization = 10.0,
cutoff = 1.7,
order_safety_factor = 0.95,
)
fact_cfg = FactorizationConfig(
tol = 1e-15,
eigentol = 1e-10,
which = :LR,
)
dos_cfg = DoSConfig(
N = 200,
R = 300,
)
vals, vecs, report = polfed(
mat,
x0,
howmany,
target;
produce_report = true,
mapping = mapping_cfg,
transform = transform_cfg,
fact = fact_cfg,
dos = dos_cfg,
)This page was generated using Literate.jl.