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