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