slepc-main 2024-11-09
Report Typos and Errors

DSNEPSetRefine

Sets the tolerance and the number of iterations of Newton iterative refinement for eigenpairs.

Synopsis

#include "slepcds.h" 
PetscErrorCode DSNEPSetRefine(DS ds,PetscReal tol,PetscInt its)
Logically Collective

Input Parameters

ds  - the direct solver context
tol  - the tolerance
its  - the number of iterations

Options Database Key

-ds_nep_refine_tol <tol>  - sets the tolerance
-ds_nep_refine_its <its>  - sets the number of Newton iterations

Notes

Iterative refinement of eigenpairs is currently used only in the contour integral method.

Use PETSC_CURRENT to retain the current value of any of the parameters. Use PETSC_DETERMINE for either argument to assign a default value computed internally.

See Also

DSNEPGetRefine()

Level

advanced

Location

src/sys/classes/ds/impls/nep/dsnep.c

Index of all DS routines
Table of Contents for all manual pages
Index of all manual pages