KSPLSQRConvergedDefault#

Determines convergence of the KSPLSQR Krylov method.

Synopsis#

#include "petscksp.h" 
PetscErrorCode KSPLSQRConvergedDefault(KSP ksp, PetscInt n, PetscReal rnorm, KSPConvergedReason *reason, void *ctx)

Collective

Input Parameters#

  • ksp - iterative context

  • n - iteration number

  • rnorm - 2-norm residual value (may be estimated)

  • ctx - convergence context which must have been created by KSPConvergedDefaultCreate()

Output Parameter#

  • reason - the convergence reason

Notes#

This is not called directly but rather is passed to KSPSetConvergenceTest(). It is used automatically by KSPLSQR

KSPConvergedDefault() is called first to check for convergence in \(A*x=b\). If that does not determine convergence then checks convergence for the least squares problem, i.e. in min{|b-A*x|}. Possible convergence for the least squares problem (which is based on the residual of the normal equations) are KSP_CONVERGED_RTOL_NORMAL norm and KSP_CONVERGED_ATOL_NORMAL.

KSP_CONVERGED_RTOL_NORMAL is returned if \(||A^T*r|| < rtol * ||A|| * ||r||\). Matrix norm \(||A||\) is iteratively refined estimate, see KSPLSQRGetNorms(). This criterion is largely compatible with that in MATLAB lsqr().

See Also#

KSP: Linear System Solvers, KSPLSQR, KSPSetConvergenceTest(), KSPSetTolerances(), KSPConvergedSkip(), KSPConvergedReason, KSPGetConvergedReason(), KSPConvergedDefaultSetUIRNorm(), KSPConvergedDefaultSetUMIRNorm(), KSPConvergedDefaultCreate(), KSPConvergedDefaultDestroy(), KSPConvergedDefault(), KSPLSQRGetNorms(), KSPLSQRSetExactMatNorm()

Level#

advanced

Location#

src/ksp/ksp/impls/lsqr/lsqr.c


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