xref: /petsc/src/snes/impls/ls/ls.c (revision eafb4bcbebb1c0621d2677f4c5d850d03cbd72eb)
1 #ifndef lint
2 static char vcid[] = "$Id: ls.c,v 1.9 1995/04/18 16:26:03 bsmith Exp bsmith $";
3 #endif
4 
5 #include <math.h>
6 #include "ls.h"
7 
8 /*
9      Implements a line search variant of Newton's Method
10     for solving systems of nonlinear equations.
11 
12     Input parameters:
13 .   snes - nonlinear context obtained from SNESCreate()
14 
15     Output Parameters:
16 .   its  - Number of global iterations until termination.
17 
18     Notes:
19     See SNESCreate() and SNESSetUp() for information on the definition and
20     initialization of the nonlinear solver context.
21 
22     This implements essentially a truncated Newton method with a
23     line search.  By default a cubic backtracking line search
24     is employed, as described in the text "Numerical Methods for
25     Unconstrained Optimization and Nonlinear Equations" by Dennis
26     and Schnabel.  See the examples in src/snes/examples.
27 */
28 
29 int SNESSolve_LS( SNES snes, int *outits )
30 {
31   SNES_LS *neP = (SNES_LS *) snes->data;
32   int     maxits, i, history_len,ierr,lits;
33   double  fnorm, gnorm, xnorm, ynorm, *history;
34   Vec     Y, X, F, G, W, TMP;
35 
36   history	= snes->conv_hist;	/* convergence history */
37   history_len	= snes->conv_hist_len;	/* convergence history length */
38   maxits	= snes->max_its;	/* maximum number of iterations */
39   X		= snes->vec_sol;		/* solution vector */
40   F		= snes->vec_res;		/* residual vector */
41   Y		= snes->work[0];		/* work vectors */
42   G		= snes->work[1];
43   W		= snes->work[2];
44 
45   ierr = SNESComputeInitialGuess(snes,X); CHKERR(ierr);  /* X <- X_0 */
46   VecNorm( X, &xnorm );		       /* xnorm = || X || */
47   ierr = SNESComputeFunction(snes,X,F); CHKERR(ierr); /* (+/-) F(X) */
48   VecNorm(F, &fnorm );	        	/* fnorm <- || F || */
49   snes->norm = fnorm;
50   if (history && history_len > 0) history[0] = fnorm;
51   if (snes->Monitor) (*snes->Monitor)(snes,0,fnorm,snes->monP);
52 
53   for ( i=0; i<maxits; i++ ) {
54        snes->iter = i+1;
55 
56        /* Solve J Y = -F, where J is Jacobian matrix */
57        (*snes->ComputeJacobian)(X,&snes->jacobian,snes->jacP);
58        ierr = SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian,0);
59        ierr = SLESSolve(snes->sles,F,Y,&lits); CHKERR(ierr);
60        ierr = (*neP->LineSearch)(snes, X, F, G, Y, W, fnorm, &ynorm, &gnorm );
61 
62        TMP = F; F = G; G = TMP;
63        TMP = X; X = Y; Y = TMP;
64        fnorm = gnorm;
65 
66        snes->norm = fnorm;
67        if (history && history_len > i+1) history[i+1] = fnorm;
68        VecNorm( X, &xnorm );		/* xnorm = || X || */
69        if (snes->Monitor) (*snes->Monitor)(snes,i+1,fnorm,snes->monP);
70 
71        /* Test for convergence */
72        if ((*snes->Converged)( snes, xnorm, ynorm, fnorm,snes->cnvP )) {
73            if (X != snes->vec_sol) VecCopy( X, snes->vec_sol );
74            break;
75        }
76   }
77   if (i == maxits) i--;
78   *outits = i+1;
79   return 0;
80 }
81 /* ------------------------------------------------------------ */
82 /*ARGSUSED*/
83 int SNESSetUp_LS(SNES snes )
84 {
85   int ierr;
86   snes->nwork = 3;
87   ierr = VecGetVecs( snes->vec_sol, snes->nwork,&snes->work ); CHKERR(ierr);
88   PLogObjectParents(snes,snes->nwork,snes->work );
89   return 0;
90 }
91 /* ------------------------------------------------------------ */
92 /*ARGSUSED*/
93 int SNESDestroy_LS(PetscObject obj)
94 {
95   SNES snes = (SNES) obj;
96   SLESDestroy(snes->sles);
97   VecFreeVecs(snes->work, snes->nwork );
98   PLogObjectDestroy(obj);
99   PETSCHEADERDESTROY(obj);
100   return 0;
101 }
102 /*@
103    SNESDefaultMonitor - Default SNES monitoring routine.  Prints the
104    residual norm at each iteration.
105 
106    Input Parameters:
107 .  snes - the SNES context
108 .  its - iteration number
109 .  fnorm - 2-norm residual value (may be estimated)
110 .  dummy - unused context
111 
112    Notes:
113    f is either the residual or its negative, depending on the user's
114    preference, as set with SNESSetFunction().
115 
116 .keywords: SNES, nonlinear, default, monitor, residual, norm
117 
118 .seealso: SNESSetMonitor()
119 @*/
120 int SNESDefaultMonitor(SNES snes,int its, double fnorm,void *dummy)
121 {
122   fprintf( stdout, "iter = %d, residual norm %g \n",its,fnorm);
123   return 0;
124 }
125 
126 /*@
127    SNESDefaultConverged - Default test for monitoring the convergence
128    of the SNES solvers.
129 
130    Input Parameters:
131 .  snes - the SNES context
132 .  xnorm - 2-norm of current iterate
133 .  pnorm - 2-norm of current step
134 .  fnorm - 2-norm of residual
135 .  dummy - unused context
136 
137    Returns:
138 $  2  if  ( fnorm < atol ),
139 $  3  if  ( pnorm < xtol*xnorm ),
140 $ -2  if  ( nres > max_res ),
141 $  0  otherwise,
142 
143    where
144 $    atol    - absolute residual norm tolerance,
145 $              set with SNESSetAbsoluteTolerance()
146 $    max_res - maximum number of residual evaluations,
147 $              set with SNESSetMaxResidualEvaluations()
148 $    nres    - number of residual evaluations
149 $    xtol    - relative residual norm tolerance,
150 $              set with SNESSetRelativeTolerance()
151 
152 .keywords: SNES, nonlinear, default, converged, convergence
153 
154 .seealso: SNESSetConvergenceTest()
155 @*/
156 int SNESDefaultConverged(SNES snes,double xnorm,double pnorm,double fnorm,
157                          void *dummy)
158 {
159   if (fnorm < snes->atol) {
160     PLogInfo((PetscObject)snes,
161       "Converged due to absolute residual norm %g < %g\n",fnorm,snes->atol);
162     return 2;
163   }
164   if (pnorm < snes->xtol*(xnorm)) {
165     PLogInfo((PetscObject)snes,
166       "Converged due to small update length  %g < %g*%g\n",
167        pnorm,snes->xtol,xnorm);
168     return 3;
169   }
170   if (snes->nresids > snes->max_resids) {
171     PLogInfo((PetscObject)snes,
172       "Exceeded maximum number of residual evaluations: %d > %d\n",
173        snes->nresids, snes->max_resids );
174     return -2;
175   }
176   return 0;
177 }
178 
179 /* ------------------------------------------------------------ */
180 /*ARGSUSED*/
181 /*@
182    SNESNoLineSearch - This routine is not a line search at all;
183    it simply uses the full Newton step.  Thus, this routine is intended
184    to serve as a template and is not recommended for general use.
185 
186    Input Parameters:
187 .  snes - nonlinear context
188 .  x - current iterate
189 .  f - residual evaluated at x
190 .  y - search direction (contains new iterate on output)
191 .  w - work vector
192 .  fnorm - 2-norm of f
193 
194    Output parameters:
195 .  g - residual evaluated at new iterate y
196 .  y - new iterate (contains search direction on input)
197 .  gnorm - 2-norm of g
198 .  ynorm - 2-norm of search length
199 
200    Returns:
201    1, indicating success of the step.
202 
203 .keywords: SNES, nonlinear, line search, cubic
204 
205 .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
206 .seealso: SNESSetLineSearchRoutine()
207 @*/
208 int SNESNoLineSearch(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w,
209                              double fnorm, double *ynorm, double *gnorm )
210 {
211   int    ierr;
212   Scalar one = 1.0;
213   VecNorm(y, ynorm );	/* ynorm = || y ||    */
214   VecAXPY(&one, x, y );	/* y <- x + y         */
215   ierr = SNESComputeFunction(snes,y,g); CHKERR(ierr);
216   VecNorm( g, gnorm ); 	/* gnorm = || g ||    */
217   return 1;
218 }
219 /* ------------------------------------------------------------------ */
220 /*@
221    SNESCubicLineSearch - This routine performs a cubic line search and
222    is the default line search method.
223 
224    Input Parameters:
225 .  snes - nonlinear context
226 .  x - current iterate
227 .  f - residual evaluated at x
228 .  y - search direction (contains new iterate on output)
229 .  w - work vector
230 .  fnorm - 2-norm of f
231 
232    Output parameters:
233 .  g - residual evaluated at new iterate y
234 .  y - new iterate (contains search direction on input)
235 .  gnorm - 2-norm of g
236 .  ynorm - 2-norm of search length
237 
238    Returns:
239    1 if the line search succeeds; 0 if the line search fails.
240 
241    Notes:
242    This line search is taken from "Numerical Methods for Unconstrained
243    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
244 
245 .keywords: SNES, nonlinear, line search, cubic
246 
247 .seealso: SNESNoLineSearch(), SNESNoLineSearch(), SNESSetLineSearchRoutine()
248 @*/
249 int SNESCubicLineSearch(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w,
250                               double fnorm, double *ynorm, double *gnorm )
251 {
252   double  steptol, initslope;
253   double  lambdaprev, gnormprev;
254   double  a, b, d, t1, t2;
255 #if defined(PETSC_COMPLEX)
256   Scalar  cinitslope,clambda;
257 #endif
258   int     ierr,count;
259   SNES_LS *neP = (SNES_LS *) snes->data;
260   Scalar  one = 1.0,scale;
261   double  maxstep,minlambda,alpha,lambda,lambdatemp;
262 
263   alpha   = neP->alpha;
264   maxstep = neP->maxstep;
265   steptol = neP->steptol;
266 
267   VecNorm(y, ynorm );
268   if (*ynorm > maxstep) {	/* Step too big, so scale back */
269     scale = maxstep/(*ynorm);
270 #if defined(PETSC_COMPLEX)
271     PLogInfo((PetscObject)snes,"Scaling step by %g\n",real(scale));
272 #else
273     PLogInfo((PetscObject)snes,"Scaling step by %g\n",scale);
274 #endif
275     VecScale(&scale, y );
276     *ynorm = maxstep;
277   }
278   minlambda = steptol/(*ynorm);
279 #if defined(PETSC_COMPLEX)
280   VecDot(f, y, &cinitslope );
281   initslope = real(cinitslope);
282 #else
283   VecDot(f, y, &initslope );
284 #endif
285   if (initslope > 0.0) initslope = -initslope;
286   if (initslope == 0.0) initslope = -1.0;
287 
288   VecCopy(y, w );
289   VecAXPY(&one, x, w );
290   ierr = SNESComputeFunction(snes,w,g); CHKERR(ierr);
291   VecNorm(g, gnorm );
292   if (*gnorm <= fnorm + alpha*initslope) {	/* Sufficient reduction */
293       VecCopy(w, y );
294       PLogInfo((PetscObject)snes,"Using full step\n");
295       return 0;
296   }
297 
298   /* Fit points with quadratic */
299   lambda = 1.0; count = 0;
300   lambdatemp = -initslope/(2.0*(*gnorm - fnorm - initslope));
301   lambdaprev = lambda;
302   gnormprev = *gnorm;
303   if (lambdatemp <= .1*lambda) {
304       lambda = .1*lambda;
305   } else lambda = lambdatemp;
306   VecCopy(x, w );
307 #if defined(PETSC_COMPLEX)
308   clambda = lambda; VecAXPY(&clambda, y, w );
309 #else
310   VecAXPY(&lambda, y, w );
311 #endif
312   ierr = SNESComputeFunction(snes,w,g); CHKERR(ierr);
313   VecNorm(g, gnorm );
314   if (*gnorm <= fnorm + alpha*initslope) {      /* sufficient reduction */
315       VecCopy(w, y );
316       PLogInfo((PetscObject)snes,"Quadratically determined step, lambda %g\n",lambda);
317       return 0;
318   }
319 
320   /* Fit points with cubic */
321   count = 1;
322   while (1) {
323       if (lambda <= minlambda) { /* bad luck; use full step */
324            PLogInfo((PetscObject)snes,"Unable to find good step length! %d \n",count);
325            PLogInfo((PetscObject)snes, "f %g fnew %g ynorm %g lambda %g \n",
326                    fnorm,*gnorm, *ynorm,lambda);
327            VecCopy(w, y );
328            return 0;
329       }
330       t1 = *gnorm - fnorm - lambda*initslope;
331       t2 = gnormprev  - fnorm - lambdaprev*initslope;
332       a = (t1/(lambda*lambda) -
333                 t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
334       b = (-lambdaprev*t1/(lambda*lambda) +
335                 lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
336       d = b*b - 3*a*initslope;
337       if (d < 0.0) d = 0.0;
338       if (a == 0.0) {
339          lambdatemp = -initslope/(2.0*b);
340       } else {
341          lambdatemp = (-b + sqrt(d))/(3.0*a);
342       }
343       if (lambdatemp > .5*lambda) {
344          lambdatemp = .5*lambda;
345       }
346       lambdaprev = lambda;
347       gnormprev = *gnorm;
348       if (lambdatemp <= .1*lambda) {
349          lambda = .1*lambda;
350       }
351       else lambda = lambdatemp;
352       VecCopy( x, w );
353 #if defined(PETSC_COMPLEX)
354       VecAXPY(&clambda, y, w );
355 #else
356       VecAXPY(&lambda, y, w );
357 #endif
358       ierr = SNESComputeFunction(snes,w,g); CHKERR(ierr);
359       VecNorm(g, gnorm );
360       if (*gnorm <= fnorm + alpha*initslope) {      /* is reduction enough */
361          VecCopy(w, y );
362          PLogInfo((PetscObject)snes,"Cubically determined step, lambda %g\n",lambda);
363          return 0;
364       }
365       count++;
366    }
367   return 0;
368 }
369 /*@
370    SNESQuadraticLineSearch - This routine performs a cubic line search.
371 
372    Input Parameters:
373 .  snes - the SNES context
374 .  x - current iterate
375 .  f - residual evaluated at x
376 .  y - search direction (contains new iterate on output)
377 .  w - work vector
378 .  fnorm - 2-norm of f
379 
380    Output parameters:
381 .  g - residual evaluated at new iterate y
382 .  y - new iterate (contains search direction on input)
383 .  gnorm - 2-norm of g
384 .  ynorm - 2-norm of search length
385 
386    Returns:
387    1 if the line search succeeds; 0 if the line search fails.
388 
389    Notes:
390    Use SNESSetLineSearchRoutine()
391    to set this routine within the SNES_NLS method.
392 
393 .keywords: SNES, nonlinear, quadratic, line search
394 
395 .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearchRoutine()
396 @*/
397 int SNESQuadraticLineSearch(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w,
398                               double fnorm, double *ynorm, double *gnorm )
399 {
400   double  steptol, initslope;
401   double  lambdaprev, gnormprev;
402 #if defined(PETSC_COMPLEX)
403   Scalar  cinitslope,clambda;
404 #endif
405   int     ierr,count;
406   SNES_LS *neP = (SNES_LS *) snes->data;
407   Scalar  one = 1.0,scale;
408   double  maxstep,minlambda,alpha,lambda,lambdatemp;
409 
410   alpha   = neP->alpha;
411   maxstep = neP->maxstep;
412   steptol = neP->steptol;
413 
414   VecNorm(y, ynorm );
415   if (*ynorm > maxstep) {	/* Step too big, so scale back */
416     scale = maxstep/(*ynorm);
417     VecScale(&scale, y );
418     *ynorm = maxstep;
419   }
420   minlambda = steptol/(*ynorm);
421 #if defined(PETSC_COMPLEX)
422   VecDot(f, y, &cinitslope );
423   initslope = real(cinitslope);
424 #else
425   VecDot(f, y, &initslope );
426 #endif
427   if (initslope > 0.0) initslope = -initslope;
428   if (initslope == 0.0) initslope = -1.0;
429 
430   VecCopy(y, w );
431   VecAXPY(&one, x, w );
432   ierr = SNESComputeFunction(snes,w,g); CHKERR(ierr);
433   VecNorm(g, gnorm );
434   if (*gnorm <= fnorm + alpha*initslope) {	/* Sufficient reduction */
435       VecCopy(w, y );
436       PLogInfo((PetscObject)snes,"Using full step\n");
437       return 0;
438   }
439 
440   /* Fit points with quadratic */
441   lambda = 1.0; count = 0;
442   count = 1;
443   while (1) {
444     if (lambda <= minlambda) { /* bad luck; use full step */
445       PLogInfo((PetscObject)snes,"Unable to find good step length! %d \n",count);
446       PLogInfo((PetscObject)snes, "f %g fnew %g ynorm %g lambda %g \n",
447                    fnorm,*gnorm, *ynorm,lambda);
448       VecCopy(w, y );
449       return 0;
450     }
451     lambdatemp = -initslope/(2.0*(*gnorm - fnorm - initslope));
452     lambdaprev = lambda;
453     gnormprev = *gnorm;
454     if (lambdatemp <= .1*lambda) {
455       lambda = .1*lambda;
456     } else lambda = lambdatemp;
457     VecCopy(x, w );
458 #if defined(PETSC_COMPLEX)
459     clambda = lambda; VecAXPY(&clambda, y, w );
460 #else
461     VecAXPY(&lambda, y, w );
462 #endif
463     ierr = SNESComputeFunction(snes,w,g); CHKERR(ierr);
464     VecNorm(g, gnorm );
465     if (*gnorm <= fnorm + alpha*initslope) {      /* sufficient reduction */
466       VecCopy(w, y );
467       PLogInfo((PetscObject)snes,"Quadratically determined step, lambda %g\n",lambda);
468       return 0;
469     }
470     count++;
471   }
472 
473   return 0;
474 }
475 /* ------------------------------------------------------------ */
476 /*@C
477    SNESSetLineSearchRoutine - Sets the line search routine to be used
478    by the method SNES_LS.
479 
480    Input Parameters:
481 .  snes - nonlinear context obtained from SNESCreate()
482 .  func - pointer to int function
483 
484    Possible routines:
485 .  SNESCubicLineSearch() - default line search
486 .  SNESQuadraticLineSearch() - quadratic line search
487 .  SNESNoLineSearch() - the full Newton step (actually not a line search)
488 
489    Calling sequence of func:
490    func (SNES snes, Vec x, Vec f, Vec g, Vec y,
491          Vec w, double fnorm, double *ynorm, double *gnorm)
492 
493     Input parameters for func:
494 .   snes - nonlinear context
495 .   x - current iterate
496 .   f - residual evaluated at x
497 .   y - search direction (contains new iterate on output)
498 .   w - work vector
499 .   fnorm - 2-norm of f
500 
501     Output parameters for func:
502 .   g - residual evaluated at new iterate y
503 .   y - new iterate (contains search direction on input)
504 .   gnorm - 2-norm of g
505 .   ynorm - 2-norm of search length
506 
507     Returned by func:
508     1 if the line search succeeds; 0 if the line search fails.
509 
510 .keywords: SNES, nonlinear, set, line search, routine
511 
512 .seealso: SNESNoLineSearch(), SNESQuadraticLineSearch(), SNESCubicLineSearch()
513 @*/
514 int SNESSetLineSearchRoutine(SNES snes,int (*func)(SNES,Vec,Vec,Vec,Vec,Vec,
515                              double,double *,double*) )
516 {
517   if ((snes)->type == SNES_NLS)
518     ((SNES_LS *)(snes->data))->LineSearch = func;
519   return 0;
520 }
521 
522 static int SNESPrintHelp_LS(SNES snes)
523 {
524   SNES_LS *ls = (SNES_LS *)snes->data;
525   fprintf(stderr,"-snes_line_search [basic,quadratic,cubic]\n");
526   fprintf(stderr,"-snes_line_search_alpha alpha (default %g)\n",ls->alpha);
527   fprintf(stderr,"-snes_line_search_maxstep max (default %g)\n",ls->maxstep);
528   fprintf(stderr,"-snes_line_search_steptol tol (default %g)\n",ls->steptol);
529   return 0;
530 }
531 
532 #include "options.h"
533 static int SNESSetFromOptions_LS(SNES snes)
534 {
535   SNES_LS *ls = (SNES_LS *)snes->data;
536   char    ver[16];
537   double  tmp;
538 
539   if (OptionsGetDouble(0,snes->prefix,"-snes_line_search_alpa",&tmp)) {
540     ls->alpha = tmp;
541   }
542   if (OptionsGetDouble(0,snes->prefix,"-snes_line_search_maxstep",&tmp)) {
543     ls->maxstep = tmp;
544   }
545   if (OptionsGetDouble(0,snes->prefix,"-snes_line_search_steptol",&tmp)) {
546     ls->steptol = tmp;
547   }
548   if (OptionsGetString(0,snes->prefix,"-snes_line_search",ver,16)) {
549     if (!strcmp(ver,"basic")) {
550       SNESSetLineSearchRoutine(snes,SNESNoLineSearch);
551     }
552     else if (!strcmp(ver,"quadratic")) {
553       SNESSetLineSearchRoutine(snes,SNESQuadraticLineSearch);
554     }
555     else if (!strcmp(ver,"cubic")) {
556       SNESSetLineSearchRoutine(snes,SNESCubicLineSearch);
557     }
558     else {SETERR(1,"Unknown line search?");}
559   }
560   return 0;
561 }
562 
563 /* ------------------------------------------------------------ */
564 int SNESCreate_LS(SNES  snes )
565 {
566   SNES_LS *neP;
567 
568   snes->type		= SNES_NLS;
569   snes->Setup		= SNESSetUp_LS;
570   snes->Solver		= SNESSolve_LS;
571   snes->destroy		= SNESDestroy_LS;
572   snes->Converged	= SNESDefaultConverged;
573   snes->PrintHelp       = SNESPrintHelp_LS;
574   snes->SetFromOptions  = SNESSetFromOptions_LS;
575 
576   neP			= NEW(SNES_LS);   CHKPTR(neP);
577   snes->data    	= (void *) neP;
578   neP->alpha		= 1.e-4;
579   neP->maxstep		= 1.e8;
580   neP->steptol		= 1.e-12;
581   neP->LineSearch       = SNESCubicLineSearch;
582   return 0;
583 }
584 
585 
586 
587 
588