xref: /petsc/src/snes/impls/tr/tr.c (revision fbbfd4723ecc59da2191385614f809b84d9232dc)
1 
2 #include <../src/snes/impls/tr/trimpl.h>                /*I   "petscsnes.h"   I*/
3 
4 typedef struct {
5   void *ctx;
6   SNES snes;
7 } SNES_TR_KSPConverged_Ctx;
8 
9 /*
10    This convergence test determines if the two norm of the
11    solution lies outside the trust region, if so it halts.
12 */
13 #undef __FUNCT__
14 #define __FUNCT__ "SNES_TR_KSPConverged_Private"
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;
23 
24   PetscFunctionBegin;
25   ierr = KSPConvergedDefault(ksp,n,rnorm,reason,ctx->ctx);CHKERRQ(ierr);
26   if (*reason) {
27     ierr = PetscInfo2(snes,"default convergence test KSP iterations=%D, rnorm=%g\n",n,(double)rnorm);CHKERRQ(ierr);
28   }
29   /* Determine norm of solution */
30   ierr = KSPBuildSolution(ksp,0,&x);CHKERRQ(ierr);
31   ierr = VecNorm(x,NORM_2,&nrm);CHKERRQ(ierr);
32   if (nrm >= neP->delta) {
33     ierr    = PetscInfo2(snes,"Ending linear iteration early, delta=%g, length=%g\n",(double)neP->delta,(double)nrm);CHKERRQ(ierr);
34     *reason = KSP_CONVERGED_STEP_LENGTH;
35   }
36   PetscFunctionReturn(0);
37 }
38 
39 #undef __FUNCT__
40 #define __FUNCT__ "SNES_TR_KSPConverged_Destroy"
41 PetscErrorCode SNES_TR_KSPConverged_Destroy(void *cctx)
42 {
43   SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx*)cctx;
44   PetscErrorCode           ierr;
45 
46   PetscFunctionBegin;
47   ierr = KSPConvergedDefaultDestroy(ctx->ctx);CHKERRQ(ierr);
48   ierr = PetscFree(ctx);CHKERRQ(ierr);
49   PetscFunctionReturn(0);
50 }
51 
52 /* ---------------------------------------------------------------- */
53 #undef __FUNCT__
54 #define __FUNCT__ "SNES_TR_Converged_Private"
55 /*
56    SNES_TR_Converged_Private -test convergence JUST for
57    the trust region tolerance.
58 
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;
63   PetscErrorCode ierr;
64 
65   PetscFunctionBegin;
66   *reason = SNES_CONVERGED_ITERATING;
67   if (neP->delta < xnorm * snes->deltatol) {
68     ierr    = PetscInfo3(snes,"Converged due to trust region param %g<%g*%g\n",(double)neP->delta,(double)xnorm,(double)snes->deltatol);CHKERRQ(ierr);
69     *reason = SNES_CONVERGED_TR_DELTA;
70   } else if (snes->nfuncs >= snes->max_funcs) {
71     ierr    = PetscInfo1(snes,"Exceeded maximum number of function evaluations: %D\n",snes->max_funcs);CHKERRQ(ierr);
72     *reason = SNES_DIVERGED_FUNCTION_COUNT;
73   }
74   PetscFunctionReturn(0);
75 }
76 
77 
78 /*
79    SNESSolve_NEWTONTR - Implements Newton's Method with a very simple trust
80    region approach for solving systems of nonlinear equations.
81 
82 
83 */
84 #undef __FUNCT__
85 #define __FUNCT__ "SNESSolve_NEWTONTR"
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 
98   PetscFunctionBegin;
99 
100   if (snes->xl || snes->xu || snes->ops->computevariablebounds) {
101     SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
102   }
103 
104   maxits = snes->max_its;               /* maximum number of iterations */
105   X      = snes->vec_sol;               /* solution vector */
106   F      = snes->vec_func;              /* residual vector */
107   Y      = snes->work[0];               /* work vectors */
108   G      = snes->work[1];
109   Ytmp   = snes->work[2];
110 
111   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
112   snes->iter = 0;
113   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
114 
115   if (!snes->vec_func_init_set) {
116     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);          /* F(X) */
117   } else snes->vec_func_init_set = PETSC_FALSE;
118 
119   ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);             /* fnorm <- || F || */
120   SNESCheckFunctionNorm(snes,fnorm);
121   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
122   snes->norm = fnorm;
123   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
124   delta      = neP->delta0*fnorm;
125   neP->delta = delta;
126   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
127   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
128 
129   /* test convergence */
130   ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
131   if (snes->reason) PetscFunctionReturn(0);
132 
133   /* Set the stopping criteria to use the More' trick. */
134   ierr = PetscOptionsGetBool(NULL,"-snes_tr_ksp_regular_convergence_test",&conv,NULL);CHKERRQ(ierr);
135   if (!conv) {
136     SNES_TR_KSPConverged_Ctx *ctx;
137     ierr      = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
138     ierr      = PetscNew(&ctx);CHKERRQ(ierr);
139     ctx->snes = snes;
140     ierr      = KSPConvergedDefaultCreate(&ctx->ctx);CHKERRQ(ierr);
141     ierr      = KSPSetConvergenceTest(ksp,SNES_TR_KSPConverged_Private,ctx,SNES_TR_KSPConverged_Destroy);CHKERRQ(ierr);
142     ierr      = PetscInfo(snes,"Using Krylov convergence test SNES_TR_KSPConverged_Private\n");CHKERRQ(ierr);
143   }
144 
145   for (i=0; i<maxits; i++) {
146 
147     /* Call general purpose update function */
148     if (snes->ops->update) {
149       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
150     }
151 
152     /* Solve J Y = F, where J is Jacobian matrix */
153     ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
154     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
155     ierr = KSPSolve(snes->ksp,F,Ytmp);CHKERRQ(ierr);
156     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
157 
158     snes->linear_its += lits;
159 
160     ierr  = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
161     ierr  = VecNorm(Ytmp,NORM_2,&nrm);CHKERRQ(ierr);
162     norm1 = nrm;
163     while (1) {
164       ierr = VecCopy(Ytmp,Y);CHKERRQ(ierr);
165       nrm  = norm1;
166 
167       /* Scale Y if need be and predict new value of F norm */
168       if (nrm >= delta) {
169         nrm    = delta/nrm;
170         gpnorm = (1.0 - nrm)*fnorm;
171         cnorm  = nrm;
172         ierr   = PetscInfo1(snes,"Scaling direction by %g\n",(double)nrm);CHKERRQ(ierr);
173         ierr   = VecScale(Y,cnorm);CHKERRQ(ierr);
174         nrm    = gpnorm;
175         ynorm  = delta;
176       } else {
177         gpnorm = 0.0;
178         ierr   = PetscInfo(snes,"Direction is in Trust Region\n");CHKERRQ(ierr);
179         ynorm  = nrm;
180       }
181       ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr);            /* Y <- X - Y */
182       ierr = VecCopy(X,snes->vec_sol_update);CHKERRQ(ierr);
183       ierr = SNESComputeFunction(snes,Y,G);CHKERRQ(ierr); /*  F(X) */
184       ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);      /* gnorm <- || g || */
185       if (fnorm == gpnorm) rho = 0.0;
186       else rho = (fnorm*fnorm - gnorm*gnorm)/(fnorm*fnorm - gpnorm*gpnorm);
187 
188       /* Update size of trust region */
189       if      (rho < neP->mu)  delta *= neP->delta1;
190       else if (rho < neP->eta) delta *= neP->delta2;
191       else                     delta *= neP->delta3;
192       ierr = PetscInfo3(snes,"fnorm=%g, gnorm=%g, ynorm=%g\n",(double)fnorm,(double)gnorm,(double)ynorm);CHKERRQ(ierr);
193       ierr = PetscInfo3(snes,"gpred=%g, rho=%g, delta=%g\n",(double)gpnorm,(double)rho,(double)delta);CHKERRQ(ierr);
194 
195       neP->delta = delta;
196       if (rho > neP->sigma) break;
197       ierr = PetscInfo(snes,"Trying again in smaller region\n");CHKERRQ(ierr);
198       /* check to see if progress is hopeless */
199       neP->itflag = PETSC_FALSE;
200       ierr        = SNES_TR_Converged_Private(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);
201       if (!reason) { ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); }
202       if (reason) {
203         /* We're not progressing, so return with the current iterate */
204         ierr     = SNESMonitor(snes,i+1,fnorm);CHKERRQ(ierr);
205         breakout = PETSC_TRUE;
206         break;
207       }
208       snes->numFailures++;
209     }
210     if (!breakout) {
211       /* Update function and solution vectors */
212       fnorm = gnorm;
213       ierr  = VecCopy(G,F);CHKERRQ(ierr);
214       ierr  = VecCopy(Y,X);CHKERRQ(ierr);
215       /* Monitor convergence */
216       ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
217       snes->iter = i+1;
218       snes->norm = fnorm;
219       ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
220       ierr       = SNESLogConvergenceHistory(snes,snes->norm,lits);CHKERRQ(ierr);
221       ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
222       /* Test for convergence, xnorm = || X || */
223       neP->itflag = PETSC_TRUE;
224       if (snes->ops->converged != SNESConvergedSkip) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); }
225       ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);
226       if (reason) break;
227     } else break;
228   }
229   if (i == maxits) {
230     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
231     if (!reason) reason = SNES_DIVERGED_MAX_IT;
232   }
233   ierr         = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
234   snes->reason = reason;
235   ierr         = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
236   PetscFunctionReturn(0);
237 }
238 /*------------------------------------------------------------*/
239 #undef __FUNCT__
240 #define __FUNCT__ "SNESSetUp_NEWTONTR"
241 static PetscErrorCode SNESSetUp_NEWTONTR(SNES snes)
242 {
243   PetscErrorCode ierr;
244 
245   PetscFunctionBegin;
246   ierr = SNESSetWorkVecs(snes,3);CHKERRQ(ierr);
247   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
248   PetscFunctionReturn(0);
249 }
250 
251 #undef __FUNCT__
252 #define __FUNCT__ "SNESReset_NEWTONTR"
253 PetscErrorCode SNESReset_NEWTONTR(SNES snes)
254 {
255 
256   PetscFunctionBegin;
257   PetscFunctionReturn(0);
258 }
259 
260 #undef __FUNCT__
261 #define __FUNCT__ "SNESDestroy_NEWTONTR"
262 static PetscErrorCode SNESDestroy_NEWTONTR(SNES snes)
263 {
264   PetscErrorCode ierr;
265 
266   PetscFunctionBegin;
267   ierr = SNESReset_NEWTONTR(snes);CHKERRQ(ierr);
268   ierr = PetscFree(snes->data);CHKERRQ(ierr);
269   PetscFunctionReturn(0);
270 }
271 /*------------------------------------------------------------*/
272 
273 #undef __FUNCT__
274 #define __FUNCT__ "SNESSetFromOptions_NEWTONTR"
275 static PetscErrorCode SNESSetFromOptions_NEWTONTR(PetscOptions *PetscOptionsObject,SNES snes)
276 {
277   SNES_NEWTONTR  *ctx = (SNES_NEWTONTR*)snes->data;
278   PetscErrorCode ierr;
279 
280   PetscFunctionBegin;
281   ierr = PetscOptionsHead(PetscOptionsObject,"SNES trust region options for nonlinear equations");CHKERRQ(ierr);
282   ierr = PetscOptionsReal("-snes_trtol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,NULL);CHKERRQ(ierr);
283   ierr = PetscOptionsReal("-snes_tr_mu","mu","None",ctx->mu,&ctx->mu,NULL);CHKERRQ(ierr);
284   ierr = PetscOptionsReal("-snes_tr_eta","eta","None",ctx->eta,&ctx->eta,NULL);CHKERRQ(ierr);
285   ierr = PetscOptionsReal("-snes_tr_sigma","sigma","None",ctx->sigma,&ctx->sigma,NULL);CHKERRQ(ierr);
286   ierr = PetscOptionsReal("-snes_tr_delta0","delta0","None",ctx->delta0,&ctx->delta0,NULL);CHKERRQ(ierr);
287   ierr = PetscOptionsReal("-snes_tr_delta1","delta1","None",ctx->delta1,&ctx->delta1,NULL);CHKERRQ(ierr);
288   ierr = PetscOptionsReal("-snes_tr_delta2","delta2","None",ctx->delta2,&ctx->delta2,NULL);CHKERRQ(ierr);
289   ierr = PetscOptionsReal("-snes_tr_delta3","delta3","None",ctx->delta3,&ctx->delta3,NULL);CHKERRQ(ierr);
290   ierr = PetscOptionsTail();CHKERRQ(ierr);
291   PetscFunctionReturn(0);
292 }
293 
294 #undef __FUNCT__
295 #define __FUNCT__ "SNESView_NEWTONTR"
296 static PetscErrorCode SNESView_NEWTONTR(SNES snes,PetscViewer viewer)
297 {
298   SNES_NEWTONTR  *tr = (SNES_NEWTONTR*)snes->data;
299   PetscErrorCode ierr;
300   PetscBool      iascii;
301 
302   PetscFunctionBegin;
303   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
304   if (iascii) {
305     ierr = PetscViewerASCIIPrintf(viewer,"  mu=%g, eta=%g, sigma=%g\n",(double)tr->mu,(double)tr->eta,(double)tr->sigma);CHKERRQ(ierr);
306     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);
307   }
308   PetscFunctionReturn(0);
309 }
310 /* ------------------------------------------------------------ */
311 /*MC
312       SNESNEWTONTR - Newton based nonlinear solver that uses a trust region
313 
314    Options Database:
315 +    -snes_trtol <tol> Trust region tolerance
316 .    -snes_tr_mu <mu>
317 .    -snes_tr_eta <eta>
318 .    -snes_tr_sigma <sigma>
319 .    -snes_tr_delta0 <delta0>
320 .    -snes_tr_delta1 <delta1>
321 .    -snes_tr_delta2 <delta2>
322 -    -snes_tr_delta3 <delta3>
323 
324    The basic algorithm is taken from "The Minpack Project", by More',
325    Sorensen, Garbow, Hillstrom, pages 88-111 of "Sources and Development
326    of Mathematical Software", Wayne Cowell, editor.
327 
328    This is intended as a model implementation, since it does not
329    necessarily have many of the bells and whistles of other
330    implementations.
331 
332    Level: intermediate
333 
334 .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESSetTrustRegionTolerance()
335 
336 M*/
337 #undef __FUNCT__
338 #define __FUNCT__ "SNESCreate_NEWTONTR"
339 PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTR(SNES snes)
340 {
341   SNES_NEWTONTR  *neP;
342   PetscErrorCode ierr;
343 
344   PetscFunctionBegin;
345   snes->ops->setup          = SNESSetUp_NEWTONTR;
346   snes->ops->solve          = SNESSolve_NEWTONTR;
347   snes->ops->destroy        = SNESDestroy_NEWTONTR;
348   snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTR;
349   snes->ops->view           = SNESView_NEWTONTR;
350   snes->ops->reset          = SNESReset_NEWTONTR;
351 
352   snes->usesksp = PETSC_TRUE;
353   snes->usespc  = PETSC_FALSE;
354 
355   ierr        = PetscNewLog(snes,&neP);CHKERRQ(ierr);
356   snes->data  = (void*)neP;
357   neP->mu     = 0.25;
358   neP->eta    = 0.75;
359   neP->delta  = 0.0;
360   neP->delta0 = 0.2;
361   neP->delta1 = 0.3;
362   neP->delta2 = 0.75;
363   neP->delta3 = 2.0;
364   neP->sigma  = 0.0001;
365   neP->itflag = PETSC_FALSE;
366   neP->rnorm0 = 0.0;
367   neP->ttol   = 0.0;
368   PetscFunctionReturn(0);
369 }
370 
371