Restarted Krylov

template<Number T, DataLocation L = Host, bool UseLanczos = false>
class Krylov

A restarted Krylov method for computing \( f(\mathbf{A} t) \vec{b} \) [1-3].

At each restart cycle k, it constructs a basis \( V_m^{(k)} \) from the Krylov subspace \( \mathcal{K}_m(\mathbf{A}, \vec{v}_{m + 1}^{(k - 1)})\) (starting with \( \vec{v}_{m + 1}^{(0)} = \vec{b} \)) and upper Hessenberg matrix \( H_m^{(k)} \) using the Arnoldi iteration. Then approximates \( f(\mathbf{A}t) \vec{b} \) as

\[ \vec{f}^{(k)} = \vec{f}^{(k - 1)} + \beta V_m^{(k)} (f(\mathbf{H}_{km}t) \vec{e}_1)_{(k-1)m:km}\]
where \( \beta = \| \vec{b} \|_2 \) and \( H_{km} \) is the combination of all \( H_m^{(j)}, j = 1, 2, ..., k \). The algorithms end when the l2 norm of the update is less than the target tolerance or when it reaches the maximum number of iterations.

References:

[1] M. Eiermann and O. G. Ernst, “A Restarted Krylov Subspace Method for the Evaluation of Matrix Functions”, SIAM J. Numer. Anal., vol. 44, no. 6, pp. 2481–2504, Jan. 2006, doi: 10.1137/050633846.

[2] M. Afanasjew, M. Eiermann, O. G. Ernst, and S. Güttel, “Implementation of a restarted Krylov subspace method for the evaluation of matrix functions”, Linear Algebra and its Applications, vol. 429, no. 10, pp. 2293–2314, Nov. 2008, doi: 10.1016/j.laa.2008.06.029.

[3] S. Güttel, D. Kressner, and K. Lund, “Limited-memory polynomial methods for large-scale matrix functions”, GAMM-Mitteilungen, vol. 43, no. 3, p. e202000019, 2020, doi: 10.1002/gamm.202000019.

Template Parameters:
  • T – Data type used by the solver.

  • L – Location of the data.

  • UseLanczos – Use the Lanczos instead of Arnoldi for constructing the basis.

Public Functions

inline Krylov(KrylovOptions params = KrylovOptions())

Initialize the solver with the parameters params.

Parameters:

params[in] solver parameters (optional). See KrylovOptions for all options.

inline void set_params(KrylovOptions params)

Set the parameters of the solver.

Parameters:

params[in] solver parameters. See ExpvTaylorOptions for all options.

inline int num_restarts() const
Returns:

The number of restarts of the last run of the solver.

Warning

doxygenclass: Cannot find class “HighFM::SymmetricKrylov” in doxygen xml output for project “HighFM” from directory: ../build/doxygen/xml/

struct KrylovOptions

Public Members

double atol = 1E-6

Tolerance. The algorithm terminates when \( \| (f(\mathbf{H}_{km}) \vec{e}_1)_{(k-1)m:km} \| \leq atol + rtol * \| \vec{f} \|_2 \) where \( \vec{f} \) is the approximation of \( f(\mathbf{A}) \vec{b} \) from the previous restart cycle.

int restart_length = 20

The size of the basis in each restart cycle.

int max_restarts = 40

Maximum number of restarts.

int num_ortho_vec = kAuto

If positive (k > 0), indicate the number of basis vectors to orthogonalize against (i.e., truncate the Arnoldi decomposition after k iterations). Otherwise, the method will perform the full orthogonalization. If the solver uses the Lanczos iteration, then this parameter has no effect.

int verbose = 0

Set the verbosity level. 0 - No report. 1 - Report the update norm

\[ \| \vec{f}_k \|_2 \]
of the restart cycle k.