Actual source code: tr.c

petsc-3.5.4 2015-05-23
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:   PetscReal           rho,fnorm,gnorm,gpnorm,xnorm=0,delta,nrm,ynorm,norm1;
 93:   PetscScalar         cnorm;
 94:   KSP                 ksp;
 95:   SNESConvergedReason reason = SNES_CONVERGED_ITERATING;
 96:   PetscBool           conv   = PETSC_FALSE,breakout = PETSC_FALSE;
 97:   PetscBool           domainerror;

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

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

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

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

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

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

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

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

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

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

163:     snes->linear_its += lits;

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

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

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

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

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

258: PetscErrorCode SNESReset_NEWTONTR(SNES snes)
259: {

262:   return(0);
263: }

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

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

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

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

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

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

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

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

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

337:    Level: intermediate

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

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

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

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

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