Actual source code: tr.c

petsc-3.4.5 2014-06-29
  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:   KSPDefaultConverged(ksp,n,rnorm,reason,ctx->ctx);
 26:   if (*reason) {
 27:     PetscInfo2(snes,"default convergence test KSP iterations=%D, rnorm=%G\n",n,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",neP->delta,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:   KSPDefaultConvergedDestroy(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",neP->delta,xnorm,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:   PetscObjectAMSTakeAccess((PetscObject)snes);
109:   snes->iter = 0;
110:   PetscObjectAMSGrantAccess((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:   if (!snes->norm_init_set) {
122:     VecNorm(F,NORM_2,&fnorm);             /* fnorm <- || F || */
123:     if (PetscIsInfOrNanReal(fnorm)) {
124:       snes->reason = SNES_DIVERGED_FNORM_NAN;
125:       return(0);
126:     }
127:   } else {
128:     fnorm               = snes->norm_init;
129:     snes->norm_init_set = PETSC_FALSE;
130:   }

132:   PetscObjectAMSTakeAccess((PetscObject)snes);
133:   snes->norm = fnorm;
134:   PetscObjectAMSGrantAccess((PetscObject)snes);
135:   delta      = neP->delta0*fnorm;
136:   neP->delta = delta;
137:   SNESLogConvergenceHistory(snes,fnorm,0);
138:   SNESMonitor(snes,0,fnorm);

140:   /* set parameter for default relative tolerance convergence test */
141:   snes->ttol = fnorm*snes->rtol;
142:   /* test convergence */
143:   (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
144:   if (snes->reason) return(0);

146:   /* Set the stopping criteria to use the More' trick. */
147:   PetscOptionsGetBool(NULL,"-snes_tr_ksp_regular_convergence_test",&conv,NULL);
148:   if (!conv) {
149:     SNES_TR_KSPConverged_Ctx *ctx;
150:     SNESGetKSP(snes,&ksp);
151:     PetscNew(SNES_TR_KSPConverged_Ctx,&ctx);
152:     ctx->snes = snes;
153:     KSPDefaultConvergedCreate(&ctx->ctx);
154:     KSPSetConvergenceTest(ksp,SNES_TR_KSPConverged_Private,ctx,SNES_TR_KSPConverged_Destroy);
155:     PetscInfo(snes,"Using Krylov convergence test SNES_TR_KSPConverged_Private\n");
156:   }

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

160:     /* Call general purpose update function */
161:     if (snes->ops->update) {
162:       (*snes->ops->update)(snes, snes->iter);
163:     }

165:     /* Solve J Y = F, where J is Jacobian matrix */
166:     SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);
167:     KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);
168:     KSPSolve(snes->ksp,F,Ytmp);
169:     KSPGetIterationNumber(snes->ksp,&lits);

171:     snes->linear_its += lits;

173:     PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);
174:     VecNorm(Ytmp,NORM_2,&nrm);
175:     norm1 = nrm;
176:     while (1) {
177:       VecCopy(Ytmp,Y);
178:       nrm  = norm1;

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

201:       /* Update size of trust region */
202:       if      (rho < neP->mu)  delta *= neP->delta1;
203:       else if (rho < neP->eta) delta *= neP->delta2;
204:       else                     delta *= neP->delta3;
205:       PetscInfo3(snes,"fnorm=%G, gnorm=%G, ynorm=%G\n",fnorm,gnorm,ynorm);
206:       PetscInfo3(snes,"gpred=%G, rho=%G, delta=%G\n",gpnorm,rho,delta);

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

259:   SNESSetWorkVecs(snes,3);
260:   SNESSetUpMatrices(snes);
261:   return(0);
262: }

266: PetscErrorCode SNESReset_NEWTONTR(SNES snes)
267: {

270:   return(0);
271: }

275: static PetscErrorCode SNESDestroy_NEWTONTR(SNES snes)
276: {

280:   SNESReset_NEWTONTR(snes);
281:   PetscFree(snes->data);
282:   return(0);
283: }
284: /*------------------------------------------------------------*/

288: static PetscErrorCode SNESSetFromOptions_NEWTONTR(SNES snes)
289: {
290:   SNES_NEWTONTR  *ctx = (SNES_NEWTONTR*)snes->data;

294:   PetscOptionsHead("SNES trust region options for nonlinear equations");
295:   PetscOptionsReal("-snes_trtol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,0);
296:   PetscOptionsReal("-snes_tr_mu","mu","None",ctx->mu,&ctx->mu,0);
297:   PetscOptionsReal("-snes_tr_eta","eta","None",ctx->eta,&ctx->eta,0);
298:   PetscOptionsReal("-snes_tr_sigma","sigma","None",ctx->sigma,&ctx->sigma,0);
299:   PetscOptionsReal("-snes_tr_delta0","delta0","None",ctx->delta0,&ctx->delta0,0);
300:   PetscOptionsReal("-snes_tr_delta1","delta1","None",ctx->delta1,&ctx->delta1,0);
301:   PetscOptionsReal("-snes_tr_delta2","delta2","None",ctx->delta2,&ctx->delta2,0);
302:   PetscOptionsReal("-snes_tr_delta3","delta3","None",ctx->delta3,&ctx->delta3,0);
303:   PetscOptionsTail();
304:   return(0);
305: }

309: static PetscErrorCode SNESView_NEWTONTR(SNES snes,PetscViewer viewer)
310: {
311:   SNES_NEWTONTR  *tr = (SNES_NEWTONTR*)snes->data;
313:   PetscBool      iascii;

316:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
317:   if (iascii) {
318:     PetscViewerASCIIPrintf(viewer,"  mu=%G, eta=%G, sigma=%G\n",tr->mu,tr->eta,tr->sigma);
319:     PetscViewerASCIIPrintf(viewer,"  delta0=%G, delta1=%G, delta2=%G, delta3=%G\n",tr->delta0,tr->delta1,tr->delta2,tr->delta3);
320:   }
321:   return(0);
322: }
323: /* ------------------------------------------------------------ */
324: /*MC
325:       SNESNEWTONTR - Newton based nonlinear solver that uses a trust region

327:    Options Database:
328: +    -snes_trtol <tol> Trust region tolerance
329: .    -snes_tr_mu <mu>
330: .    -snes_tr_eta <eta>
331: .    -snes_tr_sigma <sigma>
332: .    -snes_tr_delta0 <delta0>
333: .    -snes_tr_delta1 <delta1>
334: .    -snes_tr_delta2 <delta2>
335: -    -snes_tr_delta3 <delta3>

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

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

345:    Level: intermediate

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

349: M*/
352: PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTR(SNES snes)
353: {
354:   SNES_NEWTONTR  *neP;

358:   snes->ops->setup          = SNESSetUp_NEWTONTR;
359:   snes->ops->solve          = SNESSolve_NEWTONTR;
360:   snes->ops->destroy        = SNESDestroy_NEWTONTR;
361:   snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTR;
362:   snes->ops->view           = SNESView_NEWTONTR;
363:   snes->ops->reset          = SNESReset_NEWTONTR;

365:   snes->usesksp = PETSC_TRUE;
366:   snes->usespc  = PETSC_FALSE;

368:   PetscNewLog(snes,SNES_NEWTONTR,&neP);
369:   snes->data  = (void*)neP;
370:   neP->mu     = 0.25;
371:   neP->eta    = 0.75;
372:   neP->delta  = 0.0;
373:   neP->delta0 = 0.2;
374:   neP->delta1 = 0.3;
375:   neP->delta2 = 0.75;
376:   neP->delta3 = 2.0;
377:   neP->sigma  = 0.0001;
378:   neP->itflag = PETSC_FALSE;
379:   neP->rnorm0 = 0.0;
380:   neP->ttol   = 0.0;
381:   return(0);
382: }