Actual source code: itres.c
1: /*$Id: itres.c,v 1.52 2001/03/23 23:23:29 balay Exp $*/
3: #include src/sles/ksp/kspimpl.h
5: /*@C
6: KSPInitialResidual - Computes the residual.
8: Collective on KSP
10: Input Parameters:
11: + vsoln - solution to use in computing residual
12: . vt1, vt2 - temporary work vectors
13: . vres - calculated residual
14: . vbinvf - the result of binv^{-1} b. If null, don't do it.
15: - vb - right-hand-side vector
17: Notes:
18: This routine assumes that an iterative method, designed for
19: $ A x = b
20: will be used with a preconditioner, C, such that the actual problem is either
21: $ AC u = f (right preconditioning) or
22: $ CA x = Cf (left preconditioning).
24: Level: developer
26: .keywords: KSP, residual
28: .seealso: KSPMonitor()
29: @*/
30: int KSPInitialResidual(KSP ksp,Vec vsoln,Vec vt1,Vec vt2,Vec vres,Vec vbinvf,Vec vb)
31: {
32: Scalar mone = -1.0;
33: MatStructure pflag;
34: Mat Amat,Pmat;
35: int ierr;
39: PCGetOperators(ksp->B,&Amat,&Pmat,&pflag);
40: if (ksp->pc_side == PC_RIGHT) {
41: PCDiagonalScaleLeft(ksp->B,vb,vbinvf);
42: } else if (ksp->pc_side == PC_LEFT) {
43: KSP_PCApply(ksp,ksp->B,vb,vbinvf);
44: PCDiagonalScaleLeft(ksp->B,vbinvf,vbinvf);
45: } else {
46: SETERRQ(PETSC_ERR_SUP,"Only right and left preconditioning are currently supported");
47: }
48: if (!ksp->guess_zero) {
49: if (ksp->pc_side == PC_RIGHT) {
50: KSP_MatMult(ksp,Amat,vsoln,vt1);
51: PCDiagonalScaleLeft(ksp->B,vt1,vt1);
52: } else {
53: /* skip right scaling since current guess already has it */
54: if (!ksp->transpose_solve) {
55: MatMult(Amat,vsoln,vt2);
56: PCApply(ksp->B,vt2,vt1);
57: PCDiagonalScaleLeft(ksp->B,vt1,vt1);
58: } else {
59: MatMultTranspose(Amat,vsoln,vt2);
60: PCApplyTranspose(ksp->B,vt2,vt1);
61: PCDiagonalScaleLeft(ksp->B,vt1,vt1);
62: }
63: }
64: /* This is an extra copy for the right-inverse case */
65: VecCopy(vbinvf,vres);
66: VecAXPY(&mone,vt1,vres);
67: } else {
68: VecCopy(vbinvf,vres);
69: }
70: return(0);
71: }
73: /*@
74: KSPUnwindPreconditioner - Unwinds the preconditioning in the solution. That is,
75: takes solution to the preconditioned problem and gets the solution to the
76: original problem from it.
78: Collective on KSP
80: Input Parameters:
81: + ksp - iterative context
82: . vsoln - solution vector
83: - vt1 - temporary work vector
85: Output Parameter:
86: . vsoln - contains solution on output
88: Notes:
89: If preconditioning either symmetrically or on the right, this routine solves
90: for the correction to the unpreconditioned problem. If preconditioning on
91: the left, nothing is done.
93: Level: advanced
95: .keywords: KSP, unwind, preconditioner
97: .seealso: KSPSetPreconditionerSide()
98: @*/
99: int KSPUnwindPreconditioner(KSP ksp,Vec vsoln,Vec vt1)
100: {
105: if (ksp->pc_side == PC_RIGHT) {
106: KSP_PCApply(ksp,ksp->B,vsoln,vt1);
107: PCDiagonalScaleRight(ksp->B,vt1,vsoln);
108: } else if (ksp->pc_side == PC_SYMMETRIC) {
109: PCApplySymmetricRight(ksp->B,vsoln,vt1);
110: VecCopy(vt1,vsoln);
111: } else {
112: PCDiagonalScaleRight(ksp->B,vsoln,vsoln);
113: }
114: return(0);
115: }