Actual source code: rich.c

  2: /*          
  3:             This implements Richardson Iteration.       
  4: */
  5: #include <private/kspimpl.h>              /*I "petscksp.h" I*/
  6: #include <../src/ksp/ksp/impls/rich/richardsonimpl.h>

 10: PetscErrorCode KSPSetUp_Richardson(KSP ksp)
 11: {
 13:   KSP_Richardson *richardsonP = (KSP_Richardson*)ksp->data;

 16:   if (richardsonP->selfscale) {
 17:     KSPDefaultGetWork(ksp,4);
 18:   } else {
 19:     KSPDefaultGetWork(ksp,2);
 20:   }
 21:   return(0);
 22: }

 26: PetscErrorCode  KSPSolve_Richardson(KSP ksp)
 27: {
 29:   PetscInt       i,maxit;
 30:   MatStructure   pflag;
 31:   PetscReal      rnorm = 0.0;
 32:   PetscScalar    scale,abr,rdot;
 33:   Vec            x,b,r,z,w = PETSC_NULL,y = PETSC_NULL;
 34:   PetscInt       xs, ws;
 35:   Mat            Amat,Pmat;
 36:   KSP_Richardson *richardsonP = (KSP_Richardson*)ksp->data;
 37:   PetscBool      exists,diagonalscale;

 40:   PCGetDiagonalScale(ksp->pc,&diagonalscale);
 41:   if (diagonalscale) SETERRQ1(((PetscObject)ksp)->comm,PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);

 43:   ksp->its = 0;

 45:   PCGetOperators(ksp->pc,&Amat,&Pmat,&pflag);
 46:   x       = ksp->vec_sol;
 47:   b       = ksp->vec_rhs;
 48:   VecGetSize(x,&xs);
 49:   VecGetSize(ksp->work[0],&ws);
 50:   if (xs != ws) {
 51:     if (richardsonP->selfscale) {
 52:       KSPDefaultGetWork(ksp,4);
 53:     } else {
 54:       KSPDefaultGetWork(ksp,2);
 55:     }
 56:   }
 57:   r       = ksp->work[0];
 58:   z       = ksp->work[1];
 59:   if (richardsonP->selfscale) {
 60:     w     = ksp->work[2];
 61:     y     = ksp->work[3];
 62:   }
 63:   maxit   = ksp->max_it;

 65:   /* if user has provided fast Richardson code use that */
 66:   PCApplyRichardsonExists(ksp->pc,&exists);
 67:   if (exists && !ksp->numbermonitors && !ksp->transpose_solve) {
 68:     PCRichardsonConvergedReason reason;
 69:     PCApplyRichardson(ksp->pc,b,x,r,ksp->rtol,ksp->abstol,ksp->divtol,maxit,ksp->guess_zero,&ksp->its,&reason);
 70:     ksp->reason = (KSPConvergedReason)reason;
 71:     return(0);
 72:   }

 74:   scale   = richardsonP->scale;

 76:   if (!ksp->guess_zero) {                          /*   r <- b - A x     */
 77:     KSP_MatMult(ksp,Amat,x,r);
 78:     VecAYPX(r,-1.0,b);
 79:   } else {
 80:     VecCopy(b,r);
 81:   }

 83:   ksp->its = 0;
 84:   if (richardsonP->selfscale) {
 85:     KSP_PCApply(ksp,r,z);         /*   z <- B r          */
 86:     for (i=0; i<maxit; i++) {
 87: 
 88:       if (ksp->normtype == KSP_NORM_UNPRECONDITIONED) {
 89:         VecNorm(r,NORM_2,&rnorm); /*   rnorm <- r'*r     */
 90:         KSPMonitor(ksp,i,rnorm);
 91:         ksp->rnorm = rnorm;
 92:         KSPLogResidualHistory(ksp,rnorm);
 93:         (*ksp->converged)(ksp,i,rnorm,&ksp->reason,ksp->cnvP);
 94:         if (ksp->reason) break;
 95:       } else if (ksp->normtype == KSP_NORM_PRECONDITIONED) {
 96:         VecNorm(z,NORM_2,&rnorm); /*   rnorm <- z'*z     */
 97:         KSPMonitor(ksp,i,rnorm);
 98:         ksp->rnorm = rnorm;
 99:         KSPLogResidualHistory(ksp,rnorm);
100:         (*ksp->converged)(ksp,i,rnorm,&ksp->reason,ksp->cnvP);
101:         if (ksp->reason) break;
102:       }
103:       KSP_PCApplyBAorAB(ksp,z,y,w);  /* y = BAz = BABr */
104:       VecDotNorm2(z,y,&rdot,&abr);   /*   rdot = (Br)^T(BABR); abr = (BABr)^T (BABr) */
105:       scale = rdot/abr;

107:       VecAXPY(x,scale,z);    /*   x  <- x + scale z */
108:       VecAXPY(r,-scale,w);   /*  r <- r - scale*Az */
109:       VecAXPY(z,-scale,y);   /*  z <- z - scale*y */
110:       ksp->its++;
111:     }
112:   } else {
113:     for (i=0; i<maxit; i++) {

115:       if (ksp->normtype == KSP_NORM_UNPRECONDITIONED) {
116:         VecNorm(r,NORM_2,&rnorm); /*   rnorm <- r'*r     */
117:         KSPMonitor(ksp,i,rnorm);
118:         ksp->rnorm = rnorm;
119:         KSPLogResidualHistory(ksp,rnorm);
120:         (*ksp->converged)(ksp,i,rnorm,&ksp->reason,ksp->cnvP);
121:         if (ksp->reason) break;
122:       }

124:       KSP_PCApply(ksp,r,z);    /*   z <- B r          */

126:       if (ksp->normtype == KSP_NORM_PRECONDITIONED) {
127:         VecNorm(z,NORM_2,&rnorm); /*   rnorm <- z'*z     */
128:         KSPMonitor(ksp,i,rnorm);
129:         ksp->rnorm = rnorm;
130:         KSPLogResidualHistory(ksp,rnorm);
131:         (*ksp->converged)(ksp,i,rnorm,&ksp->reason,ksp->cnvP);
132:         if (ksp->reason) break;
133:       }
134: 
135:       VecAXPY(x,scale,z);    /*   x  <- x + scale z */
136:       ksp->its++;

138:       KSP_MatMult(ksp,Amat,x,r);      /*   r  <- b - Ax      */
139:       VecAYPX(r,-1.0,b);
140:     }
141:   }
142:   if (!ksp->reason) {
143:     if (ksp->normtype != KSP_NORM_NONE) {
144:       if (ksp->normtype == KSP_NORM_UNPRECONDITIONED){
145:         VecNorm(r,NORM_2,&rnorm);     /*   rnorm <- r'*r     */
146:       } else {
147:         KSP_PCApply(ksp,r,z);   /*   z <- B r          */
148:         VecNorm(z,NORM_2,&rnorm);     /*   rnorm <- z'*z     */
149:       }
150:       ksp->rnorm = rnorm;
151:       KSPLogResidualHistory(ksp,rnorm);
152:       KSPMonitor(ksp,i,rnorm);
153:     }
154:     if (ksp->its >= ksp->max_it) {
155:       if (ksp->normtype != KSP_NORM_NONE) {
156:         (*ksp->converged)(ksp,i,rnorm,&ksp->reason,ksp->cnvP);
157:         if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
158:       } else {
159:         ksp->reason = KSP_CONVERGED_ITS;
160:       }
161:     }
162:   }
163:   return(0);
164: }

