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