xref: /petsc/src/snes/impls/tr/tr.c (revision cc85fe4ded5189db5e5e073ce90ef04de0003fdb)
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       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
109   snes->iter = 0;
110   ierr       = PetscObjectSAWsGrantAccess((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   ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);             /* fnorm <- || F || */
122   if (PetscIsInfOrNanReal(fnorm)) {
123     snes->reason = SNES_DIVERGED_FNORM_NAN;
124     PetscFunctionReturn(0);
125   }
126 
127   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
128   snes->norm = fnorm;
129   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
130   delta      = neP->delta0*fnorm;
131   neP->delta = delta;
132   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
133   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
134 
135   /* set parameter for default relative tolerance convergence test */
136   snes->ttol = fnorm*snes->rtol;
137   /* test convergence */
138   ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
139   if (snes->reason) PetscFunctionReturn(0);
140 
141   /* Set the stopping criteria to use the More' trick. */
142   ierr = PetscOptionsGetBool(NULL,"-snes_tr_ksp_regular_convergence_test",&conv,NULL);CHKERRQ(ierr);
143   if (!conv) {
144     SNES_TR_KSPConverged_Ctx *ctx;
145     ierr      = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
146     ierr      = PetscNew(SNES_TR_KSPConverged_Ctx,&ctx);CHKERRQ(ierr);
147     ctx->snes = snes;
148     ierr      = KSPDefaultConvergedCreate(&ctx->ctx);CHKERRQ(ierr);
149     ierr      = KSPSetConvergenceTest(ksp,SNES_TR_KSPConverged_Private,ctx,SNES_TR_KSPConverged_Destroy);CHKERRQ(ierr);
150     ierr      = PetscInfo(snes,"Using Krylov convergence test SNES_TR_KSPConverged_Private\n");CHKERRQ(ierr);
151   }
152 
153   for (i=0; i<maxits; i++) {
154 
155     /* Call general purpose update function */
156     if (snes->ops->update) {
157       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
158     }
159 
160     /* Solve J Y = F, where J is Jacobian matrix */
161     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
162     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
163     ierr = KSPSolve(snes->ksp,F,Ytmp);CHKERRQ(ierr);
164     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
165 
166     snes->linear_its += lits;
167 
168     ierr  = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
169     ierr  = VecNorm(Ytmp,NORM_2,&nrm);CHKERRQ(ierr);
170     norm1 = nrm;
171     while (1) {
172       ierr = VecCopy(Ytmp,Y);CHKERRQ(ierr);
173       nrm  = norm1;
174 
175       /* Scale Y if need be and predict new value of F norm */
176       if (nrm >= delta) {
177         nrm    = delta/nrm;
178         gpnorm = (1.0 - nrm)*fnorm;
179         cnorm  = nrm;
180         ierr   = PetscInfo1(snes,"Scaling direction by %G\n",nrm);CHKERRQ(ierr);
181         ierr   = VecScale(Y,cnorm);CHKERRQ(ierr);
182         nrm    = gpnorm;
183         ynorm  = delta;
184       } else {
185         gpnorm = 0.0;
186         ierr   = PetscInfo(snes,"Direction is in Trust Region\n");CHKERRQ(ierr);
187         ynorm  = nrm;
188       }
189       ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr);            /* Y <- X - Y */
190       ierr = VecCopy(X,snes->vec_sol_update);CHKERRQ(ierr);
191       ierr = SNESComputeFunction(snes,Y,G);CHKERRQ(ierr); /*  F(X) */
192       ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);      /* gnorm <- || g || */
193       if (fnorm == gpnorm) rho = 0.0;
194       else rho = (fnorm*fnorm - gnorm*gnorm)/(fnorm*fnorm - gpnorm*gpnorm);
195 
196       /* Update size of trust region */
197       if      (rho < neP->mu)  delta *= neP->delta1;
198       else if (rho < neP->eta) delta *= neP->delta2;
199       else                     delta *= neP->delta3;
200       ierr = PetscInfo3(snes,"fnorm=%G, gnorm=%G, ynorm=%G\n",fnorm,gnorm,ynorm);CHKERRQ(ierr);
201       ierr = PetscInfo3(snes,"gpred=%G, rho=%G, delta=%G\n",gpnorm,rho,delta);CHKERRQ(ierr);
202 
203       neP->delta = delta;
204       if (rho > neP->sigma) break;
205       ierr = PetscInfo(snes,"Trying again in smaller region\n");CHKERRQ(ierr);
206       /* check to see if progress is hopeless */
207       neP->itflag = PETSC_FALSE;
208       ierr        = SNES_TR_Converged_Private(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);
209       if (!reason) { ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); }
210       if (reason) {
211         /* We're not progressing, so return with the current iterate */
212         ierr     = SNESMonitor(snes,i+1,fnorm);CHKERRQ(ierr);
213         breakout = PETSC_TRUE;
214         break;
215       }
216       snes->numFailures++;
217     }
218     if (!breakout) {
219       /* Update function and solution vectors */
220       fnorm = gnorm;
221       ierr  = VecCopy(G,F);CHKERRQ(ierr);
222       ierr  = VecCopy(Y,X);CHKERRQ(ierr);
223       /* Monitor convergence */
224       ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
225       snes->iter = i+1;
226       snes->norm = fnorm;
227       ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
228       ierr       = SNESLogConvergenceHistory(snes,snes->norm,lits);CHKERRQ(ierr);
229       ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
230       /* Test for convergence, xnorm = || X || */
231       neP->itflag = PETSC_TRUE;
232       if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); }
233       ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);
234       if (reason) break;
235     } else break;
236   }
237   if (i == maxits) {
238     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
239     if (!reason) reason = SNES_DIVERGED_MAX_IT;
240   }
241   ierr         = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
242   snes->reason = reason;
243   ierr         = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
244   PetscFunctionReturn(0);
245 }
246 /*------------------------------------------------------------*/
247 #undef __FUNCT__
248 #define __FUNCT__ "SNESSetUp_NEWTONTR"
249 static PetscErrorCode SNESSetUp_NEWTONTR(SNES snes)
250 {
251   PetscErrorCode ierr;
252 
253   PetscFunctionBegin;
254   ierr = SNESSetWorkVecs(snes,3);CHKERRQ(ierr);
255   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
256   PetscFunctionReturn(0);
257 }
258 
259 #undef __FUNCT__
260 #define __FUNCT__ "SNESReset_NEWTONTR"
261 PetscErrorCode SNESReset_NEWTONTR(SNES snes)
262 {
263 
264   PetscFunctionBegin;
265   PetscFunctionReturn(0);
266 }
267 
268 #undef __FUNCT__
269 #define __FUNCT__ "SNESDestroy_NEWTONTR"
270 static PetscErrorCode SNESDestroy_NEWTONTR(SNES snes)
271 {
272   PetscErrorCode ierr;
273 
274   PetscFunctionBegin;
275   ierr = SNESReset_NEWTONTR(snes);CHKERRQ(ierr);
276   ierr = PetscFree(snes->data);CHKERRQ(ierr);
277   PetscFunctionReturn(0);
278 }
279 /*------------------------------------------------------------*/
280 
281 #undef __FUNCT__
282 #define __FUNCT__ "SNESSetFromOptions_NEWTONTR"
283 static PetscErrorCode SNESSetFromOptions_NEWTONTR(SNES snes)
284 {
285   SNES_NEWTONTR  *ctx = (SNES_NEWTONTR*)snes->data;
286   PetscErrorCode ierr;
287 
288   PetscFunctionBegin;
289   ierr = PetscOptionsHead("SNES trust region options for nonlinear equations");CHKERRQ(ierr);
290   ierr = PetscOptionsReal("-snes_trtol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,0);CHKERRQ(ierr);
291   ierr = PetscOptionsReal("-snes_tr_mu","mu","None",ctx->mu,&ctx->mu,0);CHKERRQ(ierr);
292   ierr = PetscOptionsReal("-snes_tr_eta","eta","None",ctx->eta,&ctx->eta,0);CHKERRQ(ierr);
293   ierr = PetscOptionsReal("-snes_tr_sigma","sigma","None",ctx->sigma,&ctx->sigma,0);CHKERRQ(ierr);
294   ierr = PetscOptionsReal("-snes_tr_delta0","delta0","None",ctx->delta0,&ctx->delta0,0);CHKERRQ(ierr);
295   ierr = PetscOptionsReal("-snes_tr_delta1","delta1","None",ctx->delta1,&ctx->delta1,0);CHKERRQ(ierr);
296   ierr = PetscOptionsReal("-snes_tr_delta2","delta2","None",ctx->delta2,&ctx->delta2,0);CHKERRQ(ierr);
297   ierr = PetscOptionsReal("-snes_tr_delta3","delta3","None",ctx->delta3,&ctx->delta3,0);CHKERRQ(ierr);
298   ierr = PetscOptionsTail();CHKERRQ(ierr);
299   PetscFunctionReturn(0);
300 }
301 
302 #undef __FUNCT__
303 #define __FUNCT__ "SNESView_NEWTONTR"
304 static PetscErrorCode SNESView_NEWTONTR(SNES snes,PetscViewer viewer)
305 {
306   SNES_NEWTONTR  *tr = (SNES_NEWTONTR*)snes->data;
307   PetscErrorCode ierr;
308   PetscBool      iascii;
309 
310   PetscFunctionBegin;
311   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
312   if (iascii) {
313     ierr = PetscViewerASCIIPrintf(viewer,"  mu=%G, eta=%G, sigma=%G\n",tr->mu,tr->eta,tr->sigma);CHKERRQ(ierr);
314     ierr = PetscViewerASCIIPrintf(viewer,"  delta0=%G, delta1=%G, delta2=%G, delta3=%G\n",tr->delta0,tr->delta1,tr->delta2,tr->delta3);CHKERRQ(ierr);
315   }
316   PetscFunctionReturn(0);
317 }
318 /* ------------------------------------------------------------ */
319 /*MC
320       SNESNEWTONTR - Newton based nonlinear solver that uses a trust region
321 
322    Options Database:
323 +    -snes_trtol <tol> Trust region tolerance
324 .    -snes_tr_mu <mu>
325 .    -snes_tr_eta <eta>
326 .    -snes_tr_sigma <sigma>
327 .    -snes_tr_delta0 <delta0>
328 .    -snes_tr_delta1 <delta1>
329 .    -snes_tr_delta2 <delta2>
330 -    -snes_tr_delta3 <delta3>
331 
332    The basic algorithm is taken from "The Minpack Project", by More',
333    Sorensen, Garbow, Hillstrom, pages 88-111 of "Sources and Development
334    of Mathematical Software", Wayne Cowell, editor.
335 
336    This is intended as a model implementation, since it does not
337    necessarily have many of the bells and whistles of other
338    implementations.
339 
340    Level: intermediate
341 
342 .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESSetTrustRegionTolerance()
343 
344 M*/
345 #undef __FUNCT__
346 #define __FUNCT__ "SNESCreate_NEWTONTR"
347 PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTR(SNES snes)
348 {
349   SNES_NEWTONTR  *neP;
350   PetscErrorCode ierr;
351 
352   PetscFunctionBegin;
353   snes->ops->setup          = SNESSetUp_NEWTONTR;
354   snes->ops->solve          = SNESSolve_NEWTONTR;
355   snes->ops->destroy        = SNESDestroy_NEWTONTR;
356   snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTR;
357   snes->ops->view           = SNESView_NEWTONTR;
358   snes->ops->reset          = SNESReset_NEWTONTR;
359 
360   snes->usesksp = PETSC_TRUE;
361   snes->usespc  = PETSC_FALSE;
362 
363   ierr        = PetscNewLog(snes,SNES_NEWTONTR,&neP);CHKERRQ(ierr);
364   snes->data  = (void*)neP;
365   neP->mu     = 0.25;
366   neP->eta    = 0.75;
367   neP->delta  = 0.0;
368   neP->delta0 = 0.2;
369   neP->delta1 = 0.3;
370   neP->delta2 = 0.75;
371   neP->delta3 = 2.0;
372   neP->sigma  = 0.0001;
373   neP->itflag = PETSC_FALSE;
374   neP->rnorm0 = 0.0;
375   neP->ttol   = 0.0;
376   PetscFunctionReturn(0);
377 }
378 
379