xref: /petsc/src/snes/impls/tr/tr.c (revision 29bbc08cd5461c366aba6645924e8ff42acd1de0)
1 /*$Id: tr.c,v 1.119 2000/09/02 02:49:40 bsmith Exp bsmith $*/
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 __FUNC__
10 #define __FUNC__ /*<a name=""></a>*/"SNES_EQ_TR_KSPConverged_Private"
11 int SNES_EQ_TR_KSPConverged_Private(KSP ksp,int n,double 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_EQ_TR          *neP = (SNES_EQ_TR*)snes->data;
16   Vec                 x;
17   double              norm;
18   int                 ierr;
19 
20   PetscFunctionBegin;
21   if (snes->ksp_ewconv) {
22     if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Eisenstat-Walker onvergence 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     PLogInfo(snes,"SNES_EQ_TR_KSPConverged_Private: regular convergence test KSP iterations=%d, rnorm=%g\n",n,rnorm);
29   }
30 
31   /* Determine norm of solution */
32   ierr = KSPBuildSolution(ksp,0,&x);CHKERRQ(ierr);
33   ierr = VecNorm(x,NORM_2,&norm);CHKERRQ(ierr);
34   if (norm >= neP->delta) {
35     PLogInfo(snes,"SNES_EQ_TR_KSPConverged_Private: KSP iterations=%d, rnorm=%g\n",n,rnorm);
36     PLogInfo(snes,"SNES_EQ_TR_KSPConverged_Private: Ending linear iteration early, delta=%g, length=%g\n",neP->delta,norm);
37     *reason = KSP_CONVERGED_STEP_LENGTH;
38   }
39   PetscFunctionReturn(0);
40 }
41 
42 /*
43    SNESSolve_EQ_TR - Implements Newton's Method with a very simple trust
44    region approach for solving systems of nonlinear equations.
45 
46    The basic algorithm is taken from "The Minpack Project", by More',
47    Sorensen, Garbow, Hillstrom, pages 88-111 of "Sources and Development
48    of Mathematical Software", Wayne Cowell, editor.
49 
50    This is intended as a model implementation, since it does not
51    necessarily have many of the bells and whistles of other
52    implementations.
53 */
54 #undef __FUNC__
55 #define __FUNC__ /*<a name=""></a>*/"SNESSolve_EQ_TR"
56 static int SNESSolve_EQ_TR(SNES snes,int *its)
57 {
58   SNES_EQ_TR          *neP = (SNES_EQ_TR*)snes->data;
59   Vec                 X,F,Y,G,TMP,Ytmp;
60   int                 maxits,i,ierr,lits,breakout = 0;
61   MatStructure        flg = DIFFERENT_NONZERO_PATTERN;
62   double              rho,fnorm,gnorm,gpnorm,xnorm,delta,norm,ynorm,norm1;
63   Scalar              mone = -1.0,cnorm;
64   KSP                 ksp;
65   SLES                sles;
66   SNESConvergedReason reason;
67 
68   PetscFunctionBegin;
69   maxits	= snes->max_its;	/* maximum number of iterations */
70   X		= snes->vec_sol;	/* solution vector */
71   F		= snes->vec_func;	/* residual vector */
72   Y		= snes->work[0];	/* work vectors */
73   G		= snes->work[1];
74   Ytmp          = snes->work[2];
75 
76   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
77   snes->iter = 0;
78   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
79   ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);         /* xnorm = || X || */
80 
81   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);          /* F(X) */
82   ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);             /* fnorm <- || F || */
83   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
84   snes->norm = fnorm;
85   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
86   delta = neP->delta0*fnorm;
87   neP->delta = delta;
88   SNESLogConvHistory(snes,fnorm,0);
89   SNESMonitor(snes,0,fnorm);
90 
91  if (fnorm < snes->atol) {*its = 0; snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);}
92 
93   /* set parameter for default relative tolerance convergence test */
94   snes->ttol = fnorm*snes->rtol;
95 
96   /* Set the stopping criteria to use the More' trick. */
97   ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr);
98   ierr = SLESGetKSP(sles,&ksp);CHKERRQ(ierr);
99   ierr = KSPSetConvergenceTest(ksp,SNES_EQ_TR_KSPConverged_Private,(void *)snes);CHKERRQ(ierr);
100 
101   for (i=0; i<maxits; i++) {
102     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
103     ierr = SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
104 
105     /* Solve J Y = F, where J is Jacobian matrix */
106     ierr = SLESSolve(snes->sles,F,Ytmp,&lits);CHKERRQ(ierr);
107     snes->linear_its += lits;
108     PLogInfo(snes,"SNESSolve_EQ_TR: iter=%d, linear solve iterations=%d\n",snes->iter,lits);
109     ierr = VecNorm(Ytmp,NORM_2,&norm);CHKERRQ(ierr);
110     norm1 = norm;
111     while(1) {
112       ierr = VecCopy(Ytmp,Y);CHKERRQ(ierr);
113       norm = norm1;
114 
115       /* Scale Y if need be and predict new value of F norm */
116       if (norm >= delta) {
117         norm = delta/norm;
118         gpnorm = (1.0 - norm)*fnorm;
119         cnorm = norm;
120         PLogInfo(snes,"SNESSolve_EQ_TR: Scaling direction by %g\n",norm);
121         ierr = VecScale(&cnorm,Y);CHKERRQ(ierr);
122         norm = gpnorm;
123         ynorm = delta;
124       } else {
125         gpnorm = 0.0;
126         PLogInfo(snes,"SNESSolve_EQ_TR: Direction is in Trust Region\n");
127         ynorm = norm;
128       }
129       ierr = VecAYPX(&mone,X,Y);CHKERRQ(ierr);            /* Y <- X - Y */
130       ierr = VecCopy(X,snes->vec_sol_update_always);CHKERRQ(ierr);
131       ierr = SNESComputeFunction(snes,Y,G);CHKERRQ(ierr); /*  F(X) */
132       ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);      /* gnorm <- || g || */
133       if (fnorm == gpnorm) rho = 0.0;
134       else rho = (fnorm*fnorm - gnorm*gnorm)/(fnorm*fnorm - gpnorm*gpnorm);
135 
136       /* Update size of trust region */
137       if      (rho < neP->mu)  delta *= neP->delta1;
138       else if (rho < neP->eta) delta *= neP->delta2;
139       else                     delta *= neP->delta3;
140       PLogInfo(snes,"SNESSolve_EQ_TR: fnorm=%g, gnorm=%g, ynorm=%g\n",fnorm,gnorm,ynorm);
141       PLogInfo(snes,"SNESSolve_EQ_TR: gpred=%g, rho=%g, delta=%g\n",gpnorm,rho,delta);
142       neP->delta = delta;
143       if (rho > neP->sigma) break;
144       PLogInfo(snes,"SNESSolve_EQ_TR: Trying again in smaller region\n");
145       /* check to see if progress is hopeless */
146       neP->itflag = 0;
147       ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);
148       if (reason) {
149         /* We're not progressing, so return with the current iterate */
150         breakout = 1;
151         break;
152       }
153       snes->nfailures++;
154     }
155     if (!breakout) {
156       fnorm = gnorm;
157       ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
158       snes->iter = i+1;
159       snes->norm = fnorm;
160       ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
161       TMP = F; F = G; snes->vec_func_always = F; G = TMP;
162       TMP = X; X = Y; snes->vec_sol_always  = X; Y = TMP;
163       ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);		/* xnorm = || X || */
164       SNESLogConvHistory(snes,fnorm,lits);
165       SNESMonitor(snes,i+1,fnorm);
166 
167       /* Test for convergence */
168       neP->itflag = 1;
169       ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);
170       if (reason) {
171         break;
172       }
173     } else {
174       break;
175     }
176   }
177   /* Verify solution is in corect location */
178   if (X != snes->vec_sol) {
179     ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr);
180   }
181   if (F != snes->vec_func) {
182     ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr);
183   }
184   snes->vec_sol_always  = snes->vec_sol;
185   snes->vec_func_always = snes->vec_func;
186   if (i == maxits) {
187     PLogInfo(snes,"SNESSolve_EQ_TR: Maximum number of iterations has been reached: %d\n",maxits);
188     i--;
189     reason = SNES_DIVERGED_MAX_IT;
190   }
191   *its = i+1;
192   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
193   snes->reason = reason;
194   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
195   PetscFunctionReturn(0);
196 }
197 /*------------------------------------------------------------*/
198 #undef __FUNC__
199 #define __FUNC__ /*<a name=""></a>*/"SNESSetUp_EQ_TR"
200 static int SNESSetUp_EQ_TR(SNES snes)
201 {
202   int ierr;
203 
204   PetscFunctionBegin;
205   snes->nwork = 4;
206   ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
207   PLogObjectParents(snes,snes->nwork,snes->work);
208   snes->vec_sol_update_always = snes->work[3];
209   PetscFunctionReturn(0);
210 }
211 /*------------------------------------------------------------*/
212 #undef __FUNC__
213 #define __FUNC__ /*<a name=""></a>*/"SNESDestroy_EQ_TR"
214 static int SNESDestroy_EQ_TR(SNES snes)
215 {
216   int  ierr;
217 
218   PetscFunctionBegin;
219   if (snes->nwork) {
220     ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr);
221   }
222   ierr = PetscFree(snes->data);CHKERRQ(ierr);
223   PetscFunctionReturn(0);
224 }
225 /*------------------------------------------------------------*/
226 
227 #undef __FUNC__
228 #define __FUNC__ /*<a name=""></a>*/"SNESSetFromOptions_EQ_TR"
229 static int SNESSetFromOptions_EQ_TR(SNES snes)
230 {
231   SNES_EQ_TR *ctx = (SNES_EQ_TR *)snes->data;
232   int        ierr;
233 
234   PetscFunctionBegin;
235   ierr = OptionsHead("SNES trust region options for nonlinear equations");CHKERRQ(ierr);
236     ierr = OptionsDouble("-snes_trtol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,0);CHKERRQ(ierr);
237     ierr = OptionsDouble("-snes_eq_tr_mu","mu","None",ctx->mu,&ctx->mu,0);CHKERRQ(ierr);
238     ierr = OptionsDouble("-snes_eq_tr_eta","eta","None",ctx->eta,&ctx->eta,0);CHKERRQ(ierr);
239     ierr = OptionsDouble("-snes_eq_tr_sigma","sigma","None",ctx->sigma,&ctx->sigma,0);CHKERRQ(ierr);
240     ierr = OptionsDouble("-snes_eq_tr_delta0","delta0","None",ctx->delta0,&ctx->delta0,0);CHKERRQ(ierr);
241     ierr = OptionsDouble("-snes_eq_tr_delta1","delta1","None",ctx->delta1,&ctx->delta1,0);CHKERRQ(ierr);
242     ierr = OptionsDouble("-snes_eq_tr_delta2","delta2","None",ctx->delta2,&ctx->delta2,0);CHKERRQ(ierr);
243     ierr = OptionsDouble("-snes_eq_tr_delta3","delta3","None",ctx->delta3,&ctx->delta3,0);CHKERRQ(ierr);
244   ierr = OptionsTail();CHKERRQ(ierr);
245   PetscFunctionReturn(0);
246 }
247 
248 #undef __FUNC__
249 #define __FUNC__ /*<a name=""></a>*/"SNESView_EQ_TR"
250 static int SNESView_EQ_TR(SNES snes,Viewer viewer)
251 {
252   SNES_EQ_TR *tr = (SNES_EQ_TR *)snes->data;
253   int        ierr;
254   PetscTruth isascii;
255 
256   PetscFunctionBegin;
257   ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr);
258   if (isascii) {
259     ierr = ViewerASCIIPrintf(viewer,"  mu=%g, eta=%g, sigma=%g\n",tr->mu,tr->eta,tr->sigma);CHKERRQ(ierr);
260     ierr = ViewerASCIIPrintf(viewer,"  delta0=%g, delta1=%g, delta2=%g, delta3=%g\n",tr->delta0,tr->delta1,tr->delta2,tr->delta3);CHKERRQ(ierr);
261   } else {
262     SETERRQ1(1,"Viewer type %s not supported for SNES EQ TR",((PetscObject)viewer)->type_name);
263   }
264   PetscFunctionReturn(0);
265 }
266 
267 /* ---------------------------------------------------------------- */
268 #undef __FUNC__
269 #define __FUNC__ /*<a name=""></a>*/"SNESConverged_EQ_TR"
270 /*@C
271    SNESConverged_EQ_TR - Monitors the convergence of the trust region
272    method SNESEQTR for solving systems of nonlinear equations (default).
273 
274    Collective on SNES
275 
276    Input Parameters:
277 +  snes - the SNES context
278 .  xnorm - 2-norm of current iterate
279 .  pnorm - 2-norm of current step
280 .  fnorm - 2-norm of function
281 -  dummy - unused context
282 
283    Output Parameter:
284 .   reason - one of
285 $  SNES_CONVERGED_FNORM_ABS       - (fnorm < atol),
286 $  SNES_CONVERGED_PNORM_RELATIVE  - (pnorm < xtol*xnorm),
287 $  SNES_CONVERGED_FNORM_RELATIVE  - (fnorm < rtol*fnorm0),
288 $  SNES_DIVERGED_FUNCTION_COUNT   - (nfct > maxf),
289 $  SNES_DIVERGED_FNORM_NAN        - (fnorm == NaN),
290 $  SNES_CONVERGED_TR_DELTA        - (delta < xnorm*deltatol),
291 $  SNES_CONVERGED_ITERATING       - (otherwise)
292 
293    where
294 +    delta    - trust region paramenter
295 .    deltatol - trust region size tolerance,
296                 set with SNESSetTrustRegionTolerance()
297 .    maxf - maximum number of function evaluations,
298             set with SNESSetTolerances()
299 .    nfct - number of function evaluations,
300 .    atol - absolute function norm tolerance,
301             set with SNESSetTolerances()
302 -    xtol - relative function norm tolerance,
303             set with SNESSetTolerances()
304 
305    Level: intermediate
306 
307 .keywords: SNES, nonlinear, default, converged, convergence
308 
309 .seealso: SNESSetConvergenceTest(), SNESEisenstatWalkerConverged()
310 @*/
311 int SNESConverged_EQ_TR(SNES snes,double xnorm,double pnorm,double fnorm,SNESConvergedReason *reason,void *dummy)
312 {
313   SNES_EQ_TR *neP = (SNES_EQ_TR *)snes->data;
314   int        ierr;
315 
316   PetscFunctionBegin;
317   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
318     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only");
319   }
320 
321   if (fnorm != fnorm) {
322     PLogInfo(snes,"SNESConverged_EQ_TR:Failed to converged, function norm is NaN\n");
323     *reason = SNES_DIVERGED_FNORM_NAN;
324   } else if (neP->delta < xnorm * snes->deltatol) {
325     PLogInfo(snes,"SNESConverged_EQ_TR: Converged due to trust region param %g<%g*%g\n",neP->delta,xnorm,snes->deltatol);
326     *reason = SNES_CONVERGED_TR_DELTA;
327   } else if (neP->itflag) {
328     ierr = SNESConverged_EQ_LS(snes,xnorm,pnorm,fnorm,reason,dummy);CHKERRQ(ierr);
329   } else if (snes->nfuncs > snes->max_funcs) {
330     PLogInfo(snes,"SNESConverged_EQ_TR: Exceeded maximum number of function evaluations: %d > %d\n",snes->nfuncs,snes->max_funcs);
331     *reason = SNES_DIVERGED_FUNCTION_COUNT;
332   } else {
333     *reason = SNES_CONVERGED_ITERATING;
334   }
335   PetscFunctionReturn(0);
336 }
337 /* ------------------------------------------------------------ */
338 EXTERN_C_BEGIN
339 #undef __FUNC__
340 #define __FUNC__ /*<a name=""></a>*/"SNESCreate_EQ_TR"
341 int SNESCreate_EQ_TR(SNES snes)
342 {
343   SNES_EQ_TR *neP;
344 
345   PetscFunctionBegin;
346   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
347     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only");
348   }
349   snes->setup		= SNESSetUp_EQ_TR;
350   snes->solve		= SNESSolve_EQ_TR;
351   snes->destroy		= SNESDestroy_EQ_TR;
352   snes->converged	= SNESConverged_EQ_TR;
353   snes->setfromoptions  = SNESSetFromOptions_EQ_TR;
354   snes->view            = SNESView_EQ_TR;
355   snes->nwork           = 0;
356 
357   neP			= PetscNew(SNES_EQ_TR);CHKPTRQ(neP);
358   PLogObjectMemory(snes,sizeof(SNES_EQ_TR));
359   snes->data	        = (void*)neP;
360   neP->mu		= 0.25;
361   neP->eta		= 0.75;
362   neP->delta		= 0.0;
363   neP->delta0		= 0.2;
364   neP->delta1		= 0.3;
365   neP->delta2		= 0.75;
366   neP->delta3		= 2.0;
367   neP->sigma		= 0.0001;
368   neP->itflag		= 0;
369   neP->rnorm0		= 0;
370   neP->ttol		= 0;
371   PetscFunctionReturn(0);
372 }
373 EXTERN_C_END
374 
375