Actual source code: itres.c

petsc-3.10.5 2019-03-28
Report Typos and Errors

  2:  #include <petsc/private/kspimpl.h>

  4: /*@
  5:    KSPInitialResidual - Computes the residual. Either b - A*C*u = b - A*x with right
  6:      preconditioning or C*(b - A*x) with left preconditioning; that later
  7:      residual is often called the "preconditioned residual".

  9:    Collective on KSP

 11:    Input Parameters:
 12: +  vsoln    - solution to use in computing residual
 13: .  vt1, vt2 - temporary work vectors
 14: -  vb       - right-hand-side vector

 16:    Output Parameters:
 17: .  vres     - calculated residual

 19:    Notes:
 20:    This routine assumes that an iterative method, designed for
 21: $     A x = b
 22:    will be used with a preconditioner, C, such that the actual problem is either
 23: $     AC u = b (right preconditioning) or
 24: $     CA x = Cb (left preconditioning).
 25:    This means that the calculated residual will be scaled and/or preconditioned;
 26:    the true residual
 27: $     b-Ax
 28:    is returned in the vt2 temporary.

 30:    Level: developer

 32: .keywords: KSP, residual

 34: .seealso:  KSPMonitor()
 35: @*/

 37: PetscErrorCode  KSPInitialResidual(KSP ksp,Vec vsoln,Vec vt1,Vec vt2,Vec vres,Vec vb)
 38: {
 39:   Mat            Amat,Pmat;

 47:   if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
 48:   PCGetOperators(ksp->pc,&Amat,&Pmat);
 49:   if (!ksp->guess_zero) {
 50:     /* skip right scaling since current guess already has it */
 51:     KSP_MatMult(ksp,Amat,vsoln,vt1);
 52:     VecCopy(vb,vt2);
 53:     VecAXPY(vt2,-1.0,vt1);
 54:     if (ksp->pc_side == PC_RIGHT) {
 55:       PCDiagonalScaleLeft(ksp->pc,vt2,vres);
 56:     } else if (ksp->pc_side == PC_LEFT) {
 57:         KSP_PCApply(ksp,vt2,vres);
 58:         PCDiagonalScaleLeft(ksp->pc,vres,vres);
 59:     } else if (ksp->pc_side == PC_SYMMETRIC) {
 60:       PCApplySymmetricLeft(ksp->pc,vt2,vres);
 61:     } else SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP, "Invalid preconditioning side %d", (int)ksp->pc_side);
 62:   } else {
 63:     VecCopy(vb,vt2);
 64:     if (ksp->pc_side == PC_RIGHT) {
 65:       PCDiagonalScaleLeft(ksp->pc,vb,vres);
 66:     } else if (ksp->pc_side == PC_LEFT) {
 67:       KSP_PCApply(ksp,vb,vres);
 68:       PCDiagonalScaleLeft(ksp->pc,vres,vres);
 69:     } else if (ksp->pc_side == PC_SYMMETRIC) {
 70:       PCApplySymmetricLeft(ksp->pc, vb, vres);
 71:     } else SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP, "Invalid preconditioning side %d", (int)ksp->pc_side);
 72:   }
 73:   return(0);
 74: }

 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: KSPSetPCSide()
101: @*/
102: PetscErrorCode  KSPUnwindPreconditioner(KSP ksp,Vec vsoln,Vec vt1)
103: {

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