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