Krylov methods#
Implemented in yastn#
We provide a high-level, backend-agnostic implementation of some Krylov-based algorithms used in yastn.tn.mps.
They assume a linear operation acting on a generalized vector, with the vector being an instance of a class that includes methods
norm, vdot, add, expand_krylov_space, among others.
Examples of such a vector include yastn.Tensor (see methods).
- yastn.expmv(f, v, t=1.0, tol=1e-12, ncv=10, hermitian=False, normalize=False, return_info=False, **kwargs) vector[source]#
Calculate \(e^{(tF)}v\), where \(v\) is a vector, and \(F(v)\) is linear operator acting on \(v\).
The algorithm of: J. Niesen, W. M. Wright, ACM Trans. Math. Softw. 38, 22 (2012), Algorithm 919: A Krylov subspace algorithm for evaluating the phi-functions appearing in exponential integrators.
- Parameters:
f (Callable[[vector], vector]) – defines an action of a “square matrix” on vector.
v (vector) – input vector to apply exponential map onto.
t (number) – exponent amplitude.
tol (number) – targeted tolerance; it is used to update the time-step and size of Krylov space. The returned result should have better tolerance, as correction is included.
ncv (int) – Initial guess for the size of the Krylov space.
hermitian (bool) – Assume that
fis a hermitian operator, in which case Lanczos iterations are used. Otherwise Arnoldi iterations are used to span the Krylov space.normalize (bool) – Whether to normalize the result to unity using 2-norm.
return_info (bool) – if
True, returns(vector, info), whereinfo.ncv: guess of the Krylov-space size,info.error: estimate of error (likely over-estimate)info.krylov_steps: number of execution off(x),info.steps: number of steps to reacht,
kwargs (any) – Further parameters that are passed to
expand_krylov_space()andadd().
- yastn.eigs(f, v0, k=1, which='SR', ncv=10, maxiter=None, tol=1e-13, hermitian=False, **kwargs) tuple[array, Sequence[vectors]][source]#
Search for dominant eigenvalues of linear operator
fusing Arnoldi algorithm. Economic implementation (without restart) for internal use withinyastn.tn.mps.dmrg_().- Parameters:
f (function) – define an action of a ‘square matrix’ on the ‘vector’
v0.f(v0)should preserve the signature ofv0.v0 (Tensor) – Initial guess, ‘vector’ to span the Krylov space.
k (int) – Number of desired eigenvalues and eigenvectors. The default is 1.
which (str) – One of [
‘LM’,‘LR’,‘SR’] specifying which k eigenvectors and eigenvalues to find:‘LM’: largest magnitude,‘SM’: smallest magnitude,‘LR’: largest real part,‘SR’: smallest real part.ncv (int) – Dimension of the employed Krylov space. The default is 10. Must be greated than k.
maxiter (int) – Maximal number of restarts; (not implemented for now)
tol (float) – Stopping criterion for an expansion of the Krylov subspace. The default is
1e-13. (not implemented for now)hermitian (bool) – Assume that
fis a hermitian operator, in which case Lanczos iterations are used. Otherwise Arnoldi iterations are used to span the Krylov space.
- yastn.lin_solver(f, b, v0, ncv=10, tol=1e-13, pinv_tol=1e-13, hermitian=False, **kwargs) tuple[array, Sequence[vectors]][source]#
Search for solution of the linear equation
f(x) = b, wherexis estimated vector andf(x)is matrix-vector operation. Implementation based on pseudoinverse of Krylov expansion [1].Ref.[1] https://www.math.iit.edu/~fass/577_ch4_app.pdf
- Parameters:
f (function) – define an action of a ‘square matrix’ on the ‘vector’
v0.f(v0)should preserve the signature ofv0.b (Tensor) – free term for the linear equation.
v0 (Tensor) – Initial guess span the Krylov space.
ncv (int) – Dimension of the employed Krylov space. The default is 10.
tol (float) – Stopping criterion for an expansion of the Krylov subspace. The default is
1e-13. TODO Not implemented yet.pinv_tol (float) – Cutoff for pseudoinverve. Sets lower bound on inverted Schmidt values. The default is
1e-13.hermitian (bool) – Assume that
fis a hermitian operator, in which case Lanczos iterations are used. Otherwise Arnoldi iterations are used to span the Krylov space.
Results#
- vf: Tensor
Approximation of
vinf(v) = bproblem.- res: float
norm of the resudual vector
r = norm(f(vf) - b).
Other libraries#
With NumPy backend, it is possible to link to algorithms in
sparse.sparse.linalg,
employing yastn.Tensor.compress_to_1d() and yastn.decompress_from_1d().
See, an example.