slepc-3.18.3 2023-03-24
PEPSetConvergenceTestFunction
Sets a function to compute the error estimate used in the convergence test.
Synopsis
#include "slepcpep.h"
PetscErrorCode PEPSetConvergenceTestFunction(PEP pep,PetscErrorCode (*func)(PEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*),void* ctx,PetscErrorCode (*destroy)(void*))
Logically Collective on pep
Input Parameters
| pep | - eigensolver context obtained from PEPCreate()
|
| 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(PEP pep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
| pep | - eigensolver context obtained from PEPCreate()
|
| eigr | - real part of the eigenvalue
|
| eigi | - imaginary part of the eigenvalue
|
| res | - residual norm associated to the eigenpair
|
| errest | - (output) computed error estimate
|
| ctx | - optional context, as set by PEPSetConvergenceTestFunction()
|
Note
If the error estimate returned by the convergence test function is less than
the tolerance, then the eigenvalue is accepted as converged.
See Also
PEPSetConvergenceTest(), PEPSetTolerances()
Level
advanced
Location
src/pep/interface/pepopts.c
Index of all PEP routines
Table of Contents for all manual pages
Index of all manual pages