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