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