#include "slepcsvd.h" PetscErrorCode SVDComputeResidualNorms(SVD svd, PetscInt i, PetscReal *norm1, PetscReal *norm2)Collective on SVD
svd | - the singular value solver context | |
i | - the solution index |
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 |
Location: src/svd/interface/svdsolve.c
Index of all SVD routines
Table of Contents for all manual pages
Index of all manual pages