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