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