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