xref: /petsc/src/snes/impls/ls/ls.c (revision 6849ba73f22fecb8f92ef896a42e4e8bd4cd6965)
1 
2 #include "src/snes/impls/ls/ls.h"
3 
4 /*
5      Checks if J^T F = 0 which implies we've found a local minimum of the function,
6     but not a zero. In the case when one cannot compute J^T F we use the fact that
7     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
8     for this trick.
9 */
10 #undef __FUNCT__
11 #define __FUNCT__ "SNESLSCheckLocalMin_Private"
12 PetscErrorCode SNESLSCheckLocalMin_Private(Mat A,Vec F,Vec W,PetscReal fnorm,PetscTruth *ismin)
13 {
14   PetscReal  a1;
15   PetscErrorCode ierr;
16   PetscTruth hastranspose;
17 
18   PetscFunctionBegin;
19   *ismin = PETSC_FALSE;
20   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
21   if (hastranspose) {
22     /* Compute || J^T F|| */
23     ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr);
24     ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr);
25     PetscLogInfo(0,"SNESSolve_LS: || J^T F|| %g near zero implies found a local minimum\n",a1/fnorm);
26     if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
27   } else {
28     Vec       work;
29     PetscScalar    result;
30     PetscReal wnorm;
31 
32     ierr = VecSetRandom(PETSC_NULL,W);CHKERRQ(ierr);
33     ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr);
34     ierr = VecDuplicate(W,&work);CHKERRQ(ierr);
35     ierr = MatMult(A,W,work);CHKERRQ(ierr);
36     ierr = VecDot(F,work,&result);CHKERRQ(ierr);
37     ierr = VecDestroy(work);CHKERRQ(ierr);
38     a1   = PetscAbsScalar(result)/(fnorm*wnorm);
39     PetscLogInfo(0,"SNESSolve_LS: (F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n",a1);
40     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
41   }
42   PetscFunctionReturn(0);
43 }
44 
45 /*
46      Checks if J^T(F - J*X) = 0
47 */
48 #undef __FUNCT__
49 #define __FUNCT__ "SNESLSCheckResidual_Private"
50 PetscErrorCode SNESLSCheckResidual_Private(Mat A,Vec F,Vec X,Vec W1,Vec W2)
51 {
52   PetscReal     a1,a2;
53   PetscErrorCode ierr;
54   PetscTruth    hastranspose;
55   PetscScalar   mone = -1.0;
56 
57   PetscFunctionBegin;
58   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
59   if (hastranspose) {
60     ierr = MatMult(A,X,W1);CHKERRQ(ierr);
61     ierr = VecAXPY(&mone,F,W1);CHKERRQ(ierr);
62 
63     /* Compute || J^T W|| */
64     ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr);
65     ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr);
66     ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr);
67     if (a1 != 0) {
68       PetscLogInfo(0,"SNESSolve_LS: ||J^T(F-Ax)||/||F-AX|| %g near zero implies inconsistent rhs\n",a2/a1);
69     }
70   }
71   PetscFunctionReturn(0);
72 }
73 
74 /*  --------------------------------------------------------------------
75 
76      This file implements a truncated Newton method with a line search,
77      for solving a system of nonlinear equations, using the KSP, Vec,
78      and Mat interfaces for linear solvers, vectors, and matrices,
79      respectively.
80 
81      The following basic routines are required for each nonlinear solver:
82           SNESCreate_XXX()          - Creates a nonlinear solver context
83           SNESSetFromOptions_XXX()  - Sets runtime options
84           SNESSolve_XXX()           - Solves the nonlinear system
85           SNESDestroy_XXX()         - Destroys the nonlinear solver context
86      The suffix "_XXX" denotes a particular implementation, in this case
87      we use _LS (e.g., SNESCreate_LS, SNESSolve_LS) for solving
88      systems of nonlinear equations with a line search (LS) method.
89      These routines are actually called via the common user interface
90      routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and
91      SNESDestroy(), so the application code interface remains identical
92      for all nonlinear solvers.
93 
94      Another key routine is:
95           SNESSetUp_XXX()           - Prepares for the use of a nonlinear solver
96      by setting data structures and options.   The interface routine SNESSetUp()
97      is not usually called directly by the user, but instead is called by
98      SNESSolve() if necessary.
99 
100      Additional basic routines are:
101           SNESView_XXX()            - Prints details of runtime options that
102                                       have actually been used.
103      These are called by application codes via the interface routines
104      SNESView().
105 
106      The various types of solvers (preconditioners, Krylov subspace methods,
107      nonlinear solvers, timesteppers) are all organized similarly, so the
108      above description applies to these categories also.
109 
110     -------------------------------------------------------------------- */
111 /*
112    SNESSolve_LS - Solves a nonlinear system with a truncated Newton
113    method with a line search.
114 
115    Input Parameters:
116 .  snes - the SNES context
117 
118    Output Parameter:
119 .  outits - number of iterations until termination
120 
121    Application Interface Routine: SNESSolve()
122 
123    Notes:
124    This implements essentially a truncated Newton method with a
125    line search.  By default a cubic backtracking line search
126    is employed, as described in the text "Numerical Methods for
127    Unconstrained Optimization and Nonlinear Equations" by Dennis
128    and Schnabel.
129 */
130 #undef __FUNCT__
131 #define __FUNCT__ "SNESSolve_LS"
132 PetscErrorCode SNESSolve_LS(SNES snes)
133 {
134   SNES_LS      *neP = (SNES_LS*)snes->data;
135   PetscErrorCode ierr;
136   int          maxits,i,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 PetscErrorCode SNESSetUp_LS(SNES snes)
262 {
263   PetscErrorCode 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 PetscErrorCode SNESDestroy_LS(SNES snes)
285 {
286   PetscErrorCode 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 PetscErrorCode 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   PetscErrorCode 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 PetscErrorCode 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   PetscErrorCode 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 PetscErrorCode 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   PetscErrorCode ierr;
481   int count;
482   SNES_LS     *neP = (SNES_LS*)snes->data;
483   PetscScalar mone = -1.0,scale;
484   PetscTruth  change_y = PETSC_FALSE;
485 
486   PetscFunctionBegin;
487   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
488   *flag   = 0;
489   alpha   = neP->alpha;
490   maxstep = neP->maxstep;
491   steptol = neP->steptol;
492 
493   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
494   if (*ynorm == 0.0) {
495     PetscLogInfo(snes,"SNESCubicLineSearch: Search direction and size is 0\n");
496     *gnorm = fnorm;
497     ierr   = VecCopy(x,y);CHKERRQ(ierr);
498     ierr   = VecCopy(f,g);CHKERRQ(ierr);
499     *flag  = -1;
500     goto theend1;
501   }
502   if (*ynorm > maxstep) {	/* Step too big, so scale back */
503     scale = maxstep/(*ynorm);
504 #if defined(PETSC_USE_COMPLEX)
505     PetscLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",PetscRealPart(scale),*ynorm);
506 #else
507     PetscLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",scale,*ynorm);
508 #endif
509     ierr = VecScale(&scale,y);CHKERRQ(ierr);
510     *ynorm = maxstep;
511   }
512   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
513   minlambda = steptol/rellength;
514   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
515 #if defined(PETSC_USE_COMPLEX)
516   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
517   initslope = PetscRealPart(cinitslope);
518 #else
519   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
520 #endif
521   if (initslope > 0.0) initslope = -initslope;
522   if (initslope == 0.0) initslope = -1.0;
523 
524   ierr = VecCopy(y,w);CHKERRQ(ierr);
525   ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr);
526   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
527   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
528   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
529     ierr = VecCopy(w,y);CHKERRQ(ierr);
530     PetscLogInfo(snes,"SNESCubicLineSearch: Using full step\n");
531     goto theend1;
532   }
533 
534   /* Fit points with quadratic */
535   lambda = 1.0;
536   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
537   lambdaprev = lambda;
538   gnormprev = *gnorm;
539   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
540   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
541   else                         lambda = lambdatemp;
542   ierr   = VecCopy(x,w);CHKERRQ(ierr);
543   lambdaneg = -lambda;
544 #if defined(PETSC_USE_COMPLEX)
545   clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
546 #else
547   ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
548 #endif
549   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
550   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
551   if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
552     ierr = VecCopy(w,y);CHKERRQ(ierr);
553     PetscLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%18.16e\n",lambda);
554     goto theend1;
555   }
556 
557   /* Fit points with cubic */
558   count = 1;
559   while (count < 10000) {
560     if (lambda <= minlambda) { /* bad luck; use full step */
561       PetscLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count);
562       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);
563       ierr = VecCopy(x,y);CHKERRQ(ierr);
564       *flag = -1; break;
565     }
566     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
567     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
568     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
569     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
570     d  = b*b - 3*a*initslope;
571     if (d < 0.0) d = 0.0;
572     if (a == 0.0) {
573       lambdatemp = -initslope/(2.0*b);
574     } else {
575       lambdatemp = (-b + sqrt(d))/(3.0*a);
576     }
577     lambdaprev = lambda;
578     gnormprev  = *gnorm;
579     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
580     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
581     else                         lambda     = lambdatemp;
582     ierr = VecCopy(x,w);CHKERRQ(ierr);
583     lambdaneg = -lambda;
584 #if defined(PETSC_USE_COMPLEX)
585     clambda = lambdaneg;
586     ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
587 #else
588     ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
589 #endif
590     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
591     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
592     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */
593       ierr = VecCopy(w,y);CHKERRQ(ierr);
594       PetscLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%18.16e\n",lambda);
595       goto theend1;
596     } else {
597       PetscLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda,  lambda=%18.16e\n",lambda);
598     }
599     count++;
600   }
601   if (count >= 10000) {
602     SETERRQ(PETSC_ERR_LIB, "Lambda was decreased more than 10,000 times, so something is probably wrong with the function evaluation");
603   }
604   theend1:
605   /* Optional user-defined check for line search step validity */
606   if (neP->CheckStep) {
607     ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
608     if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */
609       ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr);
610       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
611       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
612       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
613       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
614     }
615   }
616   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
617   PetscFunctionReturn(0);
618 }
619 /* -------------------------------------------------------------------------- */
620 #undef __FUNCT__
621 #define __FUNCT__ "SNESQuadraticLineSearch"
622 /*@C
623    SNESQuadraticLineSearch - Performs a quadratic line search.
624 
625    Collective on SNES and Vec
626 
627    Input Parameters:
628 +  snes - the SNES context
629 .  lsctx - optional context for line search (not used here)
630 .  x - current iterate
631 .  f - residual evaluated at x
632 .  y - search direction (contains new iterate on output)
633 .  w - work vector
634 -  fnorm - 2-norm of f
635 
636    Output Parameters:
637 +  g - residual evaluated at new iterate y
638 .  y - new iterate (contains search direction on input)
639 .  gnorm - 2-norm of g
640 .  ynorm - 2-norm of search length
641 -  flag - 0 if line search succeeds; -1 on failure.
642 
643    Options Database Key:
644 .  -snes_ls quadratic - Activates SNESQuadraticLineSearch()
645 
646    Notes:
647    Use SNESSetLineSearch() to set this routine within the SNESLS method.
648 
649    Level: advanced
650 
651 .keywords: SNES, nonlinear, quadratic, line search
652 
653 .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
654 @*/
655 PetscErrorCode SNESQuadraticLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag)
656 {
657   /*
658      Note that for line search purposes we work with with the related
659      minimization problem:
660         min  z(x):  R^n -> R,
661      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
662    */
663   PetscReal   steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg,rellength;
664 #if defined(PETSC_USE_COMPLEX)
665   PetscScalar cinitslope,clambda;
666 #endif
667   PetscErrorCode ierr;
668   int count;
669   SNES_LS     *neP = (SNES_LS*)snes->data;
670   PetscScalar mone = -1.0,scale;
671   PetscTruth  change_y = PETSC_FALSE;
672 
673   PetscFunctionBegin;
674   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
675   *flag   = 0;
676   alpha   = neP->alpha;
677   maxstep = neP->maxstep;
678   steptol = neP->steptol;
679 
680   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
681   if (*ynorm == 0.0) {
682     PetscLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n");
683     *gnorm = fnorm;
684     ierr   = VecCopy(x,y);CHKERRQ(ierr);
685     ierr   = VecCopy(f,g);CHKERRQ(ierr);
686     *flag  = -1;
687     goto theend2;
688   }
689   if (*ynorm > maxstep) {	/* Step too big, so scale back */
690     scale = maxstep/(*ynorm);
691     ierr = VecScale(&scale,y);CHKERRQ(ierr);
692     *ynorm = maxstep;
693   }
694   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
695   minlambda = steptol/rellength;
696   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
697 #if defined(PETSC_USE_COMPLEX)
698   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
699   initslope = PetscRealPart(cinitslope);
700 #else
701   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
702 #endif
703   if (initslope > 0.0) initslope = -initslope;
704   if (initslope == 0.0) initslope = -1.0;
705 
706   ierr = VecCopy(y,w);CHKERRQ(ierr);
707   ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr);
708   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
709   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
710   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
711     ierr = VecCopy(w,y);CHKERRQ(ierr);
712     PetscLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n");
713     goto theend2;
714   }
715 
716   /* Fit points with quadratic */
717   lambda = 1.0;
718   count = 1;
719   while (PETSC_TRUE) {
720     if (lambda <= minlambda) { /* bad luck; use full step */
721       PetscLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count);
722       PetscLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",fnorm,*gnorm,*ynorm,lambda,initslope);
723       ierr = VecCopy(x,y);CHKERRQ(ierr);
724       *flag = -1; break;
725     }
726     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
727     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
728     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
729     else                         lambda     = lambdatemp;
730     ierr = VecCopy(x,w);CHKERRQ(ierr);
731     lambdaneg = -lambda;
732 #if defined(PETSC_USE_COMPLEX)
733     clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
734 #else
735     ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
736 #endif
737     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
738     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
739     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
740       ierr = VecCopy(w,y);CHKERRQ(ierr);
741       PetscLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda);
742       break;
743     }
744     count++;
745   }
746   theend2:
747   /* Optional user-defined check for line search step validity */
748   if (neP->CheckStep) {
749     ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
750     if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */
751       ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr);
752       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
753       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
754       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
755       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
756     }
757   }
758   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
759   PetscFunctionReturn(0);
760 }
761 
762 /* -------------------------------------------------------------------------- */
763 #undef __FUNCT__
764 #define __FUNCT__ "SNESSetLineSearch"
765 /*@C
766    SNESSetLineSearch - Sets the line search routine to be used
767    by the method SNESLS.
768 
769    Input Parameters:
770 +  snes - nonlinear context obtained from SNESCreate()
771 .  lsctx - optional user-defined context for use by line search
772 -  func - pointer to int function
773 
774    Collective on SNES
775 
776    Available Routines:
777 +  SNESCubicLineSearch() - default line search
778 .  SNESQuadraticLineSearch() - quadratic line search
779 .  SNESNoLineSearch() - the full Newton step (actually not a line search)
780 -  SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel)
781 
782     Options Database Keys:
783 +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
784 .   -snes_ls_alpha <alpha> - Sets alpha
785 .   -snes_ls_maxstep <max> - Sets maxstep
786 -   -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code
787                    will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length
788                    the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x)
789 
790    Calling sequence of func:
791 .vb
792    func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
793          PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,*flag)
794 .ve
795 
796     Input parameters for func:
797 +   snes - nonlinear context
798 .   lsctx - optional user-defined context for line search
799 .   x - current iterate
800 .   f - residual evaluated at x
801 .   y - search direction (contains new iterate on output)
802 .   w - work vector
803 -   fnorm - 2-norm of f
804 
805     Output parameters for func:
806 +   g - residual evaluated at new iterate y
807 .   y - new iterate (contains search direction on input)
808 .   gnorm - 2-norm of g
809 .   ynorm - 2-norm of search length
810 -   flag - set to 0 if the line search succeeds; a nonzero integer
811            on failure.
812 
813     Level: advanced
814 
815 .keywords: SNES, nonlinear, set, line search, routine
816 
817 .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(),
818           SNESSetLineSearchCheck(), SNESSetLineSearchParams(), SNESGetLineSearchParams()
819 @*/
820 PetscErrorCode SNESSetLineSearch(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void *lsctx)
821 {
822   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void*);
823 
824   PetscFunctionBegin;
825   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void (**)(void))&f);CHKERRQ(ierr);
826   if (f) {
827     ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr);
828   }
829   PetscFunctionReturn(0);
830 }
831 
832 typedef PetscErrorCode (*FCN2)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*); /* force argument to next function to not be extern C*/
833 /* -------------------------------------------------------------------------- */
834 EXTERN_C_BEGIN
835 #undef __FUNCT__
836 #define __FUNCT__ "SNESSetLineSearch_LS"
837 PetscErrorCode SNESSetLineSearch_LS(SNES snes,FCN2 func,void *lsctx)
838 {
839   PetscFunctionBegin;
840   ((SNES_LS *)(snes->data))->LineSearch = func;
841   ((SNES_LS *)(snes->data))->lsP        = lsctx;
842   PetscFunctionReturn(0);
843 }
844 EXTERN_C_END
845 /* -------------------------------------------------------------------------- */
846 #undef __FUNCT__
847 #define __FUNCT__ "SNESSetLineSearchCheck"
848 /*@C
849    SNESSetLineSearchCheck - Sets a routine to check the validity of new iterate computed
850    by the line search routine in the Newton-based method SNESLS.
851 
852    Input Parameters:
853 +  snes - nonlinear context obtained from SNESCreate()
854 .  func - pointer to int function
855 -  checkctx - optional user-defined context for use by step checking routine
856 
857    Collective on SNES
858 
859    Calling sequence of func:
860 .vb
861    int func (SNES snes, void *checkctx, Vec x, PetscTruth *flag)
862 .ve
863    where func returns an error code of 0 on success and a nonzero
864    on failure.
865 
866    Input parameters for func:
867 +  snes - nonlinear context
868 .  checkctx - optional user-defined context for use by step checking routine
869 -  x - current candidate iterate
870 
871    Output parameters for func:
872 +  x - current iterate (possibly modified)
873 -  flag - flag indicating whether x has been modified (either
874            PETSC_TRUE of PETSC_FALSE)
875 
876    Level: advanced
877 
878    Notes:
879    SNESNoLineSearch() and SNESNoLineSearchNoNorms() accept the new
880    iterate computed by the line search checking routine.  In particular,
881    these routines (1) compute a candidate iterate u_{i+1}, (2) pass control
882    to the checking routine, and then (3) compute the corresponding nonlinear
883    function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}.
884 
885    SNESQuadraticLineSearch() and SNESCubicLineSearch() also accept the
886    new iterate computed by the line search checking routine.  In particular,
887    these routines (1) compute a candidate iterate u_{i+1} as well as a
888    candidate nonlinear function f(u_{i+1}), (2) pass control to the checking
889    routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes
890    were made to the candidate iterate in the checking routine (as indicated
891    by flag=PETSC_TRUE).  The overhead of this function re-evaluation can be
892    very costly, so use this feature with caution!
893 
894 .keywords: SNES, nonlinear, set, line search check, step check, routine
895 
896 .seealso: SNESSetLineSearch()
897 @*/
898 PetscErrorCode SNESSetLineSearchCheck(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx)
899 {
900   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,void*,Vec,PetscTruth*),void*);
901 
902   PetscFunctionBegin;
903   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearchCheck_C",(void (**)(void))&f);CHKERRQ(ierr);
904   if (f) {
905     ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr);
906   }
907   PetscFunctionReturn(0);
908 }
909 /* -------------------------------------------------------------------------- */
910 typedef PetscErrorCode (*FCN)(SNES,void*,Vec,PetscTruth*); /* force argument to next function to not be extern C*/
911 EXTERN_C_BEGIN
912 #undef __FUNCT__
913 #define __FUNCT__ "SNESSetLineSearchCheck_LS"
914 PetscErrorCode SNESSetLineSearchCheck_LS(SNES snes,FCN func,void *checkctx)
915 {
916   PetscFunctionBegin;
917   ((SNES_LS *)(snes->data))->CheckStep = func;
918   ((SNES_LS *)(snes->data))->checkP    = checkctx;
919   PetscFunctionReturn(0);
920 }
921 EXTERN_C_END
922 /* -------------------------------------------------------------------------- */
923 /*
924    SNESPrintHelp_LS - Prints all options for the SNES_LS method.
925 
926    Input Parameter:
927 .  snes - the SNES context
928 
929    Application Interface Routine: SNESPrintHelp()
930 */
931 #undef __FUNCT__
932 #define __FUNCT__ "SNESPrintHelp_LS"
933 static PetscErrorCode SNESPrintHelp_LS(SNES snes,char *p)
934 {
935   SNES_LS *ls = (SNES_LS *)snes->data;
936 
937   PetscFunctionBegin;
938   (*PetscHelpPrintf)(snes->comm," method SNES_LS (ls) for systems of nonlinear equations:\n");
939   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls [cubic,quadratic,basic,basicnonorms,...]\n",p);
940   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls_alpha <alpha> (default %g)\n",p,ls->alpha);
941   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls_maxstep <max> (default %g)\n",p,ls->maxstep);
942   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls_steptol <tol> (default %g)\n",p,ls->steptol);
943   PetscFunctionReturn(0);
944 }
945 
946 /*
947    SNESView_LS - Prints info from the SNESLS data structure.
948 
949    Input Parameters:
950 .  SNES - the SNES context
951 .  viewer - visualization context
952 
953    Application Interface Routine: SNESView()
954 */
955 #undef __FUNCT__
956 #define __FUNCT__ "SNESView_LS"
957 static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer)
958 {
959   SNES_LS    *ls = (SNES_LS *)snes->data;
960   const char *cstr;
961   PetscErrorCode ierr;
962   PetscTruth iascii;
963 
964   PetscFunctionBegin;
965   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
966   if (iascii) {
967     if (ls->LineSearch == SNESNoLineSearch)             cstr = "SNESNoLineSearch";
968     else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch";
969     else if (ls->LineSearch == SNESCubicLineSearch)     cstr = "SNESCubicLineSearch";
970     else                                                cstr = "unknown";
971     ierr = PetscViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
972     ierr = PetscViewerASCIIPrintf(viewer,"  alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr);
973   } else {
974     SETERRQ1(1,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name);
975   }
976   PetscFunctionReturn(0);
977 }
978 /* -------------------------------------------------------------------------- */
979 /*
980    SNESSetFromOptions_LS - Sets various parameters for the SNESLS method.
981 
982    Input Parameter:
983 .  snes - the SNES context
984 
985    Application Interface Routine: SNESSetFromOptions()
986 */
987 #undef __FUNCT__
988 #define __FUNCT__ "SNESSetFromOptions_LS"
989 static PetscErrorCode SNESSetFromOptions_LS(SNES snes)
990 {
991   SNES_LS    *ls = (SNES_LS *)snes->data;
992   const char *lses[] = {"basic","basicnonorms","quadratic","cubic"};
993   PetscErrorCode ierr;
994   int indx;
995   PetscTruth flg;
996 
997   PetscFunctionBegin;
998   ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr);
999     ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr);
1000     ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr);
1001     ierr = PetscOptionsReal("-snes_ls_steptol","Step must be greater than","None",ls->steptol,&ls->steptol,0);CHKERRQ(ierr);
1002 
1003     ierr = PetscOptionsEList("-snes_ls","Line search used","SNESSetLineSearch",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr);
1004     if (flg) {
1005       switch (indx) {
1006       case 0:
1007         ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr);
1008         break;
1009       case 1:
1010         ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
1011         break;
1012       case 2:
1013         ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr);
1014         break;
1015       case 3:
1016         ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr);
1017         break;
1018       }
1019     }
1020   ierr = PetscOptionsTail();CHKERRQ(ierr);
1021   PetscFunctionReturn(0);
1022 }
1023 /* -------------------------------------------------------------------------- */
1024 /*MC
1025       SNESLS - Newton based nonlinear solver that uses a line search
1026 
1027    Options Database:
1028 +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
1029 .   -snes_ls_alpha <alpha> - Sets alpha
1030 .   -snes_ls_maxstep <max> - Sets maxstep
1031 -   -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code
1032                    will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length
1033                    the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x)
1034 
1035     Notes: This is the default nonlinear solver in SNES
1036 
1037    Level: beginner
1038 
1039 .seealso:  SNESCreate(), SNES, SNESSetType(), SNESTR, SNESSetLineSearch(),
1040            SNESSetLineSearchCheck(), SNESNoLineSearch(), SNESCubicLineSearch(), SNESQuadraticLineSearch(),
1041           SNESSetLineSearch(), SNESNoLineSearchNoNorms()
1042 
1043 M*/
1044 EXTERN_C_BEGIN
1045 #undef __FUNCT__
1046 #define __FUNCT__ "SNESCreate_LS"
1047 PetscErrorCode SNESCreate_LS(SNES snes)
1048 {
1049   PetscErrorCode ierr;
1050   SNES_LS *neP;
1051 
1052   PetscFunctionBegin;
1053   snes->setup		= SNESSetUp_LS;
1054   snes->solve		= SNESSolve_LS;
1055   snes->destroy		= SNESDestroy_LS;
1056   snes->converged	= SNESConverged_LS;
1057   snes->printhelp       = SNESPrintHelp_LS;
1058   snes->setfromoptions  = SNESSetFromOptions_LS;
1059   snes->view            = SNESView_LS;
1060   snes->nwork           = 0;
1061 
1062   ierr                  = PetscNew(SNES_LS,&neP);CHKERRQ(ierr);
1063   PetscLogObjectMemory(snes,sizeof(SNES_LS));
1064   snes->data    	= (void*)neP;
1065   neP->alpha		= 1.e-4;
1066   neP->maxstep		= 1.e8;
1067   neP->steptol		= 1.e-12;
1068   neP->LineSearch       = SNESCubicLineSearch;
1069   neP->lsP              = PETSC_NULL;
1070   neP->CheckStep        = PETSC_NULL;
1071   neP->checkP           = PETSC_NULL;
1072 
1073   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS",SNESSetLineSearch_LS);CHKERRQ(ierr);
1074   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS",SNESSetLineSearchCheck_LS);CHKERRQ(ierr);
1075 
1076   PetscFunctionReturn(0);
1077 }
1078 EXTERN_C_END
1079 
1080 
1081 
1082