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