xref: /petsc/src/snes/impls/ls/ls.c (revision dc357ad2be94c75fab182994dada4e07693f80fb)
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 = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);	/* fnorm <- ||F||  */
165   if (fnorm != fnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
166   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
167   snes->norm = fnorm;
168   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
169   SNESLogConvHistory(snes,fnorm,0);
170   SNESMonitor(snes,0,fnorm);
171 
172   /* set parameter for default relative tolerance convergence test */
173   snes->ttol = fnorm*snes->rtol;
174   /* test convergence */
175   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
176   if (snes->reason) PetscFunctionReturn(0);
177 
178   for (i=0; i<maxits; i++) {
179 
180     /* Call general purpose update function */
181     if (snes->ops->update) {
182       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
183     }
184 
185     /* Solve J Y = F, where J is Jacobian matrix */
186     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
187     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
188     ierr = SNES_KSPSolve(snes,snes->ksp,F,Y);CHKERRQ(ierr);
189     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
190     if (kspreason < 0) {
191       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
192         ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr);
193         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
194         break;
195       }
196     }
197     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
198     snes->linear_its += lits;
199     ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
200 
201     if (neP->precheckstep) {
202       PetscTruth changed_y = PETSC_FALSE;
203       ierr = (*neP->precheckstep)(snes,X,Y,neP->precheck,&changed_y);CHKERRQ(ierr);
204     }
205 
206     if (PetscLogPrintInfo){
207       ierr = SNESLSCheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
208     }
209 
210     /* Compute a (scaled) negative update in the line search routine:
211          Y <- X - lambda*Y
212        and evaluate G = function(Y) (depends on the line search).
213     */
214     ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr);
215     ynorm = 1; gnorm = fnorm;
216     ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
217     ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",fnorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr);
218     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
219     if (snes->domainerror) {
220       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
221       PetscFunctionReturn(0);
222     }
223     if (!lssucceed) {
224       if (++snes->numFailures >= snes->maxFailures) {
225 	PetscTruth ismin;
226         snes->reason = SNES_DIVERGED_LS_FAILURE;
227         ierr = SNESLSCheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr);
228         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
229         break;
230       }
231     }
232     /* Update function and solution vectors */
233     fnorm = gnorm;
234     ierr = VecCopy(G,F);CHKERRQ(ierr);
235     ierr = VecCopy(W,X);CHKERRQ(ierr);
236     /* Monitor convergence */
237     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
238     snes->iter = i+1;
239     snes->norm = fnorm;
240     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
241     SNESLogConvHistory(snes,snes->norm,lits);
242     SNESMonitor(snes,snes->iter,snes->norm);
243     /* Test for convergence, xnorm = || X || */
244     if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); }
245     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
246     if (snes->reason) break;
247   }
248   if (i == maxits) {
249     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
250     if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
251   }
252   PetscFunctionReturn(0);
253 }
254 /* -------------------------------------------------------------------------- */
255 /*
256    SNESSetUp_LS - Sets up the internal data structures for the later use
257    of the SNESLS nonlinear solver.
258 
259    Input Parameter:
260 .  snes - the SNES context
261 .  x - the solution vector
262 
263    Application Interface Routine: SNESSetUp()
264 
265    Notes:
266    For basic use of the SNES solvers, the user need not explicitly call
267    SNESSetUp(), since these actions will automatically occur during
268    the call to SNESSolve().
269  */
270 #undef __FUNCT__
271 #define __FUNCT__ "SNESSetUp_LS"
272 PetscErrorCode SNESSetUp_LS(SNES snes)
273 {
274   PetscErrorCode ierr;
275 
276   PetscFunctionBegin;
277   if (!snes->vec_sol_update) {
278     ierr = VecDuplicate(snes->vec_sol,&snes->vec_sol_update);CHKERRQ(ierr);
279     ierr = PetscLogObjectParent(snes,snes->vec_sol_update);CHKERRQ(ierr);
280   }
281   if (!snes->work) {
282     snes->nwork = 3;
283     ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
284     ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr);
285   }
286   PetscFunctionReturn(0);
287 }
288 /* -------------------------------------------------------------------------- */
289 /*
290    SNESDestroy_LS - Destroys the private SNES_LS context that was created
291    with SNESCreate_LS().
292 
293    Input Parameter:
294 .  snes - the SNES context
295 
296    Application Interface Routine: SNESDestroy()
297  */
298 #undef __FUNCT__
299 #define __FUNCT__ "SNESDestroy_LS"
300 PetscErrorCode SNESDestroy_LS(SNES snes)
301 {
302   PetscErrorCode ierr;
303 
304   PetscFunctionBegin;
305   if (snes->vec_sol_update) {
306     ierr = VecDestroy(snes->vec_sol_update);CHKERRQ(ierr);
307     snes->vec_sol_update = PETSC_NULL;
308   }
309   if (snes->nwork) {
310     ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr);
311     snes->nwork = 0;
312     snes->work  = PETSC_NULL;
313   }
314   ierr = PetscFree(snes->data);CHKERRQ(ierr);
315 
316   /* clear composed functions */
317   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);CHKERRQ(ierr);
318   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","",PETSC_NULL);CHKERRQ(ierr);
319   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C","",PETSC_NULL);CHKERRQ(ierr);
320 
321   PetscFunctionReturn(0);
322 }
323 /* -------------------------------------------------------------------------- */
324 #undef __FUNCT__
325 #define __FUNCT__ "SNESLineSearchNo"
326 
327 /*@C
328    SNESLineSearchNo - This routine is not a line search at all;
329    it simply uses the full Newton step.  Thus, this routine is intended
330    to serve as a template and is not recommended for general use.
331 
332    Collective on SNES and Vec
333 
334    Input Parameters:
335 +  snes - nonlinear context
336 .  lsctx - optional context for line search (not used here)
337 .  x - current iterate
338 .  f - residual evaluated at x
339 .  y - search direction
340 .  w - work vector
341 .  fnorm - 2-norm of f
342 -  xnorm - norm of x if known, otherwise 0
343 
344    Output Parameters:
345 +  g - residual evaluated at new iterate y
346 .  w - new iterate
347 .  gnorm - 2-norm of g
348 .  ynorm - 2-norm of search length
349 -  flag - PETSC_TRUE on success, PETSC_FALSE on failure
350 
351    Options Database Key:
352 .  -snes_ls basic - Activates SNESLineSearchNo()
353 
354    Level: advanced
355 
356 .keywords: SNES, nonlinear, line search, cubic
357 
358 .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(),
359           SNESLineSearchSet(), SNESLineSearchNoNorms()
360 @*/
361 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)
362 {
363   PetscErrorCode ierr;
364   SNES_LS        *neP = (SNES_LS*)snes->data;
365   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
366 
367   PetscFunctionBegin;
368   *flag = PETSC_TRUE;
369   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
370   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);         /* ynorm = || y || */
371   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
372   if (neP->postcheckstep) {
373    ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
374   }
375   if (changed_y) {
376     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
377   }
378   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
379   if (!snes->domainerror) {
380     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);  /* gnorm = || g || */
381     if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
382   }
383   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
384   PetscFunctionReturn(0);
385 }
386 /* -------------------------------------------------------------------------- */
387 #undef __FUNCT__
388 #define __FUNCT__ "SNESLineSearchNoNorms"
389 
390 /*@C
391    SNESLineSearchNoNorms - This routine is not a line search at
392    all; it simply uses the full Newton step. This version does not
393    even compute the norm of the function or search direction; this
394    is intended only when you know the full step is fine and are
395    not checking for convergence of the nonlinear iteration (for
396    example, you are running always for a fixed number of Newton steps).
397 
398    Collective on SNES and Vec
399 
400    Input Parameters:
401 +  snes - nonlinear context
402 .  lsctx - optional context for line search (not used here)
403 .  x - current iterate
404 .  f - residual evaluated at x
405 .  y - search direction
406 .  w - work vector
407 .  fnorm - 2-norm of f
408 -  xnorm - norm of x if known, otherwise 0
409 
410    Output Parameters:
411 +  g - residual evaluated at new iterate y
412 .  w - new iterate
413 .  gnorm - not changed
414 .  ynorm - not changed
415 -  flag - set to PETSC_TRUE indicating a successful line search
416 
417    Options Database Key:
418 .  -snes_ls basicnonorms - Activates SNESLineSearchNoNorms()
419 
420    Notes:
421    SNESLineSearchNoNorms() must be used in conjunction with
422    either the options
423 $     -snes_no_convergence_test -snes_max_it <its>
424    or alternatively a user-defined custom test set via
425    SNESSetConvergenceTest(); or a -snes_max_it of 1,
426    otherwise, the SNES solver will generate an error.
427 
428    During the final iteration this will not evaluate the function at
429    the solution point. This is to save a function evaluation while
430    using pseudo-timestepping.
431 
432    The residual norms printed by monitoring routines such as
433    SNESMonitorDefault() (as activated via -snes_monitor) will not be
434    correct, since they are not computed.
435 
436    Level: advanced
437 
438 .keywords: SNES, nonlinear, line search, cubic
439 
440 .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(),
441           SNESLineSearchSet(), SNESLineSearchNo()
442 @*/
443 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)
444 {
445   PetscErrorCode ierr;
446   SNES_LS        *neP = (SNES_LS*)snes->data;
447   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
448 
449   PetscFunctionBegin;
450   *flag = PETSC_TRUE;
451   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
452   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y      */
453   if (neP->postcheckstep) {
454    ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
455   }
456   if (changed_y) {
457     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
458   }
459 
460   /* don't evaluate function the last time through */
461   if (snes->iter < snes->max_its-1) {
462     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
463   }
464   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
465   PetscFunctionReturn(0);
466 }
467 /* -------------------------------------------------------------------------- */
468 #undef __FUNCT__
469 #define __FUNCT__ "SNESLineSearchCubic"
470 /*@C
471    SNESLineSearchCubic - Performs a cubic line search (default line search method).
472 
473    Collective on SNES
474 
475    Input Parameters:
476 +  snes - nonlinear context
477 .  lsctx - optional context for line search (not used here)
478 .  x - current iterate
479 .  f - residual evaluated at x
480 .  y - search direction
481 .  w - work vector
482 .  fnorm - 2-norm of f
483 -  xnorm - norm of x if known, otherwise 0
484 
485    Output Parameters:
486 +  g - residual evaluated at new iterate y
487 .  w - new iterate
488 .  gnorm - 2-norm of g
489 .  ynorm - 2-norm of search length
490 -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
491 
492    Options Database Key:
493 $  -snes_ls cubic - Activates SNESLineSearchCubic()
494 
495    Notes:
496    This line search is taken from "Numerical Methods for Unconstrained
497    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
498 
499    Level: advanced
500 
501 .keywords: SNES, nonlinear, line search, cubic
502 
503 .seealso: SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms()
504 @*/
505 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)
506 {
507   /*
508      Note that for line search purposes we work with with the related
509      minimization problem:
510         min  z(x):  R^n -> R,
511      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
512    */
513 
514   PetscReal      steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
515   PetscReal      maxstep,minlambda,alpha,lambda,lambdatemp;
516 #if defined(PETSC_USE_COMPLEX)
517   PetscScalar    cinitslope;
518 #endif
519   PetscErrorCode ierr;
520   PetscInt       count;
521   SNES_LS        *neP = (SNES_LS*)snes->data;
522   PetscScalar    scale;
523   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
524 
525   PetscFunctionBegin;
526   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
527   *flag   = PETSC_TRUE;
528   alpha   = neP->alpha;
529   maxstep = neP->maxstep;
530   steptol = snes->xtol;
531 
532   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
533   if (!*ynorm) {
534     ierr = PetscInfo(snes,"Search direction and size is 0\n");CHKERRQ(ierr);
535     *gnorm = fnorm;
536     ierr   = VecCopy(x,w);CHKERRQ(ierr);
537     ierr   = VecCopy(f,g);CHKERRQ(ierr);
538     *flag  = PETSC_FALSE;
539     goto theend1;
540   }
541   if (*ynorm > maxstep) {	/* Step too big, so scale back */
542     scale = maxstep/(*ynorm);
543     ierr = PetscInfo2(snes,"Scaling step by %G old ynorm %G\n",PetscRealPart(scale),*ynorm);CHKERRQ(ierr);
544     ierr = VecScale(y,scale);CHKERRQ(ierr);
545     *ynorm = maxstep;
546   }
547   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
548   minlambda = steptol/rellength;
549   ierr      = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
550 #if defined(PETSC_USE_COMPLEX)
551   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
552   initslope = PetscRealPart(cinitslope);
553 #else
554   ierr      = VecDot(f,w,&initslope);CHKERRQ(ierr);
555 #endif
556   if (initslope > 0.0)  initslope = -initslope;
557   if (initslope == 0.0) initslope = -1.0;
558 
559   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
560   if (snes->nfuncs >= snes->max_funcs) {
561     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
562     *flag = PETSC_FALSE;
563     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
564     goto theend1;
565   }
566   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
567   if (snes->domainerror) {
568     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
569     PetscFunctionReturn(0);
570   }
571   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
572   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
573   ierr = PetscInfo2(snes,"Initial fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr);
574   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
575     ierr = PetscInfo2(snes,"Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr);
576     goto theend1;
577   }
578   if (*ynorm < xnorm*snes->xtol) {
579     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);
580     *flag = PETSC_TRUE;
581     goto theend1;
582   }
583 
584 
585   /* Fit points with quadratic */
586   lambda     = 1.0;
587   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
588   lambdaprev = lambda;
589   gnormprev  = *gnorm;
590   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
591   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
592   else                         lambda = lambdatemp;
593 
594   ierr  = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
595   if (snes->nfuncs >= snes->max_funcs) {
596     ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr);
597     *flag = PETSC_FALSE;
598     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
599     goto theend1;
600   }
601   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
602   if (snes->domainerror) {
603     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
604     PetscFunctionReturn(0);
605   }
606   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
607   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
608   ierr = PetscInfo1(snes,"gnorm after quadratic fit %G\n",*gnorm);CHKERRQ(ierr);
609   if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
610     ierr = PetscInfo1(snes,"Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr);
611     goto theend1;
612   }
613 
614   /* Fit points with cubic */
615   count = 1;
616   while (count < 20) {
617     if (lambda <= minlambda) {
618       ierr = PetscInfo1(snes,"Unable to find good step length! After %D tries \n",count);CHKERRQ(ierr);
619       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);
620       *flag = PETSC_FALSE;
621       break;
622     }
623     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
624     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
625     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
626     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
627     d  = b*b - 3*a*initslope;
628     if (d < 0.0) d = 0.0;
629     if (a == 0.0) {
630       lambdatemp = -initslope/(2.0*b);
631     } else {
632       lambdatemp = (-b + sqrt(d))/(3.0*a);
633     }
634     lambdaprev = lambda;
635     gnormprev  = *gnorm;
636     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
637     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
638     else                         lambda     = lambdatemp;
639     ierr  = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
640     if (snes->nfuncs >= snes->max_funcs) {
641       ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
642       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);
643       *flag = PETSC_FALSE;
644       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
645       break;
646     }
647     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
648     if (snes->domainerror) {
649       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
650       PetscFunctionReturn(0);
651     }
652     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
653     if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
654     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */
655       ierr = PetscInfo2(snes,"Cubically determined step, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr);
656       break;
657     } else {
658       ierr = PetscInfo2(snes,"Cubic step no good, shrinking lambda, current gnorem %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr);
659     }
660     count++;
661   }
662   if (count >= 20) {
663     SETERRQ(PETSC_ERR_LIB, "Lambda was decreased more than 20 times, so something is probably wrong with the function evaluation");
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 (*gnorm != *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 (*gnorm != *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 (*gnorm != *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 (*gnorm != *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