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