Actual source code: itres.c

petsc-3.6.4 2016-04-12
Report Typos and Errors
  2: #include <petsc/private/kspimpl.h>   /*I "petscksp.h" I*/

  6: /*@
  7:    KSPInitialResidual - Computes the residual. Either b - A*C*u = b - A*x with right
  8:      preconditioning or C*(b - 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:   Mat            Amat,Pmat;

 49:   if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
 50:   PCGetOperators(ksp->pc,&Amat,&Pmat);
 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 SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP, "Invalid preconditioning side %d", (int)ksp->pc_side);
 68:   }
 69:   return(0);
 70: }

 74: /*@
 75:    KSPUnwindPreconditioner - Unwinds the preconditioning in the solution. That is,
 76:      takes solution to the preconditioned problem and gets the solution to the
 77:      original problem from it.

 79:    Collective on KSP

 81:    Input Parameters:
 82: +  ksp  - iterative context
 83: .  vsoln - solution vector
 84: -  vt1   - temporary work vector

 86:    Output Parameter:
 87: .  vsoln - contains solution on output

 89:    Notes:
 90:    If preconditioning either symmetrically or on the right, this routine solves
 91:    for the correction to the unpreconditioned problem.  If preconditioning on
 92:    the left, nothing is done.

 94:    Level: advanced

 96: .keywords: KSP, unwind, preconditioner

 98: .seealso: KSPSetPCSide()
 99: @*/
100: PetscErrorCode  KSPUnwindPreconditioner(KSP ksp,Vec vsoln,Vec vt1)
101: {

107:   if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
108:   if (ksp->pc_side == PC_RIGHT) {
109:     KSP_PCApply(ksp,vsoln,vt1);
110:     PCDiagonalScaleRight(ksp->pc,vt1,vsoln);
111:   } else if (ksp->pc_side == PC_SYMMETRIC) {
112:     PCApplySymmetricRight(ksp->pc,vsoln,vt1);
113:     VecCopy(vt1,vsoln);
114:   } else {
115:     PCDiagonalScaleRight(ksp->pc,vsoln,vsoln);
116:   }
117:   return(0);
118: }