slepc-main 2024-11-15
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
|
Output Parameter
Notes
The error can be computed in various ways, all of them based on the residual
norm obtained as sqrt(n1^2+n2^2) with n1 = ||A*v-sigma*u||_2 and
n2 = ||A^T*u-sigma*v||_2, where sigma is the singular value, u is the left
singular vector and v is the right singular vector.
In the case of the GSVD, the two components of the residual norm are
n1 = ||s^2*A'*u-c*B'*B*x||_2 and n2 = ||c^2*B'*v-s*A'*A*x||_2, where [u;v]
are the left singular vectors and x is the right singular vector, with
sigma=c/s.
See Also
SVDErrorType, SVDSolve()
Level
beginner
Location
src/svd/interface/svdsolve.c
Examples
src/svd/tutorials/ex15.c
src/svd/tutorials/ex15f.F90
src/svd/tutorials/ex45.c
Index of all SVD routines
Table of Contents for all manual pages
Index of all manual pages