xref: /petsc/src/snes/impls/ls/ls.c (revision 7c4f633dc6bb6149cca88d301ead35a99e103cbb)
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 (PETSC_TRUE) {
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   theend1:
666   /* Optional user-defined check for line search step validity */
667   if (neP->postcheckstep && *flag) {
668     ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
669     if (changed_y) {
670       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
671     }
672     if (changed_y || changed_w) { /* recompute the function if the step has changed */
673       ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
674       if (snes->domainerror) {
675         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
676         PetscFunctionReturn(0);
677       }
678       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
679       if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
680       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
681       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
682       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
683     }
684   }
685   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
686   PetscFunctionReturn(0);
687 }
688 /* -------------------------------------------------------------------------- */
689 #undef __FUNCT__
690 #define __FUNCT__ "SNESLineSearchQuadratic"
691 /*@C
692    SNESLineSearchQuadratic - Performs a quadratic line search.
693 
694    Collective on SNES and Vec
695 
696    Input Parameters:
697 +  snes - the SNES context
698 .  lsctx - optional context for line search (not used here)
699 .  x - current iterate
700 .  f - residual evaluated at x
701 .  y - search direction
702 .  w - work vector
703 .  fnorm - 2-norm of f
704 -  xnorm - norm of x if known, otherwise 0
705 
706    Output Parameters:
707 +  g - residual evaluated at new iterate w
708 .  w - new iterate (x + alpha*y)
709 .  gnorm - 2-norm of g
710 .  ynorm - 2-norm of search length
711 -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
712 
713    Options Database Keys:
714 +  -snes_ls quadratic - Activates SNESLineSearchQuadratic()
715 .   -snes_ls_alpha <alpha> - Sets alpha
716 .   -snes_ls_maxstep <max> - Sets maxstep
717 -   -snes_stol <steptol> - Sets steptol, this is the minimum step size that the line search code
718                    will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length
719                    the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x)
720    Notes:
721    Use SNESLineSearchSet() to set this routine within the SNESLS method.
722 
723    Level: advanced
724 
725 .keywords: SNES, nonlinear, quadratic, line search
726 
727 .seealso: SNESLineSearchCubic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms()
728 @*/
729 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)
730 {
731   /*
732      Note that for line search purposes we work with with the related
733      minimization problem:
734         min  z(x):  R^n -> R,
735      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
736    */
737   PetscReal      steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,rellength;
738 #if defined(PETSC_USE_COMPLEX)
739   PetscScalar    cinitslope;
740 #endif
741   PetscErrorCode ierr;
742   PetscInt       count;
743   SNES_LS        *neP = (SNES_LS*)snes->data;
744   PetscScalar    scale;
745   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
746 
747   PetscFunctionBegin;
748   ierr    = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
749   *flag   = PETSC_TRUE;
750   alpha   = neP->alpha;
751   maxstep = neP->maxstep;
752   steptol = snes->xtol;
753 
754   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
755   if (*ynorm == 0.0) {
756     ierr   = PetscInfo(snes,"Search direction and size is 0\n");CHKERRQ(ierr);
757     *gnorm = fnorm;
758     ierr   = VecCopy(x,w);CHKERRQ(ierr);
759     ierr   = VecCopy(f,g);CHKERRQ(ierr);
760     *flag  = PETSC_FALSE;
761     goto theend2;
762   }
763   if (*ynorm > maxstep) {	/* Step too big, so scale back */
764     scale  = maxstep/(*ynorm);
765     ierr   = VecScale(y,scale);CHKERRQ(ierr);
766     *ynorm = maxstep;
767   }
768   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
769   minlambda = steptol/rellength;
770   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
771 #if defined(PETSC_USE_COMPLEX)
772   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
773   initslope = PetscRealPart(cinitslope);
774 #else
775   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
776 #endif
777   if (initslope > 0.0)  initslope = -initslope;
778   if (initslope == 0.0) initslope = -1.0;
779   ierr = PetscInfo1(snes,"Initslope %G \n",initslope);CHKERRQ(ierr);
780 
781   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
782   if (snes->nfuncs >= snes->max_funcs) {
783     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
784     *flag = PETSC_FALSE;
785     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
786     goto theend2;
787   }
788   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
789   if (snes->domainerror) {
790     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
791     PetscFunctionReturn(0);
792   }
793   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
794   if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
795   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
796     ierr = PetscInfo(snes,"Using full step\n");CHKERRQ(ierr);
797     goto theend2;
798   }
799   if (*ynorm < xnorm*snes->xtol) {
800     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);
801     *flag = PETSC_TRUE;
802     goto theend2;
803   }
804 
805   /* Fit points with quadratic */
806   lambda = 1.0;
807   count = 1;
808   while (PETSC_TRUE) {
809     if (lambda <= minlambda) { /* bad luck; use full step */
810       ierr = PetscInfo1(snes,"Unable to find good step length! %D \n",count);CHKERRQ(ierr);
811       ierr = PetscInfo5(snes,"fnorm=%G, gnorm=%G, ynorm=%G, lambda=%G, initial slope=%G\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr);
812       ierr = VecCopy(x,w);CHKERRQ(ierr);
813       *flag = PETSC_FALSE;
814       break;
815     }
816     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
817     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
818     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
819     else                         lambda     = lambdatemp;
820 
821     ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
822     if (snes->nfuncs >= snes->max_funcs) {
823       ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
824       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);
825       *flag = PETSC_FALSE;
826       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
827       break;
828     }
829     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
830     if (snes->domainerror) {
831       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
832       PetscFunctionReturn(0);
833     }
834     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
835     if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
836     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
837       ierr = PetscInfo1(snes,"Quadratically determined step, lambda=%G\n",lambda);CHKERRQ(ierr);
838       break;
839     }
840     count++;
841   }
842   theend2:
843   /* Optional user-defined check for line search step validity */
844   if (neP->postcheckstep) {
845     ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
846     if (changed_y) {
847       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
848     }
849     if (changed_y || changed_w) { /* recompute the function if the step has changed */
850       ierr = SNESComputeFunction(snes,w,g);
851       if (snes->domainerror) {
852         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
853         PetscFunctionReturn(0);
854       }
855       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
856       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
857       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
858       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
859       if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
860     }
861   }
862   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
863   PetscFunctionReturn(0);
864 }
865 
866 /* -------------------------------------------------------------------------- */
867 #undef __FUNCT__
868 #define __FUNCT__ "SNESLineSearchSet"
869 /*@C
870    SNESLineSearchSet - Sets the line search routine to be used
871    by the method SNESLS.
872 
873    Input Parameters:
874 +  snes - nonlinear context obtained from SNESCreate()
875 .  lsctx - optional user-defined context for use by line search
876 -  func - pointer to int function
877 
878    Collective on SNES
879 
880    Available Routines:
881 +  SNESLineSearchCubic() - default line search
882 .  SNESLineSearchQuadratic() - quadratic line search
883 .  SNESLineSearchNo() - the full Newton step (actually not a line search)
884 -  SNESLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel)
885 
886     Options Database Keys:
887 +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
888 .   -snes_ls_alpha <alpha> - Sets alpha
889 .   -snes_ls_maxstep <max> - Sets maxstep
890 -   -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
891                    the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x)
892 
893    Calling sequence of func:
894 .vb
895    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)
896 .ve
897 
898     Input parameters for func:
899 +   snes - nonlinear context
900 .   lsctx - optional user-defined context for line search
901 .   x - current iterate
902 .   f - residual evaluated at x
903 .   y - search direction
904 .   w - work vector
905 -   fnorm - 2-norm of f
906 
907     Output parameters for func:
908 +   g - residual evaluated at new iterate y
909 .   w - new iterate
910 .   gnorm - 2-norm of g
911 .   ynorm - 2-norm of search length
912 -   flag - set to PETSC_TRUE if the line search succeeds; PETSC_FALSE on failure.
913 
914     Level: advanced
915 
916 .keywords: SNES, nonlinear, set, line search, routine
917 
918 .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchNoNorms(),
919           SNESLineSearchSetPostCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams(), SNESLineSearchSetPreCheck()
920 @*/
921 PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void *lsctx)
922 {
923   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void*);
924 
925   PetscFunctionBegin;
926   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSet_C",(void (**)(void))&f);CHKERRQ(ierr);
927   if (f) {
928     ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr);
929   }
930   PetscFunctionReturn(0);
931 }
932 
933 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*/
934 /* -------------------------------------------------------------------------- */
935 EXTERN_C_BEGIN
936 #undef __FUNCT__
937 #define __FUNCT__ "SNESLineSearchSet_LS"
938 PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet_LS(SNES snes,FCN2 func,void *lsctx)
939 {
940   PetscFunctionBegin;
941   ((SNES_LS *)(snes->data))->LineSearch = func;
942   ((SNES_LS *)(snes->data))->lsP        = lsctx;
943   PetscFunctionReturn(0);
944 }
945 EXTERN_C_END
946 /* -------------------------------------------------------------------------- */
947 #undef __FUNCT__
948 #define __FUNCT__ "SNESLineSearchSetPostCheck"
949 /*@C
950    SNESLineSearchSetPostCheck - Sets a routine to check the validity of new iterate computed
951    by the line search routine in the Newton-based method SNESLS.
952 
953    Input Parameters:
954 +  snes - nonlinear context obtained from SNESCreate()
955 .  func - pointer to function
956 -  checkctx - optional user-defined context for use by step checking routine
957 
958    Collective on SNES
959 
960    Calling sequence of func:
961 .vb
962    int func (SNES snes, Vec x,Vec y,Vec w,void *checkctx, PetscTruth *changed_y,PetscTruth *changed_w)
963 .ve
964    where func returns an error code of 0 on success and a nonzero
965    on failure.
966 
967    Input parameters for func:
968 +  snes - nonlinear context
969 .  checkctx - optional user-defined context for use by step checking routine
970 .  x - previous iterate
971 .  y - new search direction and length
972 -  w - current candidate iterate
973 
974    Output parameters for func:
975 +  y - search direction (possibly changed)
976 .  w - current iterate (possibly modified)
977 .  changed_y - indicates search direction was changed by this routine
978 -  changed_w - indicates current iterate was changed by this routine
979 
980    Level: advanced
981 
982    Notes: All line searches accept the new iterate computed by the line search checking routine.
983 
984    Only one of changed_y and changed_w can  be PETSC_TRUE
985 
986    On input w = x + y
987 
988    SNESLineSearchNo() and SNESLineSearchNoNorms() (1) compute a candidate iterate u_{i+1}, (2) pass control
989    to the checking routine, and then (3) compute the corresponding nonlinear
990    function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}.
991 
992    SNESLineSearchQuadratic() and SNESLineSearchCubic() (1) compute a candidate iterate u_{i+1} as well as a
993    candidate nonlinear function f(u_{i+1}), (2) pass control to the checking
994    routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes
995    were made to the candidate iterate in the checking routine (as indicated
996    by flag=PETSC_TRUE).  The overhead of this extra function re-evaluation can be
997    very costly, so use this feature with caution!
998 
999 .keywords: SNES, nonlinear, set, line search check, step check, routine
1000 
1001 .seealso: SNESLineSearchSet(), SNESLineSearchSetPreCheck()
1002 @*/
1003 PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void *checkctx)
1004 {
1005   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void*);
1006 
1007   PetscFunctionBegin;
1008   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPostCheck_C",(void (**)(void))&f);CHKERRQ(ierr);
1009   if (f) {
1010     ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr);
1011   }
1012   PetscFunctionReturn(0);
1013 }
1014 #undef __FUNCT__
1015 #define __FUNCT__ "SNESLineSearchSetPreCheck"
1016 /*@C
1017    SNESLineSearchSetPreCheck - Sets a routine to check the validity of a new direction given by the linear solve
1018          before the line search is called.
1019 
1020    Input Parameters:
1021 +  snes - nonlinear context obtained from SNESCreate()
1022 .  func - pointer to function
1023 -  checkctx - optional user-defined context for use by step checking routine
1024 
1025    Collective on SNES
1026 
1027    Calling sequence of func:
1028 .vb
1029    int func (SNES snes, Vec x,Vec y,void *checkctx, PetscTruth *changed_y)
1030 .ve
1031    where func returns an error code of 0 on success and a nonzero
1032    on failure.
1033 
1034    Input parameters for func:
1035 +  snes - nonlinear context
1036 .  checkctx - optional user-defined context for use by step checking routine
1037 .  x - previous iterate
1038 -  y - new search direction and length
1039 
1040    Output parameters for func:
1041 +  y - search direction (possibly changed)
1042 -  changed_y - indicates search direction was changed by this routine
1043 
1044    Level: advanced
1045 
1046    Notes: All line searches accept the new iterate computed by the line search checking routine.
1047 
1048 .keywords: SNES, nonlinear, set, line search check, step check, routine
1049 
1050 .seealso: SNESLineSearchSet(), SNESLineSearchSetPostCheck(), SNESSetUpdate()
1051 @*/
1052 PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPreCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,void*,PetscTruth*),void *checkctx)
1053 {
1054   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,void*,PetscTruth*),void*);
1055 
1056   PetscFunctionBegin;
1057   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPreCheck_C",(void (**)(void))&f);CHKERRQ(ierr);
1058   if (f) {
1059     ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr);
1060   }
1061   PetscFunctionReturn(0);
1062 }
1063 
1064 /* -------------------------------------------------------------------------- */
1065 typedef PetscErrorCode (*FCN1)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*); /* force argument to next function to not be extern C*/
1066 typedef PetscErrorCode (*FCN3)(SNES,Vec,Vec,void*,PetscTruth*);                 /* force argument to next function to not be extern C*/
1067 EXTERN_C_BEGIN
1068 #undef __FUNCT__
1069 #define __FUNCT__ "SNESLineSearchSetPostCheck_LS"
1070 PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck_LS(SNES snes,FCN1 func,void *checkctx)
1071 {
1072   PetscFunctionBegin;
1073   ((SNES_LS *)(snes->data))->postcheckstep = func;
1074   ((SNES_LS *)(snes->data))->postcheck     = checkctx;
1075   PetscFunctionReturn(0);
1076 }
1077 EXTERN_C_END
1078 
1079 EXTERN_C_BEGIN
1080 #undef __FUNCT__
1081 #define __FUNCT__ "SNESLineSearchSetPreCheck_LS"
1082 PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPreCheck_LS(SNES snes,FCN3 func,void *checkctx)
1083 {
1084   PetscFunctionBegin;
1085   ((SNES_LS *)(snes->data))->precheckstep = func;
1086   ((SNES_LS *)(snes->data))->precheck     = checkctx;
1087   PetscFunctionReturn(0);
1088 }
1089 EXTERN_C_END
1090 
1091 /*
1092    SNESView_LS - Prints info from the SNESLS data structure.
1093 
1094    Input Parameters:
1095 .  SNES - the SNES context
1096 .  viewer - visualization context
1097 
1098    Application Interface Routine: SNESView()
1099 */
1100 #undef __FUNCT__
1101 #define __FUNCT__ "SNESView_LS"
1102 static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer)
1103 {
1104   SNES_LS        *ls = (SNES_LS *)snes->data;
1105   const char     *cstr;
1106   PetscErrorCode ierr;
1107   PetscTruth     iascii;
1108 
1109   PetscFunctionBegin;
1110   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
1111   if (iascii) {
1112     if (ls->LineSearch == SNESLineSearchNo)             cstr = "SNESLineSearchNo";
1113     else if (ls->LineSearch == SNESLineSearchQuadratic) cstr = "SNESLineSearchQuadratic";
1114     else if (ls->LineSearch == SNESLineSearchCubic)     cstr = "SNESLineSearchCubic";
1115     else                                                cstr = "unknown";
1116     ierr = PetscViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
1117     ierr = PetscViewerASCIIPrintf(viewer,"  alpha=%G, maxstep=%G\n",ls->alpha,ls->maxstep);CHKERRQ(ierr);
1118   } else {
1119     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name);
1120   }
1121   PetscFunctionReturn(0);
1122 }
1123 /* -------------------------------------------------------------------------- */
1124 /*
1125    SNESSetFromOptions_LS - Sets various parameters for the SNESLS method.
1126 
1127    Input Parameter:
1128 .  snes - the SNES context
1129 
1130    Application Interface Routine: SNESSetFromOptions()
1131 */
1132 #undef __FUNCT__
1133 #define __FUNCT__ "SNESSetFromOptions_LS"
1134 static PetscErrorCode SNESSetFromOptions_LS(SNES snes)
1135 {
1136   SNES_LS        *ls = (SNES_LS *)snes->data;
1137   const char     *lses[] = {"basic","basicnonorms","quadratic","cubic"};
1138   PetscErrorCode ierr;
1139   PetscInt       indx;
1140   PetscTruth     flg;
1141 
1142   PetscFunctionBegin;
1143   ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr);
1144     ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr);
1145     ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr);
1146 
1147     ierr = PetscOptionsEList("-snes_ls","Line search used","SNESLineSearchSet",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr);
1148     if (flg) {
1149       switch (indx) {
1150       case 0:
1151         ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr);
1152         break;
1153       case 1:
1154         ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
1155         break;
1156       case 2:
1157         ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic,PETSC_NULL);CHKERRQ(ierr);
1158         break;
1159       case 3:
1160         ierr = SNESLineSearchSet(snes,SNESLineSearchCubic,PETSC_NULL);CHKERRQ(ierr);
1161         break;
1162       }
1163     }
1164   ierr = PetscOptionsTail();CHKERRQ(ierr);
1165   PetscFunctionReturn(0);
1166 }
1167 /* -------------------------------------------------------------------------- */
1168 /*MC
1169       SNESLS - Newton based nonlinear solver that uses a line search
1170 
1171    Options Database:
1172 +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
1173 .   -snes_ls_alpha <alpha> - Sets alpha
1174 .   -snes_ls_maxstep <max> - Sets maxstep
1175 -   -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code
1176                    will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length
1177                    the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x)
1178 
1179     Notes: This is the default nonlinear solver in SNES
1180 
1181    Level: beginner
1182 
1183 .seealso:  SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(),
1184            SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(),
1185           SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck()
1186 
1187 M*/
1188 EXTERN_C_BEGIN
1189 #undef __FUNCT__
1190 #define __FUNCT__ "SNESCreate_LS"
1191 PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_LS(SNES snes)
1192 {
1193   PetscErrorCode ierr;
1194   SNES_LS        *neP;
1195 
1196   PetscFunctionBegin;
1197   snes->ops->setup	     = SNESSetUp_LS;
1198   snes->ops->solve	     = SNESSolve_LS;
1199   snes->ops->destroy	     = SNESDestroy_LS;
1200   snes->ops->setfromoptions  = SNESSetFromOptions_LS;
1201   snes->ops->view            = SNESView_LS;
1202 
1203   ierr                  = PetscNewLog(snes,SNES_LS,&neP);CHKERRQ(ierr);
1204   snes->data    	= (void*)neP;
1205   neP->alpha		= 1.e-4;
1206   neP->maxstep		= 1.e8;
1207   neP->LineSearch       = SNESLineSearchCubic;
1208   neP->lsP              = PETSC_NULL;
1209   neP->postcheckstep    = PETSC_NULL;
1210   neP->postcheck        = PETSC_NULL;
1211   neP->precheckstep     = PETSC_NULL;
1212   neP->precheck         = PETSC_NULL;
1213 
1214   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C",
1215 					   "SNESLineSearchSet_LS",
1216 					   SNESLineSearchSet_LS);CHKERRQ(ierr);
1217   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C",
1218 					   "SNESLineSearchSetPostCheck_LS",
1219 					   SNESLineSearchSetPostCheck_LS);CHKERRQ(ierr);
1220   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C",
1221 					   "SNESLineSearchSetPreCheck_LS",
1222 					   SNESLineSearchSetPreCheck_LS);CHKERRQ(ierr);
1223 
1224   PetscFunctionReturn(0);
1225 }
1226 EXTERN_C_END
1227 
1228 
1229 
1230