slepc-3.13.0 2020-03-31
Report Typos and Errors

SVDSetConvergenceTestFunction

Sets a function to compute the error estimate used in the convergence test.

Synopsis

#include "slepcsvd.h" 
PetscErrorCode SVDSetConvergenceTestFunction(SVD svd,PetscErrorCode (*func)(SVD,PetscReal,PetscReal,PetscReal*,void*),void* ctx,PetscErrorCode (*destroy)(void*))
Logically Collective on svd

Input Parameters

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)

Calling Sequence of func

  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()

Note

If the error estimate returned by the convergence test function is less than the tolerance, then the singular value is accepted as converged.

See Also

SVDSetConvergenceTest(), SVDSetTolerances()

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