xref: /petsc/src/snes/impls/tr/tr.c (revision f05ece33fd483fcf0ae47f3ecf85fefd993a83b6)
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 = KSPDefaultConverged(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,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",neP->delta,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 = KSPDefaultConvergedDestroy(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",neP->delta,xnorm,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   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;
99 
100   PetscFunctionBegin;
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];
107 
108   ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
109   snes->iter = 0;
110   ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
111 
112   if (!snes->vec_func_init_set) {
113     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);          /* F(X) */
114     ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
115     if (domainerror) {
116       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
117       PetscFunctionReturn(0);
118     }
119   } else snes->vec_func_init_set = PETSC_FALSE;
120 
121   if (!snes->norm_init_set) {
122     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);             /* fnorm <- || F || */
123     if (PetscIsInfOrNanReal(fnorm)) {
124       snes->reason = SNES_DIVERGED_FNORM_NAN;
125       PetscFunctionReturn(0);
126     }
127   } else {
128     fnorm               = snes->norm_init;
129     snes->norm_init_set = PETSC_FALSE;
130   }
131 
132   ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
133   snes->norm = fnorm;
134   ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
135   delta      = neP->delta0*fnorm;
136   neP->delta = delta;
137   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
138   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
139 
140   /* set parameter for default relative tolerance convergence test */
141   snes->ttol = fnorm*snes->rtol;
142   /* test convergence */
143   ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
144   if (snes->reason) PetscFunctionReturn(0);
145 
146   /* Set the stopping criteria to use the More' trick. */
147   ierr = PetscOptionsGetBool(NULL,"-snes_tr_ksp_regular_convergence_test",&conv,NULL);CHKERRQ(ierr);
148   if (!conv) {
149     SNES_TR_KSPConverged_Ctx *ctx;
150     ierr      = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
151     ierr      = PetscNew(SNES_TR_KSPConverged_Ctx,&ctx);CHKERRQ(ierr);
152     ctx->snes = snes;
153     ierr      = KSPDefaultConvergedCreate(&ctx->ctx);CHKERRQ(ierr);
154     ierr      = KSPSetConvergenceTest(ksp,SNES_TR_KSPConverged_Private,ctx,SNES_TR_KSPConverged_Destroy);CHKERRQ(ierr);
155     ierr      = PetscInfo(snes,"Using Krylov convergence test SNES_TR_KSPConverged_Private\n");CHKERRQ(ierr);
156   }
157 
158   for (i=0; i<maxits; i++) {
159 
160     /* Call general purpose update function */
161     if (snes->ops->update) {
162       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
163     }
164 
165     /* Solve J Y = F, where J is Jacobian matrix */
166     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
167     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
168     ierr = KSPSolve(snes->ksp,F,Ytmp);CHKERRQ(ierr);
169     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
170 
171     snes->linear_its += lits;
172 
173     ierr  = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
174     ierr  = VecNorm(Ytmp,NORM_2,&nrm);CHKERRQ(ierr);
175     norm1 = nrm;
176     while (1) {
177       ierr = VecCopy(Ytmp,Y);CHKERRQ(ierr);
178       nrm  = norm1;
179 
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         ierr   = PetscInfo1(snes,"Scaling direction by %G\n",nrm);CHKERRQ(ierr);
186         ierr   = VecScale(Y,cnorm);CHKERRQ(ierr);
187         nrm    = gpnorm;
188         ynorm  = delta;
189       } else {
190         gpnorm = 0.0;
191         ierr   = PetscInfo(snes,"Direction is in Trust Region\n");CHKERRQ(ierr);
192         ynorm  = nrm;
193       }
194       ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr);            /* Y <- X - Y */
195       ierr = VecCopy(X,snes->vec_sol_update);CHKERRQ(ierr);
196       ierr = SNESComputeFunction(snes,Y,G);CHKERRQ(ierr); /*  F(X) */
197       ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);      /* gnorm <- || g || */
198       if (fnorm == gpnorm) rho = 0.0;
199       else rho = (fnorm*fnorm - gnorm*gnorm)/(fnorm*fnorm - gpnorm*gpnorm);
200 
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       ierr = PetscInfo3(snes,"fnorm=%G, gnorm=%G, ynorm=%G\n",fnorm,gnorm,ynorm);CHKERRQ(ierr);
206       ierr = PetscInfo3(snes,"gpred=%G, rho=%G, delta=%G\n",gpnorm,rho,delta);CHKERRQ(ierr);
207 
208       neP->delta = delta;
209       if (rho > neP->sigma) break;
210       ierr = PetscInfo(snes,"Trying again in smaller region\n");CHKERRQ(ierr);
211       /* check to see if progress is hopeless */
212       neP->itflag = PETSC_FALSE;
213       ierr        = SNES_TR_Converged_Private(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);
214       if (!reason) { ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); }
215       if (reason) {
216         /* We're not progressing, so return with the current iterate */
217         ierr     = SNESMonitor(snes,i+1,fnorm);CHKERRQ(ierr);
218         breakout = PETSC_TRUE;
219         break;
220       }
221       snes->numFailures++;
222     }
223     if (!breakout) {
224       /* Update function and solution vectors */
225       fnorm = gnorm;
226       ierr  = VecCopy(G,F);CHKERRQ(ierr);
227       ierr  = VecCopy(Y,X);CHKERRQ(ierr);
228       /* Monitor convergence */
229       ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
230       snes->iter = i+1;
231       snes->norm = fnorm;
232       ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
233       ierr       = SNESLogConvergenceHistory(snes,snes->norm,lits);CHKERRQ(ierr);
234       ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
235       /* Test for convergence, xnorm = || X || */
236       neP->itflag = PETSC_TRUE;
237       if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); }
238       ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);
239       if (reason) break;
240     } else break;
241   }
242   if (i == maxits) {
243     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
244     if (!reason) reason = SNES_DIVERGED_MAX_IT;
245   }
246   ierr         = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
247   snes->reason = reason;
248   ierr         = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
249   PetscFunctionReturn(0);
250 }
251 /*------------------------------------------------------------*/
252 #undef __FUNCT__
253 #define __FUNCT__ "SNESSetUp_NEWTONTR"
254 static PetscErrorCode SNESSetUp_NEWTONTR(SNES snes)
255 {
256   PetscErrorCode ierr;
257 
258   PetscFunctionBegin;
259   ierr = SNESSetWorkVecs(snes,3);CHKERRQ(ierr);
260   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
261   PetscFunctionReturn(0);
262 }
263 
264 #undef __FUNCT__
265 #define __FUNCT__ "SNESReset_NEWTONTR"
266 PetscErrorCode SNESReset_NEWTONTR(SNES snes)
267 {
268 
269   PetscFunctionBegin;
270   PetscFunctionReturn(0);
271 }
272 
273 #undef __FUNCT__
274 #define __FUNCT__ "SNESDestroy_NEWTONTR"
275 static PetscErrorCode SNESDestroy_NEWTONTR(SNES snes)
276 {
277   PetscErrorCode ierr;
278 
279   PetscFunctionBegin;
280   ierr = SNESReset_NEWTONTR(snes);CHKERRQ(ierr);
281   ierr = PetscFree(snes->data);CHKERRQ(ierr);
282   PetscFunctionReturn(0);
283 }
284 /*------------------------------------------------------------*/
285 
286 #undef __FUNCT__
287 #define __FUNCT__ "SNESSetFromOptions_NEWTONTR"
288 static PetscErrorCode SNESSetFromOptions_NEWTONTR(SNES snes)
289 {
290   SNES_NEWTONTR  *ctx = (SNES_NEWTONTR*)snes->data;
291   PetscErrorCode ierr;
292 
293   PetscFunctionBegin;
294   ierr = PetscOptionsHead("SNES trust region options for nonlinear equations");CHKERRQ(ierr);
295   ierr = PetscOptionsReal("-snes_trtol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,0);CHKERRQ(ierr);
296   ierr = PetscOptionsReal("-snes_tr_mu","mu","None",ctx->mu,&ctx->mu,0);CHKERRQ(ierr);
297   ierr = PetscOptionsReal("-snes_tr_eta","eta","None",ctx->eta,&ctx->eta,0);CHKERRQ(ierr);
298   ierr = PetscOptionsReal("-snes_tr_sigma","sigma","None",ctx->sigma,&ctx->sigma,0);CHKERRQ(ierr);
299   ierr = PetscOptionsReal("-snes_tr_delta0","delta0","None",ctx->delta0,&ctx->delta0,0);CHKERRQ(ierr);
300   ierr = PetscOptionsReal("-snes_tr_delta1","delta1","None",ctx->delta1,&ctx->delta1,0);CHKERRQ(ierr);
301   ierr = PetscOptionsReal("-snes_tr_delta2","delta2","None",ctx->delta2,&ctx->delta2,0);CHKERRQ(ierr);
302   ierr = PetscOptionsReal("-snes_tr_delta3","delta3","None",ctx->delta3,&ctx->delta3,0);CHKERRQ(ierr);
303   ierr = PetscOptionsTail();CHKERRQ(ierr);
304   PetscFunctionReturn(0);
305 }
306 
307 #undef __FUNCT__
308 #define __FUNCT__ "SNESView_NEWTONTR"
309 static PetscErrorCode SNESView_NEWTONTR(SNES snes,PetscViewer viewer)
310 {
311   SNES_NEWTONTR  *tr = (SNES_NEWTONTR*)snes->data;
312   PetscErrorCode ierr;
313   PetscBool      iascii;
314 
315   PetscFunctionBegin;
316   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
317   if (iascii) {
318     ierr = PetscViewerASCIIPrintf(viewer,"  mu=%G, eta=%G, sigma=%G\n",tr->mu,tr->eta,tr->sigma);CHKERRQ(ierr);
319     ierr = PetscViewerASCIIPrintf(viewer,"  delta0=%G, delta1=%G, delta2=%G, delta3=%G\n",tr->delta0,tr->delta1,tr->delta2,tr->delta3);CHKERRQ(ierr);
320   }
321   PetscFunctionReturn(0);
322 }
323 /* ------------------------------------------------------------ */
324 /*MC
325       SNESNEWTONTR - Newton based nonlinear solver that uses a trust region
326 
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>
336 
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.
340 
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.
344 
345    Level: intermediate
346 
347 .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESSetTrustRegionTolerance()
348 
349 M*/
350 #undef __FUNCT__
351 #define __FUNCT__ "SNESCreate_NEWTONTR"
352 PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTR(SNES snes)
353 {
354   SNES_NEWTONTR  *neP;
355   PetscErrorCode ierr;
356 
357   PetscFunctionBegin;
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;
364 
365   snes->usesksp = PETSC_TRUE;
366   snes->usespc  = PETSC_FALSE;
367 
368   ierr        = PetscNewLog(snes,SNES_NEWTONTR,&neP);CHKERRQ(ierr);
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   PetscFunctionReturn(0);
382 }
383 
384