Configs and Parallelization Types

This page documents the main configuration structs and execution strategy types used by polfed.

MappingConfig

Polfed.PolfedCore.MappingConfigType
MappingConfig(; kwargs...)

Configuration for the mapping stage of polfed, including:

  • mapping parallelization strategy,
  • optional matrix-entrypoint auto-optimization,
  • optional pre-rescaled mapping callback,
  • optional custom Clenshaw kernels, and
  • optional explicit spectral bounds (Emin, Emax).

Pass this struct as mapping=... in polfed.

Fields

  • parallel_strategy::Union{Parallelization,Nothing} Mapping parallelization strategy.

  • optimize_mapping::Bool Enables matrix-entrypoint optimization (polfed(mat, ...)) for structured matrices with repeated offdiagonal values. This can reduce memory traffic by building optimized mapping and Clenshaw kernels under the hood.

  • f!_rescaled::Union{Nothing,Function} Optional user-provided mapping for the rescaled operator H_tilde = (H - b)/a. Expected signature: f!_rescaled(Y, X) -> nothing, with in-place write to Y. It should support both vector and matrix X (for Lanczos and Block Lanczos).

  • clenshaw_recurrence::Union{Nothing,Function} Optional custom Clenshaw recurrence kernel. Accepted signatures:

    • (b1, b2, b3, c, k, X), or
    • (b1, b2, b3, c, X).

    It should update b1 in place and avoid allocations.

  • clenshaw_finalsum::Union{Nothing,Function} Optional custom Clenshaw final-sum kernel with signature: (b2, b3, c0, Y, X). It should write the final transformed output into Y.

  • Emin::Union{Real,Nothing}, Emax::Union{Real,Nothing} Optional unrescaled spectral bounds used for internal rescaling and target conversion. If either is nothing, POLFED estimates missing bounds with short Lanczos probes.

Keyword Defaults

  • parallel_strategy = PolfedDefaults.parallel_strategy
  • optimize_mapping = PolfedDefaults.optimize_mapping
  • f!_rescaled = nothing
  • clenshaw_recurrence = nothing
  • clenshaw_finalsum = nothing
  • Emin = nothing
  • Emax = nothing

Notes

  • optimize_mapping=true is intended for matrix input; for callback input (polfed(f!, ...)) it is ignored unless you provide custom mapping kernels.
  • If you provide f!_rescaled, keep it consistent with Emin/Emax and your original operator scaling.
  • If all mapping parallelization is handled externally (custom mapping or GPU kernel strategy), prefer NoParallel.

Example

mapping = MappingConfig(
    parallel_strategy = MulColsParallel(),
    optimize_mapping = true,
)
vals, vecs, report = polfed(mat, x0, howmany, target;
    produce_report = true,
    mapping = mapping,
)

TransformConfig

Polfed.PolfedCore.TransformConfigType
TransformConfig(; kwargs...)

Configuration for polynomial spectral transformation used by polfed.

This struct controls how the filter polynomial is constructed around the chosen target: coefficients, normalization/cutoff behavior, interval constraints, and order selection.

Pass this struct as transform=... in polfed.

Fields

  • coefficients::Function Coefficient generator used by Clenshaw polynomial evaluation. Expected callable form: coefficients(target_rescaled::Real, n::Int) -> Real. The default is Chebyshev-like coefficients from PolfedDefaults.coefficients.

  • normalization::Real Output normalization of the polynomial filter. POLFED rescales the polynomial so that P(target) = normalization.

  • cutoff::Real Threshold used when deriving interval width/order automatically. Lower values typically imply wider filters and/or larger polynomial orders.

  • left::Union{Real,Nothing}, right::Union{Real,Nothing} Optional unrescaled interval bounds. If provided, both should be set together. They are internally converted to rescaled coordinates.

  • order::Union{Int,Nothing} Optional fixed polynomial order K. If nothing, POLFED computes K from the chosen interval and cutoff.

  • order_safety_factor::Real Multiplicative safety factor used when order is computed automatically. Values below 1 reduce order (and POLFED recomputes a consistent interval). This is often used to improve robustness when requesting many eigenpairs.

  • polynomialtype::Symbol Polynomial family used by Clenshaw evaluation. Supported types follow internal Clenshaw kernels (for example :Chebyshev, :Legendre, :Hermite, :Taylor).

Keyword Defaults

  • coefficients = PolfedDefaults.coefficients
  • normalization = PolfedDefaults.normalization
  • cutoff = PolfedDefaults.cutoff
  • left = nothing
  • right = nothing
  • order = nothing
  • order_safety_factor = PolfedDefaults.order_safety_factor
  • polynomialtype = PolfedDefaults.polynomialtype

Notes

  • Do not set both explicit bounds (left/right) and explicit order: this overconstrains the filter and raises an error.
  • Plain numeric target or (:unrescaled, E) is converted to rescaled coordinates internally before calling coefficients.
  • If you use custom coefficients/polynomial families, consider validating order, interval, and convergence with produce_report=true.

Example

transform = TransformConfig(
    cutoff = 0.17,
    order_safety_factor = 0.95,
)
vals, vecs, report = polfed(mat, x0, howmany, target;
    produce_report = true,
    transform = transform,
)

FactorizationConfig

Polfed.PolfedCore.FactorizationConfigType
FactorizationConfig(; kwargs...)

Configuration structure for the Krylov factorization used within POLFED.

