Actual source code: tr.c

petsc-dev 2014-02-02
Report Typos and Errors
  2: #include <../src/snes/impls/tr/trimpl.h>                /*I   "petscsnes.h"   I*/

  4: typedef struct {
  5:   void *ctx;
  6:   SNES snes;
  7: } SNES_TR_KSPConverged_Ctx;

  9: /*
 10:    This convergence test determines if the two norm of the
 11:    solution lies outside the trust region, if so it halts.
 12: */
 15: PetscErrorCode SNES_TR_KSPConverged_Private(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason *reason,void *cctx)
 16: {
 17:   SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx*)cctx;
 18:   SNES                     snes = ctx->snes;
 19:   SNES_NEWTONTR            *neP = (SNES_NEWTONTR*)snes->data;
 20:   Vec                      x;
 21:   PetscReal                nrm;
 22:   PetscErrorCode           ierr;

 25:   KSPConvergedDefault(ksp,n,rnorm,reason,ctx->ctx);
 26:   if (*reason) {
 27:     PetscInfo2(snes,"default convergence test KSP iterations=%D, rnorm=%g\n",n,(double)rnorm);
 28:   }
 29:   /* Determine norm of solution */
 30:   KSPBuildSolution(ksp,0,&x);
 31:   VecNorm(x,NORM_2,&nrm);
 32:   if (nrm >= neP->delta) {
 33:     PetscInfo2(snes,"Ending linear iteration early, delta=%g, length=%g\n",(double)neP->delta,(double)nrm);
 34:     *reason = KSP_CONVERGED_STEP_LENGTH;
 35:   }
 36:   return(0);
 37: }

 41: PetscErrorCode SNES_TR_KSPConverged_Destroy(void *cctx)
 42: {
 43:   SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx*)cctx;
 44:   PetscErrorCode           ierr;

 47:   KSPConvergedDefaultDestroy(ctx->ctx);
 48:   PetscFree(ctx);
 49:   return(0);
 50: }

 52: /* ---------------------------------------------------------------- */
 55: /*
 56:    SNES_TR_Converged_Private -test convergence JUST for
 57:    the trust region tolerance.

 59: */
 60: static PetscErrorCode SNES_TR_Converged_Private(SNES snes,PetscInt it,PetscReal xnorm,PetscReal pnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
 61: {
 62:   SNES_NEWTONTR  *neP = (SNES_NEWTONTR*)snes->data;

 66:   *reason = SNES_CONVERGED_ITERATING;
 67:   if (neP->delta < xnorm * snes->deltatol) {
 68:     PetscInfo3(snes,"Converged due to trust region param %g<%g*%g\n",(double)neP->delta,(double)xnorm,(double)snes->deltatol);
 69:     *reason = SNES_CONVERGED_TR_DELTA;
 70:   } else if (snes->nfuncs >= snes->max_funcs) {
 71:     PetscInfo1(snes,"Exceeded maximum number of function evaluations: %D\n",snes->max_funcs);
 72:     *reason = SNES_DIVERGED_FUNCTION_COUNT;
 73:   }
 74:   return(0);
 75: }


 78: /*
 79:    SNESSolve_NEWTONTR - Implements Newton's Method with a very simple trust
 80:    region approach for solving systems of nonlinear equations.


 83: */
 86: static PetscErrorCode SNESSolve_NEWTONTR(SNES snes)
 87: {
 88:   SNES_NEWTONTR       *neP = (SNES_NEWTONTR*)snes->data;
 89:   Vec                 X,F,Y,G,Ytmp;
 90:   PetscErrorCode      ierr;
 91:   PetscInt            maxits,i,lits;
 92:   MatStructure        flg = DIFFERENT_NONZERO_PATTERN;
 93:   PetscReal           rho,fnorm,gnorm,gpnorm,xnorm=0,delta,nrm,ynorm,norm1;
 94:   PetscScalar         cnorm;
 95:   KSP                 ksp;
 96:   SNESConvergedReason reason = SNES_CONVERGED_ITERATING;
 97:   PetscBool           conv   = PETSC_FALSE,breakout = PETSC_FALSE;
 98:   PetscBool           domainerror;

101:   maxits = snes->max_its;               /* maximum number of iterations */
102:   X      = snes->vec_sol;               /* solution vector */
103:   F      = snes->vec_func;              /* residual vector */
104:   Y      = snes->work[0];               /* work vectors */
105:   G      = snes->work[1];
106:   Ytmp   = snes->work[2];

108:   PetscObjectSAWsTakeAccess((PetscObject)snes);
109:   snes->iter = 0;
110:   PetscObjectSAWsGrantAccess((PetscObject)snes);

112:   if (!snes->vec_func_init_set) {
113:     SNESComputeFunction(snes,X,F);          /* F(X) */
114:     SNESGetFunctionDomainError(snes, &domainerror);
115:     if (domainerror) {
116:       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
117:       return(0);
118:     }
119:   } else snes->vec_func_init_set = PETSC_FALSE;

121:   VecNorm(F,NORM_2,&fnorm);             /* fnorm <- || F || */
122:   if (PetscIsInfOrNanReal(fnorm)) {
123:     snes->reason = SNES_DIVERGED_FNORM_NAN;
124:     return(0);
125:   }

127:   PetscObjectSAWsTakeAccess((PetscObject)snes);
128:   snes->norm = fnorm;
129:   PetscObjectSAWsGrantAccess((PetscObject)snes);
130:   delta      = neP->delta0*fnorm;
131:   neP->delta = delta;
132:   SNESLogConvergenceHistory(snes,fnorm,0);
133:   SNESMonitor(snes,0,fnorm);

135:   /* test convergence */
136:   (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
137:   if (snes->reason) return(0);

139:   /* Set the stopping criteria to use the More' trick. */
140:   PetscOptionsGetBool(NULL,"-snes_tr_ksp_regular_convergence_test",&conv,NULL);
141:   if (!conv) {
142:     SNES_TR_KSPConverged_Ctx *ctx;
143:     SNESGetKSP(snes,&ksp);
144:     PetscNew(&ctx);
145:     ctx->snes = snes;
146:     KSPConvergedDefaultCreate(&ctx->ctx);
147:     KSPSetConvergenceTest(ksp,SNES_TR_KSPConverged_Private,ctx,SNES_TR_KSPConverged_Destroy);
148:     PetscInfo(snes,"Using Krylov convergence test SNES_TR_KSPConverged_Private\n");
149:   }

151:   for (i=0; i<maxits; i++) {

153:     /* Call general purpose update function */
154:     if (snes->ops->update) {
155:       (*snes->ops->update)(snes, snes->iter);
156:     }

158:     /* Solve J Y = F, where J is Jacobian matrix */
159:     SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);
160:     KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);
161:     KSPSolve(snes->ksp,F,Ytmp);
162:     KSPGetIterationNumber(snes->ksp,&lits);

164:     snes->linear_its += lits;

166:     PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);
167:     VecNorm(Ytmp,NORM_2,&nrm);
168:     norm1 = nrm;
169:     while (1) {
170:       VecCopy(Ytmp,Y);
171:       nrm  = norm1;

173:       /* Scale Y if need be and predict new value of F norm */
174:       if (nrm >= delta) {
175:         nrm    = delta/nrm;
176:         gpnorm = (1.0 - nrm)*fnorm;
177:         cnorm  = nrm;
178:         PetscInfo1(snes,"Scaling direction by %g\n",(double)nrm);
179:         VecScale(Y,cnorm);
180:         nrm    = gpnorm;
181:         ynorm  = delta;
182:       } else {
183:         gpnorm = 0.0;
184:         PetscInfo(snes,"Direction is in Trust Region\n");
185:         ynorm  = nrm;
186:       }
187:       VecAYPX(Y,-1.0,X);            /* Y <- X - Y */
188:       VecCopy(X,snes->vec_sol_update);
189:       SNESComputeFunction(snes,Y,G); /*  F(X) */
190:       VecNorm(G,NORM_2,&gnorm);      /* gnorm <- || g || */
191:       if (fnorm == gpnorm) rho = 0.0;
192:       else rho = (fnorm*fnorm - gnorm*gnorm)/(fnorm*fnorm - gpnorm*gpnorm);

194:       /* Update size of trust region */
195:       if      (rho < neP->mu)  delta *= neP->delta1;
196:       else if (rho < neP->eta) delta *= neP->delta2;
197:       else                     delta *= neP->delta3;
198:       PetscInfo3(snes,"fnorm=%g, gnorm=%g, ynorm=%g\n",(double)fnorm,(double)gnorm,(double)ynorm);
199:       PetscInfo3(snes,"gpred=%g, rho=%g, delta=%g\n",(double)gpnorm,(double)rho,(double)delta);

201:       neP->delta = delta;
202:       if (rho > neP->sigma) break;
203:       PetscInfo(snes,"Trying again in smaller region\n");
204:       /* check to see if progress is hopeless */
205:       neP->itflag = PETSC_FALSE;
206:       SNES_TR_Converged_Private(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);
207:       if (!reason) { (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP); }
208:       if (reason) {
209:         /* We're not progressing, so return with the current iterate */
210:         SNESMonitor(snes,i+1,fnorm);
211:         breakout = PETSC_TRUE;
212:         break;
213:       }
214:       snes->numFailures++;
215:     }
216:     if (!breakout) {
217:       /* Update function and solution vectors */
218:       fnorm = gnorm;
219:       VecCopy(G,F);
220:       VecCopy(Y,X);
221:       /* Monitor convergence */
222:       PetscObjectSAWsTakeAccess((PetscObject)snes);
223:       snes->iter = i+1;
224:       snes->norm = fnorm;
225:       PetscObjectSAWsGrantAccess((PetscObject)snes);
226:       SNESLogConvergenceHistory(snes,snes->norm,lits);
227:       SNESMonitor(snes,snes->iter,snes->norm);
228:       /* Test for convergence, xnorm = || X || */
229:       neP->itflag = PETSC_TRUE;
230:       if (snes->ops->converged != SNESConvergedSkip) { VecNorm(X,NORM_2,&xnorm); }
231:       (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);
232:       if (reason) break;
233:     } else break;
234:   }
235:   if (i == maxits) {
236:     PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);
237:     if (!reason) reason = SNES_DIVERGED_MAX_IT;
238:   }
239:   PetscObjectSAWsTakeAccess((PetscObject)snes);
240:   snes->reason = reason;
241:   PetscObjectSAWsGrantAccess((PetscObject)snes);
242:   return(0);
243: }
244: /*------------------------------------------------------------*/
247: static PetscErrorCode SNESSetUp_NEWTONTR(SNES snes)
248: {

252:   SNESSetWorkVecs(snes,3);
253:   SNESSetUpMatrices(snes);
254:   return(0);
255: }

