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