This struct defines options that control the iterative factorization process, including re-orthogonalization strategies, basis representation, and convergence tolerances. These parameters are passed to the underlying Lanczos or Block Lanczos solver.

Fields

  • rot::ReOrthTechnique Re-orthogonalization technique. Options:

    • FullRO() (default): Full re-orthogonalization.
    • PRO(): Partial re-orthogonalization.

    Since the polynomial spectral transformation typically dominates the runtime, the overhead of full re-orthogonalization is negligible, making FullRO() generally preferred.

  • basistype::Type{<:OrthonormalBasis} Type of orthonormal basis used. Default: PolfedDefaults.basistype. Options:

    • MatrixBasis: Standard CPU-based orthonormal basis.
    • HybridMatrixBasis: Enables hybrid CPU–GPU computation (slower due to data transfer overhead).
  • which::Symbol Specifies which eigenvalues to target (:LR or :SR). Default: PolfedDefaults.which.

  • tol::Real Convergence tolerance for the Krylov factorization. Default: PolfedDefaults.tol. Convergence is assessed based on the residual norm of the Ritz pairs.

  • eigentol::Union{Real,Nothing} Convergence tolerance for eigenvectors. Default: PolfedDefaults.eigentol. Computed using the L₂ norm of the eigenvector residual.

Example

fact_cfg = FactorizationConfig(tol=1e-12, which=:SR)

Notes

  • The tol parameter typically controls when Ritz values are accepted as converged.
  • For most applications, full re-orthogonalization is recommended for numerical stability.
  • The hybrid basis is mainly useful when the matrix–vector products are performed on the GPU.

DoSConfig

Polfed.PolfedCore.DoSConfigType
DoSConfig(; kwargs...)

Configuration structure for Density of States (DoS) calculations in POLFED.

This struct specifies parameters for stochastic estimation of the density of states using random vector sampling and Chebyshev moment expansion. The configuration affects the resolution and statistical accuracy of the DoS computation.

Fields

  • N::Int Number of Chebyshev moments (or energy points) used in the DoS calculation. Higher values improve spectral resolution at increased computational cost. Default: PolfedDefaults.N (typically 250).

  • R::Int Number of random vectors used to stochastically estimate the DoS. Increasing R reduces statistical noise in the estimated spectrum. Default: PolfedDefaults.R (typically 300).

  • kernel::Symbol Kernel used in the KPM expansion. Default: PolfedDefaults.kernel (typically :Jackson).

Example

```julia

Use default DoS settings

dos_cfg = DoSConfig()

Custom configuration

doscfgcustom = DoSConfig(N=500, R=400)

Notes

  • The statistical accuracy of the estimated DoS scales as 1/√R, while the spectral resolution improves roughly as 1/N.
  • The DoS computation is typically performed as a preparatory step before spectral filtering, to identify the spectral interval around the target energy and determine an appropriate polynomial order.
  • For large sparse systems, moderate R values (100–400) usually provide a good balance between computational cost and accuracy.

Parallelization Strategy Types

Polfed.PolfedCore.ParallelizationType

Abstract base type for all parallelization strategies.

All parallelization types should inherit from Parallelization. It is used to define interfaces and enable multiple dispatch based on the chosen strategy.

Polfed.PolfedCore.MulColsParallelType

Represents multi-column parallelization, where independent columns of a matrix or data are processed in parallel.

Each column can be handled by a separate worker or thread, providing speedup when the number of columns is large. This strategy assumes that column operations are independent.

Polfed.PolfedCore.TwoLevelParallelType

Represents a two-level parallelization strategy.

  1. Inter-column parallelism: Columns are distributed among a pool of worker processes.
  2. Intra-column parallelism: Within each column, multiple threads (or tasks) perform computations in parallel.

Fields

  • worker_pool::WorkerPool: A collection of worker processes handling column-level tasks.
  • nt_per_col::Int: Number of threads per column for intra-column parallelism.

Constructor

TwoLevelParallel(nt_per_col::Int)

Creates a TwoLevelParallel instance with a new worker pool and the specified number of threads per column.

This strategy is ideal for large-scale computations where both the number of columns and the per-column workload are significant.

Polfed.PolfedCore.NoParallelType

Represents a serial execution strategy with no parallelization.

Useful for debugging, testing, when parallelism is not required, and mostly when all of the parallelization is handled by the user (for example, when using GPUs).

Processing Unit Types

Polfed.PolfedCore.CPUType
CPU() -> CPU

Construct a CPU processing-unit adapter containing matrix/vector types and allocation/random helper functions.

GPU processing is available through GPU() in CUDA-enabled environments.

Polfed.Lanczos.PartialROType
PartialRO(skip::Integer)

Partial reorthogonalization strategy. Reorthogonalization is performed only on selected iterations controlled by skip.

Polfed.Lanczos.MatrixBasisType
MatrixBasis{E}(maxdim::Int, x0::AbstractVecOrMat, pu::ProcessingUnit)

Dense column-major basis storage for Lanczos vectors/blocks.

Polfed.Lanczos.HybridMatrixBasisType
HybridMatrixBasis(maxdim_gpu::Int, maxdim_cpu::Int, x0::AbstractVecOrMat{T}) where {T}

Two-tier basis storage backend: fills GPU matrix columns first, then spills to CPU matrix columns.

Supported target Specifications

  • :maxdos
  • :middle
  • (:offset, frac)
  • (:unrescaled, E) or plain E::Real
  • (:rescaled, e)

This page was generated using Literate.jl.