petsc-3.12.5 2020-03-29
PCFieldSplitSetGKBNu
Sets the scalar value nu >= 0 in the transformation H = A00 + nu*A01*A01' of the (1,1) block in the Golub-Kahan bidiagonalization preconditioner.
Synopsis
#include "petscpc.h"
PetscErrorCode PCFieldSplitSetGKBNu(PC pc,PetscReal nu)
Collective on PC
Notes
This shift is in general done to obtain better convergence properties for the outer loop of the algorithm. This is often achieved by chosing nu sufficiently big. However,
if nu is chosen too big, the matrix H might be badly conditioned and the solution of the linear system Hx = b in the inner loop gets difficult. It is therefore
necessary to find a good balance in between the convergence of the inner and outer loop.
For nu = 0, no shift is done. In this case A00 has to be positive definite. The matrix N in [Ar13] is then chosen as identity.
[Ar13] Generalized Golub-Kahan bidiagonalization and stopping criteria, SIAM J. Matrix Anal. Appl., Vol. 34, No. 2, pp. 571-592, 2013.
Input Parameters
| pc | - the preconditioner context
|
| nu | - the shift parameter
|
Options Database
-pc_fieldsplit_gkb_nu -default is 1
See Also
PCFIELDSPLIT, PCFieldSplitSetGKBDelay(), PCFieldSplitSetGKBTol(), PCFieldSplitSetGKBMaxit()
Level
intermediate
Location
src/ksp/pc/impls/fieldsplit/fieldsplit.c
Implementations
PCFieldSplitSetGKBNu_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