Actual source code: itres.c
2: #include <private/kspimpl.h> /*I "petscksp.h" I*/
6: /*@
7: KSPInitialResidual - Computes the residual. Either b - A*C*x with right
8: preconditioning or C*b - C*A*x with left preconditioning; that later
9: residual is often called the "preconditioned residual".
11: Collective on KSP
13: Input Parameters:
14: + vsoln - solution to use in computing residual
15: . vt1, vt2 - temporary work vectors
16: - vb - right-hand-side vector
18: Output Parameters:
19: . vres - calculated residual
21: Notes:
22: This routine assumes that an iterative method, designed for
23: $ A x = b
24: will be used with a preconditioner, C, such that the actual problem is either
25: $ AC u = b (right preconditioning) or
26: $ CA x = Cb (left preconditioning).
27: This means that the calculated residual will be scaled and/or preconditioned;
28: the true residual
29: $ b-Ax
30: is returned in the vt2 temporary.
32: Level: developer
34: .keywords: KSP, residual
36: .seealso: KSPMonitor()
37: @*/
39: PetscErrorCode KSPInitialResidual(KSP ksp,Vec vsoln,Vec vt1,Vec vt2,Vec vres,Vec vb)
40: {
41: MatStructure pflag;
42: Mat Amat,Pmat;
50: if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
51: PCGetOperators(ksp->pc,&Amat,&Pmat,&pflag);
52: if (!ksp->guess_zero) {
53: /* skip right scaling since current guess already has it */
54: KSP_MatMult(ksp,Amat,vsoln,vt1);
55: VecCopy(vb,vt2);
56: VecAXPY(vt2,-1.0,vt1);
57: (ksp->pc_side == PC_RIGHT)?(VecCopy(vt2,vres)):(KSP_PCApply(ksp,vt2,vres));
58: PCDiagonalScaleLeft(ksp->pc,vres,vres);
59: } else {
60: VecCopy(vb,vt2);
61: if (ksp->pc_side == PC_RIGHT) {
62: PCDiagonalScaleLeft(ksp->pc,vb,vres);
63: } else if (ksp->pc_side == PC_LEFT) {
64: KSP_PCApply(ksp,vb,vres);
65: PCDiagonalScaleLeft(ksp->pc,vres,vres);
66: } else if (ksp->pc_side == PC_SYMMETRIC) {
67: PCApplySymmetricLeft(ksp->pc, vb, vres);
68: } else {
69: SETERRQ1(((PetscObject)ksp)->comm,PETSC_ERR_SUP, "Invalid preconditioning side %d", (int)ksp->pc_side);
70: }
71: }
72: return(0);
73: }
77: /*@
78: KSPUnwindPreconditioner - Unwinds the preconditioning in the solution. That is,
79: takes solution to the preconditioned problem and gets the solution to the
80: original problem from it.
82: Collective on KSP
84: Input Parameters:
85: + ksp - iterative context
86: . vsoln - solution vector
87: - vt1 - temporary work vector
89: Output Parameter:
90: . vsoln - contains solution on output
92: Notes:
93: If preconditioning either symmetrically or on the right, this routine solves
94: for the correction to the unpreconditioned problem. If preconditioning on
95: the left, nothing is done.
97: Level: advanced
99: .keywords: KSP, unwind, preconditioner
101: .seealso: KSPSetPCSide()
102: @*/
103: PetscErrorCode KSPUnwindPreconditioner(KSP ksp,Vec vsoln,Vec vt1)
104: {
110: if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
111: if (ksp->pc_side == PC_RIGHT) {
112: KSP_PCApply(ksp,vsoln,vt1);
113: PCDiagonalScaleRight(ksp->pc,vt1,vsoln);
114: } else if (ksp->pc_side == PC_SYMMETRIC) {
115: PCApplySymmetricRight(ksp->pc,vsoln,vt1);
116: VecCopy(vt1,vsoln);
117: } else {
118: PCDiagonalScaleRight(ksp->pc,vsoln,vsoln);
119: }
120: return(0);
121: }