Actual source code: rich.c

petsc-master 2015-01-29
Report Typos and Errors
  2: /*
  3:             This implements Richardson Iteration.
  4: */
  5: #include <petsc-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:     KSPSetWorkVecs(ksp,4);
 18:   } else {
 19:     KSPSetWorkVecs(ksp,2);
 20:   }
 21:   return(0);
 22: }

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

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

 42:   ksp->its = 0;

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

 64:   /* if user has provided fast Richardson code use that */
 65:   PCApplyRichardsonExists(ksp->pc,&exists);
 66:   if (exists && richardsonP->scale == 1.0 && !ksp->numbermonitors && !ksp->transpose_solve & !ksp->nullsp) {
 67:     PCRichardsonConvergedReason reason;
 68:     PCApplyRichardson(ksp->pc,b,x,r,ksp->rtol,ksp->abstol,ksp->divtol,maxit,ksp->guess_zero,&ksp->its,&reason);
 69:     ksp->reason = (KSPConvergedReason)reason;
 70:     return(0);
 71:   } else if (exists && !ksp->numbermonitors && !ksp->transpose_solve & !ksp->nullsp) {
 72:     PetscInfo(ksp,"KSPSolve_Richardson: Warning, skipping optimized PCApplyRichardson() because scale factor is not 1.0\n");
 73:   }

 75:   scale = richardsonP->scale;

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

 84:   ksp->its = 0;
 85:   if (richardsonP->selfscale) {
 86:     KSP_PCApply(ksp,r,z);         /*   z <- B r          */
 87:     for (i=0; i<maxit; i++) {

 89:       if (ksp->normtype == KSP_NORM_UNPRECONDITIONED) {
 90:         VecNorm(r,NORM_2,&rnorm); /*   rnorm <- r'*r     */
 91:         KSPMonitor(ksp,i,rnorm);
 92:         ksp->rnorm = rnorm;
 93:         KSPLogResidualHistory(ksp,rnorm);
 94:         (*ksp->converged)(ksp,i,rnorm,&ksp->reason,ksp->cnvP);
 95:         if (ksp->reason) break;
 96:       } else if (ksp->normtype == KSP_NORM_PRECONDITIONED) {
 97:         VecNorm(z,NORM_2,&rnorm); /*   rnorm <- z'*z     */
 98:         KSPMonitor(ksp,i,rnorm);
 99:         ksp->rnorm = rnorm;
100:         KSPLogResidualHistory(ksp,rnorm);
101:         (*ksp->converged)(ksp,i,rnorm,&ksp->reason,ksp->cnvP);
102:         if (ksp->reason) break;
103:       }
104:       KSP_PCApplyBAorAB(ksp,z,y,w); /* y = BAz = BABr */
105:       VecDotNorm2(z,y,&rdot,&abr);   /*   rdot = (Br)^T(BABR); abr = (BABr)^T (BABr) */
106:       scale = rdot/abr;
107:       PetscInfo1(ksp,"Self-scale factor %g\n",(double)PetscRealPart(scale));
108:       VecAXPY(x,scale,z);   /*   x  <- x + scale z */
109:       VecAXPY(r,-scale,w);  /*  r <- r - scale*Az */
110:       VecAXPY(z,-scale,y);  /*  z <- z - scale*y */
111:       ksp->its++;
112:     }
113:   } else {
114:     for (i=0; i<maxit; i++) {

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

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

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

136:       VecAXPY(x,scale,z);    /*   x  <- x + scale z */
137:       ksp->its++;

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

171: PetscErrorCode KSPView_Richardson(KSP ksp,PetscViewer viewer)
172: {
173:   KSP_Richardson *richardsonP = (KSP_Richardson*)ksp->data;
175:   PetscBool      iascii;

178:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
179:   if (iascii) {
180:     if (richardsonP->selfscale) {
181:       PetscViewerASCIIPrintf(viewer,"  Richardson: using self-scale best computed damping factor\n");
182:     } else {
183:       PetscViewerASCIIPrintf(viewer,"  Richardson: damping factor=%g\n",(double)richardsonP->scale);
184:     }
185:   }
186:   return(0);
187: }

191: PetscErrorCode KSPSetFromOptions_Richardson(PetscOptions *PetscOptionsObject,KSP ksp)
192: {
193:   KSP_Richardson *rich = (KSP_Richardson*)ksp->data;
195:   PetscReal      tmp;
196:   PetscBool      flg,flg2;

199:   PetscOptionsHead(PetscOptionsObject,"KSP Richardson Options");
200:   PetscOptionsReal("-ksp_richardson_scale","damping factor","KSPRichardsonSetScale",rich->scale,&tmp,&flg);
201:   if (flg) { KSPRichardsonSetScale(ksp,tmp); }
202:   PetscOptionsBool("-ksp_richardson_self_scale","dynamically determine optimal damping factor","KSPRichardsonSetSelfScale",rich->selfscale,&flg2,&flg);
203:   if (flg) { KSPRichardsonSetSelfScale(ksp,flg2); }
204:   PetscOptionsTail();
205:   return(0);
206: }

210: PetscErrorCode KSPDestroy_Richardson(KSP ksp)
211: {

215:   PetscObjectComposeFunction((PetscObject)ksp,"KSPRichardsonSetScale_C",NULL);
216:   KSPDestroyDefault(ksp);
217:   return(0);
218: }

222: static PetscErrorCode  KSPRichardsonSetScale_Richardson(KSP ksp,PetscReal scale)
223: {
224:   KSP_Richardson *richardsonP;

227:   richardsonP        = (KSP_Richardson*)ksp->data;
228:   richardsonP->scale = scale;
229:   return(0);
230: }

234: static PetscErrorCode  KSPRichardsonSetSelfScale_Richardson(KSP ksp,PetscBool selfscale)
235: {
236:   KSP_Richardson *richardsonP;

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

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

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

250:    Level: beginner

252:    Notes: x^{n+1} = x^{n} + scale*B(b - A x^{n})

254:           Here B is the application of the preconditioner

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

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

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

267:     Supports only left preconditioning

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

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

278: M*/

282: PETSC_EXTERN PetscErrorCode KSPCreate_Richardson(KSP ksp)
283: {
285:   KSP_Richardson *richardsonP;

288:   PetscNewLog(ksp,&richardsonP);
289:   ksp->data = (void*)richardsonP;

291:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);
292:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,2);

294:   ksp->ops->setup          = KSPSetUp_Richardson;
295:   ksp->ops->solve          = KSPSolve_Richardson;
296:   ksp->ops->destroy        = KSPDestroy_Richardson;
297:   ksp->ops->buildsolution  = KSPBuildSolutionDefault;
298:   ksp->ops->buildresidual  = KSPBuildResidualDefault;
299:   ksp->ops->view           = KSPView_Richardson;
300:   ksp->ops->setfromoptions = KSPSetFromOptions_Richardson;

302:   PetscObjectComposeFunction((PetscObject)ksp,"KSPRichardsonSetScale_C",KSPRichardsonSetScale_Richardson);
303:   PetscObjectComposeFunction((PetscObject)ksp,"KSPRichardsonSetSelfScale_C",KSPRichardsonSetSelfScale_Richardson);

305:   richardsonP->scale = 1.0;
306:   return(0);
307: }