PCFieldSplitSetGKBNu#
Sets the scalar value nu >= 0 in the transformation H = A00 + nuA01A01’ of the (1,1) block in the Golub-Kahan bidiagonalization preconditioner [Ari13] in PCFIELDSPLIT
Synopsis#
#include "petscpc.h"
PetscErrorCode PCFieldSplitSetGKBNu(PC pc, PetscReal nu)
Collective
Input Parameters#
pc - the preconditioner context
nu - the shift parameter
Options Database Key#
-pc_fieldsplit_gkb_nu
- default is 1
Notes#
This shift is in general done to obtain better convergence properties for the outer loop of the algorithm. This is often achieved by choosing nu
sufficiently large. However,
if nu
is chosen too large, the matrix H might be badly conditioned and the solution of the linear system \(Hx = b\) in the inner loop becomes 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 [Ari13] is then chosen as identity.
References#
See Also#
Solving Block Matrices, PC
, PCFIELDSPLIT
, PCFieldSplitSetGKBDelay()
, PCFieldSplitSetGKBTol()
, PCFieldSplitSetGKBMaxit()
Level#
intermediate
Location#
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