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