Reports, Logging, and Defaults
Reporting is enabled through produce_report=true in polfed. The returned Report provides spectral-transform, factorization, and benchmark diagnostics.
Report Types and Display Functions
Polfed.PolfedCore.display_report — Functiondisplay_report(report::Report;
use_colors::Bool=true,
include_spectral_transform::Bool=true,
include_factorization::Bool=true,
show_convergence_details::Bool=false,
include_benchmark::Bool=true)Displays a comprehensive report for a Polfed.PolfedCore.Report object.
Keyword Arguments
use_colors::Bool=true: Enable or disable ANSI color formatting.include_spectral_transform::Bool=true: Include the spectral transform sectiondisplay_spectral_report.include_factorization::Bool=true: Include the factorization sectiondisplay_factorization_report.show_convergence_details::Bool=false: Iftrue, display detailed convergence information in the factorization report.include_benchmark::Bool=true: Include the benchmark sectiondisplay_benchmark_report.
Notes
The show_timings flag for the factorization section is always set to false in this combined display function.
Polfed.PolfedCore.Report — Typestruct ReportA container for aggregating the main results of a simulation or computation. One can deisplay it with display_report function.
Fields
SpectralTransformReport- Tells you about spectral transformation parameter like rescaled targeted interval, order of polynomial, number of matrix-vector multiplications etc.FactorizationReport- Provides insights into the factorization process, including convergence behavior, iteration counts and errors of obtained eigenpairs.BenchmarkReport- Summarizes performance such as total walltime and CPU time, as well as the times spend in different parts of the factorization algorithm.
Polfed.Lanczos.FactorizationReport — TypeFactorizationReport(convergenceinfo::ConvergenceInfo,
factorization::KrylovFactorization,
iterator::LanczosIterator,
walltimes::Vector{<:Real},
cputimes::Vector{<:Real})Stores detailed information about a single factorization run, including basis statistics, convergence results, and timing data.
This report is typically constructed internally after a Lanczos or block-Lanczos factorization completes, and it provides the core diagnostic information used by higher-level reporting tools.
Arguments
convergenceinfo::ConvergenceInfo: Summary of convergence behavior, including tolerances, number of converged eigenpairs, iteration counts, and convergence history.factorization::KrylovFactorization: The factorization object (e.g.,LanczosFactorizationorBlockLanczosFactorization) containing the computed Krylov basis and internal parameters.iterator::LanczosIterator: The iterator used to build the Krylov subspace, providing information about reorthogonalization techniques and iteration limits.walltimes::Vector{<:Real}: Wall-clock times (in seconds) for each major code segment of the factorization routine.cputimes::Vector{<:Real}: CPU times (in seconds) corresponding to the same code segments as inwalltimes.
Fields
factorizationtype::String: Type of the factorization method used (e.g.,"LanczosFactorization").blocksize::Integer: Block size used in the factorization (1 for standard Lanczos).basistype::Type{<:OrthonormalBasis}: Type of orthonormal basis used (e.g.,MatrixBasis,HybridMatrixBasis).ncols_cpu_reserved,ncols_gpu_reserved::Integer: Number of basis vectors reserved on CPU and GPU, respectively.ncols_cpu_needed,ncols_gpu_needed::Integer: Number of basis vectors actually used during the run.howmany::Integer: Number of eigenpairs requested.itersneeded::Integer: Number of iterations performed before convergence.itersreserved::Integer: Number of iterations preallocated or reserved.rot::ReOrthTechnique: Reorthogonalization method used (e.g.,FullReorth,PartialReorth).tol,residual::Real: Lanczos convergence tolerance and final residual norm.converged_tol::Integer: Number of eigenpairs satisfying the Lanczos convergence tolerance.eigentol,eigenresidual::Real: Eigenvalue convergence tolerance and corresponding residual.converged_eigentol::Integer: Number of eigenpairs satisfying the eigenvalue convergence tolerance.numofchecks::Integer: Number of convergence checks performed during the run.converged_history::Vector{<:Integer}: Number of converged eigenpairs at each check.krylovbasisdim_history::Vector{<:Integer}: Krylov subspace dimension recorded at each check.maxresidual_history::Vector: Maximum residual per convergence check.walltimes,cputimes::Vector{<:Real}: Timings for individual stages (e.g., mapping, reorthogonalization, convergence check).
Notes
- This report is primarily used by
display_factorization_reportto summarize performance and convergence behavior. - Supports both CPU-only and hybrid CPU/GPU basis types.
Polfed.Lanczos.display_factorization_report — Functiondisplay_factorization_report(report::FactorizationReport;
use_colors::Bool=true,
show_convergence_details::Bool=false,
show_timings::Bool=true)Displays a detailed summary of a FactorizationReport object, including convergence, iteration, and timing information. Optionally includes a convergence history table for fine-grained inspection.
Arguments
report::FactorizationReport: The report object to display.use_colors::Bool=true: Iftrue, enables ANSI color formatting for terminal readability.show_convergence_details::Bool=false: Iftrue, prints a formatted table showing convergence evolution across checks.show_timings::Bool=true: Iftrue, includes timing summaries for walltime and CPU time, including stage-wise percentages.
Output
The function prints:
- Convergence summary – how many eigenpairs reached the desired tolerance.
- Iteration details – total iterations used vs. reserved.
- Basis information – type of basis and reorthogonalization strategy.
- Timing summary – breakdown of walltime and CPU time.
- (Optional) Convergence table – evolution of convergence metrics if
show_convergence_details=true.
Notes
- Uses color highlighting to indicate convergence success:
- ✅ Green → all requested eigenpairs converged.
- ⚠️ Yellow → moderate margin to reserved iterations.
- ❌ Red → insufficient convergence or heavy overestimation.
- Timing percentages are grouped as:
(Mapping, Reorthogonalization, Convergence check, others) - To suppress colors for non-TTY environments (e.g., when redirecting output to file), set
use_colors=false.
Polfed.PolfedCore.SpectralTransformReport — Typestruct SpectralTransformReportHolds parameters and statistics for the spectral transformation phase of the computation.
Constructor
SpectralTransformReport(transform_plan,
mapping_plan,
fact::FactorizationReport)Fields
target::Real: The target eigenvalue (if any).target_strategy::String: The user-facing target strategy used to choose the target.left::Real: Left endpoint of the spectral interval.right::Real: Right endpoint of the spectral interval.polynomialtype::String: The polynomial type used (e.g.,"Chebyshev").order::Integer: Order of the polynomial.order_safety_factor::Real: Safety multiplier applied to the order.parallel_strategy::Parallelization: Parallelization strategy used.howmany::Integer: Number of requested eigenpairs.howmany_ininterval::Integer: Number of eigenvalues found in the target interval.matrixvec_muls::Integer: Number of matrix–vector multiplications performed.optimize_mapping::Bool: Whether automatic mapping optimization was enabled.
Notes
The field howmany_ininterval may be updated later, after results are known.
Polfed.PolfedCore.display_spectral_report — Functiondisplay_spectral_report(report::SpectralTransformReport; use_colors::Bool=true)Pretty-prints a SpectralTransformReport summarizing key parameters and statistics.
Keyword Arguments
use_colors::Bool=true: Enable or disable ANSI color formatting.
Output
Displays:
- Target strategy, target energy, and number of eigenpairs requested.
- Spectral interval and width.
- Polynomial type and order.
- Number of matrix–vector multiplications performed.
- Parallelization strategy and automatic optimization status.
Polfed.PolfedCore.BenchmarkReport — Typestruct BenchmarkReportAggregates performance statistics from a single POLFED run, such as total walltime and CPU time, as well as the times spend in different parts of the factorization algorithm. It serves as a high-level summary of the computational efficiency of a polfed.
Fields
timings::TimingReport: Stores walltime and CPU time statistics across the full simulation and individual factorization stages.memory_cpu::Union{MemoryReport, Nothing}: Memory usage report for CPU computations (ornothingif not applicable).memory_gpu::Union{MemoryReport, Nothing}: Memory usage report for GPU computations (ornothingif not applicable).
Constructors
```julia BenchmarkReport( factreport::FactorizationReport, totalwt::Real, total_ct::Real, x0::AbstractVecOrMat{T}, pu::ProcessingUnit ) where {T<:Number}
Constructs a new BenchmarkReport from a FactorizationReport, recording total walltime and CPU time along with detailed memory usage estimates.
- The TimingReport is created automatically from factreport.walltimes and factreport.cputimes.
- The MemoryReport entries are generated based on the basistype in fact_report and the given processing unit pu (CPU() or GPU()).
- If a HybridMatrixBasis is detected, both CPU and GPU memory usage are recorded.
Polfed.PolfedCore.display_benchmark_report — Functiondisplay_benchmark_report(report::BenchmarkReport;
show_timings::Bool=true,
show_memory::Bool=false,
use_colors::Bool=true)Pretty-prints a BenchmarkReport with optional colorized output.
Keyword Arguments
show_timings::Bool=true: Display timing information (walltime and CPU time).show_memory::Bool=false: Display memory usage details.use_colors::Bool=true: Use ANSI colors in the printed output.
Output
- Summarizes total walltime and CPU time.
- Shows distribution of time across factorization stages.
- Optionally displays memory usage for CPU and/or GPU computations.
Logging Levels and Helpers
Polfed.Common.POLFED_SILENT_LEVEL — ConstantPOLFED_SILENT_LEVEL::Int = 0Disable POLFED logging output.
Polfed.Common.POLFED_WARN_LEVEL — ConstantPOLFED_WARN_LEVEL::Int = 1Emit warning-level POLFED logs.
Polfed.Common.POLFED_INFO_LEVEL — ConstantPOLFED_INFO_LEVEL::Int = 2Emit informational and warning POLFED logs.
Polfed.Common.POLFED_DEBUG_LEVEL — ConstantPOLFED_DEBUG_LEVEL::Int = 3Emit debug, info, and warning POLFED logs.
Polfed.Common.verbosity — Constantverbosity::Base.RefValue{Int}Runtime POLFED logging level. Typical values:
0(POLFED_SILENT_LEVEL)1(POLFED_WARN_LEVEL)2(POLFED_INFO_LEVEL)3(POLFED_DEBUG_LEVEL)
Polfed.Common.should_log — Functionshould_log(level::Int) -> BoolReturn true if current verbosity allows logging at level.
Polfed.Common.polfed_log — Functionpolfed_log(level::Int, msg; kwargs...) -> nothingEmit a structured POLFED log message gated by verbosity.
Log routing:
- debug level:
@infowith[DEBUG]prefix, - info level:
@info, - warn level:
@warn.
PolfedDefaults Module
Polfed.PolfedCore.PolfedDefaults — ModulePolfedDefaults — Default settings for POLFED.jl
This module contains all default constants and functions used across the POLFED framework, including:
- Krylov subspace dimension prediction
- Reporting and mapping options
- Spectral transformation defaults (polynomial type, coefficients, cutoff, etc.)
- Lanczos algorithm defaults (rotation type, basis, tolerances)
- Density-of-states calculation defaults (kernel, number of points, random vectors)
All defaults are intended to provide reasonable starting points for typical simulations, but can be overridden by user-specified options.
Default values:
produce_report = falseoptimize_mapping = falsepolynomialtype = :Chebyshevcutoff = 0.17normalization = 1.0order_safety_factor = 0.97parallel_strategy = MulColsParallel()overestimate_iters = 1.25rot = FullRO()basistype = MatrixBasistol = 1e-14eigentol = 1e-9which = :LRkernel = :JacksonN = 250R = 300
Krylov Dimension Prediction
Polfed.PolfedCore.PolfedDefaults.expectedkrylovdim — Functionexpectedkrylovdim(howmany::Int, blocksize::Int, η::Real)
Predicts the expected Krylov subspace dimension for POLFED given the number of vectors (howmany), block size (blocksize), and a safety factor η.
Returns an integer, rounded up using ceil(Int64, ...)
expectedkrylovdim(howmany, blocksize, η) = ceil(Int64, (20.427*blocksize + 1.696*howmany)*η)Example:
PolfedDefaults.expectedkrylovdim(10, 5, 1.1) # => 240 (example)Reporting and Mapping Defaults
Polfed.PolfedCore.PolfedDefaults.produce_report — Constantproduce_report::Bool = false
Whether to produce a summary report at the end of a POLFED run.
Polfed.PolfedCore.PolfedDefaults.optimize_mapping — Constantoptimize_mapping::Bool = false
Whether to optimize the Hamiltonian mapping before computation.
Spectral Transformation Defaults
Polfed.PolfedCore.PolfedDefaults.coefficients — Functioncoefficients(λ::T, n::Int) where {T<:Real}
Function that produces polynomial coefficients for spectral transformations. Default is Chebyshev polynomial:
- For n = 0: coefficient is 2 * cos(0 * acos(λ)) = 2
- For n > 0: coefficient is cos(n * acos(λ))
const coefficients(λ::T, n::Int) where {T<:Real} = T((2 - ==(n,0)) * cos(n * acos(λ)))Polfed.PolfedCore.PolfedDefaults.polynomialtype — Constantpolynomialtype::Symbol = :Chebyshev
Default polynomial type for spectral transformations. Currently only :Chebyshev is supported.
Polfed.PolfedCore.PolfedDefaults.cutoff — Constantcutoff::Real = 0.17
Default cutoff value for the polynomial expansion.
Polfed.PolfedCore.PolfedDefaults.normalization — Constantnormalization::Real = 1.0
Spectral transformation is normalized so that P_K(target) = normalization.
Polfed.PolfedCore.PolfedDefaults.order_safety_factor — Constantorder_safety_factor::Real=0.97
Safety factor used to reduce the polynomial order to ensure, all of the howmany target eigenvalues are captured within the filter.
Polfed.PolfedCore.PolfedDefaults.parallel_strategy — Constantparallel_strategy::Parallelization = MulColsParallel()
Default parallelization strategy for spectral transformations. Default is MulColsParallel().
Polfed.PolfedCore.PolfedDefaults.overestimate_iters — Constantoverestimate_iters::Real = 1.25
Factor used to overestimate the number of iterations needed for convergence of the algorithm. In particular, the predicted Krylov dimension expectedkrylovdim is multiplied by this factor to ensure sufficient iterations.
Note
Consider the case where not all of the requested eigenpairs are converged, and you check the Polfed.PolfedCore.Report at the end of the run (with Polfed.PolfedCore.display_report). There you can see that not all eigenpairs are converged, and under Polfed.Lanczos.FactorizationReport you can see that all iterations were used. This is usually a clear sign that one should increase overestimate_iters factor.
Factorization Defaults
Polfed.PolfedCore.PolfedDefaults.rot — Constantrot::FullRO = FullRO()
Default reorthogonalization type for the Lanczos algorithm, another option is to use PartialRO(skip) but it does not bring any benefit in most cases.
Polfed.PolfedCore.PolfedDefaults.basistype — Typebasistype::Type = MatrixBasis
Default basis type for the Lanczos algorithm, and also the only one.
Polfed.PolfedCore.PolfedDefaults.tol — Constanttol::Real = 1e-14
Default convergence tolerance for Lanczos iterations.
Polfed.PolfedCore.PolfedDefaults.eigentol — Constanteigentol::Real = 1e-9
Default eigenvalue tolerance for Lanczos eigenvectors.
Polfed.PolfedCore.PolfedDefaults.which — Constantwhich::Symbol = :LR
Default selection of eigenvalues: :LR or :SR.
DoS Defaults
Polfed.PolfedCore.PolfedDefaults.kernel — Constantkernel::Symbol = :Jackson
Default kernel used in density-of-states calculations. Available options are also :Lorentz, :Fejer, :LanczosK, :WangZunger, and :Dirichlet.
Polfed.PolfedCore.PolfedDefaults.N — ConstantN::Int = 250
Default order of Chebyshev polynomial expansion in density-of-states calculations.
Polfed.PolfedCore.PolfedDefaults.R — ConstantR::Int = 300
Default number of random vectors used in density-of-states calculations.
This page was generated using Literate.jl.