Actual source code: itres.c
1: /*$Id: itres.c,v 1.54 2001/08/07 03:03:45 balay Exp $*/
3: #include src/sles/ksp/kspimpl.h
5: #undef __FUNCT__
7: /*@C
8: KSPInitialResidual - Computes the residual.
10: Collective on KSP
12: Input Parameters:
13: + vsoln - solution to use in computing residual
14: . vt1, vt2 - temporary work vectors
15: . vres - calculated residual
16: - vb - right-hand-side vector
18: Notes:
19: This routine assumes that an iterative method, designed for
20: $ A x = b
21: will be used with a preconditioner, C, such that the actual problem is either
22: $ AC u = f (right preconditioning) or
23: $ CA x = Cf (left preconditioning).
25: Level: developer
27: .keywords: KSP, residual
29: .seealso: KSPMonitor()
30: @*/
31: int KSPInitialResidual(KSP ksp,Vec vsoln,Vec vt1,Vec vt2,Vec vres,Vec vb)
32: {
33: PetscScalar mone = -1.0;
34: MatStructure pflag;
35: Mat Amat,Pmat;
36: int ierr;
40: PCGetOperators(ksp->B,&Amat,&Pmat,&pflag);
41: if (!ksp->guess_zero) {
42: /* skip right scaling since current guess already has it */
43: KSP_MatMult(ksp,Amat,vsoln,vt1);
44: VecCopy(vb,vt2);
45: VecAXPY(&mone,vt1,vt2);
46: (ksp->pc_side == PC_RIGHT)?(VecCopy(vt2,vres)):(KSP_PCApply(ksp,ksp->B,vt2,vres));
47: PCDiagonalScaleLeft(ksp->B,vres,vres);
48: } else {
49: if (ksp->pc_side == PC_RIGHT) {
50: PCDiagonalScaleLeft(ksp->B,vb,vres);
51: } else if (ksp->pc_side == PC_LEFT) {
52: KSP_PCApply(ksp,ksp->B,vb,vres);
53: PCDiagonalScaleLeft(ksp->B,vres,vres);
54: } else if (ksp->pc_side == PC_SYMMETRIC) {
55: PCApplySymmetricLeft(ksp->B, vb, vres);
56: } else {
57: SETERRQ1(PETSC_ERR_SUP, "Invalid preconditioning side %d", ksp->pc_side);
58: }
60: }
61: return(0);
62: }
64: #undef __FUNCT__
66: /*@
67: KSPUnwindPreconditioner - Unwinds the preconditioning in the solution. That is,
68: takes solution to the preconditioned problem and gets the solution to the
69: original problem from it.
71: Collective on KSP
73: Input Parameters:
74: + ksp - iterative context
75: . vsoln - solution vector
76: - vt1 - temporary work vector
78: Output Parameter:
79: . vsoln - contains solution on output
81: Notes:
82: If preconditioning either symmetrically or on the right, this routine solves
83: for the correction to the unpreconditioned problem. If preconditioning on
84: the left, nothing is done.
86: Level: advanced
88: .keywords: KSP, unwind, preconditioner
90: .seealso: KSPSetPreconditionerSide()
91: @*/
92: int KSPUnwindPreconditioner(KSP ksp,Vec vsoln,Vec vt1)
93: {
98: if (ksp->pc_side == PC_RIGHT) {
99: KSP_PCApply(ksp,ksp->B,vsoln,vt1);
100: PCDiagonalScaleRight(ksp->B,vt1,vsoln);
101: } else if (ksp->pc_side == PC_SYMMETRIC) {
102: PCApplySymmetricRight(ksp->B,vsoln,vt1);
103: VecCopy(vt1,vsoln);
104: } else {
105: PCDiagonalScaleRight(ksp->B,vsoln,vsoln);
106: }
107: return(0);
108: }