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