xref: /petsc/src/snes/impls/ls/ls.c (revision be1d678a52e6eff2808b2fa31ae986cdbf03c9fe)
1 
2 #include "src/snes/impls/ls/ls.h"
3 
4 /*
5      Checks if J^T F = 0 which implies we've found a local minimum of the function,
6     but not a zero. In the case when one cannot compute J^T F we use the fact that
7     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
8     for this trick.
9 */
10 #undef __FUNCT__
11 #define __FUNCT__ "SNESLSCheckLocalMin_Private"
12 PetscErrorCode SNESLSCheckLocalMin_Private(Mat A,Vec F,Vec W,PetscReal fnorm,PetscTruth *ismin)
13 {
14   PetscReal      a1;
15   PetscErrorCode ierr;
16   PetscTruth     hastranspose;
17 
18   PetscFunctionBegin;
19   *ismin = PETSC_FALSE;
20   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
21   if (hastranspose) {
22     /* Compute || J^T F|| */
23     ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr);
24     ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr);
25     ierr = PetscLogInfo((0,"SNESSolve_LS: || J^T F|| %g near zero implies found a local minimum\n",a1/fnorm));CHKERRQ(ierr);
26     if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
27   } else {
28     Vec         work;
29     PetscScalar result;
30     PetscReal   wnorm;
31 
32     ierr = VecSetRandom(PETSC_NULL,W);CHKERRQ(ierr);
33     ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr);
34     ierr = VecDuplicate(W,&work);CHKERRQ(ierr);
35     ierr = MatMult(A,W,work);CHKERRQ(ierr);
36     ierr = VecDot(F,work,&result);CHKERRQ(ierr);
37     ierr = VecDestroy(work);CHKERRQ(ierr);
38     a1   = PetscAbsScalar(result)/(fnorm*wnorm);
39     ierr = PetscLogInfo((0,"SNESSolve_LS: (F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n",a1));CHKERRQ(ierr);
40     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
41   }
42   PetscFunctionReturn(0);
43 }
44 
45 /*
46      Checks if J^T(F - J*X) = 0
47 */
48 #undef __FUNCT__
49 #define __FUNCT__ "SNESLSCheckResidual_Private"
50 PetscErrorCode SNESLSCheckResidual_Private(Mat A,Vec F,Vec X,Vec W1,Vec W2)
51 {
52   PetscReal     a1,a2;
53   PetscErrorCode ierr;
54   PetscTruth    hastranspose;
55   PetscScalar   mone = -1.0;
56 
57   PetscFunctionBegin;
58   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
59   if (hastranspose) {
60     ierr = MatMult(A,X,W1);CHKERRQ(ierr);
61     ierr = VecAXPY(&mone,F,W1);CHKERRQ(ierr);
62 
63     /* Compute || J^T W|| */
64     ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr);
65     ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr);
66     ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr);
67     if (a1 != 0) {
68       ierr = PetscLogInfo((0,"SNESSolve_LS: ||J^T(F-Ax)||/||F-AX|| %g near zero implies inconsistent rhs\n",a2/a1));CHKERRQ(ierr);
69     }
70   }
71   PetscFunctionReturn(0);
72 }
73 
74 /*  --------------------------------------------------------------------
75 
76      This file implements a truncated Newton method with a line search,
77      for solving a system of nonlinear equations, using the KSP, Vec,
78      and Mat interfaces for linear solvers, vectors, and matrices,
79      respectively.
80 
81      The following basic routines are required for each nonlinear solver:
82           SNESCreate_XXX()          - Creates a nonlinear solver context
83           SNESSetFromOptions_XXX()  - Sets runtime options
84           SNESSolve_XXX()           - Solves the nonlinear system
85           SNESDestroy_XXX()         - Destroys the nonlinear solver context
86      The suffix "_XXX" denotes a particular implementation, in this case
87      we use _LS (e.g., SNESCreate_LS, SNESSolve_LS) for solving
88      systems of nonlinear equations with a line search (LS) method.
89      These routines are actually called via the common user interface
90      routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and
91      SNESDestroy(), so the application code interface remains identical
92      for all nonlinear solvers.
93 
94      Another key routine is:
95           SNESSetUp_XXX()           - Prepares for the use of a nonlinear solver
96      by setting data structures and options.   The interface routine SNESSetUp()
97      is not usually called directly by the user, but instead is called by
98      SNESSolve() if necessary.
99 
100      Additional basic routines are:
101           SNESView_XXX()            - Prints details of runtime options that
102                                       have actually been used.
103      These are called by application codes via the interface routines
104      SNESView().
105 
106      The various types of solvers (preconditioners, Krylov subspace methods,
107      nonlinear solvers, timesteppers) are all organized similarly, so the
108      above description applies to these categories also.
109 
110     -------------------------------------------------------------------- */
111 /*
112    SNESSolve_LS - Solves a nonlinear system with a truncated Newton
113    method with a line search.
114 
115    Input Parameters:
116 .  snes - the SNES context
117 
118    Output Parameter:
119 .  outits - number of iterations until termination
120 
121    Application Interface Routine: SNESSolve()
122 
123    Notes:
124    This implements essentially a truncated Newton method with a
125    line search.  By default a cubic backtracking line search
126    is employed, as described in the text "Numerical Methods for
127    Unconstrained Optimization and Nonlinear Equations" by Dennis
128    and Schnabel.
129 */
130 #undef __FUNCT__
131 #define __FUNCT__ "SNESSolve_LS"
132 PetscErrorCode SNESSolve_LS(SNES snes)
133 {
134   SNES_LS        *neP = (SNES_LS*)snes->data;
135   PetscErrorCode ierr;
136   PetscInt       maxits,i,lits;
137   PetscTruth     lssucceed;
138   MatStructure   flg = DIFFERENT_NONZERO_PATTERN;
139   PetscReal      fnorm,gnorm,xnorm,ynorm;
140   Vec            Y,X,F,G,W,TMP;
141   KSP            ksp;
142 
143   PetscFunctionBegin;
144   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
145   snes->reason  = SNES_CONVERGED_ITERATING;
146 
147   maxits	= snes->max_its;	/* maximum number of iterations */
148   X		= snes->vec_sol;	/* solution vector */
149   F		= snes->vec_func;	/* residual vector */
150   Y		= snes->work[0];	/* work vectors */
151   G		= snes->work[1];
152   W		= snes->work[2];
153 
154   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
155   snes->iter = 0;
156   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
157   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);  /*  F(X)      */
158   ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);	/* fnorm <- ||F||  */
159   if (fnorm != fnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
160   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
161   snes->norm = fnorm;
162   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
163   SNESLogConvHistory(snes,fnorm,0);
164   SNESMonitor(snes,0,fnorm);
165 
166   if (fnorm < snes->abstol) {snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);}
167 
168   /* set parameter for default relative tolerance convergence test */
169   snes->ttol = fnorm*snes->rtol;
170 
171   for (i=0; i<maxits; i++) {
172 
173     /* Call general purpose update function */
174     if (snes->update) {
175       ierr = (*snes->update)(snes, snes->iter);CHKERRQ(ierr);
176     }
177 
178     /* Solve J Y = F, where J is Jacobian matrix */
179     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
180     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
181     ierr = KSPSolve(snes->ksp,F,Y);CHKERRQ(ierr);
182     ierr = KSPGetIterationNumber(ksp,&lits);CHKERRQ(ierr);
183 
184     if (PetscLogPrintInfo){
185       ierr = SNESLSCheckResidual_Private(snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
186     }
187 
188     /* should check what happened to the linear solve? */
189     snes->linear_its += lits;
190     ierr = PetscLogInfo((snes,"SNESSolve_LS: iter=%D, linear solve iterations=%D\n",snes->iter,lits));CHKERRQ(ierr);
191 
192     /* Compute a (scaled) negative update in the line search routine:
193          Y <- X - lambda*Y
194        and evaluate G(Y) = function(Y))
195     */
196     ierr = VecCopy(Y,snes->vec_sol_update_always);CHKERRQ(ierr);
197     ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
198     ierr = PetscLogInfo((snes,"SNESSolve_LS: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",fnorm,gnorm,ynorm,(int)lssucceed));CHKERRQ(ierr);
199     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
200 
201     TMP = F; F = G; snes->vec_func_always = F; G = TMP;
202     TMP = X; X = Y; snes->vec_sol_always = X;  Y = TMP;
203     fnorm = gnorm;
204 
205     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
206     snes->iter = i+1;
207     snes->norm = fnorm;
208     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
209     SNESLogConvHistory(snes,fnorm,lits);
210     SNESMonitor(snes,i+1,fnorm);
211 
212     if (!lssucceed) {
213       PetscTruth ismin;
214 
215       if (++snes->numFailures >= snes->maxFailures) {
216         snes->reason = SNES_DIVERGED_LS_FAILURE;
217         ierr = SNESLSCheckLocalMin_Private(snes->jacobian,F,W,fnorm,&ismin);CHKERRQ(ierr);
218         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
219         break;
220       }
221     }
222 
223     /* Test for convergence */
224     if (snes->converged) {
225       ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm = || X || */
226       ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
227       if (snes->reason) {
228         break;
229       }
230     }
231   }
232   if (X != snes->vec_sol) {
233     ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr);
234   }
235   if (F != snes->vec_func) {
236     ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr);
237   }
238   snes->vec_sol_always  = snes->vec_sol;
239   snes->vec_func_always = snes->vec_func;
240   if (i == maxits) {
241     ierr = PetscLogInfo((snes,"SNESSolve_LS: Maximum number of iterations has been reached: %D\n",maxits));CHKERRQ(ierr);
242     snes->reason = SNES_DIVERGED_MAX_IT;
243   }
244   PetscFunctionReturn(0);
245 }
246 /* -------------------------------------------------------------------------- */
247 /*
248    SNESSetUp_LS - Sets up the internal data structures for the later use
249    of the SNESLS nonlinear solver.
250 
251    Input Parameter:
252 .  snes - the SNES context
253 .  x - the solution vector
254 
255    Application Interface Routine: SNESSetUp()
256 
257    Notes:
258    For basic use of the SNES solvers, the user need not explicitly call
259    SNESSetUp(), since these actions will automatically occur during
260    the call to SNESSolve().
261  */
262 #undef __FUNCT__
263 #define __FUNCT__ "SNESSetUp_LS"
264 PetscErrorCode SNESSetUp_LS(SNES snes)
265 {
266   PetscErrorCode ierr;
267 
268   PetscFunctionBegin;
269   snes->nwork = 4;
270   ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
271   ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr);
272   snes->vec_sol_update_always = snes->work[3];
273   PetscFunctionReturn(0);
274 }
275 /* -------------------------------------------------------------------------- */
276 /*
277    SNESDestroy_LS - Destroys the private SNES_LS context that was created
278    with SNESCreate_LS().
279 
280    Input Parameter:
281 .  snes - the SNES context
282 
283    Application Interface Routine: SNESDestroy()
284  */
285 #undef __FUNCT__
286 #define __FUNCT__ "SNESDestroy_LS"
287 PetscErrorCode SNESDestroy_LS(SNES snes)
288 {
289   PetscErrorCode ierr;
290 
291   PetscFunctionBegin;
292   if (snes->nwork) {
293     ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr);
294   }
295   ierr = PetscFree(snes->data);CHKERRQ(ierr);
296   PetscFunctionReturn(0);
297 }
298 /* -------------------------------------------------------------------------- */
299 #undef __FUNCT__
300 #define __FUNCT__ "SNESNoLineSearch"
301 
302 /*@C
303    SNESNoLineSearch - This routine is not a line search at all;
304    it simply uses the full Newton step.  Thus, this routine is intended
305    to serve as a template and is not recommended for general use.
306 
307    Collective on SNES and Vec
308 
309    Input Parameters:
310 +  snes - nonlinear context
311 .  lsctx - optional context for line search (not used here)
312 .  x - current iterate
313 .  f - residual evaluated at x
314 .  y - search direction (contains new iterate on output)
315 .  w - work vector
316 -  fnorm - 2-norm of f
317 
318    Output Parameters:
319 +  g - residual evaluated at new iterate y
320 .  y - new iterate (contains search direction on input)
321 .  gnorm - 2-norm of g
322 .  ynorm - 2-norm of search length
323 -  flag - PETSC_TRUE on success, PETSC_FALSE on failure
324 
325    Options Database Key:
326 .  -snes_ls basic - Activates SNESNoLineSearch()
327 
328    Level: advanced
329 
330 .keywords: SNES, nonlinear, line search, cubic
331 
332 .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
333           SNESSetLineSearch(), SNESNoLineSearchNoNorms()
334 @*/
335 PetscErrorCode SNESNoLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
336 {
337   PetscErrorCode ierr;
338   PetscScalar    mone = -1.0;
339   SNES_LS        *neP = (SNES_LS*)snes->data;
340   PetscTruth     change_y = PETSC_FALSE;
341 
342   PetscFunctionBegin;
343   *flag = PETSC_TRUE;
344   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
345   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);  /* ynorm = || y || */
346   ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr);            /* y <- y - x      */
347   if (neP->CheckStep) {
348    ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
349   }
350   ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y)    */
351   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);  /* gnorm = || g || */
352   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
353   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
354   PetscFunctionReturn(0);
355 }
356 /* -------------------------------------------------------------------------- */
357 #undef __FUNCT__
358 #define __FUNCT__ "SNESNoLineSearchNoNorms"
359 
360 /*@C
361    SNESNoLineSearchNoNorms - This routine is not a line search at
362    all; it simply uses the full Newton step. This version does not
363    even compute the norm of the function or search direction; this
364    is intended only when you know the full step is fine and are
365    not checking for convergence of the nonlinear iteration (for
366    example, you are running always for a fixed number of Newton steps).
367 
368    Collective on SNES and Vec
369 
370    Input Parameters:
371 +  snes - nonlinear context
372 .  lsctx - optional context for line search (not used here)
373 .  x - current iterate
374 .  f - residual evaluated at x
375 .  y - search direction (contains new iterate on output)
376 .  w - work vector
377 -  fnorm - 2-norm of f
378 
379    Output Parameters:
380 +  g - residual evaluated at new iterate y
381 .  gnorm - not changed
382 .  ynorm - not changed
383 -  flag - set to PETSC_TRUE indicating a successful line search
384 
385    Options Database Key:
386 .  -snes_ls basicnonorms - Activates SNESNoLineSearchNoNorms()
387 
388    Notes:
389    SNESNoLineSearchNoNorms() must be used in conjunction with
390    either the options
391 $     -snes_no_convergence_test -snes_max_it <its>
392    or alternatively a user-defined custom test set via
393    SNESSetConvergenceTest(); or a -snes_max_it of 1,
394    otherwise, the SNES solver will generate an error.
395 
396    During the final iteration this will not evaluate the function at
397    the solution point. This is to save a function evaluation while
398    using pseudo-timestepping.
399 
400    The residual norms printed by monitoring routines such as
401    SNESDefaultMonitor() (as activated via -snes_monitor) will not be
402    correct, since they are not computed.
403 
404    Level: advanced
405 
406 .keywords: SNES, nonlinear, line search, cubic
407 
408 .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
409           SNESSetLineSearch(), SNESNoLineSearch()
410 @*/
411 PetscErrorCode SNESNoLineSearchNoNorms(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
412 {
413   PetscErrorCode ierr;
414   PetscScalar mone = -1.0;
415   SNES_LS     *neP = (SNES_LS*)snes->data;
416   PetscTruth  change_y = PETSC_FALSE;
417 
418   PetscFunctionBegin;
419   *flag = PETSC_TRUE;
420   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
421   ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr);            /* y <- y - x      */
422   if (neP->CheckStep) {
423    ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
424   }
425 
426   /* don't evaluate function the last time through */
427   if (snes->iter < snes->max_its-1) {
428     ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y)    */
429   }
430   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
431   PetscFunctionReturn(0);
432 }
433 /* -------------------------------------------------------------------------- */
434 #undef __FUNCT__
435 #define __FUNCT__ "SNESCubicLineSearch"
436 /*@C
437    SNESCubicLineSearch - Performs a cubic line search (default line search method).
438 
439    Collective on SNES
440 
441    Input Parameters:
442 +  snes - nonlinear context
443 .  lsctx - optional context for line search (not used here)
444 .  x - current iterate
445 .  f - residual evaluated at x
446 .  y - search direction (contains new iterate on output)
447 .  w - work vector
448 -  fnorm - 2-norm of f
449 
450    Output Parameters:
451 +  g - residual evaluated at new iterate y
452 .  y - new iterate (contains search direction on input)
453 .  gnorm - 2-norm of g
454 .  ynorm - 2-norm of search length
455 -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
456 
457    Options Database Key:
458 $  -snes_ls cubic - Activates SNESCubicLineSearch()
459 
460    Notes:
461    This line search is taken from "Numerical Methods for Unconstrained
462    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
463 
464    Level: advanced
465 
466 .keywords: SNES, nonlinear, line search, cubic
467 
468 .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
469 @*/
470 PetscErrorCode SNESCubicLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
471 {
472   /*
473      Note that for line search purposes we work with with the related
474      minimization problem:
475         min  z(x):  R^n -> R,
476      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
477    */
478 
479   PetscReal   steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
480   PetscReal   maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg;
481 #if defined(PETSC_USE_COMPLEX)
482   PetscScalar cinitslope,clambda;
483 #endif
484   PetscErrorCode ierr;
485   PetscInt count;
486   SNES_LS     *neP = (SNES_LS*)snes->data;
487   PetscScalar mone = -1.0,scale;
488   PetscTruth  change_y = PETSC_FALSE;
489 
490   PetscFunctionBegin;
491   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
492   *flag   = PETSC_TRUE;
493   alpha   = neP->alpha;
494   maxstep = neP->maxstep;
495   steptol = neP->steptol;
496 
497   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
498   if (*ynorm == 0.0) {
499     ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Search direction and size is 0\n"));CHKERRQ(ierr);
500     *gnorm = fnorm;
501     ierr   = VecCopy(x,y);CHKERRQ(ierr);
502     ierr   = VecCopy(f,g);CHKERRQ(ierr);
503     *flag  = PETSC_FALSE;
504     goto theend1;
505   }
506   if (*ynorm > maxstep) {	/* Step too big, so scale back */
507     scale = maxstep/(*ynorm);
508 #if defined(PETSC_USE_COMPLEX)
509     ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",PetscRealPart(scale),*ynorm));CHKERRQ(ierr);
510 #else
511     ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",scale,*ynorm));CHKERRQ(ierr);
512 #endif
513     ierr = VecScale(&scale,y);CHKERRQ(ierr);
514     *ynorm = maxstep;
515   }
516   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
517   minlambda = steptol/rellength;
518   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
519 #if defined(PETSC_USE_COMPLEX)
520   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
521   initslope = PetscRealPart(cinitslope);
522 #else
523   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
524 #endif
525   if (initslope > 0.0) initslope = -initslope;
526   if (initslope == 0.0) initslope = -1.0;
527 
528   ierr = VecCopy(y,w);CHKERRQ(ierr);
529   ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr);
530   if (snes->nfuncs >= snes->max_funcs) {
531     ierr = PetscLogInfo((snes,"SNESCubicLineSearch:Exceeded maximum function evaluations, while checking full step length!\n"));CHKERRQ(ierr);
532     ierr = VecCopy(w,y);CHKERRQ(ierr);
533     *flag = PETSC_FALSE;
534     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
535     goto theend1;
536   }
537   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
538   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
539   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
540   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
541     ierr = VecCopy(w,y);CHKERRQ(ierr);
542     ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Using full step\n"));CHKERRQ(ierr);
543     goto theend1;
544   }
545 
546   /* Fit points with quadratic */
547   lambda = 1.0;
548   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
549   lambdaprev = lambda;
550   printf("tmp %g ddsd %g %g %g %g\n",lambdatemp,((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope),*gnorm,fnorm,initslope);
551   gnormprev = *gnorm;
552   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
553   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
554   else                         lambda = lambdatemp;
555   ierr   = VecCopy(x,w);CHKERRQ(ierr);
556   lambdaneg = -lambda;
557 #if defined(PETSC_USE_COMPLEX)
558   clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
559 #else
560   ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
561 #endif
562   if (snes->nfuncs >= snes->max_funcs) {
563     ierr = PetscLogInfo((snes,"SNESCubicLineSearch:Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n"));CHKERRQ(ierr);
564     ierr = VecCopy(w,y);CHKERRQ(ierr);
565     *flag = PETSC_FALSE;
566     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
567     goto theend1;
568   }
569   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
570   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
571   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
572   if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
573     ierr = VecCopy(w,y);CHKERRQ(ierr);
574     ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%18.16e\n",lambda));CHKERRQ(ierr);
575     goto theend1;
576   }
577 
578   /* Fit points with cubic */
579   count = 1;
580   while (count < 10000) {
581     if (lambda <= minlambda) { /* bad luck; use full step */
582       ierr = PetscLogInfo((snes,"SNESCubicLineSearch:Unable to find good step length! %D \n",count));CHKERRQ(ierr);
583       ierr = PetscLogInfo((snes,"SNESCubicLineSearch:fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope));CHKERRQ(ierr);
584       ierr = VecCopy(x,y);CHKERRQ(ierr);
585       *flag = PETSC_FALSE;
586       break;
587     }
588     printf("lamd %g %g\n",lambda,lambdaprev);
589     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
590     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
591     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
592     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
593     d  = b*b - 3*a*initslope;
594     if (d < 0.0) d = 0.0;
595     if (a == 0.0) {
596       lambdatemp = -initslope/(2.0*b);
597     } else {
598       lambdatemp = (-b + sqrt(d))/(3.0*a);
599     }
600     lambdaprev = lambda;
601     gnormprev  = *gnorm;
602     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
603     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
604     else                         lambda     = lambdatemp;
605     ierr = VecCopy(x,w);CHKERRQ(ierr);
606     lambdaneg = -lambda;
607 #if defined(PETSC_USE_COMPLEX)
608     clambda = lambdaneg;
609     ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
610 #else
611     ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
612 #endif
613     if (snes->nfuncs >= snes->max_funcs) {
614       ierr = PetscLogInfo((snes,"SNESCubicLineSearch:Exceeded maximum function evaluations, while looking for good step length! %D \n",count));CHKERRQ(ierr);
615       ierr = PetscLogInfo((snes,"SNESCubicLineSearch:fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope));CHKERRQ(ierr);
616       ierr = VecCopy(w,y);CHKERRQ(ierr);
617       *flag = PETSC_FALSE;
618       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
619       break;
620     }
621     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
622     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
623     if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
624     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */
625       ierr = VecCopy(w,y);CHKERRQ(ierr);
626       ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Cubically determined step, lambda=%18.16e\n",lambda));CHKERRQ(ierr);
627       break;
628     } else {
629       ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda,  lambda=%18.16e\n",lambda));CHKERRQ(ierr);
630     }
631     count++;
632   }
633   if (count >= 10000) {
634     SETERRQ(PETSC_ERR_LIB, "Lambda was decreased more than 10,000 times, so something is probably wrong with the function evaluation");
635   }
636   theend1:
637   /* Optional user-defined check for line search step validity */
638   if (neP->CheckStep && *flag) {
639     ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
640     if (change_y) { /* recompute the function if the step has changed */
641       ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr);
642       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
643       if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
644       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
645       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
646       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
647     }
648   }
649   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
650   PetscFunctionReturn(0);
651 }
652 /* -------------------------------------------------------------------------- */
653 #undef __FUNCT__
654 #define __FUNCT__ "SNESQuadraticLineSearch"
655 /*@C
656    SNESQuadraticLineSearch - Performs a quadratic line search.
657 
658    Collective on SNES and Vec
659 
660    Input Parameters:
661 +  snes - the SNES context
662 .  lsctx - optional context for line search (not used here)
663 .  x - current iterate
664 .  f - residual evaluated at x
665 .  y - search direction (contains new iterate on output)
666 .  w - work vector
667 -  fnorm - 2-norm of f
668 
669    Output Parameters:
670 +  g - residual evaluated at new iterate y
671 .  y - new iterate (contains search direction on input)
672 .  gnorm - 2-norm of g
673 .  ynorm - 2-norm of search length
674 -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
675 
676    Options Database Key:
677 .  -snes_ls quadratic - Activates SNESQuadraticLineSearch()
678 
679    Notes:
680    Use SNESSetLineSearch() to set this routine within the SNESLS method.
681 
682    Level: advanced
683 
684 .keywords: SNES, nonlinear, quadratic, line search
685 
686 .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
687 @*/
688 PetscErrorCode SNESQuadraticLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
689 {
690   /*
691      Note that for line search purposes we work with with the related
692      minimization problem:
693         min  z(x):  R^n -> R,
694      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
695    */
696   PetscReal   steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg,rellength;
697 #if defined(PETSC_USE_COMPLEX)
698   PetscScalar cinitslope,clambda;
699 #endif
700   PetscErrorCode ierr;
701   PetscInt count;
702   SNES_LS     *neP = (SNES_LS*)snes->data;
703   PetscScalar mone = -1.0,scale;
704   PetscTruth  change_y = PETSC_FALSE;
705 
706   PetscFunctionBegin;
707   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
708   *flag   = PETSC_TRUE;
709   alpha   = neP->alpha;
710   maxstep = neP->maxstep;
711   steptol = neP->steptol;
712 
713   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
714   if (*ynorm == 0.0) {
715     ierr = PetscLogInfo((snes,"SNESQuadraticLineSearch: Search direction and size is 0\n"));CHKERRQ(ierr);
716     *gnorm = fnorm;
717     ierr   = VecCopy(x,y);CHKERRQ(ierr);
718     ierr   = VecCopy(f,g);CHKERRQ(ierr);
719     *flag  = PETSC_FALSE;
720     goto theend2;
721   }
722   if (*ynorm > maxstep) {	/* Step too big, so scale back */
723     scale = maxstep/(*ynorm);
724     ierr = VecScale(&scale,y);CHKERRQ(ierr);
725     *ynorm = maxstep;
726   }
727   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
728   minlambda = steptol/rellength;
729   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
730 #if defined(PETSC_USE_COMPLEX)
731   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
732   initslope = PetscRealPart(cinitslope);
733 #else
734   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
735 #endif
736   if (initslope > 0.0) initslope = -initslope;
737   if (initslope == 0.0) initslope = -1.0;
738 
739   ierr = VecCopy(y,w);CHKERRQ(ierr);
740   ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr);
741   if (snes->nfuncs >= snes->max_funcs) {
742     ierr = PetscLogInfo((snes,"SNESQuadraticLineSearch:Exceeded maximum function evaluations, while checking full step length!\n"));CHKERRQ(ierr);
743     ierr = VecCopy(w,y);CHKERRQ(ierr);
744     *flag = PETSC_FALSE;
745     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
746     goto theend2;
747   }
748   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
749   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
750   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
751   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
752     ierr = VecCopy(w,y);CHKERRQ(ierr);
753     ierr = PetscLogInfo((snes,"SNESQuadraticLineSearch: Using full step\n"));CHKERRQ(ierr);
754     goto theend2;
755   }
756 
757   /* Fit points with quadratic */
758   lambda = 1.0;
759   count = 1;
760   while (PETSC_TRUE) {
761     if (lambda <= minlambda) { /* bad luck; use full step */
762       ierr = PetscLogInfo((snes,"SNESQuadraticLineSearch:Unable to find good step length! %D \n",count));CHKERRQ(ierr);
763       ierr = PetscLogInfo((snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",fnorm,*gnorm,*ynorm,lambda,initslope));CHKERRQ(ierr);
764       ierr = VecCopy(x,y);CHKERRQ(ierr);
765       *flag = PETSC_FALSE;
766       break;
767     }
768     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
769     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
770     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
771     else                         lambda     = lambdatemp;
772     ierr = VecCopy(x,w);CHKERRQ(ierr);
773     lambdaneg = -lambda;
774 #if defined(PETSC_USE_COMPLEX)
775     clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
776 #else
777     ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
778 #endif
779     if (snes->nfuncs >= snes->max_funcs) {
780       ierr = PetscLogInfo((snes,"SNESQuadraticLineSearch:Exceeded maximum function evaluations, while looking for good step length! %D \n",count));CHKERRQ(ierr);
781       ierr = PetscLogInfo((snes,"SNESQuadraticLineSearch:fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope));CHKERRQ(ierr);
782       ierr = VecCopy(w,y);CHKERRQ(ierr);
783       *flag = PETSC_FALSE;
784       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
785       break;
786     }
787     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
788     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
789     if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
790     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
791       ierr = VecCopy(w,y);CHKERRQ(ierr);
792       ierr = PetscLogInfo((snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda));CHKERRQ(ierr);
793       break;
794     }
795     count++;
796   }
797   theend2:
798   /* Optional user-defined check for line search step validity */
799   if (neP->CheckStep) {
800     ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
801     if (change_y) { /* recompute the function if the step has changed */
802       ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr);
803       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
804       if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
805       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
806       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
807       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
808     }
809   }
810   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
811   PetscFunctionReturn(0);
812 }
813 
814 /* -------------------------------------------------------------------------- */
815 #undef __FUNCT__
816 #define __FUNCT__ "SNESSetLineSearch"
817 /*@C
818    SNESSetLineSearch - Sets the line search routine to be used
819    by the method SNESLS.
820 
821    Input Parameters:
822 +  snes - nonlinear context obtained from SNESCreate()
823 .  lsctx - optional user-defined context for use by line search
824 -  func - pointer to int function
825 
826    Collective on SNES
827 
828    Available Routines:
829 +  SNESCubicLineSearch() - default line search
830 .  SNESQuadraticLineSearch() - quadratic line search
831 .  SNESNoLineSearch() - the full Newton step (actually not a line search)
832 -  SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel)
833 
834     Options Database Keys:
835 +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
836 .   -snes_ls_alpha <alpha> - Sets alpha
837 .   -snes_ls_maxstep <max> - Sets maxstep
838 -   -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code
839                    will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length
840                    the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x)
841 
842    Calling sequence of func:
843 .vb
844    func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
845          PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
846 .ve
847 
848     Input parameters for func:
849 +   snes - nonlinear context
850 .   lsctx - optional user-defined context for line search
851 .   x - current iterate
852 .   f - residual evaluated at x
853 .   y - search direction (contains new iterate on output)
854 .   w - work vector
855 -   fnorm - 2-norm of f
856 
857     Output parameters for func:
858 +   g - residual evaluated at new iterate y
859 .   y - new iterate (contains search direction on input)
860 .   gnorm - 2-norm of g
861 .   ynorm - 2-norm of search length
862 -   flag - set to PETSC_TRUE if the line search succeeds; PETSC_FALSE on failure.
863 
864     Level: advanced
865 
866 .keywords: SNES, nonlinear, set, line search, routine
867 
868 .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(),
869           SNESSetLineSearchCheck(), SNESSetLineSearchParams(), SNESGetLineSearchParams()
870 @*/
871 PetscErrorCode SNESSetLineSearch(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void *lsctx)
872 {
873   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void*);
874 
875   PetscFunctionBegin;
876   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void (**)(void))&f);CHKERRQ(ierr);
877   if (f) {
878     ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr);
879   }
880   PetscFunctionReturn(0);
881 }
882 
883 typedef PetscErrorCode (*FCN2)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*); /* force argument to next function to not be extern C*/
884 /* -------------------------------------------------------------------------- */
885 EXTERN_C_BEGIN
886 #undef __FUNCT__
887 #define __FUNCT__ "SNESSetLineSearch_LS"
888 PetscErrorCode SNESSetLineSearch_LS(SNES snes,FCN2 func,void *lsctx)
889 {
890   PetscFunctionBegin;
891   ((SNES_LS *)(snes->data))->LineSearch = func;
892   ((SNES_LS *)(snes->data))->lsP        = lsctx;
893   PetscFunctionReturn(0);
894 }
895 EXTERN_C_END
896 /* -------------------------------------------------------------------------- */
897 #undef __FUNCT__
898 #define __FUNCT__ "SNESSetLineSearchCheck"
899 /*@C
900    SNESSetLineSearchCheck - Sets a routine to check the validity of new iterate computed
901    by the line search routine in the Newton-based method SNESLS.
902 
903    Input Parameters:
904 +  snes - nonlinear context obtained from SNESCreate()
905 .  func - pointer to int function
906 -  checkctx - optional user-defined context for use by step checking routine
907 
908    Collective on SNES
909 
910    Calling sequence of func:
911 .vb
912    int func (SNES snes, void *checkctx, Vec x, PetscTruth *flag)
913 .ve
914    where func returns an error code of 0 on success and a nonzero
915    on failure.
916 
917    Input parameters for func:
918 +  snes - nonlinear context
919 .  checkctx - optional user-defined context for use by step checking routine
920 -  x - current candidate iterate
921 
922    Output parameters for func:
923 +  x - current iterate (possibly modified)
924 -  flag - flag indicating whether x has been modified (either
925            PETSC_TRUE of PETSC_FALSE)
926 
927    Level: advanced
928 
929    Notes:
930    SNESNoLineSearch() and SNESNoLineSearchNoNorms() accept the new
931    iterate computed by the line search checking routine.  In particular,
932    these routines (1) compute a candidate iterate u_{i+1}, (2) pass control
933    to the checking routine, and then (3) compute the corresponding nonlinear
934    function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}.
935 
936    SNESQuadraticLineSearch() and SNESCubicLineSearch() also accept the
937    new iterate computed by the line search checking routine.  In particular,
938    these routines (1) compute a candidate iterate u_{i+1} as well as a
939    candidate nonlinear function f(u_{i+1}), (2) pass control to the checking
940    routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes
941    were made to the candidate iterate in the checking routine (as indicated
942    by flag=PETSC_TRUE).  The overhead of this function re-evaluation can be
943    very costly, so use this feature with caution!
944 
945 .keywords: SNES, nonlinear, set, line search check, step check, routine
946 
947 .seealso: SNESSetLineSearch()
948 @*/
949 PetscErrorCode SNESSetLineSearchCheck(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx)
950 {
951   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,void*,Vec,PetscTruth*),void*);
952 
953   PetscFunctionBegin;
954   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearchCheck_C",(void (**)(void))&f);CHKERRQ(ierr);
955   if (f) {
956     ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr);
957   }
958   PetscFunctionReturn(0);
959 }
960 /* -------------------------------------------------------------------------- */
961 typedef PetscErrorCode (*FCN)(SNES,void*,Vec,PetscTruth*); /* force argument to next function to not be extern C*/
962 EXTERN_C_BEGIN
963 #undef __FUNCT__
964 #define __FUNCT__ "SNESSetLineSearchCheck_LS"
965 PetscErrorCode SNESSetLineSearchCheck_LS(SNES snes,FCN func,void *checkctx)
966 {
967   PetscFunctionBegin;
968   ((SNES_LS *)(snes->data))->CheckStep = func;
969   ((SNES_LS *)(snes->data))->checkP    = checkctx;
970   PetscFunctionReturn(0);
971 }
972 EXTERN_C_END
973 /* -------------------------------------------------------------------------- */
974 /*
975    SNESPrintHelp_LS - Prints all options for the SNES_LS method.
976 
977    Input Parameter:
978 .  snes - the SNES context
979 
980    Application Interface Routine: SNESPrintHelp()
981 */
982 #undef __FUNCT__
983 #define __FUNCT__ "SNESPrintHelp_LS"
984 static PetscErrorCode SNESPrintHelp_LS(SNES snes,char *p)
985 {
986   SNES_LS *ls = (SNES_LS *)snes->data;
987 
988   PetscFunctionBegin;
989   (*PetscHelpPrintf)(snes->comm," method SNES_LS (ls) for systems of nonlinear equations:\n");
990   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls [cubic,quadratic,basic,basicnonorms,...]\n",p);
991   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls_alpha <alpha> (default %g)\n",p,ls->alpha);
992   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls_maxstep <max> (default %g)\n",p,ls->maxstep);
993   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls_steptol <tol> (default %g)\n",p,ls->steptol);
994   PetscFunctionReturn(0);
995 }
996 
997 /*
998    SNESView_LS - Prints info from the SNESLS data structure.
999 
1000    Input Parameters:
1001 .  SNES - the SNES context
1002 .  viewer - visualization context
1003 
1004    Application Interface Routine: SNESView()
1005 */
1006 #undef __FUNCT__
1007 #define __FUNCT__ "SNESView_LS"
1008 static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer)
1009 {
1010   SNES_LS    *ls = (SNES_LS *)snes->data;
1011   const char *cstr;
1012   PetscErrorCode ierr;
1013   PetscTruth iascii;
1014 
1015   PetscFunctionBegin;
1016   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
1017   if (iascii) {
1018     if (ls->LineSearch == SNESNoLineSearch)             cstr = "SNESNoLineSearch";
1019     else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch";
1020     else if (ls->LineSearch == SNESCubicLineSearch)     cstr = "SNESCubicLineSearch";
1021     else                                                cstr = "unknown";
1022     ierr = PetscViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
1023     ierr = PetscViewerASCIIPrintf(viewer,"  alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr);
1024   } else {
1025     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name);
1026   }
1027   PetscFunctionReturn(0);
1028 }
1029 /* -------------------------------------------------------------------------- */
1030 /*
1031    SNESSetFromOptions_LS - Sets various parameters for the SNESLS method.
1032 
1033    Input Parameter:
1034 .  snes - the SNES context
1035 
1036    Application Interface Routine: SNESSetFromOptions()
1037 */
1038 #undef __FUNCT__
1039 #define __FUNCT__ "SNESSetFromOptions_LS"
1040 static PetscErrorCode SNESSetFromOptions_LS(SNES snes)
1041 {
1042   SNES_LS    *ls = (SNES_LS *)snes->data;
1043   const char *lses[] = {"basic","basicnonorms","quadratic","cubic"};
1044   PetscErrorCode ierr;
1045   PetscInt indx;
1046   PetscTruth flg;
1047 
1048   PetscFunctionBegin;
1049   ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr);
1050     ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr);
1051     ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr);
1052     ierr = PetscOptionsReal("-snes_ls_steptol","Step must be greater than","None",ls->steptol,&ls->steptol,0);CHKERRQ(ierr);
1053 
1054     ierr = PetscOptionsEList("-snes_ls","Line search used","SNESSetLineSearch",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr);
1055     if (flg) {
1056       switch (indx) {
1057       case 0:
1058         ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr);
1059         break;
1060       case 1:
1061         ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
1062         break;
1063       case 2:
1064         ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr);
1065         break;
1066       case 3:
1067         ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr);
1068         break;
1069       }
1070     }
1071   ierr = PetscOptionsTail();CHKERRQ(ierr);
1072   PetscFunctionReturn(0);
1073 }
1074 /* -------------------------------------------------------------------------- */
1075 /*MC
1076       SNESLS - Newton based nonlinear solver that uses a line search
1077 
1078    Options Database:
1079 +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
1080 .   -snes_ls_alpha <alpha> - Sets alpha
1081 .   -snes_ls_maxstep <max> - Sets maxstep
1082 -   -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code
1083                    will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length
1084                    the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x)
1085 
1086     Notes: This is the default nonlinear solver in SNES
1087 
1088    Level: beginner
1089 
1090 .seealso:  SNESCreate(), SNES, SNESSetType(), SNESTR, SNESSetLineSearch(),
1091            SNESSetLineSearchCheck(), SNESNoLineSearch(), SNESCubicLineSearch(), SNESQuadraticLineSearch(),
1092           SNESSetLineSearch(), SNESNoLineSearchNoNorms()
1093 
1094 M*/
1095 EXTERN_C_BEGIN
1096 #undef __FUNCT__
1097 #define __FUNCT__ "SNESCreate_LS"
1098 PetscErrorCode SNESCreate_LS(SNES snes)
1099 {
1100   PetscErrorCode ierr;
1101   SNES_LS *neP;
1102 
1103   PetscFunctionBegin;
1104   snes->setup		= SNESSetUp_LS;
1105   snes->solve		= SNESSolve_LS;
1106   snes->destroy		= SNESDestroy_LS;
1107   snes->converged	= SNESConverged_LS;
1108   snes->printhelp       = SNESPrintHelp_LS;
1109   snes->setfromoptions  = SNESSetFromOptions_LS;
1110   snes->view            = SNESView_LS;
1111   snes->nwork           = 0;
1112 
1113   ierr                  = PetscNew(SNES_LS,&neP);CHKERRQ(ierr);
1114   ierr = PetscLogObjectMemory(snes,sizeof(SNES_LS));CHKERRQ(ierr);
1115   snes->data    	= (void*)neP;
1116   neP->alpha		= 1.e-4;
1117   neP->maxstep		= 1.e8;
1118   neP->steptol		= 1.e-12;
1119   neP->LineSearch       = SNESCubicLineSearch;
1120   neP->lsP              = PETSC_NULL;
1121   neP->CheckStep        = PETSC_NULL;
1122   neP->checkP           = PETSC_NULL;
1123 
1124   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS",SNESSetLineSearch_LS);CHKERRQ(ierr);
1125   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS",SNESSetLineSearchCheck_LS);CHKERRQ(ierr);
1126 
1127   PetscFunctionReturn(0);
1128 }
1129 EXTERN_C_END
1130 
1131 
1132 
1133