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