KSPCGSetObjectiveTarget#

Sets the target value for the CG quadratic model

Synopsis#

#include "petscksp.h" 
PetscErrorCode KSPCGSetObjectiveTarget(KSP ksp, PetscReal obj)

Logically Collective

Input Parameters#

  • ksp - the iterative context

  • obj - the objective value (0 is the default)

Notes#

The KSPSolve() will stop when the current value of the objective function \(\frac{1}{2} x^H_k A x_k - b^H x_k\) is smaller than obj if obj is negative. Otherwise the test is ignored. This is used for example inside SNESNEWTONTR.

See Also#

KSP: Linear System Solvers, KSP, KSPCG, KSPNASH, KSPSTCG, KSPGLTR, SNESNEWTONTR

Level#

advanced

Location#

src/ksp/ksp/impls/cg/cgtype.c

Implementations#

KSPCGSetObjectiveTarget_CG() in src/ksp/ksp/impls/cg/cg.c


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