SVDComputeError#

Computes the error (based on the residual norm) associated with the i-th singular triplet.

Synopsis#

#include "slepcsvd.h" 
PetscErrorCode SVDComputeError(SVD svd,PetscInt i,SVDErrorType type,PetscReal *error)

Collective

Input Parameters#

  • svd - the singular value solver context

  • i - the solution index

  • type - the type of error to compute, see SVDErrorType

Output Parameter#

  • error - the error

Notes#

The error can be computed in various ways, all of them based on the residual norm obtained as \(\sqrt{\eta_1^2+\eta_2^2}\) with \(\eta_1 = \|Av-\sigma u\|_2\) and \(\eta_2 = \|A^*u-\sigma v\|_2\), where \((\sigma,u,v)\) is the approximate singular triplet.

In the case of the GSVD, the two components of the residual norm are \(\eta_1 = \|s^2 A^*u-cB^*Bx\|_2\) and \(\eta_2 = ||c^2 B^*v-sA^*Ax||_2\), where \((\sigma,u,v,x)\) is the approximate generalized singular quadruple, with \(\sigma=c/s\).

See Also#

SVD: Singular Value Decomposition, SVDErrorType, SVDSolve()

Level#

beginner

Location#

src/svd/interface/svdsolve.c

Examples#

src/svd/tutorials/ex45.c
src/svd/tutorials/ex15.c
src/svd/tutorials/ex15f.F90


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