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