#include "slepcsvd.h" PetscErrorCode SVDSetConvergenceTestFunction(SVD svd,PetscErrorCode (*func)(SVD,PetscReal,PetscReal,PetscReal*,void*),void* ctx,PetscErrorCode (*destroy)(void*))Logically Collective on svd
svd | - singular value solver context obtained from SVDCreate() | |
func | - a pointer to the convergence test function | |
ctx | - context for private data for the convergence routine (may be null) | |
destroy | - a routine for destroying the context (may be null) |
func(SVD svd,PetscReal sigma,PetscReal res,PetscReal *errest,void *ctx)
svd | - singular value solver context obtained from SVDCreate() | |
sigma | - computed singular value | |
res | - residual norm associated to the singular triplet | |
errest | - (output) computed error estimate | |
ctx | - optional context, as set by SVDSetConvergenceTestFunction() |
Location: src/svd/interface/svdopts.c
Index of all SVD routines
Table of Contents for all manual pages
Index of all manual pages