Actual source code: rich.c

petsc-master 2016-06-26
Report Typos and Errors
  2: /*
  3:             This implements Richardson Iteration.
  4: */
  5: #include <../src/ksp/ksp/impls/rich/richardsonimpl.h>     /*I "petscksp.h" I*/

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

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

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

 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:   MatGetNullSpace(Pmat,&nullsp);
 67:   if (exists && richardsonP->scale == 1.0 && !ksp->numbermonitors && !ksp->transpose_solve & !nullsp) {
 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:   } else if (exists && !ksp->numbermonitors && !ksp->transpose_solve & !nullsp) {
 73:     PetscInfo(ksp,"KSPSolve_Richardson: Warning, skipping optimized PCApplyRichardson() because scale factor is not 1.0\n");
 74:   }

 76:   scale = richardsonP->scale;

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

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

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

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

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

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

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

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

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

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

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

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

211: PetscErrorCode KSPDestroy_Richardson(KSP ksp)
212: {

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

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

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

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

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

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

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

251:    Level: beginner

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

255:           Here B is the application of the preconditioner

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

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: .  1. - L. F. Richardson, "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:   Philosophical Transactions of the Royal Society of London. Series A,
274:   Containing Papers of a Mathematical or Physical Character, Vol. 210, 1911 (1911).

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

279: M*/

283: PETSC_EXTERN PetscErrorCode KSPCreate_Richardson(KSP ksp)
284: {
286:   KSP_Richardson *richardsonP;

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

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

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

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

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