SVDComputeResidualNorms

Computes the norms of the residual vectors associated with the i-th computed singular triplet.

Synopsis

#include "slepcsvd.h" 
PetscErrorCode SVDComputeResidualNorms(SVD svd, PetscInt i, PetscReal *norm1, PetscReal *norm2)
Collective on SVD

Input Parameters

svd - the singular value solver context
i - the solution index

Output Parameters

norm1 - the norm ||A*v-sigma*u||_2 where sigma is the singular value, u and v are the left and right singular vectors.
norm2 - the norm ||A^T*u-sigma*v||_2 with the same sigma, u and v

Note

The index i should be a value between 0 and nconv-1 (see SVDGetConverged()). Both output parameters can be PETSC_NULL on input if not needed.

See Also

SVDSolve(), SVDGetConverged(), SVDComputeRelativeError()

Location: src/svd/interface/svdsolve.c
Index of all SVD routines
Table of Contents for all manual pages
Index of all manual pages