xref: /petsc/src/snes/impls/ls/ls.c (revision 1184f6bfc0267870dcbedb6aa2d30e108ed61be2)
1 #ifdef PETSC_RCS_HEADER
2 static char vcid[] = "$Id: ls.c,v 1.127 1999/03/14 22:00:47 curfman Exp bsmith $";
3 #endif
4 
5 #include "src/snes/impls/ls/ls.h"
6 
7 /*  --------------------------------------------------------------------
8 
9      This file implements a truncated Newton method with a line search,
10      for solving a system of nonlinear equations, using the SLES, Vec,
11      and Mat interfaces for linear solvers, vectors, and matrices,
12      respectively.
13 
14      The following basic routines are required for each nonlinear solver:
15           SNESCreate_XXX()          - Creates a nonlinear solver context
16           SNESSetFromOptions_XXX()  - Sets runtime options
17           SNESSolve_XXX()           - Solves the nonlinear system
18           SNESDestroy_XXX()         - Destroys the nonlinear solver context
19      The suffix "_XXX" denotes a particular implementation, in this case
20      we use _EQ_LS (e.g., SNESCreate_EQ_LS, SNESSolve_EQ_LS) for solving
21      systems of nonlinear equations (EQ) with a line search (LS) method.
22      These routines are actually called via the common user interface
23      routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and
24      SNESDestroy(), so the application code interface remains identical
25      for all nonlinear solvers.
26 
27      Another key routine is:
28           SNESSetUp_XXX()           - Prepares for the use of a nonlinear solver
29      by setting data structures and options.   The interface routine SNESSetUp()
30      is not usually called directly by the user, but instead is called by
31      SNESSolve() if necessary.
32 
33      Additional basic routines are:
34           SNESPrintHelp_XXX()       - Prints nonlinear solver runtime options
35           SNESView_XXX()            - Prints details of runtime options that
36                                       have actually been used.
37      These are called by application codes via the interface routines
38      SNESPrintHelp() and SNESView().
39 
40      The various types of solvers (preconditioners, Krylov subspace methods,
41      nonlinear solvers, timesteppers) are all organized similarly, so the
42      above description applies to these categories also.
43 
44     -------------------------------------------------------------------- */
45 /*
46    SNESSolve_EQ_LS - Solves a nonlinear system with a truncated Newton
47    method with a line search.
48 
49    Input Parameters:
50 .  snes - the SNES context
51 
52    Output Parameter:
53 .  outits - number of iterations until termination
54 
55    Application Interface Routine: SNESSolve()
56 
57    Notes:
58    This implements essentially a truncated Newton method with a
59    line search.  By default a cubic backtracking line search
60    is employed, as described in the text "Numerical Methods for
61    Unconstrained Optimization and Nonlinear Equations" by Dennis
62    and Schnabel.
63 */
64 #undef __FUNC__
65 #define __FUNC__ "SNESSolve_EQ_LS"
66 int SNESSolve_EQ_LS(SNES snes,int *outits)
67 {
68   SNES_LS       *neP = (SNES_LS *) snes->data;
69   int           maxits, i, ierr, lits, lsfail;
70   MatStructure  flg = DIFFERENT_NONZERO_PATTERN;
71   double        fnorm, gnorm, xnorm, ynorm;
72   Vec           Y, X, F, G, W, TMP;
73 
74   PetscFunctionBegin;
75   maxits	= snes->max_its;	/* maximum number of iterations */
76   X		= snes->vec_sol;	/* solution vector */
77   F		= snes->vec_func;	/* residual vector */
78   Y		= snes->work[0];	/* work vectors */
79   G		= snes->work[1];
80   W		= snes->work[2];
81 
82   ierr = SNESComputeFunction(snes,X,F); CHKERRQ(ierr);  /*  F(X)      */
83   ierr = VecNorm(F,NORM_2,&fnorm); CHKERRQ(ierr);	/* fnorm <- ||F||  */
84   PetscAMSTakeAccess(snes);
85   snes->iter = 0;
86   snes->norm = fnorm;
87   PetscAMSGrantAccess(snes);
88   SNESLogConvHistory(snes,fnorm,0);
89   SNESMonitor(snes,0,fnorm);
90 
91   if (fnorm < snes->atol) {*outits = 0; PetscFunctionReturn(0);}
92 
93   /* set parameter for default relative tolerance convergence test */
94   snes->ttol = fnorm*snes->rtol;
95 
96   for ( i=0; i<maxits; i++ ) {
97 
98     /* Solve J Y = F, where J is Jacobian matrix */
99     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg); CHKERRQ(ierr);
100     ierr = SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,flg); CHKERRQ(ierr);
101     ierr = SLESSolve(snes->sles,F,Y,&lits); CHKERRQ(ierr);
102     snes->linear_its += PetscAbsInt(lits);
103     PLogInfo(snes,"SNESSolve_EQ_LS: iter=%d, linear solve iterations=%d\n",snes->iter,lits);
104 
105     /* Compute a (scaled) negative update in the line search routine:
106          Y <- X - lambda*Y
107        and evaluate G(Y) = function(Y))
108     */
109     ierr = VecCopy(Y,snes->vec_sol_update_always); CHKERRQ(ierr);
110     ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lsfail); CHKERRQ(ierr);
111     PLogInfo(snes,"SNESSolve_EQ_LS: fnorm=%g, gnorm=%g, ynorm=%g, lsfail=%d\n",fnorm,gnorm,ynorm,lsfail);
112     if (lsfail) snes->nfailures++;
113 
114     TMP = F; F = G; snes->vec_func_always = F; G = TMP;
115     TMP = X; X = Y; snes->vec_sol_always = X;  Y = TMP;
116     fnorm = gnorm;
117 
118     PetscAMSTakeAccess(snes);
119     snes->iter = i+1;
120     snes->norm = fnorm;
121     PetscAMSGrantAccess(snes);
122     SNESLogConvHistory(snes,fnorm,lits);
123     SNESMonitor(snes,i+1,fnorm);
124 
125     /* Test for convergence */
126     if (snes->converged) {
127       ierr = VecNorm(X,NORM_2,&xnorm); CHKERRQ(ierr);	/* xnorm = || X || */
128       if ((*snes->converged)(snes,xnorm,ynorm,fnorm,snes->cnvP)) {
129         break;
130       }
131     }
132   }
133   if (X != snes->vec_sol) {
134     ierr = VecCopy(X,snes->vec_sol); CHKERRQ(ierr);
135     snes->vec_sol_always  = snes->vec_sol;
136     snes->vec_func_always = snes->vec_func;
137   }
138   if (i == maxits) {
139     PLogInfo(snes,"SNESSolve_EQ_LS: Maximum number of iterations has been reached: %d\n",maxits);
140     i--;
141   }
142   *outits = i+1;
143   PetscFunctionReturn(0);
144 }
145 /* -------------------------------------------------------------------------- */
146 /*
147    SNESSetUp_EQ_LS - Sets up the internal data structures for the later use
148    of the SNES_EQ_LS nonlinear solver.
149 
150    Input Parameter:
151 .  snes - the SNES context
152 .  x - the solution vector
153 
154    Application Interface Routine: SNESSetUp()
155 
156    Notes:
157    For basic use of the SNES solvers the user need not explicitly call
158    SNESSetUp(), since these actions will automatically occur during
159    the call to SNESSolve().
160  */
161 #undef __FUNC__
162 #define __FUNC__ "SNESSetUp_EQ_LS"
163 int SNESSetUp_EQ_LS(SNES snes)
164 {
165   int ierr;
166 
167   PetscFunctionBegin;
168   snes->nwork = 4;
169   ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
170   PLogObjectParents(snes,snes->nwork,snes->work);
171   snes->vec_sol_update_always = snes->work[3];
172   PetscFunctionReturn(0);
173 }
174 /* -------------------------------------------------------------------------- */
175 /*
176    SNESDestroy_EQ_LS - Destroys the private SNES_LS context that was created
177    with SNESCreate_EQ_LS().
178 
179    Input Parameter:
180 .  snes - the SNES context
181 
182    Application Interface Routine: SNESDestroy()
183  */
184 #undef __FUNC__
185 #define __FUNC__ "SNESDestroy_EQ_LS"
186 int SNESDestroy_EQ_LS(SNES snes)
187 {
188   int  ierr;
189 
190   PetscFunctionBegin;
191   if (snes->nwork) {
192     ierr = VecDestroyVecs(snes->work,snes->nwork); CHKERRQ(ierr);
193   }
194   PetscFree(snes->data);
195   PetscFunctionReturn(0);
196 }
197 /* -------------------------------------------------------------------------- */
198 #undef __FUNC__
199 #define __FUNC__ "SNESNoLineSearch"
200 
201 /*@C
202    SNESNoLineSearch - This routine is not a line search at all;
203    it simply uses the full Newton step.  Thus, this routine is intended
204    to serve as a template and is not recommended for general use.
205 
206    Collective on SNES and Vec
207 
208    Input Parameters:
209 +  snes - nonlinear context
210 .  lsctx - optional context for line search (not used here)
211 .  x - current iterate
212 .  f - residual evaluated at x
213 .  y - search direction (contains new iterate on output)
214 .  w - work vector
215 -  fnorm - 2-norm of f
216 
217    Output Parameters:
218 +  g - residual evaluated at new iterate y
219 .  y - new iterate (contains search direction on input)
220 .  gnorm - 2-norm of g
221 .  ynorm - 2-norm of search length
222 -  flag - set to 0, indicating a successful line search
223 
224    Options Database Key:
225 .  -snes_eq_ls basic - Activates SNESNoLineSearch()
226 
227    Level: advanced
228 
229 .keywords: SNES, nonlinear, line search, cubic
230 
231 .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
232           SNESSetLineSearch(), SNESNoLineSearchNoNorms()
233 @*/
234 int SNESNoLineSearch(SNES snes, void *lsctx, Vec x, Vec f, Vec g, Vec y, Vec w,
235                      double fnorm, double *ynorm, double *gnorm,int *flag )
236 {
237   int        ierr;
238   Scalar     mone = -1.0;
239   SNES_LS    *neP = (SNES_LS *) snes->data;
240   PetscTruth change_y = PETSC_FALSE;
241 
242   PetscFunctionBegin;
243   *flag = 0;
244   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
245   ierr = VecAYPX(&mone,x,y); CHKERRQ(ierr);            /* y <- y - x      */
246   if (neP->CheckStep) {
247    ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y); CHKERRQ(ierr);
248   }
249   ierr = SNESComputeFunction(snes,y,g); CHKERRQ(ierr); /* Compute F(y)    */
250   ierr = VecNormBegin(y,NORM_2,ynorm); CHKERRQ(ierr);  /* ynorm = || y || */
251   ierr = VecNormBegin(g,NORM_2,gnorm); CHKERRQ(ierr);  /* gnorm = || g || */
252   ierr = VecNormEnd(y,NORM_2,ynorm); CHKERRQ(ierr);
253   ierr = VecNormEnd(g,NORM_2,gnorm); CHKERRQ(ierr);
254   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
255   PetscFunctionReturn(0);
256 }
257 /* -------------------------------------------------------------------------- */
258 #undef __FUNC__
259 #define __FUNC__ "SNESNoLineSearchNoNorms"
260 
261 /*@C
262    SNESNoLineSearchNoNorms - This routine is not a line search at
263    all; it simply uses the full Newton step. This version does not
264    even compute the norm of the function or search direction; this
265    is intended only when you know the full step is fine and are
266    not checking for convergence of the nonlinear iteration (for
267    example, you are running always for a fixed number of Newton steps).
268 
269    Collective on SNES and Vec
270 
271    Input Parameters:
272 +  snes - nonlinear context
273 .  lsctx - optional context for line search (not used here)
274 .  x - current iterate
275 .  f - residual evaluated at x
276 .  y - search direction (contains new iterate on output)
277 .  w - work vector
278 -  fnorm - 2-norm of f
279 
280    Output Parameters:
281 +  g - residual evaluated at new iterate y
282 .  gnorm - not changed
283 .  ynorm - not changed
284 -  flag - set to 0, indicating a successful line search
285 
286    Options Database Key:
287 .  -snes_eq_ls basicnonorms - Activates SNESNoLineSearchNoNorms()
288 
289    Level: advanced
290 
291 .keywords: SNES, nonlinear, line search, cubic
292 
293 .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
294           SNESSetLineSearch(), SNESNoLineSearch()
295 @*/
296 int SNESNoLineSearchNoNorms(SNES snes, void *lsctx, Vec x, Vec f, Vec g, Vec y, Vec w,
297                      double fnorm, double *ynorm, double *gnorm,int *flag )
298 {
299   int        ierr;
300   Scalar     mone = -1.0;
301   SNES_LS    *neP = (SNES_LS *) snes->data;
302   PetscTruth change_y = PETSC_FALSE;
303 
304   PetscFunctionBegin;
305   *flag = 1;
306   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
307   ierr = VecAYPX(&mone,x,y); CHKERRQ(ierr);            /* y <- y - x      */
308   if (neP->CheckStep) {
309    ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y); CHKERRQ(ierr);
310   }
311   ierr = SNESComputeFunction(snes,y,g); CHKERRQ(ierr); /* Compute F(y)    */
312   *gnorm = 0.0;
313   *ynorm = 0.0;
314   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
315   PetscFunctionReturn(0);
316 }
317 /* -------------------------------------------------------------------------- */
318 #undef __FUNC__
319 #define __FUNC__ "SNESCubicLineSearch"
320 /*@C
321    SNESCubicLineSearch - Performs a cubic line search (default line search method).
322 
323    Collective on SNES
324 
325    Input Parameters:
326 +  snes - nonlinear context
327 .  lsctx - optional context for line search (not used here)
328 .  x - current iterate
329 .  f - residual evaluated at x
330 .  y - search direction (contains new iterate on output)
331 .  w - work vector
332 -  fnorm - 2-norm of f
333 
334    Output Parameters:
335 +  g - residual evaluated at new iterate y
336 .  y - new iterate (contains search direction on input)
337 .  gnorm - 2-norm of g
338 .  ynorm - 2-norm of search length
339 -  flag - 0 if line search succeeds; -1 on failure.
340 
341    Options Database Key:
342 $  -snes_eq_ls cubic - Activates SNESCubicLineSearch()
343 
344    Notes:
345    This line search is taken from "Numerical Methods for Unconstrained
346    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
347 
348    Level: advanced
349 
350 .keywords: SNES, nonlinear, line search, cubic
351 
352 .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
353 @*/
354 int SNESCubicLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
355                         double fnorm,double *ynorm,double *gnorm,int *flag)
356 {
357   /*
358      Note that for line search purposes we work with with the related
359      minimization problem:
360         min  z(x):  R^n -> R,
361      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
362    */
363 
364   double     steptol, initslope, lambdaprev, gnormprev, a, b, d, t1, t2;
365   double     maxstep, minlambda, alpha, lambda, lambdatemp, lambdaneg;
366 #if defined(USE_PETSC_COMPLEX)
367   Scalar     cinitslope, clambda;
368 #endif
369   int        ierr, count;
370   SNES_LS    *neP = (SNES_LS *) snes->data;
371   Scalar     mone = -1.0, scale;
372   PetscTruth change_y = PETSC_FALSE;
373 
374   PetscFunctionBegin;
375   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
376   *flag   = 0;
377   alpha   = neP->alpha;
378   maxstep = neP->maxstep;
379   steptol = neP->steptol;
380 
381   ierr = VecNorm(y,NORM_2,ynorm); CHKERRQ(ierr);
382   if (*ynorm < snes->atol) {
383     PLogInfo(snes,"SNESCubicLineSearch: Search direction and size are nearly 0\n");
384     *gnorm = fnorm;
385     ierr = VecCopy(x,y); CHKERRQ(ierr);
386     ierr = VecCopy(f,g); CHKERRQ(ierr);
387     goto theend1;
388   }
389   if (*ynorm > maxstep) {	/* Step too big, so scale back */
390     scale = maxstep/(*ynorm);
391 #if defined(USE_PETSC_COMPLEX)
392     PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",PetscReal(scale));
393 #else
394     PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",scale);
395 #endif
396     ierr = VecScale(&scale,y); CHKERRQ(ierr);
397     *ynorm = maxstep;
398   }
399   minlambda = steptol/(*ynorm);
400   ierr = MatMult(snes->jacobian,y,w); CHKERRQ(ierr);
401 #if defined(USE_PETSC_COMPLEX)
402   ierr = VecDot(f,w,&cinitslope); CHKERRQ(ierr);
403   initslope = PetscReal(cinitslope);
404 #else
405   ierr = VecDot(f,w,&initslope); CHKERRQ(ierr);
406 #endif
407   if (initslope > 0.0) initslope = -initslope;
408   if (initslope == 0.0) initslope = -1.0;
409 
410   ierr = VecCopy(y,w); CHKERRQ(ierr);
411   ierr = VecAYPX(&mone,x,w); CHKERRQ(ierr);
412   ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
413   ierr = VecNorm(g,NORM_2,gnorm);
414   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */
415     ierr = VecCopy(w,y); CHKERRQ(ierr);
416     PLogInfo(snes,"SNESCubicLineSearch: Using full step\n");
417     goto theend1;
418   }
419 
420   /* Fit points with quadratic */
421   lambda = 1.0; count = 0;
422   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
423   lambdaprev = lambda;
424   gnormprev = *gnorm;
425   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
426   else lambda = lambdatemp;
427   ierr   = VecCopy(x,w); CHKERRQ(ierr);
428   lambdaneg = -lambda;
429 #if defined(USE_PETSC_COMPLEX)
430   clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr);
431 #else
432   ierr = VecAXPY(&lambdaneg,y,w); CHKERRQ(ierr);
433 #endif
434   ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
435   ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);
436   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */
437     ierr = VecCopy(w,y); CHKERRQ(ierr);
438     PLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%g\n",lambda);
439     goto theend1;
440   }
441 
442   /* Fit points with cubic */
443   count = 1;
444   while (1) {
445     if (lambda <= minlambda) { /* bad luck; use full step */
446       PLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count);
447       PLogInfo(snes,"SNESCubicLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",
448                fnorm,*gnorm,*ynorm,lambda,initslope);
449       ierr = VecCopy(w,y); CHKERRQ(ierr);
450       *flag = -1; break;
451     }
452     t1 = ((*gnorm)*(*gnorm) - fnorm*fnorm)*0.5 - lambda*initslope;
453     t2 = (gnormprev*gnormprev  - fnorm*fnorm)*0.5 - lambdaprev*initslope;
454     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
455     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
456     d  = b*b - 3*a*initslope;
457     if (d < 0.0) d = 0.0;
458     if (a == 0.0) {
459       lambdatemp = -initslope/(2.0*b);
460     } else {
461       lambdatemp = (-b + sqrt(d))/(3.0*a);
462     }
463     if (lambdatemp > .5*lambda) {
464       lambdatemp = .5*lambda;
465     }
466     lambdaprev = lambda;
467     gnormprev = *gnorm;
468     if (lambdatemp <= .1*lambda) {
469       lambda = .1*lambda;
470     }
471     else lambda = lambdatemp;
472     ierr = VecCopy( x, w ); CHKERRQ(ierr);
473     lambdaneg = -lambda;
474 #if defined(USE_PETSC_COMPLEX)
475     clambda = lambdaneg;
476     ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr);
477 #else
478     ierr = VecAXPY(&lambdaneg,y,w); CHKERRQ(ierr);
479 #endif
480     ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
481     ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);
482     if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* is reduction enough? */
483       ierr = VecCopy(w,y); CHKERRQ(ierr);
484       PLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%g\n",lambda);
485       goto theend1;
486     } else {
487       PLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda,  lambda=%g\n",lambda);
488     }
489     count++;
490   }
491   theend1:
492   /* Optional user-defined check for line search step validity */
493   if (neP->CheckStep) {
494     ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y); CHKERRQ(ierr);
495     if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */
496       ierr = SNESComputeFunction(snes,y,g); CHKERRQ(ierr);
497       ierr = VecNormBegin(y,NORM_2,ynorm); CHKERRQ(ierr);
498       ierr = VecNormBegin(g,NORM_2,gnorm); CHKERRQ(ierr);
499       ierr = VecNormEnd(y,NORM_2,ynorm); CHKERRQ(ierr);
500       ierr = VecNormEnd(g,NORM_2,gnorm); CHKERRQ(ierr);
501     }
502   }
503   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
504   PetscFunctionReturn(0);
505 }
506 /* -------------------------------------------------------------------------- */
507 #undef __FUNC__
508 #define __FUNC__ "SNESQuadraticLineSearch"
509 /*@C
510    SNESQuadraticLineSearch - Performs a quadratic line search.
511 
512    Collective on SNES and Vec
513 
514    Input Parameters:
515 +  snes - the SNES context
516 .  lsctx - optional context for line search (not used here)
517 .  x - current iterate
518 .  f - residual evaluated at x
519 .  y - search direction (contains new iterate on output)
520 .  w - work vector
521 -  fnorm - 2-norm of f
522 
523    Output Parameters:
524 +  g - residual evaluated at new iterate y
525 .  y - new iterate (contains search direction on input)
526 .  gnorm - 2-norm of g
527 .  ynorm - 2-norm of search length
528 -  flag - 0 if line search succeeds; -1 on failure.
529 
530    Options Database Key:
531 .  -snes_eq_ls quadratic - Activates SNESQuadraticLineSearch()
532 
533    Notes:
534    Use SNESSetLineSearch() to set this routine within the SNES_EQ_LS method.
535 
536    Level: advanced
537 
538 .keywords: SNES, nonlinear, quadratic, line search
539 
540 .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
541 @*/
542 int SNESQuadraticLineSearch(SNES snes, void *lsctx, Vec x, Vec f, Vec g, Vec y, Vec w,
543                            double fnorm, double *ynorm, double *gnorm,int *flag)
544 {
545   /*
546      Note that for line search purposes we work with with the related
547      minimization problem:
548         min  z(x):  R^n -> R,
549      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
550    */
551   double     steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp;
552 #if defined(USE_PETSC_COMPLEX)
553   Scalar     cinitslope,clambda;
554 #endif
555   int        ierr, count;
556   SNES_LS    *neP = (SNES_LS *) snes->data;
557   Scalar     mone = -1.0,scale;
558   PetscTruth change_y = PETSC_FALSE;
559 
560   PetscFunctionBegin;
561   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
562   *flag   = 0;
563   alpha   = neP->alpha;
564   maxstep = neP->maxstep;
565   steptol = neP->steptol;
566 
567   VecNorm(y, NORM_2,ynorm );
568   if (*ynorm < snes->atol) {
569     PLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n");
570     *gnorm = fnorm;
571     ierr = VecCopy(x,y); CHKERRQ(ierr);
572     ierr = VecCopy(f,g); CHKERRQ(ierr);
573     goto theend2;
574   }
575   if (*ynorm > maxstep) {	/* Step too big, so scale back */
576     scale = maxstep/(*ynorm);
577     ierr = VecScale(&scale,y); CHKERRQ(ierr);
578     *ynorm = maxstep;
579   }
580   minlambda = steptol/(*ynorm);
581   ierr = MatMult(snes->jacobian,y,w); CHKERRQ(ierr);
582 #if defined(USE_PETSC_COMPLEX)
583   ierr = VecDot(f,w,&cinitslope); CHKERRQ(ierr);
584   initslope = PetscReal(cinitslope);
585 #else
586   ierr = VecDot(f,w,&initslope); CHKERRQ(ierr);
587 #endif
588   if (initslope > 0.0) initslope = -initslope;
589   if (initslope == 0.0) initslope = -1.0;
590 
591   ierr = VecCopy(y,w); CHKERRQ(ierr);
592   ierr = VecAYPX(&mone,x,w); CHKERRQ(ierr);
593   ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
594   ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);
595   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */
596     ierr = VecCopy(w,y); CHKERRQ(ierr);
597     PLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n");
598     goto theend2;
599   }
600 
601   /* Fit points with quadratic */
602   lambda = 1.0; count = 0;
603   count = 1;
604   while (1) {
605     if (lambda <= minlambda) { /* bad luck; use full step */
606       PLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count);
607       PLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",
608                fnorm,*gnorm,*ynorm,lambda,initslope);
609       ierr = VecCopy(w,y); CHKERRQ(ierr);
610       *flag = -1; break;
611     }
612     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
613     if (lambdatemp <= .1*lambda) {
614       lambda = .1*lambda;
615     } else lambda = lambdatemp;
616     ierr = VecCopy(x,w); CHKERRQ(ierr);
617     lambda = -lambda;
618 #if defined(USE_PETSC_COMPLEX)
619     clambda = lambda; ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr);
620 #else
621     ierr = VecAXPY(&lambda,y,w); CHKERRQ(ierr);
622 #endif
623     ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
624     ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);
625     if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */
626       ierr = VecCopy(w,y); CHKERRQ(ierr);
627       PLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda);
628       break;
629     }
630     count++;
631   }
632   theend2:
633   /* Optional user-defined check for line search step validity */
634   if (neP->CheckStep) {
635     ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y); CHKERRQ(ierr);
636     if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */
637       ierr = SNESComputeFunction(snes,y,g); CHKERRQ(ierr);
638       ierr = VecNormBegin(y,NORM_2,ynorm); CHKERRQ(ierr);
639       ierr = VecNormBegin(g,NORM_2,gnorm); CHKERRQ(ierr);
640       ierr = VecNormEnd(y,NORM_2,ynorm); CHKERRQ(ierr);
641       ierr = VecNormEnd(g,NORM_2,gnorm); CHKERRQ(ierr);
642     }
643   }
644   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
645   PetscFunctionReturn(0);
646 }
647 /* -------------------------------------------------------------------------- */
648 #undef __FUNC__
649 #define __FUNC__ "SNESSetLineSearch"
650 /*@C
651    SNESSetLineSearch - Sets the line search routine to be used
652    by the method SNES_EQ_LS.
653 
654    Input Parameters:
655 +  snes - nonlinear context obtained from SNESCreate()
656 .  lsctx - optional user-defined context for use by line search
657 -  func - pointer to int function
658 
659    Collective on SNES
660 
661    Available Routines:
662 +  SNESCubicLineSearch() - default line search
663 .  SNESQuadraticLineSearch() - quadratic line search
664 .  SNESNoLineSearch() - the full Newton step (actually not a line search)
665 -  SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel)
666 
667     Options Database Keys:
668 +   -snes_eq_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
669 .   -snes_eq_ls_alpha <alpha> - Sets alpha
670 .   -snes_eq_ls_maxstep <max> - Sets maxstep
671 -   -snes_eq_ls_steptol <steptol> - Sets steptol
672 
673    Calling sequence of func:
674 .vb
675    func (SNES snes, void *lsctx, Vec x, Vec f, Vec g, Vec y, Vec w,
676          double fnorm, double *ynorm, double *gnorm, *flag)
677 .ve
678 
679     Input parameters for func:
680 +   snes - nonlinear context
681 .   lsctx - optional user-defined context for line search
682 .   x - current iterate
683 .   f - residual evaluated at x
684 .   y - search direction (contains new iterate on output)
685 .   w - work vector
686 -   fnorm - 2-norm of f
687 
688     Output parameters for func:
689 +   g - residual evaluated at new iterate y
690 .   y - new iterate (contains search direction on input)
691 .   gnorm - 2-norm of g
692 .   ynorm - 2-norm of search length
693 -   flag - set to 0 if the line search succeeds; a nonzero integer
694            on failure.
695 
696     Level: advanced
697 
698 .keywords: SNES, nonlinear, set, line search, routine
699 
700 .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(), SNESSetLineSearchCheck()
701 @*/
702 int SNESSetLineSearch(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,double,double*,double*,int*),void *lsctx)
703 {
704   int ierr, (*f)(SNES,int (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,double,double*,double*,int*),void*);
705 
706   PetscFunctionBegin;
707   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void **)&f);CHKERRQ(ierr);
708   if (f) {
709     ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr);
710   }
711   PetscFunctionReturn(0);
712 }
713 /* -------------------------------------------------------------------------- */
714 EXTERN_C_BEGIN
715 #undef __FUNC__
716 #define __FUNC__ "SNESSetLineSearch_LS"
717 int SNESSetLineSearch_LS(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,
718                          double,double*,double*,int*),void *lsctx)
719 {
720   PetscFunctionBegin;
721   ((SNES_LS *)(snes->data))->LineSearch = func;
722   ((SNES_LS *)(snes->data))->lsP = lsctx;
723   PetscFunctionReturn(0);
724 }
725 EXTERN_C_END
726 /* -------------------------------------------------------------------------- */
727 #undef __FUNC__
728 #define __FUNC__ "SNESSetLineSearchCheck"
729 /*@C
730    SNESSetLineSearchCheck - Sets a routine to check the new iterate computed
731    by the line search routine in the Newton-based method SNES_EQ_LS.
732 
733    Input Parameters:
734 +  snes - nonlinear context obtained from SNESCreate()
735 .  func - pointer to int function
736 -  checkctx - optional user-defined context for use by step checking routine
737 
738    Collective on SNES
739 
740    Calling sequence of func:
741 .vb
742    func (SNES snes, void *lsctx, Vec x, PetscTruth *flag)
743 .ve
744 
745    Input parameters for func:
746 +  snes - nonlinear context
747 .  checkctx - optional user-defined context for use by step checking routine
748 -  x - current candidate iterate
749 
750    Output parameters for func:
751 +  x - current iterate (possibly modified)
752 -  flag - flag indicating whether x has been modified (either
753            PETSC_TRUE of PETSC_FALSE)
754 
755    Level: advanced
756 
757    Notes:
758    The user-defined line search checking routine is available for
759    use in conjunction with SNESNoLineSearch(), SNESNoLineSearchNoNorms(),
760    which (1) compute a candidate iterate u_{i+1}, (2) pass control to the
761    checking routine, and then (3) compute the corresponding function f(u_{i+1})
762    with the (possibly altered) iterate u_{i+1}.
763 
764    The user-defined line search checking routine is also available for
765    use in conjunction with SNESQuadraticLineSearch() and SNESCubicLineSearch().
766    These routines (1) compute a candidate iterate u_{i+1} as well as a
767    candidate function f(u_{i+1}), (2) pass control to the checking routine,
768    and then (3) force a re-evaluation of f(u_{i+1}) if any changes were made
769    to the candidate iterate in the checking routine (as indicated by
770    flag=PETSC_TRUE).  The overhead of this re-evaluation can be costly, so
771    this feature with caution.
772 
773 .keywords: SNES, nonlinear, set, line search check, step check, routine
774 
775 .seealso: SNESSetLineSearch()
776 @*/
777 int SNESSetLineSearchCheck(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx)
778 {
779   int ierr, (*f)(SNES,int (*)(SNES,void*,Vec,PetscTruth*),void*);
780 
781   PetscFunctionBegin;
782   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearchCheck_C",(void **)&f);CHKERRQ(ierr);
783   if (f) {
784     ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr);
785   }
786   PetscFunctionReturn(0);
787 }
788 /* -------------------------------------------------------------------------- */
789 EXTERN_C_BEGIN
790 #undef __FUNC__
791 #define __FUNC__ "SNESSetLineSearchCheck_LS"
792 int SNESSetLineSearchCheck_LS(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx)
793 {
794   PetscFunctionBegin;
795   ((SNES_LS *)(snes->data))->CheckStep = func;
796   ((SNES_LS *)(snes->data))->checkP = checkctx;
797   PetscFunctionReturn(0);
798 }
799 EXTERN_C_END
800 /* -------------------------------------------------------------------------- */
801 /*
802    SNESPrintHelp_EQ_LS - Prints all options for the SNES_EQ_LS method.
803 
804    Input Parameter:
805 .  snes - the SNES context
806 
807    Application Interface Routine: SNESPrintHelp()
808 */
809 #undef __FUNC__
810 #define __FUNC__ "SNESPrintHelp_EQ_LS"
811 static int SNESPrintHelp_EQ_LS(SNES snes,char *p)
812 {
813   SNES_LS *ls = (SNES_LS *)snes->data;
814 
815   PetscFunctionBegin;
816   (*PetscHelpPrintf)(snes->comm," method SNES_EQ_LS (ls) for systems of nonlinear equations:\n");
817   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls [cubic,quadratic,basic,basicnonorms,...]\n",p);
818   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_alpha <alpha> (default %g)\n",p,ls->alpha);
819   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_maxstep <max> (default %g)\n",p,ls->maxstep);
820   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_steptol <tol> (default %g)\n",p,ls->steptol);
821   PetscFunctionReturn(0);
822 }
823 /* -------------------------------------------------------------------------- */
824 /*
825    SNESView_EQ_LS - Prints info from the SNES_EQ_LS data structure.
826 
827    Input Parameters:
828 .  SNES - the SNES context
829 .  viewer - visualization context
830 
831    Application Interface Routine: SNESView()
832 */
833 #undef __FUNC__
834 #define __FUNC__ "SNESView_EQ_LS"
835 static int SNESView_EQ_LS(SNES snes,Viewer viewer)
836 {
837   SNES_LS    *ls = (SNES_LS *)snes->data;
838   char       *cstr;
839   int        ierr;
840   ViewerType vtype;
841 
842   PetscFunctionBegin;
843   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
844   if (PetscTypeCompare(vtype,ASCII_VIEWER)) {
845     if (ls->LineSearch == SNESNoLineSearch)             cstr = "SNESNoLineSearch";
846     else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch";
847     else if (ls->LineSearch == SNESCubicLineSearch)     cstr = "SNESCubicLineSearch";
848     else                                                cstr = "unknown";
849     ViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);
850     ViewerASCIIPrintf(viewer,"  alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);
851   } else {
852     SETERRQ(1,1,"Viewer type not supported for this object");
853   }
854   PetscFunctionReturn(0);
855 }
856 /* -------------------------------------------------------------------------- */
857 /*
858    SNESSetFromOptions_EQ_LS - Sets various parameters for the SNES_EQ_LS method.
859 
860    Input Parameter:
861 .  snes - the SNES context
862 
863    Application Interface Routine: SNESSetFromOptions()
864 */
865 #undef __FUNC__
866 #define __FUNC__ "SNESSetFromOptions_EQ_LS"
867 static int SNESSetFromOptions_EQ_LS(SNES snes)
868 {
869   SNES_LS *ls = (SNES_LS *)snes->data;
870   char    ver[16];
871   double  tmp;
872   int     ierr,flg;
873 
874   PetscFunctionBegin;
875   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_alpha",&tmp, &flg);CHKERRQ(ierr);
876   if (flg) {
877     ls->alpha = tmp;
878   }
879   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_maxstep",&tmp, &flg);CHKERRQ(ierr);
880   if (flg) {
881     ls->maxstep = tmp;
882   }
883   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_steptol",&tmp, &flg);CHKERRQ(ierr);
884   if (flg) {
885     ls->steptol = tmp;
886   }
887   ierr = OptionsGetString(snes->prefix,"-snes_eq_ls",ver,16, &flg); CHKERRQ(ierr);
888   if (flg) {
889     if (!PetscStrcmp(ver,"basic")) {
890       ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr);
891     } else if (!PetscStrcmp(ver,"basicnonorms")) {
892       ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
893     } else if (!PetscStrcmp(ver,"quadratic")) {
894       ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr);
895     } else if (!PetscStrcmp(ver,"cubic")) {
896       ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr);
897     }
898     else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown line search");}
899   }
900   PetscFunctionReturn(0);
901 }
902 /* -------------------------------------------------------------------------- */
903 /*
904    SNESCreate_EQ_LS - Creates a nonlinear solver context for the SNES_EQ_LS method,
905    SNES_LS, and sets this as the private data within the generic nonlinear solver
906    context, SNES, that was created within SNESCreate().
907 
908 
909    Input Parameter:
910 .  snes - the SNES context
911 
912    Application Interface Routine: SNESCreate()
913  */
914 EXTERN_C_BEGIN
915 #undef __FUNC__
916 #define __FUNC__ "SNESCreate_EQ_LS"
917 int SNESCreate_EQ_LS(SNES snes)
918 {
919   int     ierr;
920   SNES_LS *neP;
921 
922   PetscFunctionBegin;
923   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
924     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
925   }
926 
927   snes->setup		= SNESSetUp_EQ_LS;
928   snes->solve		= SNESSolve_EQ_LS;
929   snes->destroy		= SNESDestroy_EQ_LS;
930   snes->converged	= SNESConverged_EQ_LS;
931   snes->printhelp       = SNESPrintHelp_EQ_LS;
932   snes->setfromoptions  = SNESSetFromOptions_EQ_LS;
933   snes->view            = SNESView_EQ_LS;
934   snes->nwork           = 0;
935 
936   neP			= PetscNew(SNES_LS);   CHKPTRQ(neP);
937   PLogObjectMemory(snes,sizeof(SNES_LS));
938   snes->data    	= (void *) neP;
939   neP->alpha		= 1.e-4;
940   neP->maxstep		= 1.e8;
941   neP->steptol		= 1.e-12;
942   neP->LineSearch       = SNESCubicLineSearch;
943   neP->lsP              = PETSC_NULL;
944   neP->CheckStep        = PETSC_NULL;
945   neP->checkP           = PETSC_NULL;
946 
947   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS",
948                     (void*)SNESSetLineSearch_LS);CHKERRQ(ierr);
949   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS",
950                     (void*)SNESSetLineSearchCheck_LS);CHKERRQ(ierr);
951 
952   PetscFunctionReturn(0);
953 }
954 EXTERN_C_END
955 
956 
957 
958