KSPInitialResidual#
Computes the residual. Either b - ACu = b - Ax with right preconditioning or C(b - A*x) with left preconditioning; the latter residual is often called the “preconditioned residual”.
Synopsis#
#include "petscksp.h"
PetscErrorCode KSPInitialResidual(KSP ksp, Vec vsoln, Vec vt1, Vec vt2, Vec vres, Vec vb)
Collective
Input Parameters#
vsoln - solution to use in computing residual
vt1, vt2 - temporary work vectors
vb - right-hand-side vector
Output Parameter#
vres - calculated residual
Note#
This routine assumes that an iterative method, designed for
A x = b
will be used with a preconditioner, C, such that the actual problem is either
AC u = b (right preconditioning) or
CA x = Cb (left preconditioning).
This means that the calculated residual will be scaled and/or preconditioned; the true residual
b-Ax
is returned in the vt2
temporary.
See Also#
Level#
developer
Location#
Index of all KSP routines
Table of Contents for all manual pages
Index of all manual pages