PCFieldSplitSetGKBTol#

Sets the solver tolerance for the generalized Golub-Kahan bidiagonalization preconditioner [Ari13] in PCFIELDSPLIT

Synopsis#

#include "petscpc.h" 
PetscErrorCode PCFieldSplitSetGKBTol(PC pc, PetscReal tolerance)

Collective

Input Parameters#

  • pc - the preconditioner context

  • tolerance - the solver tolerance

Options Database Key#

  • -pc_fieldsplit_gkb_tol - default is 1e-5

Note#

The generalized GKB algorithm [Ari13] uses a lower bound estimate of the error in energy norm as stopping criterion. It stops once the lower bound estimate undershoots the required solver tolerance. Although the actual error might be bigger than this estimate, the stopping criterion is satisfactory in practical cases.

References#

[Ari13] (1,2)

Mario Arioli. Generalized Golub–Kahan bidiagonalization and stopping criteria. SIAM Journal on Matrix Analysis and Applications, 34(2):571–592, 2013.

See Also#

Solving Block Matrices with PCFIELDSPLIT, PC, PCFIELDSPLIT, PCFieldSplitSetGKBDelay(), PCFieldSplitSetGKBNu(), PCFieldSplitSetGKBMaxit()

Level#

intermediate

Location#

src/ksp/pc/impls/fieldsplit/fieldsplit.c

Implementations#

PCFieldSplitSetGKBTol_FieldSplit() in src/ksp/pc/impls/fieldsplit/fieldsplit.c


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