PCFieldSplitSetGKBDelay#

Sets the delay in the lower bound error estimate in the generalized Golub-Kahan bidiagonalization [Ari13] in PCFIELDSPLIT preconditioner.

Synopsis#

#include "petscpc.h" 
PetscErrorCode PCFieldSplitSetGKBDelay(PC pc, PetscInt delay)

Collective

Input Parameters#

  • pc - the preconditioner context

  • delay - the delay window in the lower bound estimate

Options Database Key#

  • -pc_fieldsplit_gkb_delay - default is 5

Notes#

The algorithm uses a lower bound estimate of the error in energy norm as stopping criterion. The lower bound of the error \( ||u-u^k||_H \) is expressed as a truncated sum. The error at iteration k can only be measured at iteration (k + delay), and thus the algorithm needs at least (delay + 1) iterations to stop.

For more details on the generalized Golub-Kahan bidiagonalization method and its lower bound stopping criterion, please refer to [Ari13]

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, PC, PCFIELDSPLIT, PCFieldSplitSetGKBNu(), PCFieldSplitSetGKBTol(), PCFieldSplitSetGKBMaxit()

Level#

intermediate

Location#

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

Implementations#

PCFieldSplitSetGKBDelay_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