Actual source code: tr.c

petsc-3.7.3 2016-08-01
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;

 99:   if (snes->xl || snes->xu || snes->ops->computevariablebounds) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);

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:   } else snes->vec_func_init_set = PETSC_FALSE;

116:   VecNorm(F,NORM_2,&fnorm);             /* fnorm <- || F || */
117:   SNESCheckFunctionNorm(snes,fnorm);
118:   VecNorm(X,NORM_2,&xnorm);             /* fnorm <- || F || */
119:   PetscObjectSAWsTakeAccess((PetscObject)snes);
120:   snes->norm = fnorm;
121:   PetscObjectSAWsGrantAccess((PetscObject)snes);
122:   delta      = xnorm ? neP->delta0*xnorm : neP->delta0;
123:   neP->delta = delta;
124:   SNESLogConvergenceHistory(snes,fnorm,0);
125:   SNESMonitor(snes,0,fnorm);

127:   /* test convergence */
128:   (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
129:   if (snes->reason) return(0);

131:   /* Set the stopping criteria to use the More' trick. */
132:   PetscOptionsGetBool(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_tr_ksp_regular_convergence_test",&conv,NULL);
133:   if (!conv) {
134:     SNES_TR_KSPConverged_Ctx *ctx;
135:     SNESGetKSP(snes,&ksp);
136:     PetscNew(&ctx);
137:     ctx->snes = snes;
138:     KSPConvergedDefaultCreate(&ctx->ctx);
139:     KSPSetConvergenceTest(ksp,SNES_TR_KSPConverged_Private,ctx,SNES_TR_KSPConverged_Destroy);
140:     PetscInfo(snes,"Using Krylov convergence test SNES_TR_KSPConverged_Private\n");
141:   }

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

145:     /* Call general purpose update function */
146:     if (snes->ops->update) {
147:       (*snes->ops->update)(snes, snes->iter);
148:     }

150:     /* Solve J Y = F, where J is Jacobian matrix */
151:     SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);
152:     KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);
153:     KSPSolve(snes->ksp,F,Ytmp);
154:     KSPGetIterationNumber(snes->ksp,&lits);

156:     snes->linear_its += lits;

158:     PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);
159:     VecNorm(Ytmp,NORM_2,&nrm);
160:     norm1 = nrm;
161:     while (1) {
162:       VecCopy(Ytmp,Y);
163:       nrm  = norm1;

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

186:       /* Update size of trust region */
187:       if      (rho < neP->mu)  delta *= neP->delta1;
188:       else if (rho < neP->eta) delta *= neP->delta2;
189:       else                     delta *= neP->delta3;
190:       PetscInfo3(snes,"fnorm=%g, gnorm=%g, ynorm=%g\n",(double)fnorm,(double)gnorm,(double)ynorm);
191:       PetscInfo3(snes,"gpred=%g, rho=%g, delta=%g\n",(double)gpnorm,(double)rho,(double)delta);

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

244:   SNESSetWorkVecs(snes,3);
245:   SNESSetUpMatrices(snes);
246:   return(0);
247: }

251: PetscErrorCode SNESReset_NEWTONTR(SNES snes)
252: {

255:   return(0);
256: }

260: static PetscErrorCode SNESDestroy_NEWTONTR(SNES snes)
261: {

265:   SNESReset_NEWTONTR(snes);
266:   PetscFree(snes->data);
267:   return(0);
268: }
269: /*------------------------------------------------------------*/

273: static PetscErrorCode SNESSetFromOptions_NEWTONTR(PetscOptionItems *PetscOptionsObject,SNES snes)
274: {
275:   SNES_NEWTONTR  *ctx = (SNES_NEWTONTR*)snes->data;

279:   PetscOptionsHead(PetscOptionsObject,"SNES trust region options for nonlinear equations");
280:   PetscOptionsReal("-snes_trtol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,NULL);
281:   PetscOptionsReal("-snes_tr_mu","mu","None",ctx->mu,&ctx->mu,NULL);
282:   PetscOptionsReal("-snes_tr_eta","eta","None",ctx->eta,&ctx->eta,NULL);
283:   PetscOptionsReal("-snes_tr_sigma","sigma","None",ctx->sigma,&ctx->sigma,NULL);
284:   PetscOptionsReal("-snes_tr_delta0","delta0","None",ctx->delta0,&ctx->delta0,NULL);
285:   PetscOptionsReal("-snes_tr_delta1","delta1","None",ctx->delta1,&ctx->delta1,NULL);
286:   PetscOptionsReal("-snes_tr_delta2","delta2","None",ctx->delta2,&ctx->delta2,NULL);
287:   PetscOptionsReal("-snes_tr_delta3","delta3","None",ctx->delta3,&ctx->delta3,NULL);
288:   PetscOptionsTail();
289:   return(0);
290: }

294: static PetscErrorCode SNESView_NEWTONTR(SNES snes,PetscViewer viewer)
295: {
296:   SNES_NEWTONTR  *tr = (SNES_NEWTONTR*)snes->data;
298:   PetscBool      iascii;

301:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
302:   if (iascii) {
303:     PetscViewerASCIIPrintf(viewer,"  mu=%g, eta=%g, sigma=%g\n",(double)tr->mu,(double)tr->eta,(double)tr->sigma);
304:     PetscViewerASCIIPrintf(viewer,"  delta0=%g, delta1=%g, delta2=%g, delta3=%g\n",(double)tr->delta0,(double)tr->delta1,(double)tr->delta2,(double)tr->delta3);
305:   }
306:   return(0);
307: }
308: /* ------------------------------------------------------------ */
309: /*MC
310:       SNESNEWTONTR - Newton based nonlinear solver that uses a trust region

312:    Options Database:
313: +    -snes_trtol <tol> Trust region tolerance
314: .    -snes_tr_mu <mu>
315: .    -snes_tr_eta <eta>
316: .    -snes_tr_sigma <sigma>
317: .    -snes_tr_delta0 <delta0> -  initial size of the trust region is delta0*norm2(x)
318: .    -snes_tr_delta1 <delta1>
319: .    -snes_tr_delta2 <delta2>
320: -    -snes_tr_delta3 <delta3>

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

326:    Level: intermediate

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

330: M*/
333: PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTR(SNES snes)
334: {
335:   SNES_NEWTONTR  *neP;

339:   snes->ops->setup          = SNESSetUp_NEWTONTR;
340:   snes->ops->solve          = SNESSolve_NEWTONTR;
341:   snes->ops->destroy        = SNESDestroy_NEWTONTR;
342:   snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTR;
343:   snes->ops->view           = SNESView_NEWTONTR;
344:   snes->ops->reset          = SNESReset_NEWTONTR;

346:   snes->usesksp = PETSC_TRUE;
347:   snes->usespc  = PETSC_FALSE;

349:   PetscNewLog(snes,&neP);
350:   snes->data  = (void*)neP;
351:   neP->mu     = 0.25;
352:   neP->eta    = 0.75;
353:   neP->delta  = 0.0;
354:   neP->delta0 = 0.2;
355:   neP->delta1 = 0.3;
356:   neP->delta2 = 0.75;
357:   neP->delta3 = 2.0;
358:   neP->sigma  = 0.0001;
359:   neP->itflag = PETSC_FALSE;
360:   neP->rnorm0 = 0.0;
361:   neP->ttol   = 0.0;
362:   return(0);
363: }