xref: /petsc/src/snes/impls/ls/ls.c (revision 606d414c6e034e02e67059b83ebaefc3ebe99698)
1 #ifdef PETSC_RCS_HEADER
2 static char vcid[] = "$Id: ls.c,v 1.137 1999/06/08 22:53:15 bsmith 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    Notes:
290    SNESNoLineSearchNoNorms() must be used in conjunction with
291    either the options
292 $     -snes_no_convergence_test -snes_max_it <its>
293    or alternatively a user-defined custom test set via
294    SNESSetConvergenceTest(); otherwise, the SNES solver will
295    generate an error.
296 
297    The residual norms printed by monitoring routines such as
298    SNESDefaultMonitor() (as activated via -snes_monitor) will not be
299    correct, since they are not computed.
300 
301    Level: advanced
302 
303 .keywords: SNES, nonlinear, line search, cubic
304 
305 .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
306           SNESSetLineSearch(), SNESNoLineSearch()
307 @*/
308 int SNESNoLineSearchNoNorms(SNES snes, void *lsctx, Vec x, Vec f, Vec g, Vec y, Vec w,
309                      double fnorm, double *ynorm, double *gnorm,int *flag )
310 {
311   int        ierr;
312   Scalar     mone = -1.0;
313   SNES_LS    *neP = (SNES_LS *) snes->data;
314   PetscTruth change_y = PETSC_FALSE;
315 
316   PetscFunctionBegin;
317   *flag = 0;
318   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
319   ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr);            /* y <- y - x      */
320   if (neP->CheckStep) {
321    ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
322   }
323   ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y)    */
324   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
325   PetscFunctionReturn(0);
326 }
327 /* -------------------------------------------------------------------------- */
328 #undef __FUNC__
329 #define __FUNC__ "SNESCubicLineSearch"
330 /*@C
331    SNESCubicLineSearch - Performs a cubic line search (default line search method).
332 
333    Collective on SNES
334 
335    Input Parameters:
336 +  snes - nonlinear context
337 .  lsctx - optional context for line search (not used here)
338 .  x - current iterate
339 .  f - residual evaluated at x
340 .  y - search direction (contains new iterate on output)
341 .  w - work vector
342 -  fnorm - 2-norm of f
343 
344    Output Parameters:
345 +  g - residual evaluated at new iterate y
346 .  y - new iterate (contains search direction on input)
347 .  gnorm - 2-norm of g
348 .  ynorm - 2-norm of search length
349 -  flag - 0 if line search succeeds; -1 on failure.
350 
351    Options Database Key:
352 $  -snes_eq_ls cubic - Activates SNESCubicLineSearch()
353 
354    Notes:
355    This line search is taken from "Numerical Methods for Unconstrained
356    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
357 
358    Level: advanced
359 
360 .keywords: SNES, nonlinear, line search, cubic
361 
362 .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
363 @*/
364 int SNESCubicLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
365                         double fnorm,double *ynorm,double *gnorm,int *flag)
366 {
367   /*
368      Note that for line search purposes we work with with the related
369      minimization problem:
370         min  z(x):  R^n -> R,
371      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
372    */
373 
374   double     steptol, initslope, lambdaprev, gnormprev, a, b, d, t1, t2;
375   double     maxstep, minlambda, alpha, lambda, lambdatemp, lambdaneg;
376 #if defined(PETSC_USE_COMPLEX)
377   Scalar     cinitslope, clambda;
378 #endif
379   int        ierr, count;
380   SNES_LS    *neP = (SNES_LS *) snes->data;
381   Scalar     mone = -1.0, scale;
382   PetscTruth change_y = PETSC_FALSE;
383 
384   PetscFunctionBegin;
385   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
386   *flag   = 0;
387   alpha   = neP->alpha;
388   maxstep = neP->maxstep;
389   steptol = neP->steptol;
390 
391   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
392   if (*ynorm < snes->atol) {
393     PLogInfo(snes,"SNESCubicLineSearch: Search direction and size are nearly 0\n");
394     *gnorm = fnorm;
395     ierr = VecCopy(x,y);CHKERRQ(ierr);
396     ierr = VecCopy(f,g);CHKERRQ(ierr);
397     goto theend1;
398   }
399   if (*ynorm > maxstep) {	/* Step too big, so scale back */
400     scale = maxstep/(*ynorm);
401 #if defined(PETSC_USE_COMPLEX)
402     PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",PetscReal(scale));
403 #else
404     PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",scale);
405 #endif
406     ierr = VecScale(&scale,y);CHKERRQ(ierr);
407     *ynorm = maxstep;
408   }
409   minlambda = steptol/(*ynorm);
410   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
411 #if defined(PETSC_USE_COMPLEX)
412   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
413   initslope = PetscReal(cinitslope);
414 #else
415   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
416 #endif
417   if (initslope > 0.0) initslope = -initslope;
418   if (initslope == 0.0) initslope = -1.0;
419 
420   ierr = VecCopy(y,w);CHKERRQ(ierr);
421   ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr);
422   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
423   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
424   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */
425     ierr = VecCopy(w,y);CHKERRQ(ierr);
426     PLogInfo(snes,"SNESCubicLineSearch: Using full step\n");
427     goto theend1;
428   }
429 
430   /* Fit points with quadratic */
431   lambda = 1.0;
432   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
433   lambdaprev = lambda;
434   gnormprev = *gnorm;
435   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
436   else lambda = lambdatemp;
437   ierr   = VecCopy(x,w);CHKERRQ(ierr);
438   lambdaneg = -lambda;
439 #if defined(PETSC_USE_COMPLEX)
440   clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
441 #else
442   ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
443 #endif
444   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
445   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
446   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */
447     ierr = VecCopy(w,y);CHKERRQ(ierr);
448     PLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%g\n",lambda);
449     goto theend1;
450   }
451 
452   /* Fit points with cubic */
453   count = 1;
454   while (1) {
455     if (lambda <= minlambda) { /* bad luck; use full step */
456       PLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count);
457       PLogInfo(snes,"SNESCubicLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",
458                fnorm,*gnorm,*ynorm,lambda,initslope);
459       ierr = VecCopy(w,y);CHKERRQ(ierr);
460       *flag = -1; break;
461     }
462     t1 = ((*gnorm)*(*gnorm) - fnorm*fnorm)*0.5 - lambda*initslope;
463     t2 = (gnormprev*gnormprev  - fnorm*fnorm)*0.5 - lambdaprev*initslope;
464     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
465     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
466     d  = b*b - 3*a*initslope;
467     if (d < 0.0) d = 0.0;
468     if (a == 0.0) {
469       lambdatemp = -initslope/(2.0*b);
470     } else {
471       lambdatemp = (-b + sqrt(d))/(3.0*a);
472     }
473     if (lambdatemp > .5*lambda) {
474       lambdatemp = .5*lambda;
475     }
476     lambdaprev = lambda;
477     gnormprev = *gnorm;
478     if (lambdatemp <= .1*lambda) {
479       lambda = .1*lambda;
480     }
481     else lambda = lambdatemp;
482     ierr = VecCopy( x, w );CHKERRQ(ierr);
483     lambdaneg = -lambda;
484 #if defined(PETSC_USE_COMPLEX)
485     clambda = lambdaneg;
486     ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
487 #else
488     ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
489 #endif
490     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
491     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
492     if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* is reduction enough? */
493       ierr = VecCopy(w,y);CHKERRQ(ierr);
494       PLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%g\n",lambda);
495       goto theend1;
496     } else {
497       PLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda,  lambda=%g\n",lambda);
498     }
499     count++;
500   }
501   theend1:
502   /* Optional user-defined check for line search step validity */
503   if (neP->CheckStep) {
504     ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
505     if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */
506       ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr);
507       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
508       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
509       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
510       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
511     }
512   }
513   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
514   PetscFunctionReturn(0);
515 }
516 /* -------------------------------------------------------------------------- */
517 #undef __FUNC__
518 #define __FUNC__ "SNESQuadraticLineSearch"
519 /*@C
520    SNESQuadraticLineSearch - Performs a quadratic line search.
521 
522    Collective on SNES and Vec
523 
524    Input Parameters:
525 +  snes - the SNES context
526 .  lsctx - optional context for line search (not used here)
527 .  x - current iterate
528 .  f - residual evaluated at x
529 .  y - search direction (contains new iterate on output)
530 .  w - work vector
531 -  fnorm - 2-norm of f
532 
533    Output Parameters:
534 +  g - residual evaluated at new iterate y
535 .  y - new iterate (contains search direction on input)
536 .  gnorm - 2-norm of g
537 .  ynorm - 2-norm of search length
538 -  flag - 0 if line search succeeds; -1 on failure.
539 
540    Options Database Key:
541 .  -snes_eq_ls quadratic - Activates SNESQuadraticLineSearch()
542 
543    Notes:
544    Use SNESSetLineSearch() to set this routine within the SNES_EQ_LS method.
545 
546    Level: advanced
547 
548 .keywords: SNES, nonlinear, quadratic, line search
549 
550 .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
551 @*/
552 int SNESQuadraticLineSearch(SNES snes, void *lsctx, Vec x, Vec f, Vec g, Vec y, Vec w,
553                            double fnorm, double *ynorm, double *gnorm,int *flag)
554 {
555   /*
556      Note that for line search purposes we work with with the related
557      minimization problem:
558         min  z(x):  R^n -> R,
559      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
560    */
561   double     steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp;
562 #if defined(PETSC_USE_COMPLEX)
563   Scalar     cinitslope,clambda;
564 #endif
565   int        ierr, count;
566   SNES_LS    *neP = (SNES_LS *) snes->data;
567   Scalar     mone = -1.0,scale;
568   PetscTruth change_y = PETSC_FALSE;
569 
570   PetscFunctionBegin;
571   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
572   *flag   = 0;
573   alpha   = neP->alpha;
574   maxstep = neP->maxstep;
575   steptol = neP->steptol;
576 
577   VecNorm(y, NORM_2,ynorm );
578   if (*ynorm < snes->atol) {
579     PLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n");
580     *gnorm = fnorm;
581     ierr = VecCopy(x,y);CHKERRQ(ierr);
582     ierr = VecCopy(f,g);CHKERRQ(ierr);
583     goto theend2;
584   }
585   if (*ynorm > maxstep) {	/* Step too big, so scale back */
586     scale = maxstep/(*ynorm);
587     ierr = VecScale(&scale,y);CHKERRQ(ierr);
588     *ynorm = maxstep;
589   }
590   minlambda = steptol/(*ynorm);
591   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
592 #if defined(PETSC_USE_COMPLEX)
593   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
594   initslope = PetscReal(cinitslope);
595 #else
596   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
597 #endif
598   if (initslope > 0.0) initslope = -initslope;
599   if (initslope == 0.0) initslope = -1.0;
600 
601   ierr = VecCopy(y,w);CHKERRQ(ierr);
602   ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr);
603   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
604   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
605   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */
606     ierr = VecCopy(w,y);CHKERRQ(ierr);
607     PLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n");
608     goto theend2;
609   }
610 
611   /* Fit points with quadratic */
612   lambda = 1.0;
613   count = 1;
614   while (1) {
615     if (lambda <= minlambda) { /* bad luck; use full step */
616       PLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count);
617       PLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",
618                fnorm,*gnorm,*ynorm,lambda,initslope);
619       ierr = VecCopy(w,y);CHKERRQ(ierr);
620       *flag = -1; break;
621     }
622     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
623     if (lambdatemp <= .1*lambda) {
624       lambda = .1*lambda;
625     } else lambda = lambdatemp;
626     ierr = VecCopy(x,w);CHKERRQ(ierr);
627     lambda = -lambda;
628 #if defined(PETSC_USE_COMPLEX)
629     clambda = lambda; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
630 #else
631     ierr = VecAXPY(&lambda,y,w);CHKERRQ(ierr);
632 #endif
633     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
634     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
635     if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */
636       ierr = VecCopy(w,y);CHKERRQ(ierr);
637       PLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda);
638       break;
639     }
640     count++;
641   }
642   theend2:
643   /* Optional user-defined check for line search step validity */
644   if (neP->CheckStep) {
645     ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
646     if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */
647       ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr);
648       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
649       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
650       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
651       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
652     }
653   }
654   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
655   PetscFunctionReturn(0);
656 }
657 /* -------------------------------------------------------------------------- */
658 #undef __FUNC__
659 #define __FUNC__ "SNESSetLineSearch"
660 /*@C
661    SNESSetLineSearch - Sets the line search routine to be used
662    by the method SNES_EQ_LS.
663 
664    Input Parameters:
665 +  snes - nonlinear context obtained from SNESCreate()
666 .  lsctx - optional user-defined context for use by line search
667 -  func - pointer to int function
668 
669    Collective on SNES
670 
671    Available Routines:
672 +  SNESCubicLineSearch() - default line search
673 .  SNESQuadraticLineSearch() - quadratic line search
674 .  SNESNoLineSearch() - the full Newton step (actually not a line search)
675 -  SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel)
676 
677     Options Database Keys:
678 +   -snes_eq_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
679 .   -snes_eq_ls_alpha <alpha> - Sets alpha
680 .   -snes_eq_ls_maxstep <max> - Sets maxstep
681 -   -snes_eq_ls_steptol <steptol> - Sets steptol
682 
683    Calling sequence of func:
684 .vb
685    func (SNES snes, void *lsctx, Vec x, Vec f, Vec g, Vec y, Vec w,
686          double fnorm, double *ynorm, double *gnorm, *flag)
687 .ve
688 
689     Input parameters for func:
690 +   snes - nonlinear context
691 .   lsctx - optional user-defined context for line search
692 .   x - current iterate
693 .   f - residual evaluated at x
694 .   y - search direction (contains new iterate on output)
695 .   w - work vector
696 -   fnorm - 2-norm of f
697 
698     Output parameters for func:
699 +   g - residual evaluated at new iterate y
700 .   y - new iterate (contains search direction on input)
701 .   gnorm - 2-norm of g
702 .   ynorm - 2-norm of search length
703 -   flag - set to 0 if the line search succeeds; a nonzero integer
704            on failure.
705 
706     Level: advanced
707 
708 .keywords: SNES, nonlinear, set, line search, routine
709 
710 .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(), SNESSetLineSearchCheck(), SNESSetLineSearchParams(), SNESGetLineSearchParams()
711 @*/
712 int SNESSetLineSearch(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,double,double*,double*,int*),void *lsctx)
713 {
714   int ierr, (*f)(SNES,int (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,double,double*,double*,int*),void*);
715 
716   PetscFunctionBegin;
717   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void **)&f);CHKERRQ(ierr);
718   if (f) {
719     ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr);
720   }
721   PetscFunctionReturn(0);
722 }
723 /* -------------------------------------------------------------------------- */
724 EXTERN_C_BEGIN
725 #undef __FUNC__
726 #define __FUNC__ "SNESSetLineSearch_LS"
727 int SNESSetLineSearch_LS(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,
728                          double,double*,double*,int*),void *lsctx)
729 {
730   PetscFunctionBegin;
731   ((SNES_LS *)(snes->data))->LineSearch = func;
732   ((SNES_LS *)(snes->data))->lsP = lsctx;
733   PetscFunctionReturn(0);
734 }
735 EXTERN_C_END
736 /* -------------------------------------------------------------------------- */
737 #undef __FUNC__
738 #define __FUNC__ "SNESSetLineSearchCheck"
739 /*@C
740    SNESSetLineSearchCheck - Sets a routine to check the validity of new iterate computed
741    by the line search routine in the Newton-based method SNES_EQ_LS.
742 
743    Input Parameters:
744 +  snes - nonlinear context obtained from SNESCreate()
745 .  func - pointer to int function
746 -  checkctx - optional user-defined context for use by step checking routine
747 
748    Collective on SNES
749 
750    Calling sequence of func:
751 .vb
752    int func (SNES snes, void *checkctx, Vec x, PetscTruth *flag)
753 .ve
754    where func returns an error code of 0 on success and a nonzero
755    on failure.
756 
757    Input parameters for func:
758 +  snes - nonlinear context
759 .  checkctx - optional user-defined context for use by step checking routine
760 -  x - current candidate iterate
761 
762    Output parameters for func:
763 +  x - current iterate (possibly modified)
764 -  flag - flag indicating whether x has been modified (either
765            PETSC_TRUE of PETSC_FALSE)
766 
767    Level: advanced
768 
769    Notes:
770    SNESNoLineSearch() and SNESNoLineSearchNoNorms() accept the new
771    iterate computed by the line search checking routine.  In particular,
772    these routines (1) compute a candidate iterate u_{i+1}, (2) pass control
773    to the checking routine, and then (3) compute the corresponding nonlinear
774    function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}.
775 
776    SNESQuadraticLineSearch() and SNESCubicLineSearch() also accept the
777    new iterate computed by the line search checking routine.  In particular,
778    these routines (1) compute a candidate iterate u_{i+1} as well as a
779    candidate nonlinear function f(u_{i+1}), (2) pass control to the checking
780    routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes
781    were made to the candidate iterate in the checking routine (as indicated
782    by flag=PETSC_TRUE).  The overhead of this function re-evaluation can be
783    very costly, so use this feature with caution!
784 
785 .keywords: SNES, nonlinear, set, line search check, step check, routine
786 
787 .seealso: SNESSetLineSearch()
788 @*/
789 int SNESSetLineSearchCheck(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx)
790 {
791   int ierr, (*f)(SNES,int (*)(SNES,void*,Vec,PetscTruth*),void*);
792 
793   PetscFunctionBegin;
794   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearchCheck_C",(void **)&f);CHKERRQ(ierr);
795   if (f) {
796     ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr);
797   }
798   PetscFunctionReturn(0);
799 }
800 /* -------------------------------------------------------------------------- */
801 EXTERN_C_BEGIN
802 #undef __FUNC__
803 #define __FUNC__ "SNESSetLineSearchCheck_LS"
804 int SNESSetLineSearchCheck_LS(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx)
805 {
806   PetscFunctionBegin;
807   ((SNES_LS *)(snes->data))->CheckStep = func;
808   ((SNES_LS *)(snes->data))->checkP = checkctx;
809   PetscFunctionReturn(0);
810 }
811 EXTERN_C_END
812 /* -------------------------------------------------------------------------- */
813 /*
814    SNESPrintHelp_EQ_LS - Prints all options for the SNES_EQ_LS method.
815 
816    Input Parameter:
817 .  snes - the SNES context
818 
819    Application Interface Routine: SNESPrintHelp()
820 */
821 #undef __FUNC__
822 #define __FUNC__ "SNESPrintHelp_EQ_LS"
823 static int SNESPrintHelp_EQ_LS(SNES snes,char *p)
824 {
825   SNES_LS *ls = (SNES_LS *)snes->data;
826   int     ierr;
827 
828   PetscFunctionBegin;
829   ierr = (*PetscHelpPrintf)(snes->comm," method SNES_EQ_LS (ls) for systems of nonlinear equations:\n");CHKERRQ(ierr);
830   ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls [cubic,quadratic,basic,basicnonorms,...]\n",p);CHKERRQ(ierr);
831   ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_alpha <alpha> (default %g)\n",p,ls->alpha);CHKERRQ(ierr);
832   ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_maxstep <max> (default %g)\n",p,ls->maxstep);CHKERRQ(ierr);
833   ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_steptol <tol> (default %g)\n",p,ls->steptol);CHKERRQ(ierr);
834   PetscFunctionReturn(0);
835 }
836 /* -------------------------------------------------------------------------- */
837 /*
838    SNESView_EQ_LS - Prints info from the SNES_EQ_LS data structure.
839 
840    Input Parameters:
841 .  SNES - the SNES context
842 .  viewer - visualization context
843 
844    Application Interface Routine: SNESView()
845 */
846 #undef __FUNC__
847 #define __FUNC__ "SNESView_EQ_LS"
848 static int SNESView_EQ_LS(SNES snes,Viewer viewer)
849 {
850   SNES_LS    *ls = (SNES_LS *)snes->data;
851   char       *cstr;
852   int        ierr;
853   ViewerType vtype;
854 
855   PetscFunctionBegin;
856   ierr = ViewerGetType(viewer,&vtype);CHKERRQ(ierr);
857   if (PetscTypeCompare(vtype,ASCII_VIEWER)) {
858     if (ls->LineSearch == SNESNoLineSearch)             cstr = "SNESNoLineSearch";
859     else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch";
860     else if (ls->LineSearch == SNESCubicLineSearch)     cstr = "SNESCubicLineSearch";
861     else                                                cstr = "unknown";
862     ierr = ViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
863     ierr = ViewerASCIIPrintf(viewer,"  alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr);
864   } else {
865     SETERRQ(1,1,"Viewer type not supported for this object");
866   }
867   PetscFunctionReturn(0);
868 }
869 /* -------------------------------------------------------------------------- */
870 /*
871    SNESSetFromOptions_EQ_LS - Sets various parameters for the SNES_EQ_LS method.
872 
873    Input Parameter:
874 .  snes - the SNES context
875 
876    Application Interface Routine: SNESSetFromOptions()
877 */
878 #undef __FUNC__
879 #define __FUNC__ "SNESSetFromOptions_EQ_LS"
880 static int SNESSetFromOptions_EQ_LS(SNES snes)
881 {
882   SNES_LS *ls = (SNES_LS *)snes->data;
883   char    ver[16];
884   double  tmp;
885   int     ierr,flg;
886 
887   PetscFunctionBegin;
888   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_alpha",&tmp, &flg);CHKERRQ(ierr);
889   if (flg) {
890     ls->alpha = tmp;
891   }
892   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_maxstep",&tmp, &flg);CHKERRQ(ierr);
893   if (flg) {
894     ls->maxstep = tmp;
895   }
896   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_steptol",&tmp, &flg);CHKERRQ(ierr);
897   if (flg) {
898     ls->steptol = tmp;
899   }
900   ierr = OptionsGetString(snes->prefix,"-snes_eq_ls",ver,16, &flg);CHKERRQ(ierr);
901   if (flg) {
902     if (!PetscStrcmp(ver,"basic")) {
903       ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr);
904     } else if (!PetscStrcmp(ver,"basicnonorms")) {
905       ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
906     } else if (!PetscStrcmp(ver,"quadratic")) {
907       ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr);
908     } else if (!PetscStrcmp(ver,"cubic")) {
909       ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr);
910     }
911     else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown line search");}
912   }
913   PetscFunctionReturn(0);
914 }
915 /* -------------------------------------------------------------------------- */
916 /*
917    SNESCreate_EQ_LS - Creates a nonlinear solver context for the SNES_EQ_LS method,
918    SNES_LS, and sets this as the private data within the generic nonlinear solver
919    context, SNES, that was created within SNESCreate().
920 
921 
922    Input Parameter:
923 .  snes - the SNES context
924 
925    Application Interface Routine: SNESCreate()
926  */
927 EXTERN_C_BEGIN
928 #undef __FUNC__
929 #define __FUNC__ "SNESCreate_EQ_LS"
930 int SNESCreate_EQ_LS(SNES snes)
931 {
932   int     ierr;
933   SNES_LS *neP;
934 
935   PetscFunctionBegin;
936   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
937     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
938   }
939 
940   snes->setup		= SNESSetUp_EQ_LS;
941   snes->solve		= SNESSolve_EQ_LS;
942   snes->destroy		= SNESDestroy_EQ_LS;
943   snes->converged	= SNESConverged_EQ_LS;
944   snes->printhelp       = SNESPrintHelp_EQ_LS;
945   snes->setfromoptions  = SNESSetFromOptions_EQ_LS;
946   snes->view            = SNESView_EQ_LS;
947   snes->nwork           = 0;
948 
949   neP			= PetscNew(SNES_LS);CHKPTRQ(neP);
950   PLogObjectMemory(snes,sizeof(SNES_LS));
951   snes->data    	= (void *) neP;
952   neP->alpha		= 1.e-4;
953   neP->maxstep		= 1.e8;
954   neP->steptol		= 1.e-12;
955   neP->LineSearch       = SNESCubicLineSearch;
956   neP->lsP              = PETSC_NULL;
957   neP->CheckStep        = PETSC_NULL;
958   neP->checkP           = PETSC_NULL;
959 
960   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS",
961                     (void*)SNESSetLineSearch_LS);CHKERRQ(ierr);
962   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS",
963                     (void*)SNESSetLineSearchCheck_LS);CHKERRQ(ierr);
964 
965   PetscFunctionReturn(0);
966 }
967 EXTERN_C_END
968 
969 
970 
971