168: PetscErrorCode KSPView_Richardson(KSP ksp,PetscViewer viewer)
169: {
170:   KSP_Richardson *richardsonP = (KSP_Richardson*)ksp->data;
172:   PetscBool      iascii;

175:   PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
176:   if (iascii) {
177:     PetscViewerASCIIPrintf(viewer,"  Richardson: damping factor=%G\n",richardsonP->scale);
178:   } else {
179:     SETERRQ1(((PetscObject)ksp)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for KSP Richardson",((PetscObject)viewer)->type_name);
180:   }
181:   return(0);
182: }

186: PetscErrorCode KSPSetFromOptions_Richardson(KSP ksp)
187: {
188:   KSP_Richardson *rich = (KSP_Richardson*)ksp->data;
190:   PetscReal      tmp;
191:   PetscBool      flg,flg2;

194:   PetscOptionsHead("KSP Richardson Options");
195:     PetscOptionsReal("-ksp_richardson_scale","damping factor","KSPRichardsonSetScale",rich->scale,&tmp,&flg);
196:     if (flg) { KSPRichardsonSetScale(ksp,tmp); }
197:     PetscOptionsBool("-ksp_richardson_self_scale","dynamically determine optimal damping factor","KSPRichardsonSetSelfScale",rich->selfscale,&flg2,&flg);
198:     if (flg) { KSPRichardsonSetSelfScale(ksp,flg2); }
199:   PetscOptionsTail();
200:   return(0);
201: }

205: PetscErrorCode KSPDestroy_Richardson(KSP ksp)
206: {

210:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPRichardsonSetScale_C","",PETSC_NULL);
211:   KSPDefaultDestroy(ksp);
212:   return(0);
213: }

218: PetscErrorCode  KSPRichardsonSetScale_Richardson(KSP ksp,PetscReal scale)
219: {
220:   KSP_Richardson *richardsonP;

223:   richardsonP = (KSP_Richardson*)ksp->data;
224:   richardsonP->scale = scale;
225:   return(0);
226: }

232: PetscErrorCode  KSPRichardsonSetSelfScale_Richardson(KSP ksp,PetscBool  selfscale)
233: {
234:   KSP_Richardson *richardsonP;

237:   richardsonP = (KSP_Richardson*)ksp->data;
238:   richardsonP->selfscale = selfscale;
239:   return(0);
240: }

243: /*MC
244:      KSPRICHARDSON - The preconditioned Richardson iterative method

246:    Options Database Keys:
247: .   -ksp_richardson_scale - damping factor on the correction (defaults to 1.0)

249:    Level: beginner

251:    Notes: x^{n+1} = x^{n} + scale*B(b - A x^{n})
252:  
253:           Here B is the application of the preconditioner

255:           This method often (usually) will not converge unless scale is very small. It
256: is described in


259:    Notes: For some preconditioners, currently SOR, the convergence test is skipped to improve speed,
260:     thus it always iterates the maximum number of iterations you've selected. When -ksp_monitor 
261:     (or any other monitor) is turned on, the norm is computed at each iteration and so the convergence test is run unless
262:     you specifically call KSPSetNormType(ksp,KSP_NORM_NONE);

264:          For some preconditioners, currently PCMG and PCHYPRE with BoomerAMG if -ksp_monitor (and also
265:     any other monitor) is not turned on then the convergence test is done by the preconditioner itself and
266:     so the solver may run more or fewer iterations then if -ksp_monitor is selected.

268:     Supports only left preconditioning

270:   References:
271:   "The Approximate Arithmetical Solution by Finite Differences of Physical Problems Involving
272:    Differential Equations, with an Application to the Stresses in a Masonry Dam",
273:   L. F. Richardson, Philosophical Transactions of the Royal Society of London. Series A,
274:   Containing Papers of a Mathematical or Physical Character, Vol. 210, 1911 (1911), pp. 307-357.

276: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP,
277:            KSPRichardsonSetScale()

279: M*/

284: PetscErrorCode  KSPCreate_Richardson(KSP ksp)
285: {
287:   KSP_Richardson *richardsonP;

290:   PetscNewLog(ksp,KSP_Richardson,&richardsonP);
291:   ksp->data                        = (void*)richardsonP;

293:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,2);
294:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,1);

296:   ksp->ops->setup                  = KSPSetUp_Richardson;
297:   ksp->ops->solve                  = KSPSolve_Richardson;
298:   ksp->ops->destroy                = KSPDestroy_Richardson;
299:   ksp->ops->buildsolution          = KSPDefaultBuildSolution;
300:   ksp->ops->buildresidual          = KSPDefaultBuildResidual;
301:   ksp->ops->view                   = KSPView_Richardson;
302:   ksp->ops->setfromoptions         = KSPSetFromOptions_Richardson;

304:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPRichardsonSetScale_C","KSPRichardsonSetScale_Richardson",KSPRichardsonSetScale_Richardson);
305:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPRichardsonSetSelfScale_C","KSPRichardsonSetSelfScale_Richardson",KSPRichardsonSetSelfScale_Richardson);
306:   richardsonP->scale               = 1.0;
307:   return(0);
308: }