259: PetscErrorCode SNESReset_NEWTONTR(SNES snes)
260: {

263:   return(0);
264: }

268: static PetscErrorCode SNESDestroy_NEWTONTR(SNES snes)
269: {

273:   SNESReset_NEWTONTR(snes);
274:   PetscFree(snes->data);
275:   return(0);
276: }
277: /*------------------------------------------------------------*/

281: static PetscErrorCode SNESSetFromOptions_NEWTONTR(SNES snes)
282: {
283:   SNES_NEWTONTR  *ctx = (SNES_NEWTONTR*)snes->data;

287:   PetscOptionsHead("SNES trust region options for nonlinear equations");
288:   PetscOptionsReal("-snes_trtol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,0);
289:   PetscOptionsReal("-snes_tr_mu","mu","None",ctx->mu,&ctx->mu,0);
290:   PetscOptionsReal("-snes_tr_eta","eta","None",ctx->eta,&ctx->eta,0);
291:   PetscOptionsReal("-snes_tr_sigma","sigma","None",ctx->sigma,&ctx->sigma,0);
292:   PetscOptionsReal("-snes_tr_delta0","delta0","None",ctx->delta0,&ctx->delta0,0);
293:   PetscOptionsReal("-snes_tr_delta1","delta1","None",ctx->delta1,&ctx->delta1,0);
294:   PetscOptionsReal("-snes_tr_delta2","delta2","None",ctx->delta2,&ctx->delta2,0);
295:   PetscOptionsReal("-snes_tr_delta3","delta3","None",ctx->delta3,&ctx->delta3,0);
296:   PetscOptionsTail();
297:   return(0);
298: }

302: static PetscErrorCode SNESView_NEWTONTR(SNES snes,PetscViewer viewer)
303: {
304:   SNES_NEWTONTR  *tr = (SNES_NEWTONTR*)snes->data;
306:   PetscBool      iascii;

309:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
310:   if (iascii) {
311:     PetscViewerASCIIPrintf(viewer,"  mu=%g, eta=%g, sigma=%g\n",(double)tr->mu,(double)tr->eta,(double)tr->sigma);
312:     PetscViewerASCIIPrintf(viewer,"  delta0=%g, delta1=%g, delta2=%g, delta3=%g\n",(double)tr->delta0,(double)tr->delta1,(double)tr->delta2,(double)tr->delta3);
313:   }
314:   return(0);
315: }
316: /* ------------------------------------------------------------ */
317: /*MC
318:       SNESNEWTONTR - Newton based nonlinear solver that uses a trust region

320:    Options Database:
321: +    -snes_trtol <tol> Trust region tolerance
322: .    -snes_tr_mu <mu>
323: .    -snes_tr_eta <eta>
324: .    -snes_tr_sigma <sigma>
325: .    -snes_tr_delta0 <delta0>
326: .    -snes_tr_delta1 <delta1>
327: .    -snes_tr_delta2 <delta2>
328: -    -snes_tr_delta3 <delta3>

330:    The basic algorithm is taken from "The Minpack Project", by More',
331:    Sorensen, Garbow, Hillstrom, pages 88-111 of "Sources and Development
332:    of Mathematical Software", Wayne Cowell, editor.

334:    This is intended as a model implementation, since it does not
335:    necessarily have many of the bells and whistles of other
336:    implementations.

338:    Level: intermediate

340: .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESSetTrustRegionTolerance()

342: M*/
345: PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTR(SNES snes)
346: {
347:   SNES_NEWTONTR  *neP;

351:   snes->ops->setup          = SNESSetUp_NEWTONTR;
352:   snes->ops->solve          = SNESSolve_NEWTONTR;
353:   snes->ops->destroy        = SNESDestroy_NEWTONTR;
354:   snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTR;
355:   snes->ops->view           = SNESView_NEWTONTR;
356:   snes->ops->reset          = SNESReset_NEWTONTR;

358:   snes->usesksp = PETSC_TRUE;
359:   snes->usespc  = PETSC_FALSE;

361:   PetscNewLog(snes,&neP);
362:   snes->data  = (void*)neP;
363:   neP->mu     = 0.25;
364:   neP->eta    = 0.75;
365:   neP->delta  = 0.0;
366:   neP->delta0 = 0.2;
367:   neP->delta1 = 0.3;
368:   neP->delta2 = 0.75;
369:   neP->delta3 = 2.0;
370:   neP->sigma  = 0.0001;
371:   neP->itflag = PETSC_FALSE;
372:   neP->rnorm0 = 0.0;
373:   neP->ttol   = 0.0;
374:   return(0);
375: }