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