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