xref: /petsc/src/snes/impls/tr/tr.c (revision c93683566a9f093a797b09ebeb4ca4f8a209c4a3)
1 
2 #include <../src/snes/impls/tr/trimpl.h>                /*I   "petscsnes.h"   I*/
3 
4 typedef struct {
5   SNES           snes;
6   /*  Information on the regular SNES convergence test; which may have been user provided */
7   PetscErrorCode (*convtest)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
8   PetscErrorCode (*convdestroy)(void*);
9   void           *convctx;
10 } SNES_TR_KSPConverged_Ctx;
11 
12 static PetscErrorCode SNESTR_KSPConverged_Private(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason *reason,void *cctx)
13 {
14   SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx*)cctx;
15   SNES                     snes = ctx->snes;
16   SNES_NEWTONTR            *neP = (SNES_NEWTONTR*)snes->data;
17   Vec                      x;
18   PetscReal                nrm;
19   PetscErrorCode           ierr;
20 
21   PetscFunctionBegin;
22   ierr = (*ctx->convtest)(ksp,n,rnorm,reason,ctx->convctx);CHKERRQ(ierr);
23   if (*reason) {
24     ierr = PetscInfo2(snes,"Default or user provided convergence test KSP iterations=%D, rnorm=%g\n",n,(double)rnorm);CHKERRQ(ierr);
25   }
26   /* Determine norm of solution */
27   ierr = KSPBuildSolution(ksp,0,&x);CHKERRQ(ierr);
28   ierr = VecNorm(x,NORM_2,&nrm);CHKERRQ(ierr);
29   if (nrm >= neP->delta) {
30     ierr    = PetscInfo2(snes,"Ending linear iteration early, delta=%g, length=%g\n",(double)neP->delta,(double)nrm);CHKERRQ(ierr);
31     *reason = KSP_CONVERGED_STEP_LENGTH;
32   }
33   PetscFunctionReturn(0);
34 }
35 
36 static PetscErrorCode SNESTR_KSPConverged_Destroy(void *cctx)
37 {
38   SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx*)cctx;
39   PetscErrorCode           ierr;
40 
41   PetscFunctionBegin;
42   ierr = (*ctx->convdestroy)(ctx->convctx);CHKERRQ(ierr);
43   ierr = PetscFree(ctx);CHKERRQ(ierr);
44   PetscFunctionReturn(0);
45 }
46 
47 /* ---------------------------------------------------------------- */
48 /*
49    SNESTR_Converged_Private -test convergence JUST for
50    the trust region tolerance.
51 
52 */
53 static PetscErrorCode SNESTR_Converged_Private(SNES snes,PetscInt it,PetscReal xnorm,PetscReal pnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
54 {
55   SNES_NEWTONTR  *neP = (SNES_NEWTONTR*)snes->data;
56   PetscErrorCode ierr;
57 
58   PetscFunctionBegin;
59   *reason = SNES_CONVERGED_ITERATING;
60   if (neP->delta < xnorm * snes->deltatol) {
61     ierr    = PetscInfo3(snes,"Converged due to trust region param %g<%g*%g\n",(double)neP->delta,(double)xnorm,(double)snes->deltatol);CHKERRQ(ierr);
62     *reason = SNES_DIVERGED_TR_DELTA;
63   } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
64     ierr    = PetscInfo1(snes,"Exceeded maximum number of function evaluations: %D\n",snes->max_funcs);CHKERRQ(ierr);
65     *reason = SNES_DIVERGED_FUNCTION_COUNT;
66   }
67   PetscFunctionReturn(0);
68 }
69 
70 /*@C
71    SNESNewtonTRSetPreCheck - Sets a user function that is called before the search step has been determined.
72        Allows the user a chance to change or override the decision of the line search routine.
73 
74    Logically Collective on snes
75 
76    Input Parameters:
77 +  snes - the nonlinear solver object
78 .  func - [optional] function evaluation routine, see SNESNewtonTRPreCheck()  for the calling sequence
79 -  ctx  - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
80 
81    Level: intermediate
82 
83    Note: This function is called BEFORE the function evaluation within the SNESNEWTONTR solver.
84 
85 .seealso: SNESNewtonTRPreCheck(), SNESNewtonTRGetPreCheck()
86 @*/
87 PetscErrorCode  SNESNewtonTRSetPreCheck(SNES snes, PetscErrorCode (*func)(SNES,Vec,Vec,PetscBool*,void*),void *ctx)
88 {
89   SNES_NEWTONTR  *tr = (SNES_NEWTONTR*)snes->data;
90 
91   PetscFunctionBegin;
92   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
93   if (func) tr->precheck    = func;
94   if (ctx)  tr->precheckctx = ctx;
95   PetscFunctionReturn(0);
96 }
97 
98 /*@C
99    SNESNewtonTRGetPreCheck - Gets the pre-check function
100 
101    Not collective
102 
103    Input Parameter:
104 .  snes - the nonlinear solver context
105 
106    Output Parameters:
107 +  func - [optional] function evaluation routine, see for the calling sequence SNESNewtonTRPreCheck()
108 -  ctx  - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
109 
110    Level: intermediate
111 
112 .seealso: SNESNewtonTRSetPreCheck(), SNESNewtonTRPreCheck()
113 @*/
114 PetscErrorCode  SNESNewtonTRGetPreCheck(SNES snes, PetscErrorCode (**func)(SNES,Vec,Vec,PetscBool*,void*),void **ctx)
115 {
116   SNES_NEWTONTR  *tr = (SNES_NEWTONTR*)snes->data;
117 
118   PetscFunctionBegin;
119   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
120   if (func) *func = tr->precheck;
121   if (ctx)  *ctx  = tr->precheckctx;
122   PetscFunctionReturn(0);
123 }
124 
125 /*@C
126    SNESNewtonTRSetPostCheck - Sets a user function that is called after the search step has been determined but before the next
127        function evaluation. Allows the user a chance to change or override the decision of the line search routine
128 
129    Logically Collective on snes
130 
131    Input Parameters:
132 +  snes - the nonlinear solver object
133 .  func - [optional] function evaluation routine, see SNESNewtonTRPostCheck()  for the calling sequence
134 -  ctx  - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
135 
136    Level: intermediate
137 
138    Note: This function is called BEFORE the function evaluation within the SNESNEWTONTR solver while the function set in
139    SNESLineSearchSetPostCheck() is called AFTER the function evaluation.
140 
141 .seealso: SNESNewtonTRPostCheck(), SNESNewtonTRGetPostCheck()
142 @*/
143 PetscErrorCode  SNESNewtonTRSetPostCheck(SNES snes, PetscErrorCode (*func)(SNES,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void *ctx)
144 {
145   SNES_NEWTONTR  *tr = (SNES_NEWTONTR*)snes->data;
146 
147   PetscFunctionBegin;
148   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
149   if (func) tr->postcheck    = func;
150   if (ctx)  tr->postcheckctx = ctx;
151   PetscFunctionReturn(0);
152 }
153 
154 /*@C
155    SNESNewtonTRGetPostCheck - Gets the post-check function
156 
157    Not collective
158 
159    Input Parameter:
160 .  snes - the nonlinear solver context
161 
162    Output Parameters:
163 +  func - [optional] function evaluation routine, see for the calling sequence SNESNewtonTRPostCheck()
164 -  ctx  - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
165 
166    Level: intermediate
167 
168 .seealso: SNESNewtonTRSetPostCheck(), SNESNewtonTRPostCheck()
169 @*/
170 PetscErrorCode  SNESNewtonTRGetPostCheck(SNES snes, PetscErrorCode (**func)(SNES,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void **ctx)
171 {
172   SNES_NEWTONTR  *tr = (SNES_NEWTONTR*)snes->data;
173 
174   PetscFunctionBegin;
175   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
176   if (func) *func = tr->postcheck;
177   if (ctx)  *ctx  = tr->postcheckctx;
178   PetscFunctionReturn(0);
179 }
180 
181 /*@C
182    SNESNewtonTRPreCheck - Called before the step has been determined in SNESNEWTONTR
183 
184    Logically Collective on snes
185 
186    Input Parameters:
187 +  snes - the solver
188 .  X - The last solution
189 -  Y - The step direction
190 
191    Output Parameters:
192 .  changed_Y - Indicator that the step direction Y has been changed.
193 
194    Level: developer
195 
196 .seealso: SNESNewtonTRSetPreCheck(), SNESNewtonTRGetPreCheck()
197 @*/
198 static PetscErrorCode SNESNewtonTRPreCheck(SNES snes,Vec X,Vec Y,PetscBool *changed_Y)
199 {
200   SNES_NEWTONTR  *tr = (SNES_NEWTONTR*)snes->data;
201   PetscErrorCode ierr;
202 
203   PetscFunctionBegin;
204   *changed_Y = PETSC_FALSE;
205   if (tr->precheck) {
206     ierr = (*tr->precheck)(snes,X,Y,changed_Y,tr->precheckctx);CHKERRQ(ierr);
207     PetscValidLogicalCollectiveBool(snes,*changed_Y,4);
208   }
209   PetscFunctionReturn(0);
210 }
211 
212 /*@C
213    SNESNewtonTRPostCheck - Called after the step has been determined in SNESNEWTONTR but before the function evaluation
214 
215    Logically Collective on snes
216 
217    Input Parameters:
218 +  snes - the solver.  X - The last solution
219 .  Y - The full step direction
220 -  W - The updated solution, W = X + tao*Y for some tao
221 
222    Output Parameters:
223 .  changed_W - Indicator if the new candidate solution W has been changed.
224 
225    Notes:
226      If Y is changed then W is recomputed as X + tao*Y where tao is the current update value in SNESNEWTONTR
227 
228    Level: developer
229 
230 .seealso: SNESNewtonTRSetPostCheck(), SNESNewtonTRGetPostCheck()
231 @*/
232 static PetscErrorCode SNESNewtonTRPostCheck(SNES snes,Vec X,Vec Y,Vec W,PetscBool *changed_Y,PetscBool *changed_W)
233 {
234   SNES_NEWTONTR  *tr = (SNES_NEWTONTR*)snes->data;
235   PetscErrorCode ierr;
236 
237   PetscFunctionBegin;
238   *changed_Y = PETSC_FALSE;
239   *changed_W = PETSC_FALSE;
240   if (tr->postcheck) {
241     ierr = (*tr->postcheck)(snes,X,Y,W,changed_Y,changed_W,tr->postcheckctx);CHKERRQ(ierr);
242     PetscValidLogicalCollectiveBool(snes,*changed_Y,5);
243     PetscValidLogicalCollectiveBool(snes,*changed_W,6);
244   }
245   PetscFunctionReturn(0);
246 }
247 
248 /*
249    SNESSolve_NEWTONTR - Implements Newton's Method with a very simple trust
250    region approach for solving systems of nonlinear equations.
251 
252 
253 */
254 static PetscErrorCode SNESSolve_NEWTONTR(SNES snes)
255 {
256   SNES_NEWTONTR            *neP = (SNES_NEWTONTR*)snes->data;
257   Vec                      X,F,Y,G,Ytmp,W;
258   PetscErrorCode           ierr;
259   PetscInt                 maxits,i,lits;
260   PetscReal                rho,fnorm,gnorm,gpnorm,xnorm=0,delta,nrm,ynorm,norm1;
261   PetscScalar              cnorm;
262   KSP                      ksp;
263   SNESConvergedReason      reason = SNES_CONVERGED_ITERATING;
264   PetscBool                breakout = PETSC_FALSE;
265   SNES_TR_KSPConverged_Ctx *ctx;
266   PetscErrorCode           (*convtest)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
267 
268   PetscFunctionBegin;
269   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);
270 
271   maxits = snes->max_its;               /* maximum number of iterations */
272   X      = snes->vec_sol;               /* solution vector */
273   F      = snes->vec_func;              /* residual vector */
274   Y      = snes->work[0];               /* work vectors */
275   G      = snes->work[1];
276   Ytmp   = snes->work[2];
277   W      = snes->work[3];
278 
279   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
280   snes->iter = 0;
281   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
282 
283   /* Set the linear stopping criteria to use the More' trick. */
284   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
285   ierr = KSPGetConvergenceTest(ksp,&convtest,NULL,NULL);CHKERRQ(ierr);
286   if (convtest != SNESTR_KSPConverged_Private) {
287     ierr                  = PetscNew(&ctx);CHKERRQ(ierr);
288     ctx->snes             = snes;
289     ierr                  = KSPGetAndClearConvergenceTest(ksp,&ctx->convtest,&ctx->convctx,&ctx->convdestroy);CHKERRQ(ierr);
290     ierr                  = KSPSetConvergenceTest(ksp,SNESTR_KSPConverged_Private,ctx,SNESTR_KSPConverged_Destroy);CHKERRQ(ierr);
291     ierr                  = PetscInfo(snes,"Using Krylov convergence test SNESTR_KSPConverged_Private\n");CHKERRQ(ierr);
292   }
293 
294   if (!snes->vec_func_init_set) {
295     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);          /* F(X) */
296   } else snes->vec_func_init_set = PETSC_FALSE;
297 
298   ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);             /* fnorm <- || F || */
299   SNESCheckFunctionNorm(snes,fnorm);
300   ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);             /* fnorm <- || F || */
301   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
302   snes->norm = fnorm;
303   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
304   delta      = xnorm ? neP->delta0*xnorm : neP->delta0;
305   neP->delta = delta;
306   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
307   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
308 
309   /* test convergence */
310   ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
311   if (snes->reason) PetscFunctionReturn(0);
312 
313 
314   for (i=0; i<maxits; i++) {
315 
316     /* Call general purpose update function */
317     if (snes->ops->update) {
318       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
319     }
320 
321     /* Solve J Y = F, where J is Jacobian matrix */
322     ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
323     SNESCheckJacobianDomainerror(snes);
324     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
325     ierr = KSPSolve(snes->ksp,F,Ytmp);CHKERRQ(ierr);
326     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
327 
328     snes->linear_its += lits;
329 
330     ierr  = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
331     ierr  = VecNorm(Ytmp,NORM_2,&nrm);CHKERRQ(ierr);
332     norm1 = nrm;
333     while (1) {
334       PetscBool changed_y;
335       PetscBool changed_w;
336       ierr = VecCopy(Ytmp,Y);CHKERRQ(ierr);
337       nrm  = norm1;
338 
339       /* Scale Y if need be and predict new value of F norm */
340       if (nrm >= delta) {
341         nrm    = delta/nrm;
342         gpnorm = (1.0 - nrm)*fnorm;
343         cnorm  = nrm;
344         ierr   = PetscInfo1(snes,"Scaling direction by %g\n",(double)nrm);CHKERRQ(ierr);
345         ierr   = VecScale(Y,cnorm);CHKERRQ(ierr);
346         nrm    = gpnorm;
347         ynorm  = delta;
348       } else {
349         gpnorm = 0.0;
350         ierr   = PetscInfo(snes,"Direction is in Trust Region\n");CHKERRQ(ierr);
351         ynorm  = nrm;
352       }
353       ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr);
354       /* PreCheck() allows for updates to Y prior to W <- X - Y */
355       ierr = SNESNewtonTRPreCheck(snes,X,Y,&changed_y);CHKERRQ(ierr);
356       ierr = VecWAXPY(W,-1.0,Y,X);CHKERRQ(ierr);         /* W <- X - Y */
357       ierr = SNESNewtonTRPostCheck(snes,X,Y,W,&changed_y,&changed_w);CHKERRQ(ierr);
358       if (changed_y) ierr = VecWAXPY(W,-1.0,Y,X);CHKERRQ(ierr);
359       ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); /*  F(X) */
360       ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);      /* gnorm <- || g || */
361       SNESCheckFunctionNorm(snes,gnorm);
362       if (fnorm == gpnorm) rho = 0.0;
363       else rho = (fnorm*fnorm - gnorm*gnorm)/(fnorm*fnorm - gpnorm*gpnorm);
364 
365       /* Update size of trust region */
366       if      (rho < neP->mu)  delta *= neP->delta1;
367       else if (rho < neP->eta) delta *= neP->delta2;
368       else                     delta *= neP->delta3;
369       ierr = PetscInfo3(snes,"fnorm=%g, gnorm=%g, ynorm=%g\n",(double)fnorm,(double)gnorm,(double)ynorm);CHKERRQ(ierr);
370       ierr = PetscInfo3(snes,"gpred=%g, rho=%g, delta=%g\n",(double)gpnorm,(double)rho,(double)delta);CHKERRQ(ierr);
371 
372       neP->delta = delta;
373       if (rho > neP->sigma) break;
374       ierr = PetscInfo(snes,"Trying again in smaller region\n");CHKERRQ(ierr);
375       /* check to see if progress is hopeless */
376       neP->itflag = PETSC_FALSE;
377       ierr        = SNESTR_Converged_Private(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);
378       if (!reason) {ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);}
379       if (reason == SNES_CONVERGED_SNORM_RELATIVE) reason = SNES_DIVERGED_INNER;
380       if (reason) {
381         /* We're not progressing, so return with the current iterate */
382         ierr     = SNESMonitor(snes,i+1,fnorm);CHKERRQ(ierr);
383         breakout = PETSC_TRUE;
384         break;
385       }
386       snes->numFailures++;
387     }
388     if (!breakout) {
389       /* Update function and solution vectors */
390       fnorm = gnorm;
391       ierr  = VecCopy(G,F);CHKERRQ(ierr);
392       ierr  = VecCopy(W,X);CHKERRQ(ierr);
393       /* Monitor convergence */
394       ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
395       snes->iter = i+1;
396       snes->norm = fnorm;
397       snes->xnorm = xnorm;
398       snes->ynorm = ynorm;
399       ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
400       ierr       = SNESLogConvergenceHistory(snes,snes->norm,lits);CHKERRQ(ierr);
401       ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
402       /* Test for convergence, xnorm = || X || */
403       neP->itflag = PETSC_TRUE;
404       if (snes->ops->converged != SNESConvergedSkip) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); }
405       ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);
406       if (reason) break;
407     } else break;
408   }
409   if (i == maxits) {
410     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
411     if (!reason) reason = SNES_DIVERGED_MAX_IT;
412   }
413   ierr         = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
414   snes->reason = reason;
415   ierr         = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
416   PetscFunctionReturn(0);
417 }
418 /*------------------------------------------------------------*/
419 static PetscErrorCode SNESSetUp_NEWTONTR(SNES snes)
420 {
421   PetscErrorCode ierr;
422 
423   PetscFunctionBegin;
424   ierr = SNESSetWorkVecs(snes,4);CHKERRQ(ierr);
425   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
426   PetscFunctionReturn(0);
427 }
428 
429 PetscErrorCode SNESReset_NEWTONTR(SNES snes)
430 {
431 
432   PetscFunctionBegin;
433   PetscFunctionReturn(0);
434 }
435 
436 static PetscErrorCode SNESDestroy_NEWTONTR(SNES snes)
437 {
438   PetscErrorCode ierr;
439 
440   PetscFunctionBegin;
441   ierr = SNESReset_NEWTONTR(snes);CHKERRQ(ierr);
442   ierr = PetscFree(snes->data);CHKERRQ(ierr);
443   PetscFunctionReturn(0);
444 }
445 /*------------------------------------------------------------*/
446 
447 static PetscErrorCode SNESSetFromOptions_NEWTONTR(PetscOptionItems *PetscOptionsObject,SNES snes)
448 {
449   SNES_NEWTONTR  *ctx = (SNES_NEWTONTR*)snes->data;
450   PetscErrorCode ierr;
451 
452   PetscFunctionBegin;
453   ierr = PetscOptionsHead(PetscOptionsObject,"SNES trust region options for nonlinear equations");CHKERRQ(ierr);
454   ierr = PetscOptionsReal("-snes_trtol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,NULL);CHKERRQ(ierr);
455   ierr = PetscOptionsReal("-snes_tr_mu","mu","None",ctx->mu,&ctx->mu,NULL);CHKERRQ(ierr);
456   ierr = PetscOptionsReal("-snes_tr_eta","eta","None",ctx->eta,&ctx->eta,NULL);CHKERRQ(ierr);
457   ierr = PetscOptionsReal("-snes_tr_sigma","sigma","None",ctx->sigma,&ctx->sigma,NULL);CHKERRQ(ierr);
458   ierr = PetscOptionsReal("-snes_tr_delta0","delta0","None",ctx->delta0,&ctx->delta0,NULL);CHKERRQ(ierr);
459   ierr = PetscOptionsReal("-snes_tr_delta1","delta1","None",ctx->delta1,&ctx->delta1,NULL);CHKERRQ(ierr);
460   ierr = PetscOptionsReal("-snes_tr_delta2","delta2","None",ctx->delta2,&ctx->delta2,NULL);CHKERRQ(ierr);
461   ierr = PetscOptionsReal("-snes_tr_delta3","delta3","None",ctx->delta3,&ctx->delta3,NULL);CHKERRQ(ierr);
462   ierr = PetscOptionsTail();CHKERRQ(ierr);
463   PetscFunctionReturn(0);
464 }
465 
466 static PetscErrorCode SNESView_NEWTONTR(SNES snes,PetscViewer viewer)
467 {
468   SNES_NEWTONTR  *tr = (SNES_NEWTONTR*)snes->data;
469   PetscErrorCode ierr;
470   PetscBool      iascii;
471 
472   PetscFunctionBegin;
473   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
474   if (iascii) {
475     ierr = PetscViewerASCIIPrintf(viewer,"  Trust region tolerance (-snes_trtol)\n",(double)snes->deltatol);CHKERRQ(ierr);
476     ierr = PetscViewerASCIIPrintf(viewer,"  mu=%g, eta=%g, sigma=%g\n",(double)tr->mu,(double)tr->eta,(double)tr->sigma);CHKERRQ(ierr);
477     ierr = PetscViewerASCIIPrintf(viewer,"  delta0=%g, delta1=%g, delta2=%g, delta3=%g\n",(double)tr->delta0,(double)tr->delta1,(double)tr->delta2,(double)tr->delta3);CHKERRQ(ierr);
478   }
479   PetscFunctionReturn(0);
480 }
481 /* ------------------------------------------------------------ */
482 /*MC
483       SNESNEWTONTR - Newton based nonlinear solver that uses a trust region
484 
485    Options Database:
486 +    -snes_trtol <tol> - trust region tolerance
487 .    -snes_tr_mu <mu> - trust region parameter
488 .    -snes_tr_eta <eta> - trust region parameter
489 .    -snes_tr_sigma <sigma> - trust region parameter
490 .    -snes_tr_delta0 <delta0> -  initial size of the trust region is delta0*norm2(x)
491 .    -snes_tr_delta1 <delta1> - trust region parameter
492 .    -snes_tr_delta2 <delta2> - trust region parameter
493 -    -snes_tr_delta3 <delta3> - trust region parameter
494 
495    The basic algorithm is taken from "The Minpack Project", by More',
496    Sorensen, Garbow, Hillstrom, pages 88-111 of "Sources and Development
497    of Mathematical Software", Wayne Cowell, editor.
498 
499    Level: intermediate
500 
501 .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESSetTrustRegionTolerance()
502 
503 M*/
504 PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTR(SNES snes)
505 {
506   SNES_NEWTONTR  *neP;
507   PetscErrorCode ierr;
508 
509   PetscFunctionBegin;
510   snes->ops->setup          = SNESSetUp_NEWTONTR;
511   snes->ops->solve          = SNESSolve_NEWTONTR;
512   snes->ops->destroy        = SNESDestroy_NEWTONTR;
513   snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTR;
514   snes->ops->view           = SNESView_NEWTONTR;
515   snes->ops->reset          = SNESReset_NEWTONTR;
516 
517   snes->usesksp = PETSC_TRUE;
518   snes->usesnpc = PETSC_FALSE;
519 
520   snes->alwayscomputesfinalresidual = PETSC_TRUE;
521 
522   ierr        = PetscNewLog(snes,&neP);CHKERRQ(ierr);
523   snes->data  = (void*)neP;
524   neP->mu     = 0.25;
525   neP->eta    = 0.75;
526   neP->delta  = 0.0;
527   neP->delta0 = 0.2;
528   neP->delta1 = 0.3;
529   neP->delta2 = 0.75;
530   neP->delta3 = 2.0;
531   neP->sigma  = 0.0001;
532   neP->itflag = PETSC_FALSE;
533   neP->rnorm0 = 0.0;
534   neP->ttol   = 0.0;
535   PetscFunctionReturn(0);
536 }
537 
538