Statistics

Krylov.SimpleStatsType

Type for storing statistics returned by the majority of (block) Krylov solvers.

The fields are as follows:

  • niter: The total number of iterations completed by the solver;
  • solved: Indicates whether the solver successfully reached convergence (true if solved, false otherwise);
  • inconsistent: Flags whether the system was detected as inconsistent (i.e., when b is not in the range of A);
  • indefinite: Flags whether the system was detected as indefinite (i.e., when A is not positive definite);
  • npcCount: The number of nonpositive curvature directions encountered during the solve;
  • residuals: A vector containing the residual norms at each iteration;
  • Aresiduals: A vector of A'-residual norms at each iteration;
  • Acond: An estimate of the condition number of matrix A.
  • allocation_timer: The elapsed time (in seconds) spent on allocations;
  • timer: The elapsed time (in seconds) taken by the solver to complete all iterations;
  • status: A string indicating the outcome of the solve, providing additional details beyond solved.
source
Krylov.LanczosStatsType

Type for storing statistics returned by CG-LANCZOS. The fields are as follows:

  • niter
  • solved
  • residuals
  • indefinite
  • Anorm
  • Acond
  • allocation_timer
  • timer
  • status
source
Krylov.LanczosShiftStatsType

Type for storing statistics returned by CG-LANCZOS-SHIFT and CGLS-LANCZOS-SHIFT. The fields are as follows:

  • niter
  • solved
  • residuals
  • indefinite
  • Anorm
  • Acond
  • allocation_timer
  • timer
  • status
source
Krylov.SymmlqStatsType

Type for storing statistics returned by SYMMLQ. The fields are as follows:

  • niter
  • solved
  • residuals
  • residualscg
  • errors
  • errorscg
  • Anorm
  • Acond
  • allocation_timer
  • timer
  • status
source
Krylov.AdjointStatsType

Type for storing statistics returned by adjoint systems solvers BiLQR and TriLQR. The fields are as follows:

  • niter
  • solved_primal
  • solved_dual
  • residuals_primal
  • residuals_dual
  • allocation_timer
  • timer
  • status
source
Krylov.LNLQStatsType

Type for storing statistics returned by the LNLQ method. The fields are as follows:

  • niter
  • solved
  • residuals
  • errorwithbnd
  • errorbndx
  • errorbndy
  • allocation_timer
  • timer
  • status
source
Krylov.LSLQStatsType

Type for storing statistics returned by the LSLQ method. The fields are as follows:

  • niter
  • solved
  • inconsistent
  • residuals
  • Aresiduals
  • err_lbnds
  • errorwithbnd
  • errubndslq
  • errubndscg
  • allocation_timer
  • timer
  • status
source
Krylov.LsmrStatsType

Type for storing statistics returned by LSMR. The fields are as follows:

  • niter
  • solved
  • inconsistent
  • residuals
  • Aresiduals
  • Acond
  • Anorm
  • xNorm
  • allocation_timer
  • timer
  • status
source

Workspaces

Utilities

Krylov.roots_quadraticFunction
roots = roots_quadratic(q₂, q₁, q₀; nitref)

Find the real roots of the quadratic

q(x) = q₂ x² + q₁ x + q₀,

where q₂, q₁ and q₀ are real. Care is taken to avoid numerical cancellation. Optionally, nitref steps of iterative refinement may be performed to improve accuracy. By default, nitref=1.

source
Krylov.sym_givensFunction
(c, s, ρ) = sym_givens(a, b)

Numerically stable symmetric Givens reflection. Given a and b reals, return (c, s, ρ) such that

[ c  s ] [ a ] = [ ρ ]
[ s -c ] [ b ] = [ 0 ].
source

Numerically stable symmetric Givens reflection. Given a and b complexes, return (c, s, ρ) with c real and (s, ρ) complexes such that

[ c   s ] [ a ] = [ ρ ]
[ s̅  -c ] [ b ] = [ 0 ].
source
Krylov.to_boundaryFunction
roots = to_boundary(n, x, d, radius; flip, xNorm2, dNorm2)

Given a trust-region radius radius, a vector x lying inside the trust-region and a direction d, return σ1 and σ2 such that

‖x + σi d‖ = radius, i = 1, 2

in the Euclidean norm. n is the length of vectors x and d. If known, ‖x‖² and ‖d‖² may be supplied with xNorm2 and dNorm2.

If flip is set to true, σ1 and σ2 are computed such that

‖x - σi d‖ = radius, i = 1, 2.
source
Krylov.vec2strFunction
s = vec2str(x; ndisp)

Display an array in the form

[ -3.0e-01 -5.1e-01  1.9e-01 ... -2.3e-01 -4.4e-01  2.4e-01 ]

with (ndisp - 1)/2 elements on each side.

source
Krylov.ktypeofFunction
S = ktypeof(v)

Return the most relevant storage type S based on the type of v.

source
Krylov.vector_to_matrixFunction
M = vector_to_matrix(S)

Return the dense matrix storage type M related to the dense vector storage type S.

source
Krylov.matrix_to_vectorFunction
S = matrix_to_vector(M)

Return the dense vector storage type S related to the dense matrix storage type M.

source