slepc-main 2024-11-09
Report Typos and Errors

BVMatLanczos

Computes a Lanczos factorization associated with a matrix.

Synopsis

#include "slepcbv.h"   
PetscErrorCode BVMatLanczos(BV V,Mat A,Mat T,PetscInt k,PetscInt *m,PetscReal *beta,PetscBool *breakdown)
Collective

Input Parameters

V  - basis vectors context
A  - the matrix
T  - (optional) the tridiagonal matrix
k  - number of locked columns
m  - dimension of the Lanczos basis, may be modified

Output Parameters

beta  - (optional) norm of last vector before normalization
breakdown  - (optional) flag indicating that breakdown occurred

Notes

Computes an m-step Lanczos factorization for matrix A, with full reorthogonalization. At each Lanczos step, the corresponding Lanczos vector is orthogonalized with respect to all previous Lanczos vectors. This is equivalent to computing an m-step Arnoldi factorization and exploting symmetry of the operator.

The first k columns are assumed to be locked and therefore they are not modified. On exit, the following relation is satisfied

                   A * V - V * T = beta*v_m * e_m^T

where the columns of V are the Lanczos vectors (which are B-orthonormal), T is a real symmetric tridiagonal matrix, and e_m is the m-th vector of the canonical basis. On exit, beta contains the B-norm of V[m] before normalization. The T matrix is stored in a special way, its first column contains the diagonal elements, and its second column the off-diagonal ones. In complex scalars, the elements are stored as PetscReal and thus occupy only the first column of the Mat object. This is the same storage scheme used in matrix DS_MAT_T obtained with DSGetMat().

The breakdown flag indicates that orthogonalization failed, see BVOrthonormalizeColumn(). In that case, on exit m contains the index of the column that failed.

The values of k and m are not restricted to the active columns of V.

To create a Lanczos factorization from scratch, set k=0 and make sure the first column contains the normalized initial vector.

See Also

BVMatArnoldi(), BVSetActiveColumns(), BVOrthonormalizeColumn(), DSGetMat()

Level

advanced

Location

src/sys/classes/bv/interface/bvkrylov.c

Index of all BV routines
Table of Contents for all manual pages
Index of all manual pages