xref: /petsc/src/snes/impls/tr/tr.c (revision 2a6744eb01855f5aa328eb8fdf4b0d01e72ad151)
1 #define PETSCSNES_DLL
2 
3 #include "src/snes/impls/tr/tr.h"                /*I   "petscsnes.h"   I*/
4 
5 /*
6    This convergence test determines if the two norm of the
7    solution lies outside the trust region, if so it halts.
8 */
9 #undef __FUNCT__
10 #define __FUNCT__ "SNES_TR_KSPConverged_Private"
11 PetscErrorCode SNES_TR_KSPConverged_Private(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason *reason,void *ctx)
12 {
13   SNES                snes = (SNES) ctx;
14   SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx*)snes->kspconvctx;
15   SNES_TR             *neP = (SNES_TR*)snes->data;
16   Vec                 x;
17   PetscReal           nrm;
18   PetscErrorCode      ierr;
19 
20   PetscFunctionBegin;
21   if (snes->ksp_ewconv) {
22     if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Eisenstat-Walker convergence context not created");
23     if (!n) {ierr = SNES_KSP_EW_ComputeRelativeTolerance_Private(snes,ksp);CHKERRQ(ierr);}
24     kctx->lresid_last = rnorm;
25   }
26   ierr = KSPDefaultConverged(ksp,n,rnorm,reason,ctx);CHKERRQ(ierr);
27   if (*reason) {
28     ierr = PetscInfo2(snes,"regular convergence test KSP iterations=%D, rnorm=%G\n",n,rnorm);CHKERRQ(ierr);
29   }
30 
31   /* Determine norm of solution */
32   ierr = KSPBuildSolution(ksp,0,&x);CHKERRQ(ierr);
33   ierr = VecNorm(x,NORM_2,&nrm);CHKERRQ(ierr);
34   if (nrm >= neP->delta) {
35     ierr = PetscInfo2(snes,"KSP iterations=%D, rnorm=%G\n",n,rnorm);CHKERRQ(ierr);
36     ierr = PetscInfo2(snes,"Ending linear iteration early, delta=%G, length=%G\n",neP->delta,nrm);CHKERRQ(ierr);
37     *reason = KSP_CONVERGED_STEP_LENGTH;
38   }
39   PetscFunctionReturn(0);
40 }
41 
42 /*
43    SNESSolve_TR - Implements Newton's Method with a very simple trust
44    region approach for solving systems of nonlinear equations.
45 
46 
47 */
48 #undef __FUNCT__
49 #define __FUNCT__ "SNESSolve_TR"
50 static PetscErrorCode SNESSolve_TR(SNES snes)
51 {
52   SNES_TR             *neP = (SNES_TR*)snes->data;
53   Vec                 X,F,Y,G,TMP,Ytmp;
54   PetscErrorCode      ierr;
55   PetscInt            maxits,i,lits;
56   MatStructure        flg = DIFFERENT_NONZERO_PATTERN;
57   PetscReal           rho,fnorm,gnorm,gpnorm,xnorm,delta,nrm,ynorm,norm1;
58   PetscScalar         cnorm;
59   KSP                 ksp;
60   SNESConvergedReason reason;
61   PetscTruth          conv,breakout = PETSC_FALSE;
62 
63   PetscFunctionBegin;
64   maxits	= snes->max_its;	/* maximum number of iterations */
65   X		= snes->vec_sol;	/* solution vector */
66   F		= snes->vec_func;	/* residual vector */
67   Y		= snes->work[0];	/* work vectors */
68   G		= snes->work[1];
69   Ytmp          = snes->work[2];
70 
71   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
72   snes->iter = 0;
73   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
74   ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);         /* xnorm = || X || */
75 
76   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);          /* F(X) */
77   ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);             /* fnorm <- || F || */
78   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
79   snes->norm = fnorm;
80   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
81   delta = neP->delta0*fnorm;
82   neP->delta = delta;
83   SNESLogConvHistory(snes,fnorm,0);
84   SNESMonitor(snes,0,fnorm);
85   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
86 
87  if (fnorm < snes->abstol) {snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);}
88 
89   /* set parameter for default relative tolerance convergence test */
90   snes->ttol = fnorm*snes->rtol;
91 
92   /* Set the stopping criteria to use the More' trick. */
93   ierr = PetscOptionsHasName(PETSC_NULL,"-snes_tr_ksp_regular_convergence_test",&conv);CHKERRQ(ierr);
94   if (!conv) {
95     ierr = KSPSetConvergenceTest(ksp,SNES_TR_KSPConverged_Private,(void*)snes);CHKERRQ(ierr);
96     ierr = PetscInfo(snes,"Using Krylov convergence test SNES_TR_KSPConverged_Private\n");CHKERRQ(ierr);
97   }
98 
99   for (i=0; i<maxits; i++) {
100 
101     /* Call general purpose update function */
102     if (snes->update) {
103       ierr = (*snes->update)(snes, snes->iter);CHKERRQ(ierr);
104     }
105 
106     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
107     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
108 
109     /* Solve J Y = F, where J is Jacobian matrix */
110     ierr = KSPSolve(snes->ksp,F,Ytmp);CHKERRQ(ierr);
111     ierr = KSPGetIterationNumber(ksp,&lits);CHKERRQ(ierr);
112     snes->linear_its += lits;
113     ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
114     ierr = VecNorm(Ytmp,NORM_2,&nrm);CHKERRQ(ierr);
115     norm1 = nrm;
116     while(1) {
117       ierr = VecCopy(Ytmp,Y);CHKERRQ(ierr);
118       nrm = norm1;
119 
120       /* Scale Y if need be and predict new value of F norm */
121       if (nrm >= delta) {
122         nrm = delta/nrm;
123         gpnorm = (1.0 - nrm)*fnorm;
124         cnorm = nrm;
125         ierr = PetscInfo1(snes,"Scaling direction by %G\n",nrm);CHKERRQ(ierr);
126         ierr = VecScale(Y,cnorm);CHKERRQ(ierr);
127         nrm = gpnorm;
128         ynorm = delta;
129       } else {
130         gpnorm = 0.0;
131         ierr = PetscInfo(snes,"Direction is in Trust Region\n");CHKERRQ(ierr);
132         ynorm = nrm;
133       }
134       ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr);            /* Y <- X - Y */
135       ierr = VecCopy(X,snes->vec_sol_update_always);CHKERRQ(ierr);
136       ierr = SNESComputeFunction(snes,Y,G);CHKERRQ(ierr); /*  F(X) */
137       ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);      /* gnorm <- || g || */
138       if (fnorm == gpnorm) rho = 0.0;
139       else rho = (fnorm*fnorm - gnorm*gnorm)/(fnorm*fnorm - gpnorm*gpnorm);
140 
141       /* Update size of trust region */
142       if      (rho < neP->mu)  delta *= neP->delta1;
143       else if (rho < neP->eta) delta *= neP->delta2;
144       else                     delta *= neP->delta3;
145       ierr = PetscInfo3(snes,"fnorm=%G, gnorm=%G, ynorm=%G\n",fnorm,gnorm,ynorm);CHKERRQ(ierr);
146       ierr = PetscInfo3(snes,"gpred=%G, rho=%G, delta=%G\n",gpnorm,rho,delta);CHKERRQ(ierr);
147       neP->delta = delta;
148       if (rho > neP->sigma) break;
149       ierr = PetscInfo(snes,"Trying again in smaller region\n");CHKERRQ(ierr);
150       /* check to see if progress is hopeless */
151       neP->itflag = PETSC_FALSE;
152       ierr = (*snes->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);
153       if (reason) {
154         /* We're not progressing, so return with the current iterate */
155         SNESMonitor(snes,i+1,fnorm);
156         breakout = PETSC_TRUE;
157         break;
158       }
159       snes->numFailures++;
160     }
161     if (!breakout) {
162       fnorm = gnorm;
163       ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
164       snes->iter = i+1;
165       snes->norm = fnorm;
166       ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
167       TMP = F; F = G; snes->vec_func_always = F; G = TMP;
168       TMP = X; X = Y; snes->vec_sol_always  = X; Y = TMP;
169       ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);		/* xnorm = || X || */
170       SNESLogConvHistory(snes,fnorm,lits);
171       SNESMonitor(snes,i+1,fnorm);
172 
173       /* Test for convergence */
174       neP->itflag = PETSC_TRUE;
175       ierr = (*snes->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);
176       if (reason) {
177         break;
178       }
179     } else {
180       break;
181     }
182   }
183   /* Verify solution is in corect location */
184   if (X != snes->vec_sol) {
185     ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr);
186   }
187   if (F != snes->vec_func) {
188     ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr);
189   }
190   snes->vec_sol_always  = snes->vec_sol;
191   snes->vec_func_always = snes->vec_func;
192   if (i == maxits) {
193     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
194     reason = SNES_DIVERGED_MAX_IT;
195   }
196   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
197   snes->reason = reason;
198   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
199   PetscFunctionReturn(0);
200 }
201 /*------------------------------------------------------------*/
202 #undef __FUNCT__
203 #define __FUNCT__ "SNESSetUp_TR"
204 static PetscErrorCode SNESSetUp_TR(SNES snes)
205 {
206   PetscErrorCode ierr;
207 
208   PetscFunctionBegin;
209   snes->nwork = 4;
210   ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
211   ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr);
212   snes->vec_sol_update_always = snes->work[3];
213   PetscFunctionReturn(0);
214 }
215 /*------------------------------------------------------------*/
216 #undef __FUNCT__
217 #define __FUNCT__ "SNESDestroy_TR"
218 static PetscErrorCode SNESDestroy_TR(SNES snes)
219 {
220   PetscErrorCode ierr;
221 
222   PetscFunctionBegin;
223   if (snes->nwork) {
224     ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr);
225   }
226   ierr = PetscFree(snes->data);CHKERRQ(ierr);
227   PetscFunctionReturn(0);
228 }
229 /*------------------------------------------------------------*/
230 
231 #undef __FUNCT__
232 #define __FUNCT__ "SNESSetFromOptions_TR"
233 static PetscErrorCode SNESSetFromOptions_TR(SNES snes)
234 {
235   SNES_TR *ctx = (SNES_TR *)snes->data;
236   PetscErrorCode ierr;
237 
238   PetscFunctionBegin;
239   ierr = PetscOptionsHead("SNES trust region options for nonlinear equations");CHKERRQ(ierr);
240     ierr = PetscOptionsReal("-snes_trtol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,0);CHKERRQ(ierr);
241     ierr = PetscOptionsReal("-snes_tr_mu","mu","None",ctx->mu,&ctx->mu,0);CHKERRQ(ierr);
242     ierr = PetscOptionsReal("-snes_tr_eta","eta","None",ctx->eta,&ctx->eta,0);CHKERRQ(ierr);
243     ierr = PetscOptionsReal("-snes_tr_sigma","sigma","None",ctx->sigma,&ctx->sigma,0);CHKERRQ(ierr);
244     ierr = PetscOptionsReal("-snes_tr_delta0","delta0","None",ctx->delta0,&ctx->delta0,0);CHKERRQ(ierr);
245     ierr = PetscOptionsReal("-snes_tr_delta1","delta1","None",ctx->delta1,&ctx->delta1,0);CHKERRQ(ierr);
246     ierr = PetscOptionsReal("-snes_tr_delta2","delta2","None",ctx->delta2,&ctx->delta2,0);CHKERRQ(ierr);
247     ierr = PetscOptionsReal("-snes_tr_delta3","delta3","None",ctx->delta3,&ctx->delta3,0);CHKERRQ(ierr);
248   ierr = PetscOptionsTail();CHKERRQ(ierr);
249   PetscFunctionReturn(0);
250 }
251 
252 #undef __FUNCT__
253 #define __FUNCT__ "SNESView_TR"
254 static PetscErrorCode SNESView_TR(SNES snes,PetscViewer viewer)
255 {
256   SNES_TR *tr = (SNES_TR *)snes->data;
257   PetscErrorCode ierr;
258   PetscTruth iascii;
259 
260   PetscFunctionBegin;
261   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
262   if (iascii) {
263     ierr = PetscViewerASCIIPrintf(viewer,"  mu=%G, eta=%G, sigma=%G\n",tr->mu,tr->eta,tr->sigma);CHKERRQ(ierr);
264     ierr = PetscViewerASCIIPrintf(viewer,"  delta0=%G, delta1=%G, delta2=%G, delta3=%G\n",tr->delta0,tr->delta1,tr->delta2,tr->delta3);CHKERRQ(ierr);
265   } else {
266     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ TR",((PetscObject)viewer)->type_name);
267   }
268   PetscFunctionReturn(0);
269 }
270 
271 /* ---------------------------------------------------------------- */
272 #undef __FUNCT__
273 #define __FUNCT__ "SNESConverged_TR"
274 /*@C
275    SNESConverged_TR - Monitors the convergence of the trust region
276    method SNESTR for solving systems of nonlinear equations (default).
277 
278    Collective on SNES
279 
280    Input Parameters:
281 +  snes - the SNES context
282 .  xnorm - 2-norm of current iterate
283 .  pnorm - 2-norm of current step
284 .  fnorm - 2-norm of function
285 -  dummy - unused context
286 
287    Output Parameter:
288 .   reason - one of
289 $  SNES_CONVERGED_FNORM_ABS       - (fnorm < abstol),
290 $  SNES_CONVERGED_PNORM_RELATIVE  - (pnorm < xtol*xnorm),
291 $  SNES_CONVERGED_FNORM_RELATIVE  - (fnorm < rtol*fnorm0),
292 $  SNES_DIVERGED_FUNCTION_COUNT   - (nfct > maxf),
293 $  SNES_DIVERGED_FNORM_NAN        - (fnorm == NaN),
294 $  SNES_CONVERGED_TR_DELTA        - (delta < xnorm*deltatol),
295 $  SNES_CONVERGED_ITERATING       - (otherwise)
296 
297    where
298 +    delta    - trust region paramenter
299 .    deltatol - trust region size tolerance,
300                 set with SNESSetTrustRegionTolerance()
301 .    maxf - maximum number of function evaluations,
302             set with SNESSetTolerances()
303 .    nfct - number of function evaluations,
304 .    abstol - absolute function norm tolerance,
305             set with SNESSetTolerances()
306 -    xtol - relative function norm tolerance,
307             set with SNESSetTolerances()
308 
309    Level: intermediate
310 
311 .keywords: SNES, nonlinear, default, converged, convergence
312 
313 .seealso: SNESSetConvergenceTest(), SNESEisenstatWalkerConverged()
314 @*/
315 PetscErrorCode PETSCSNES_DLLEXPORT SNESConverged_TR(SNES snes,PetscInt it,PetscReal xnorm,PetscReal pnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
316 {
317   SNES_TR *neP = (SNES_TR *)snes->data;
318   PetscErrorCode ierr;
319 
320   PetscFunctionBegin;
321   if (fnorm != fnorm) {
322     ierr = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr);
323     *reason = SNES_DIVERGED_FNORM_NAN;
324   } else if (neP->delta < xnorm * snes->deltatol) {
325     ierr = PetscInfo3(snes,"Converged due to trust region param %G<%G*%G\n",neP->delta,xnorm,snes->deltatol);CHKERRQ(ierr);
326     *reason = SNES_CONVERGED_TR_DELTA;
327   } else if (neP->itflag) {
328     ierr = SNESConverged_LS(snes,it,xnorm,pnorm,fnorm,reason,dummy);CHKERRQ(ierr);
329   } else if (snes->nfuncs >= snes->max_funcs) {
330     ierr = PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);CHKERRQ(ierr);
331     *reason = SNES_DIVERGED_FUNCTION_COUNT;
332   } else {
333     *reason = SNES_CONVERGED_ITERATING;
334   }
335   PetscFunctionReturn(0);
336 }
337 /* ------------------------------------------------------------ */
338 /*MC
339       SNESTR - Newton based nonlinear solver that uses a trust region
340 
341    Options Database:
342 +    -snes_trtol <tol> Trust region tolerance
343 .    -snes_tr_mu <mu>
344 .    -snes_tr_eta <eta>
345 .    -snes_tr_sigma <sigma>
346 .    -snes_tr_delta0 <delta0>
347 .    -snes_tr_delta1 <delta1>
348 .    -snes_tr_delta2 <delta2>
349 -    -snes_tr_delta3 <delta3>
350 
351    The basic algorithm is taken from "The Minpack Project", by More',
352    Sorensen, Garbow, Hillstrom, pages 88-111 of "Sources and Development
353    of Mathematical Software", Wayne Cowell, editor.
354 
355    This is intended as a model implementation, since it does not
356    necessarily have many of the bells and whistles of other
357    implementations.
358 
359    Level: intermediate
360 
361 .seealso:  SNESCreate(), SNES, SNESSetType(), SNESLS, SNESSetTrustRegionTolerance()
362 
363 M*/
364 EXTERN_C_BEGIN
365 #undef __FUNCT__
366 #define __FUNCT__ "SNESCreate_TR"
367 PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_TR(SNES snes)
368 {
369   SNES_TR        *neP;
370   PetscErrorCode ierr;
371 
372   PetscFunctionBegin;
373   snes->setup		= SNESSetUp_TR;
374   snes->solve		= SNESSolve_TR;
375   snes->destroy		= SNESDestroy_TR;
376   snes->converged	= SNESConverged_TR;
377   snes->setfromoptions  = SNESSetFromOptions_TR;
378   snes->view            = SNESView_TR;
379   snes->nwork           = 0;
380 
381   ierr			= PetscNew(SNES_TR,&neP);CHKERRQ(ierr);
382   ierr = PetscLogObjectMemory(snes,sizeof(SNES_TR));CHKERRQ(ierr);
383   snes->data	        = (void*)neP;
384   neP->mu		= 0.25;
385   neP->eta		= 0.75;
386   neP->delta		= 0.0;
387   neP->delta0		= 0.2;
388   neP->delta1		= 0.3;
389   neP->delta2		= 0.75;
390   neP->delta3		= 2.0;
391   neP->sigma		= 0.0001;
392   neP->itflag		= PETSC_FALSE;
393   neP->rnorm0		= 0;
394   neP->ttol		= 0;
395   PetscFunctionReturn(0);
396 }
397 EXTERN_C_END
398 
399