xref: /petsc/src/snes/impls/ls/ls.c (revision 76be9ce4a233aaa47cda2bc7f5a27cd7faabecaa)
1 #ifdef PETSC_RCS_HEADER
2 static char vcid[] = "$Id: ls.c,v 1.98 1998/01/06 20:12:26 bsmith Exp bsmith $";
3 #endif
4 
5 #include <math.h>
6 #include "src/snes/impls/ls/ls.h"
7 #include "pinclude/pviewer.h"
8 
9 /*
10      Implements a line search variant of Newton's Method
11     for solving systems of nonlinear equations.
12 
13     Input parameters:
14 .   snes - nonlinear context obtained from SNESCreate()
15 
16     Output Parameters:
17 .   outits  - Number of global iterations until termination.
18 
19     Notes:
20     This implements essentially a truncated Newton method with a
21     line search.  By default a cubic backtracking line search
22     is employed, as described in the text "Numerical Methods for
23     Unconstrained Optimization and Nonlinear Equations" by Dennis
24     and Schnabel.
25 */
26 
27 #undef __FUNC__
28 #define __FUNC__ "SNESSolve_EQ_LS"
29 int SNESSolve_EQ_LS(SNES snes,int *outits)
30 {
31   SNES_LS       *neP = (SNES_LS *) snes->data;
32   int           maxits, i, history_len, ierr, lits, lsfail;
33   MatStructure  flg = DIFFERENT_NONZERO_PATTERN;
34   double        fnorm, gnorm, xnorm, ynorm, *history;
35   Vec           Y, X, F, G, W, TMP;
36 
37   PetscFunctionBegin;
38   history	= snes->conv_hist;	/* convergence history */
39   history_len	= snes->conv_hist_size;	/* convergence history length */
40   maxits	= snes->max_its;	/* maximum number of iterations */
41   X		= snes->vec_sol;	/* solution vector */
42   F		= snes->vec_func;	/* residual vector */
43   Y		= snes->work[0];	/* work vectors */
44   G		= snes->work[1];
45   W		= snes->work[2];
46 
47   snes->iter = 0;
48   ierr = SNESComputeFunction(snes,X,F); CHKERRQ(ierr);  /*  F(X)      */
49   ierr = VecNorm(F,NORM_2,&fnorm); CHKERRQ(ierr);	/* fnorm <- ||F||  */
50   snes->norm = fnorm;
51   if (history) history[0] = fnorm;
52   SNESMonitor(snes,0,fnorm);
53 
54   if (fnorm < snes->atol) {*outits = 0; PetscFunctionReturn(0);}
55 
56   /* set parameter for default relative tolerance convergence test */
57   snes->ttol = fnorm*snes->rtol;
58 
59   for ( i=0; i<maxits; i++ ) {
60     snes->iter = i+1;
61 
62     /* Solve J Y = F, where J is Jacobian matrix */
63     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg); CHKERRQ(ierr);
64     ierr = SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,flg); CHKERRQ(ierr);
65     ierr = SLESSolve(snes->sles,F,Y,&lits); CHKERRQ(ierr);
66     snes->linear_its += PetscAbsInt(lits);
67     PLogInfo(snes,"SNESSolve_EQ_LS: iter=%d, linear solve iterations=%d\n",snes->iter,lits);
68 
69     /* Compute a (scaled) negative update in the line search routine:
70          Y <- X - lambda*Y
71        and evaluate G(Y) = function(Y))
72     */
73     ierr = VecCopy(Y,snes->vec_sol_update_always); CHKERRQ(ierr);
74     ierr = (*neP->LineSearch)(snes,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lsfail); CHKERRQ(ierr);
75     PLogInfo(snes,"SNESSolve_EQ_LS: fnorm=%g, gnorm=%g, ynorm=%g, lsfail=%d\n",fnorm,gnorm,ynorm,lsfail);
76     if (lsfail) snes->nfailures++;
77 
78     TMP = F; F = G; snes->vec_func_always = F; G = TMP;
79     TMP = X; X = Y; snes->vec_sol_always = X;  Y = TMP;
80     fnorm = gnorm;
81 
82     snes->norm = fnorm;
83     if (history && history_len > i+1) history[i+1] = fnorm;
84     SNESMonitor(snes,i+1,fnorm);
85 
86     /* Test for convergence */
87     if (snes->converged) {
88       ierr = VecNorm(X,NORM_2,&xnorm); CHKERRQ(ierr);	/* xnorm = || X || */
89       if ((*snes->converged)(snes,xnorm,ynorm,fnorm,snes->cnvP)) {
90         break;
91       }
92     }
93   }
94   if (X != snes->vec_sol) {
95     ierr = VecCopy(X,snes->vec_sol); CHKERRQ(ierr);
96     snes->vec_sol_always  = snes->vec_sol;
97     snes->vec_func_always = snes->vec_func;
98   }
99   if (i == maxits) {
100     PLogInfo(snes,"SNESSolve_EQ_LS: Maximum number of iterations has been reached: %d\n",maxits);
101     i--;
102   }
103   if (history) snes->conv_act_size = (history_len < i+1) ? history_len : i+1;
104   *outits = i+1;
105   PetscFunctionReturn(0);
106 }
107 /* ------------------------------------------------------------ */
108 #undef __FUNC__
109 #define __FUNC__ "SNESSetUp_EQ_LS"
110 int SNESSetUp_EQ_LS(SNES snes )
111 {
112   int ierr;
113 
114   PetscFunctionBegin;
115   snes->nwork = 4;
116   ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
117   PLogObjectParents(snes,snes->nwork,snes->work);
118   snes->vec_sol_update_always = snes->work[3];
119   PetscFunctionReturn(0);
120 }
121 /* ------------------------------------------------------------ */
122 #undef __FUNC__
123 #define __FUNC__ "SNESDestroy_EQ_LS"
124 int SNESDestroy_EQ_LS(PetscObject obj)
125 {
126   SNES snes = (SNES) obj;
127   int  ierr;
128 
129   PetscFunctionBegin;
130   if (snes->nwork) {
131     ierr = VecDestroyVecs(snes->work,snes->nwork); CHKERRQ(ierr);
132   }
133   PetscFree(snes->data);
134   PetscFunctionReturn(0);
135 }
136 /* ------------------------------------------------------------ */
137 #undef __FUNC__
138 #define __FUNC__ "SNESNoLineSearch"
139 /*ARGSUSED*/
140 /*@C
141    SNESNoLineSearch - This routine is not a line search at all;
142    it simply uses the full Newton step.  Thus, this routine is intended
143    to serve as a template and is not recommended for general use.
144 
145    Input Parameters:
146 .  snes - nonlinear context
147 .  x - current iterate
148 .  f - residual evaluated at x
149 .  y - search direction (contains new iterate on output)
150 .  w - work vector
151 .  fnorm - 2-norm of f
152 
153    Output Parameters:
154 .  g - residual evaluated at new iterate y
155 .  y - new iterate (contains search direction on input)
156 .  gnorm - 2-norm of g
157 .  ynorm - 2-norm of search length
158 .  flag - set to 0, indicating a successful line search
159 
160    Options Database Key:
161 $  -snes_eq_ls basic
162 
163 .keywords: SNES, nonlinear, line search, cubic
164 
165 .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
166           SNESSetLineSearch()
167 @*/
168 int SNESNoLineSearch(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w,
169                      double fnorm, double *ynorm, double *gnorm,int *flag )
170 {
171   int    ierr;
172   Scalar mone = -1.0;
173 
174   PetscFunctionBegin;
175   *flag = 0;
176   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
177   ierr = VecNorm(y,NORM_2,ynorm); CHKERRQ(ierr);       /* ynorm = || y || */
178   ierr = VecAYPX(&mone,x,y); CHKERRQ(ierr);            /* y <- y - x      */
179   ierr = SNESComputeFunction(snes,y,g); CHKERRQ(ierr); /* Compute F(y)    */
180   ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);       /* gnorm = || g || */
181   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
182   PetscFunctionReturn(0);
183 }
184 /* ------------------------------------------------------------------ */
185 #undef __FUNC__
186 #define __FUNC__ "SNESNoLineSearchNoNorms"
187 /*ARGSUSED*/
188 /*@C
189    SNESNoLineSearchNoNorms - This routine is not a line search at
190    all; it simply uses the full Newton step. This version does not
191    even compute the norm of the function or search direction; this
192    is intended only when you know the full step is fine and are
193    not checking for convergence of the nonlinear iteration (for
194    example, you are running always for a fixed number of Newton
195    steps).
196 
197    Input Parameters:
198 .  snes - nonlinear context
199 .  x - current iterate
200 .  f - residual evaluated at x
201 .  y - search direction (contains new iterate on output)
202 .  w - work vector
203 .  fnorm - 2-norm of f
204 
205    Output Parameters:
206 .  g - residual evaluated at new iterate y
207 .  gnorm - not changed
208 .  ynorm - not changed
209 .  flag - set to 0, indicating a successful line search
210 
211    Options Database Key:
212 $  -snes_eq_ls basicnonorms
213 
214 .keywords: SNES, nonlinear, line search, cubic
215 
216 .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
217           SNESSetLineSearch(), SNESNoLineSearch()
218 @*/
219 int SNESNoLineSearchNoNorms(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w,
220                      double fnorm, double *ynorm, double *gnorm,int *flag )
221 {
222   int    ierr;
223   Scalar mone = -1.0;
224 
225   PetscFunctionBegin;
226   *flag = 0;
227   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
228   ierr = VecAYPX(&mone,x,y); CHKERRQ(ierr);            /* y <- y - x      */
229   ierr = SNESComputeFunction(snes,y,g); CHKERRQ(ierr); /* Compute F(y)    */
230   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
231   PetscFunctionReturn(0);
232 }
233 /* ------------------------------------------------------------------ */
234 #undef __FUNC__
235 #define __FUNC__ "SNESCubicLineSearch"
236 /*@C
237    SNESCubicLineSearch - Performs a cubic line search (default line search method).
238 
239    Input Parameters:
240 .  snes - nonlinear context
241 .  x - current iterate
242 .  f - residual evaluated at x
243 .  y - search direction (contains new iterate on output)
244 .  w - work vector
245 .  fnorm - 2-norm of f
246 
247    Output Parameters:
248 .  g - residual evaluated at new iterate y
249 .  y - new iterate (contains search direction on input)
250 .  gnorm - 2-norm of g
251 .  ynorm - 2-norm of search length
252 .  flag - 0 if line search succeeds; -1 on failure.
253 
254    Options Database Key:
255 $  -snes_eq_ls cubic
256 
257    Notes:
258    This line search is taken from "Numerical Methods for Unconstrained
259    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
260 
261 .keywords: SNES, nonlinear, line search, cubic
262 
263 .seealso: SNESNoLineSearch(), SNESNoLineSearch(), SNESSetLineSearch()
264 @*/
265 int SNESCubicLineSearch(SNES snes,Vec x,Vec f,Vec g,Vec y,Vec w,
266                         double fnorm,double *ynorm,double *gnorm,int *flag)
267 {
268   /*
269      Note that for line search purposes we work with with the related
270      minimization problem:
271         min  z(x):  R^n -> R,
272      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
273    */
274 
275   double  steptol, initslope, lambdaprev, gnormprev, a, b, d, t1, t2;
276   double  maxstep, minlambda, alpha, lambda, lambdatemp, lambdaneg;
277 #if defined(USE_PETSC_COMPLEX)
278   Scalar  cinitslope, clambda;
279 #endif
280   int     ierr, count;
281   SNES_LS *neP = (SNES_LS *) snes->data;
282   Scalar  mone = -1.0,scale;
283 
284   PetscFunctionBegin;
285   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
286   *flag   = 0;
287   alpha   = neP->alpha;
288   maxstep = neP->maxstep;
289   steptol = neP->steptol;
290 
291   ierr = VecNorm(y,NORM_2,ynorm); CHKERRQ(ierr);
292   if (*ynorm < snes->atol) {
293     PLogInfo(snes,"SNESCubicLineSearch: Search direction and size are nearly 0\n");
294     *gnorm = fnorm;
295     ierr = VecCopy(x,y); CHKERRQ(ierr);
296     ierr = VecCopy(f,g); CHKERRQ(ierr);
297     goto theend1;
298   }
299   if (*ynorm > maxstep) {	/* Step too big, so scale back */
300     scale = maxstep/(*ynorm);
301 #if defined(USE_PETSC_COMPLEX)
302     PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",real(scale));
303 #else
304     PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",scale);
305 #endif
306     ierr = VecScale(&scale,y); CHKERRQ(ierr);
307     *ynorm = maxstep;
308   }
309   minlambda = steptol/(*ynorm);
310   ierr = MatMult(snes->jacobian,y,w); CHKERRQ(ierr);
311 #if defined(USE_PETSC_COMPLEX)
312   ierr = VecDot(f,w,&cinitslope); CHKERRQ(ierr);
313   initslope = real(cinitslope);
314 #else
315   ierr = VecDot(f,w,&initslope); CHKERRQ(ierr);
316 #endif
317   if (initslope > 0.0) initslope = -initslope;
318   if (initslope == 0.0) initslope = -1.0;
319 
320   ierr = VecCopy(y,w); CHKERRQ(ierr);
321   ierr = VecAYPX(&mone,x,w); CHKERRQ(ierr);
322   ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
323   ierr = VecNorm(g,NORM_2,gnorm);
324   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */
325     ierr = VecCopy(w,y); CHKERRQ(ierr);
326     PLogInfo(snes,"SNESCubicLineSearch: Using full step\n");
327     goto theend1;
328   }
329 
330   /* Fit points with quadratic */
331   lambda = 1.0; count = 0;
332   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
333   lambdaprev = lambda;
334   gnormprev = *gnorm;
335   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
336   else lambda = lambdatemp;
337   ierr   = VecCopy(x,w); CHKERRQ(ierr);
338   lambdaneg = -lambda;
339 #if defined(USE_PETSC_COMPLEX)
340   clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr);
341 #else
342   ierr = VecAXPY(&lambdaneg,y,w); CHKERRQ(ierr);
343 #endif
344   ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
345   ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);
346   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */
347     ierr = VecCopy(w,y); CHKERRQ(ierr);
348     PLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%g\n",lambda);
349     goto theend1;
350   }
351 
352   /* Fit points with cubic */
353   count = 1;
354   while (1) {
355     if (lambda <= minlambda) { /* bad luck; use full step */
356       PLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count);
357       PLogInfo(snes,"SNESCubicLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",
358                fnorm,*gnorm,*ynorm,lambda,initslope);
359       ierr = VecCopy(w,y); CHKERRQ(ierr);
360       *flag = -1; break;
361     }
362     t1 = ((*gnorm)*(*gnorm) - fnorm*fnorm)*0.5 - lambda*initslope;
363     t2 = (gnormprev*gnormprev  - fnorm*fnorm)*0.5 - lambdaprev*initslope;
364     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
365     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
366     d  = b*b - 3*a*initslope;
367     if (d < 0.0) d = 0.0;
368     if (a == 0.0) {
369       lambdatemp = -initslope/(2.0*b);
370     } else {
371       lambdatemp = (-b + sqrt(d))/(3.0*a);
372     }
373     if (lambdatemp > .5*lambda) {
374       lambdatemp = .5*lambda;
375     }
376     lambdaprev = lambda;
377     gnormprev = *gnorm;
378     if (lambdatemp <= .1*lambda) {
379       lambda = .1*lambda;
380     }
381     else lambda = lambdatemp;
382     ierr = VecCopy( x, w ); CHKERRQ(ierr);
383     lambdaneg = -lambda;
384 #if defined(USE_PETSC_COMPLEX)
385     clambda = lambdaneg;
386     ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr);
387 #else
388     ierr = VecAXPY(&lambdaneg,y,w); CHKERRQ(ierr);
389 #endif
390     ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
391     ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);
392     if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* is reduction enough? */
393       ierr = VecCopy(w,y); CHKERRQ(ierr);
394       PLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%g\n",lambda);
395       goto theend1;
396     } else {
397       PLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda,  lambda=%g\n",lambda);
398     }
399     count++;
400   }
401   theend1:
402   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
403   PetscFunctionReturn(0);
404 }
405 /* ------------------------------------------------------------------ */
406 #undef __FUNC__
407 #define __FUNC__ "SNESQuadraticLineSearch"
408 /*@C
409    SNESQuadraticLineSearch - Performs a quadratic line search.
410 
411    Input Parameters:
412 .  snes - the SNES context
413 .  x - current iterate
414 .  f - residual evaluated at x
415 .  y - search direction (contains new iterate on output)
416 .  w - work vector
417 .  fnorm - 2-norm of f
418 
419    Output Parameters:
420 .  g - residual evaluated at new iterate y
421 .  y - new iterate (contains search direction on input)
422 .  gnorm - 2-norm of g
423 .  ynorm - 2-norm of search length
424 .  flag - 0 if line search succeeds; -1 on failure.
425 
426    Options Database Key:
427 $  -snes_eq_ls quadratic
428 
429    Notes:
430    Use SNESSetLineSearch()
431    to set this routine within the SNES_EQ_LS method.
432 
433 .keywords: SNES, nonlinear, quadratic, line search
434 
435 .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch()
436 @*/
437 int SNESQuadraticLineSearch(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w,
438                            double fnorm, double *ynorm, double *gnorm,int *flag)
439 {
440   /*
441      Note that for line search purposes we work with with the related
442      minimization problem:
443         min  z(x):  R^n -> R,
444      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
445    */
446   double  steptol,initslope,lambdaprev,gnormprev,maxstep,minlambda,alpha,lambda,lambdatemp;
447 #if defined(USE_PETSC_COMPLEX)
448   Scalar  cinitslope,clambda;
449 #endif
450   int     ierr,count;
451   SNES_LS *neP = (SNES_LS *) snes->data;
452   Scalar  mone = -1.0,scale;
453 
454   PetscFunctionBegin;
455   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
456   *flag   = 0;
457   alpha   = neP->alpha;
458   maxstep = neP->maxstep;
459   steptol = neP->steptol;
460 
461   VecNorm(y, NORM_2,ynorm );
462   if (*ynorm < snes->atol) {
463     PLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n");
464     *gnorm = fnorm;
465     ierr = VecCopy(x,y); CHKERRQ(ierr);
466     ierr = VecCopy(f,g); CHKERRQ(ierr);
467     goto theend2;
468   }
469   if (*ynorm > maxstep) {	/* Step too big, so scale back */
470     scale = maxstep/(*ynorm);
471     ierr = VecScale(&scale,y); CHKERRQ(ierr);
472     *ynorm = maxstep;
473   }
474   minlambda = steptol/(*ynorm);
475   ierr = MatMult(snes->jacobian,y,w); CHKERRQ(ierr);
476 #if defined(USE_PETSC_COMPLEX)
477   ierr = VecDot(f,w,&cinitslope); CHKERRQ(ierr);
478   initslope = real(cinitslope);
479 #else
480   ierr = VecDot(f,w,&initslope); CHKERRQ(ierr);
481 #endif
482   if (initslope > 0.0) initslope = -initslope;
483   if (initslope == 0.0) initslope = -1.0;
484 
485   ierr = VecCopy(y,w); CHKERRQ(ierr);
486   ierr = VecAYPX(&mone,x,w); CHKERRQ(ierr);
487   ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
488   ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);
489   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */
490     ierr = VecCopy(w,y); CHKERRQ(ierr);
491     PLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n");
492     goto theend2;
493   }
494 
495   /* Fit points with quadratic */
496   lambda = 1.0; count = 0;
497   count = 1;
498   while (1) {
499     if (lambda <= minlambda) { /* bad luck; use full step */
500       PLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count);
501       PLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",
502                fnorm,*gnorm,*ynorm,lambda,initslope);
503       ierr = VecCopy(w,y); CHKERRQ(ierr);
504       *flag = -1; break;
505     }
506     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
507     lambdaprev = lambda;
508     gnormprev = *gnorm;
509     if (lambdatemp <= .1*lambda) {
510       lambda = .1*lambda;
511     } else lambda = lambdatemp;
512     ierr = VecCopy(x,w); CHKERRQ(ierr);
513     lambda = -lambda;
514 #if defined(USE_PETSC_COMPLEX)
515     clambda = lambda; ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr);
516 #else
517     ierr = VecAXPY(&lambda,y,w); CHKERRQ(ierr);
518 #endif
519     ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
520     ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);
521     if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */
522       ierr = VecCopy(w,y); CHKERRQ(ierr);
523       PLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda);
524       break;
525     }
526     count++;
527   }
528   theend2:
529   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
530   PetscFunctionReturn(0);
531 }
532 /* ------------------------------------------------------------ */
533 #undef __FUNC__
534 #define __FUNC__ "SNESSetLineSearch"
535 /*@C
536    SNESSetLineSearch - Sets the line search routine to be used
537    by the method SNES_EQ_LS.
538 
539    Input Parameters:
540 .  snes - nonlinear context obtained from SNESCreate()
541 .  func - pointer to int function
542 
543    Available Routines:
544 .  SNESCubicLineSearch() - default line search
545 .  SNESQuadraticLineSearch() - quadratic line search
546 .  SNESNoLineSearch() - the full Newton step (actually not a line search)
547 
548     Options Database Keys:
549 $   -snes_eq_ls [basic,quadratic,cubic]
550 $   -snes_eq_ls_alpha <alpha>
551 $   -snes_eq_ls_maxstep <max>
552 $   -snes_eq_ls_steptol <steptol>
553 
554    Calling sequence of func:
555    func (SNES snes, Vec x, Vec f, Vec g, Vec y,
556          Vec w, double fnorm, double *ynorm,
557          double *gnorm, *flag)
558 
559     Input parameters for func:
560 .   snes - nonlinear context
561 .   x - current iterate
562 .   f - residual evaluated at x
563 .   y - search direction (contains new iterate on output)
564 .   w - work vector
565 .   fnorm - 2-norm of f
566 
567     Output parameters for func:
568 .   g - residual evaluated at new iterate y
569 .   y - new iterate (contains search direction on input)
570 .   gnorm - 2-norm of g
571 .   ynorm - 2-norm of search length
572 .   flag - set to 0 if the line search succeeds; a nonzero integer
573            on failure.
574 
575 .keywords: SNES, nonlinear, set, line search, routine
576 
577 .seealso: SNESNoLineSearch(), SNESQuadraticLineSearch(), SNESCubicLineSearch()
578 @*/
579 int SNESSetLineSearch(SNES snes,int (*func)(SNES,Vec,Vec,Vec,Vec,Vec,
580                              double,double*,double*,int*))
581 {
582   PetscFunctionBegin;
583   if ((snes)->type == SNES_EQ_LS) ((SNES_LS *)(snes->data))->LineSearch = func;
584   PetscFunctionReturn(0);
585 }
586 /* ------------------------------------------------------------------ */
587 #undef __FUNC__
588 #define __FUNC__ "SNESPrintHelp_EQ_LS"
589 static int SNESPrintHelp_EQ_LS(SNES snes,char *p)
590 {
591   SNES_LS *ls = (SNES_LS *)snes->data;
592 
593   PetscFunctionBegin;
594   (*PetscHelpPrintf)(snes->comm," method SNES_EQ_LS (ls) for systems of nonlinear equations:\n");
595   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls [basic,quadratic,cubic]\n",p);
596   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_alpha <alpha> (default %g)\n",p,ls->alpha);
597   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_maxstep <max> (default %g)\n",p,ls->maxstep);
598   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_steptol <tol> (default %g)\n",p,ls->steptol);
599   PetscFunctionReturn(0);
600 }
601 /* ------------------------------------------------------------------ */
602 #undef __FUNC__
603 #define __FUNC__ "SNESView_EQ_LS"
604 static int SNESView_EQ_LS(PetscObject obj,Viewer viewer)
605 {
606   SNES       snes = (SNES)obj;
607   SNES_LS    *ls = (SNES_LS *)snes->data;
608   FILE       *fd;
609   char       *cstr;
610   int        ierr;
611   ViewerType vtype;
612 
613   PetscFunctionBegin;
614   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
615   if (vtype  == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
616     ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
617     if (ls->LineSearch == SNESNoLineSearch) cstr = "SNESNoLineSearch";
618     else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch";
619     else if (ls->LineSearch == SNESCubicLineSearch) cstr = "SNESCubicLineSearch";
620     else cstr = "unknown";
621     PetscFPrintf(snes->comm,fd,"    line search variant: %s\n",cstr);
622     PetscFPrintf(snes->comm,fd,"    alpha=%g, maxstep=%g, steptol=%g\n",
623                  ls->alpha,ls->maxstep,ls->steptol);
624   }
625   PetscFunctionReturn(0);
626 }
627 /* ------------------------------------------------------------------ */
628 #undef __FUNC__
629 #define __FUNC__ "SNESSetFromOptions_EQ_LS"
630 static int SNESSetFromOptions_EQ_LS(SNES snes)
631 {
632   SNES_LS *ls = (SNES_LS *)snes->data;
633   char    ver[16];
634   double  tmp;
635   int     ierr,flg;
636 
637   PetscFunctionBegin;
638   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_alpha",&tmp, &flg);CHKERRQ(ierr);
639   if (flg) {
640     ls->alpha = tmp;
641   }
642   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_maxstep",&tmp, &flg);CHKERRQ(ierr);
643   if (flg) {
644     ls->maxstep = tmp;
645   }
646   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_steptol",&tmp, &flg);CHKERRQ(ierr);
647   if (flg) {
648     ls->steptol = tmp;
649   }
650   ierr = OptionsGetString(snes->prefix,"-snes_eq_ls",ver,16, &flg); CHKERRQ(ierr);
651   if (flg) {
652     if (!PetscStrcmp(ver,"basic")) {
653       SNESSetLineSearch(snes,SNESNoLineSearch);
654     }
655     else if (!PetscStrcmp(ver,"basicnonorms")) {
656       SNESSetLineSearch(snes,SNESNoLineSearchNoNorms);
657     }
658     else if (!PetscStrcmp(ver,"quadratic")) {
659       SNESSetLineSearch(snes,SNESQuadraticLineSearch);
660     }
661     else if (!PetscStrcmp(ver,"cubic")) {
662       SNESSetLineSearch(snes,SNESCubicLineSearch);
663     }
664     else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown line search");}
665   }
666   PetscFunctionReturn(0);
667 }
668 /* ------------------------------------------------------------ */
669 #undef __FUNC__
670 #define __FUNC__ "SNESCreate_EQ_LS"
671 int SNESCreate_EQ_LS(SNES  snes )
672 {
673   SNES_LS *neP;
674 
675   PetscFunctionBegin;
676   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
677     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
678   }
679   snes->type		= SNES_EQ_LS;
680   snes->setup		= SNESSetUp_EQ_LS;
681   snes->solve		= SNESSolve_EQ_LS;
682   snes->destroy		= SNESDestroy_EQ_LS;
683   snes->converged	= SNESConverged_EQ_LS;
684   snes->printhelp       = SNESPrintHelp_EQ_LS;
685   snes->setfromoptions  = SNESSetFromOptions_EQ_LS;
686   snes->view            = SNESView_EQ_LS;
687   snes->nwork           = 0;
688 
689   neP			= PetscNew(SNES_LS);   CHKPTRQ(neP);
690   PLogObjectMemory(snes,sizeof(SNES_LS));
691   snes->data    	= (void *) neP;
692   neP->alpha		= 1.e-4;
693   neP->maxstep		= 1.e8;
694   neP->steptol		= 1.e-12;
695   neP->LineSearch       = SNESCubicLineSearch;
696   PetscFunctionReturn(0);
697 }
698 
699