Actual source code: rich.c

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

  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:   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++) {

 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;
106:       PetscInfo1(ksp,"Self-scale factor %g\n",(double)PetscRealPart(scale));
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:       }

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

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

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

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

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

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

209: PetscErrorCode KSPDestroy_Richardson(KSP ksp)
210: {

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

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

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

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

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

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})

253:           Here B is the application of the preconditioner

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

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

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

266:     Supports only left preconditioning

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

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

277: M*/

281: PETSC_EXTERN PetscErrorCode KSPCreate_Richardson(KSP ksp)
282: {
284:   KSP_Richardson *richardsonP;

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

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

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

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

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