xref: /petsc/src/snes/impls/ls/ls.c (revision 8e3fc8c070edf75bc3ff0904c3a8a3a272dc5aef)
1 
2 #include <../src/snes/impls/ls/lsimpl.h>
3 
4 /*
5      Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function,
6     || 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
7     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
8     for this trick. One assumes that the probability that W is in the null space of J is very, very small.
9 */
10 #undef __FUNCT__
11 #define __FUNCT__ "SNESLSCheckLocalMin_Private"
12 PetscErrorCode SNESLSCheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscBool  *ismin)
13 {
14   PetscReal      a1;
15   PetscErrorCode ierr;
16   PetscBool      hastranspose;
17 
18   PetscFunctionBegin;
19   *ismin = PETSC_FALSE;
20   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
21   if (hastranspose) {
22     /* Compute || J^T F|| */
23     ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr);
24     ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr);
25     ierr = PetscInfo1(snes,"|| J^T F|| %14.12e near zero implies found a local minimum\n",(double)(a1/fnorm));CHKERRQ(ierr);
26     if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
27   } else {
28     Vec         work;
29     PetscScalar result;
30     PetscReal   wnorm;
31 
32     ierr = VecSetRandom(W,PETSC_NULL);CHKERRQ(ierr);
33     ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr);
34     ierr = VecDuplicate(W,&work);CHKERRQ(ierr);
35     ierr = MatMult(A,W,work);CHKERRQ(ierr);
36     ierr = VecDot(F,work,&result);CHKERRQ(ierr);
37     ierr = VecDestroy(&work);CHKERRQ(ierr);
38     a1   = PetscAbsScalar(result)/(fnorm*wnorm);
39     ierr = PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %14.12e near zero implies found a local minimum\n",(double)a1);CHKERRQ(ierr);
40     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
41   }
42   PetscFunctionReturn(0);
43 }
44 
45 /*
46      Checks if J^T(F - J*X) = 0
47 */
48 #undef __FUNCT__
49 #define __FUNCT__ "SNESLSCheckResidual_Private"
50 PetscErrorCode SNESLSCheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2)
51 {
52   PetscReal      a1,a2;
53   PetscErrorCode ierr;
54   PetscBool      hastranspose;
55 
56   PetscFunctionBegin;
57   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
58   if (hastranspose) {
59     ierr = MatMult(A,X,W1);CHKERRQ(ierr);
60     ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr);
61 
62     /* Compute || J^T W|| */
63     ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr);
64     ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr);
65     ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr);
66     if (a1 != 0.0) {
67       ierr = PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %14.12e near zero implies inconsistent rhs\n",(double)(a2/a1));CHKERRQ(ierr);
68     }
69   }
70   PetscFunctionReturn(0);
71 }
72 
73 /*  --------------------------------------------------------------------
74 
75      This file implements a truncated Newton method with a line search,
76      for solving a system of nonlinear equations, using the KSP, Vec,
77      and Mat interfaces for linear solvers, vectors, and matrices,
78      respectively.
79 
80      The following basic routines are required for each nonlinear solver:
81           SNESCreate_XXX()          - Creates a nonlinear solver context
82           SNESSetFromOptions_XXX()  - Sets runtime options
83           SNESSolve_XXX()           - Solves the nonlinear system
84           SNESDestroy_XXX()         - Destroys the nonlinear solver context
85      The suffix "_XXX" denotes a particular implementation, in this case
86      we use _LS (e.g., SNESCreate_LS, SNESSolve_LS) for solving
87      systems of nonlinear equations with a line search (LS) method.
88      These routines are actually called via the common user interface
89      routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and
90      SNESDestroy(), so the application code interface remains identical
91      for all nonlinear solvers.
92 
93      Another key routine is:
94           SNESSetUp_XXX()           - Prepares for the use of a nonlinear solver
95      by setting data structures and options.   The interface routine SNESSetUp()
96      is not usually called directly by the user, but instead is called by
97      SNESSolve() if necessary.
98 
99      Additional basic routines are:
100           SNESView_XXX()            - Prints details of runtime options that
101                                       have actually been used.
102      These are called by application codes via the interface routines
103      SNESView().
104 
105      The various types of solvers (preconditioners, Krylov subspace methods,
106      nonlinear solvers, timesteppers) are all organized similarly, so the
107      above description applies to these categories also.
108 
109     -------------------------------------------------------------------- */
110 /*
111    SNESSolve_LS - Solves a nonlinear system with a truncated Newton
112    method with a line search.
113 
114    Input Parameters:
115 .  snes - the SNES context
116 
117    Output Parameter:
118 .  outits - number of iterations until termination
119 
120    Application Interface Routine: SNESSolve()
121 
122    Notes:
123    This implements essentially a truncated Newton method with a
124    line search.  By default a cubic backtracking line search
125    is employed, as described in the text "Numerical Methods for
126    Unconstrained Optimization and Nonlinear Equations" by Dennis
127    and Schnabel.
128 */
129 #undef __FUNCT__
130 #define __FUNCT__ "SNESSolve_LS"
131 PetscErrorCode SNESSolve_LS(SNES snes)
132 {
133   PetscErrorCode     ierr;
134   PetscInt           maxits,i,lits;
135   PetscBool          lssucceed;
136   MatStructure       flg = DIFFERENT_NONZERO_PATTERN;
137   PetscReal          fnorm,gnorm,xnorm=0,ynorm;
138   Vec                Y,X,F,G,W;
139   KSPConvergedReason kspreason;
140 
141   PetscFunctionBegin;
142   snes->numFailures            = 0;
143   snes->numLinearSolveFailures = 0;
144   snes->reason                 = SNES_CONVERGED_ITERATING;
145 
146   maxits	= snes->max_its;	/* maximum number of iterations */
147   X		= snes->vec_sol;	/* solution vector */
148   F		= snes->vec_func;	/* residual vector */
149   Y		= snes->work[0];	/* work vectors */
150   G		= snes->work[1];
151   W		= snes->work[2];
152 
153   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
154   snes->iter = 0;
155   snes->norm = 0.0;
156   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
157   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
158   if (snes->domainerror) {
159     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
160     PetscFunctionReturn(0);
161   }
162   ierr = VecNormBegin(F,NORM_2,&fnorm);CHKERRQ(ierr);	/* fnorm <- ||F||  */
163   ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm <- ||x||  */
164   ierr = VecNormEnd(F,NORM_2,&fnorm);CHKERRQ(ierr);
165   ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr);
166   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
167   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
168   snes->norm = fnorm;
169   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
170   SNESLogConvHistory(snes,fnorm,0);
171   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
172 
173   /* set parameter for default relative tolerance convergence test */
174   snes->ttol = fnorm*snes->rtol;
175   /* test convergence */
176   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
177   if (snes->reason) PetscFunctionReturn(0);
178 
179   for (i=0; i<maxits; i++) {
180 
181     /* Call general purpose update function */
182     if (snes->ops->update) {
183       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
184     }
185 
186     /* Solve J Y = F, where J is Jacobian matrix */
187     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
188     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
189     ierr = SNES_KSPSolve(snes,snes->ksp,F,Y);CHKERRQ(ierr);
190     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
191     if (kspreason < 0) {
192       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
193         ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr);
194         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
195         break;
196       }
197     }
198     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
199     snes->linear_its += lits;
200     ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
201 
202     if (snes->ops->precheckstep) {
203       PetscBool  changed_y = PETSC_FALSE;
204       ierr = (*snes->ops->precheckstep)(snes,X,Y,snes->precheck,&changed_y);CHKERRQ(ierr);
205     }
206 
207     if (PetscLogPrintInfo){
208       ierr = SNESLSCheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
209     }
210 
211     /* Compute a (scaled) negative update in the line search routine:
212          Y <- X - lambda*Y
213        and evaluate G = function(Y) (depends on the line search).
214     */
215     ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr);
216     ynorm = 1; gnorm = fnorm;
217     ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Y,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
218     ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",(double)fnorm,(double)gnorm,(double)ynorm,(int)lssucceed);CHKERRQ(ierr);
219     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
220     if (snes->domainerror) {
221       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
222       PetscFunctionReturn(0);
223     }
224     if (!lssucceed) {
225       if (++snes->numFailures >= snes->maxFailures) {
226 	PetscBool  ismin;
227         snes->reason = SNES_DIVERGED_LINE_SEARCH;
228         ierr = SNESLSCheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr);
229         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
230         break;
231       }
232     }
233     /* Update function and solution vectors */
234     fnorm = gnorm;
235     ierr = VecCopy(G,F);CHKERRQ(ierr);
236     ierr = VecCopy(W,X);CHKERRQ(ierr);
237     /* Monitor convergence */
238     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
239     snes->iter = i+1;
240     snes->norm = fnorm;
241     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
242     SNESLogConvHistory(snes,snes->norm,lits);
243     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
244     /* Test for convergence, xnorm = || X || */
245     if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); }
246     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
247     if (snes->reason) break;
248   }
249   if (i == maxits) {
250     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
251     if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
252   }
253   PetscFunctionReturn(0);
254 }
255 /* -------------------------------------------------------------------------- */
256 /*
257    SNESSetUp_LS - Sets up the internal data structures for the later use
258    of the SNESLS nonlinear solver.
259 
260    Input Parameter:
261 .  snes - the SNES context
262 .  x - the solution vector
263 
264    Application Interface Routine: SNESSetUp()
265 
266    Notes:
267    For basic use of the SNES solvers, the user need not explicitly call
268    SNESSetUp(), since these actions will automatically occur during
269    the call to SNESSolve().
270  */
271 #undef __FUNCT__
272 #define __FUNCT__ "SNESSetUp_LS"
273 PetscErrorCode SNESSetUp_LS(SNES snes)
274 {
275   PetscErrorCode ierr;
276 
277   PetscFunctionBegin;
278   ierr = SNESDefaultGetWork(snes,3);CHKERRQ(ierr);
279   PetscFunctionReturn(0);
280 }
281 /* -------------------------------------------------------------------------- */
282 
283 #undef __FUNCT__
284 #define __FUNCT__ "SNESReset_LS"
285 PetscErrorCode SNESReset_LS(SNES snes)
286 {
287   PetscErrorCode ierr;
288 
289   PetscFunctionBegin;
290   if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);}
291   PetscFunctionReturn(0);
292 }
293 
294 /*
295    SNESDestroy_LS - Destroys the private SNES_LS context that was created
296    with SNESCreate_LS().
297 
298    Input Parameter:
299 .  snes - the SNES context
300 
301    Application Interface Routine: SNESDestroy()
302  */
303 #undef __FUNCT__
304 #define __FUNCT__ "SNESDestroy_LS"
305 PetscErrorCode SNESDestroy_LS(SNES snes)
306 {
307   PetscErrorCode ierr;
308 
309   PetscFunctionBegin;
310   ierr = SNESReset_LS(snes);CHKERRQ(ierr);
311   ierr = PetscFree(snes->data);CHKERRQ(ierr);
312 
313   PetscFunctionReturn(0);
314 }
315 /* -------------------------------------------------------------------------- */
316 #undef __FUNCT__
317 #define __FUNCT__ "SNESLineSearchNo_LS"
318 
319 /*@C
320    SNESLineSearchNo - This routine is not a line search at all;
321    it simply uses the full Newton step.  Thus, this routine is intended
322    to serve as a template and is not recommended for general use.
323 
324    Logically Collective on SNES and Vec
325 
326    Input Parameters:
327 +  snes - nonlinear context
328 .  lsctx - optional context for line search (not used here)
329 .  x - current iterate
330 .  f - residual evaluated at x
331 .  y - search direction
332 .  fnorm - 2-norm of f
333 -  xnorm - norm of x if known, otherwise 0
334 
335    Output Parameters:
336 +  g - residual evaluated at new iterate y
337 .  w - new iterate
338 .  gnorm - 2-norm of g
339 .  ynorm - 2-norm of search length
340 -  flag - PETSC_TRUE on success, PETSC_FALSE on failure
341 
342    Options Database Key:
343 .  -snes_ls basic - Activates SNESLineSearchNo()
344 
345    Level: advanced
346 
347 .keywords: SNES, nonlinear, line search, cubic
348 
349 .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(),
350           SNESLineSearchSet(), SNESLineSearchNoNorms()
351 @*/
352 PetscErrorCode  SNESLineSearchNo_LS(SNES snes,void *lsctx,Vec x,Vec f,Vec y,PetscReal fnorm,PetscReal xnorm,Vec g,Vec w,PetscReal *ynorm,PetscReal *gnorm,PetscBool  *flag)
353 {
354   PetscErrorCode ierr;
355   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
356 
357   PetscFunctionBegin;
358   *flag = PETSC_TRUE;
359   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
360   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);         /* ynorm = || y || */
361   ierr = VecWAXPY(w,-snes->damping,y,x);CHKERRQ(ierr);            /* w <- x - y   */
362   if (snes->ops->postcheckstep) {
363    ierr = (*snes->ops->postcheckstep)(snes,x,y,w,snes->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
364   }
365   if (changed_y) {
366     ierr = VecWAXPY(w,-snes->damping,y,x);CHKERRQ(ierr);            /* w <- x - y   */
367   }
368   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
369   if (!snes->domainerror) {
370     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);  /* gnorm = || g || */
371     if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
372   }
373   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
374   PetscFunctionReturn(0);
375 }
376 /* -------------------------------------------------------------------------- */
377 #undef __FUNCT__
378 #define __FUNCT__ "SNESLineSearchNoNorms_LS"
379 
380 /*@C
381    SNESLineSearchNoNorms - This routine is not a line search at
382    all; it simply uses the full Newton step. This version does not
383    even compute the norm of the function or search direction; this
384    is intended only when you know the full step is fine and are
385    not checking for convergence of the nonlinear iteration (for
386    example, you are running always for a fixed number of Newton steps).
387 
388    Logically Collective on SNES and Vec
389 
390    Input Parameters:
391 +  snes - nonlinear context
392 .  lsctx - optional context for line search (not used here)
393 .  x - current iterate
394 .  f - residual evaluated at x
395 .  y - search direction
396 .  fnorm - 2-norm of f
397 -  xnorm - norm of x if known, otherwise 0
398 
399    Output Parameters:
400 +  g - residual evaluated at new iterate y
401 .  w - new iterate
402 .  gnorm - not changed
403 .  ynorm - not changed
404 -  flag - set to PETSC_TRUE indicating a successful line search
405 
406    Options Database Key:
407 .  -snes_ls basicnonorms - Activates SNESLineSearchNoNorms()
408 
409    Notes:
410    SNESLineSearchNoNorms() must be used in conjunction with
411    either the options
412 $     -snes_no_convergence_test -snes_max_it <its>
413    or alternatively a user-defined custom test set via
414    SNESSetConvergenceTest(); or a -snes_max_it of 1,
415    otherwise, the SNES solver will generate an error.
416 
417    During the final iteration this will not evaluate the function at
418    the solution point. This is to save a function evaluation while
419    using pseudo-timestepping.
420 
421    The residual norms printed by monitoring routines such as
422    SNESMonitorDefault() (as activated via -snes_monitor) will not be
423    correct, since they are not computed.
424 
425    Level: advanced
426 
427 .keywords: SNES, nonlinear, line search, cubic
428 
429 .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(),
430           SNESLineSearchSet(), SNESLineSearchNo()
431 @*/
432 PetscErrorCode  SNESLineSearchNoNorms_LS(SNES snes,void *lsctx,Vec x,Vec f,Vec y,PetscReal fnorm,PetscReal xnorm,Vec g,Vec w,PetscReal *ynorm,PetscReal *gnorm,PetscBool  *flag)
433 {
434   PetscErrorCode ierr;
435   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
436 
437   PetscFunctionBegin;
438   *flag = PETSC_TRUE;
439   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
440   ierr = VecWAXPY(w,-snes->damping,y,x);CHKERRQ(ierr);            /* w <- x - y      */
441   if (snes->ops->postcheckstep) {
442    ierr = (*snes->ops->postcheckstep)(snes,x,y,w,snes->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
443   }
444   if (changed_y) {
445     ierr = VecWAXPY(w,-snes->damping,y,x);CHKERRQ(ierr);            /* w <- x - y   */
446   }
447 
448   /* don't evaluate function the last time through */
449   if (snes->iter < snes->max_its-1) {
450     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
451   }
452   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
453   PetscFunctionReturn(0);
454 }
455 /* -------------------------------------------------------------------------- */
456 #undef __FUNCT__
457 #define __FUNCT__ "SNESLineSearchCubic_LS"
458 /*@C
459    SNESLineSearchCubic - Performs a cubic line search (default line search method).
460 
461    Collective on SNES
462 
463    Input Parameters:
464 +  snes - nonlinear context
465 .  lsctx - optional context for line search (not used here)
466 .  x - current iterate
467 .  f - residual evaluated at x
468 .  y - search direction
469 .  fnorm - 2-norm of f
470 -  xnorm - norm of x if known, otherwise 0
471 
472    Output Parameters:
473 +  g - residual evaluated at new iterate y
474 .  w - new iterate
475 .  gnorm - 2-norm of g
476 .  ynorm - 2-norm of search length
477 -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
478 
479    Options Database Key:
480 +  -snes_ls cubic - Activates SNESLineSearchCubic()
481 .   -snes_ls_alpha <alpha> - Sets alpha
482 .   -snes_ls_maxstep <maxstep> - Sets the maximum stepsize the line search will use (if the 2-norm(y) > maxstep then scale y to be y = (maxstep/2-norm(y)) *y)
483 -   -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda/ max_i ( y[i]/x[i] )
484 
485 
486    Notes:
487    This line search is taken from "Numerical Methods for Unconstrained
488    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
489 
490    Level: advanced
491 
492 .keywords: SNES, nonlinear, line search, cubic
493 
494 .seealso: SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms()
495 @*/
496 PetscErrorCode  SNESLineSearchCubic_LS(SNES snes,void *lsctx,Vec x,Vec f,Vec y,PetscReal fnorm,PetscReal xnorm,Vec g,Vec w,PetscReal *ynorm,PetscReal *gnorm,PetscBool  *flag)
497 {
498   /*
499      Note that for line search purposes we work with with the related
500      minimization problem:
501         min  z(x):  R^n -> R,
502      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
503    */
504 
505   PetscReal      initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
506   PetscReal      minlambda,lambda,lambdatemp;
507 #if defined(PETSC_USE_COMPLEX)
508   PetscScalar    cinitslope;
509 #endif
510   PetscErrorCode ierr;
511   PetscInt       count;
512   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
513   MPI_Comm       comm;
514 
515   PetscFunctionBegin;
516   ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
517   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
518   *flag   = PETSC_TRUE;
519 
520   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
521   if (*ynorm == 0.0) {
522     if (snes->ls_monitor) {
523       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
524       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Initial direction and size is 0\n");CHKERRQ(ierr);
525       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
526     }
527     *gnorm = fnorm;
528     ierr   = VecCopy(x,w);CHKERRQ(ierr);
529     ierr   = VecCopy(f,g);CHKERRQ(ierr);
530     *flag  = PETSC_FALSE;
531     goto theend1;
532   }
533   if (*ynorm > snes->maxstep) {	/* Step too big, so scale back */
534     if (snes->ls_monitor) {
535       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
536       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Scaling step by %14.12e old ynorm %14.12e\n",(double)(snes->maxstep/(*ynorm)),(double)(*ynorm));CHKERRQ(ierr);
537       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
538     }
539     ierr = VecScale(y,snes->maxstep/(*ynorm));CHKERRQ(ierr);
540     *ynorm = snes->maxstep;
541   }
542   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
543   minlambda = snes->steptol/rellength;
544   ierr      = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
545 #if defined(PETSC_USE_COMPLEX)
546   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
547   initslope = PetscRealPart(cinitslope);
548 #else
549   ierr      = VecDot(f,w,&initslope);CHKERRQ(ierr);
550 #endif
551   if (initslope > 0.0)  initslope = -initslope;
552   if (initslope == 0.0) initslope = -1.0;
553 
554   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
555   if (snes->nfuncs >= snes->max_funcs) {
556     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
557     *flag = PETSC_FALSE;
558     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
559     goto theend1;
560   }
561   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
562   if (snes->domainerror) {
563     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
564     PetscFunctionReturn(0);
565   }
566   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
567   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
568   ierr = PetscInfo2(snes,"Initial fnorm %14.12e gnorm %14.12e\n",(double)fnorm,(double)*gnorm);CHKERRQ(ierr);
569   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + snes->ls_alpha*initslope) { /* Sufficient reduction */
570     if (snes->ls_monitor) {
571       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
572       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Using full step: fnorm %14.12e gnorm %14.12e\n",(double)fnorm,(double)*gnorm);CHKERRQ(ierr);
573       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
574     }
575     goto theend1;
576   }
577 
578   /* Fit points with quadratic */
579   lambda     = 1.0;
580   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
581   lambdaprev = lambda;
582   gnormprev  = *gnorm;
583   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
584   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
585   else                         lambda = lambdatemp;
586 
587   ierr  = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
588   if (snes->nfuncs >= snes->max_funcs) {
589     ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr);
590     *flag = PETSC_FALSE;
591     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
592     goto theend1;
593   }
594   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
595   if (snes->domainerror) {
596     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
597     PetscFunctionReturn(0);
598   }
599   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
600   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
601   if (snes->ls_monitor) {
602     ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
603     ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: gnorm after quadratic fit %14.12e\n",(double)*gnorm);CHKERRQ(ierr);
604     ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
605   }
606   if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*snes->ls_alpha*initslope) { /* sufficient reduction */
607     if (snes->ls_monitor) {
608       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
609       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr);
610       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
611     }
612     goto theend1;
613   }
614 
615   /* Fit points with cubic */
616   count = 1;
617   while (PETSC_TRUE) {
618     if (lambda <= minlambda) {
619       if (snes->ls_monitor) {
620         ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
621 	ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr);
622 	ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",(double)fnorm,(double)*gnorm,(double)*ynorm,(double)minlambda,(double)lambda,(double)initslope);CHKERRQ(ierr);
623         ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
624       }
625       *flag = PETSC_FALSE;
626       break;
627     }
628     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
629     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
630     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
631     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
632     d  = b*b - 3*a*initslope;
633     if (d < 0.0) d = 0.0;
634     if (a == 0.0) {
635       lambdatemp = -initslope/(2.0*b);
636     } else {
637       lambdatemp = (-b + sqrt(d))/(3.0*a);
638     }
639     lambdaprev = lambda;
640     gnormprev  = *gnorm;
641     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
642     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
643     else                         lambda     = lambdatemp;
644     ierr  = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
645     if (snes->nfuncs >= snes->max_funcs) {
646       ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
647       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);
648       *flag = PETSC_FALSE;
649       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
650       break;
651     }
652     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
653     if (snes->domainerror) {
654       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
655       PetscFunctionReturn(0);
656     }
657     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
658     if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
659     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*snes->ls_alpha*initslope) { /* is reduction enough? */
660       if (snes->ls_monitor) {
661 	ierr = PetscPrintf(comm,"    Line search: Cubically determined step, current gnorm %14.12e lambda=%18.16e\n",(double)*gnorm,(double)lambda);CHKERRQ(ierr);
662       }
663       break;
664     } else {
665       if (snes->ls_monitor) {
666         ierr = PetscPrintf(comm,"    Line search: Cubic step no good, shrinking lambda, current gnorm %12.12e lambda=%18.16e\n",(double)*gnorm,(double)lambda);CHKERRQ(ierr);
667       }
668     }
669     count++;
670   }
671   theend1:
672   /* Optional user-defined check for line search step validity */
673   if (snes->ops->postcheckstep && *flag) {
674     ierr = (*snes->ops->postcheckstep)(snes,x,y,w,snes->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
675     if (changed_y) {
676       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
677     }
678     if (changed_y || changed_w) { /* recompute the function if the step has changed */
679       ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
680       if (snes->domainerror) {
681         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
682         PetscFunctionReturn(0);
683       }
684       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
685       if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
686       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
687       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
688       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
689     }
690   }
691   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
692   PetscFunctionReturn(0);
693 }
694 /* -------------------------------------------------------------------------- */
695 #undef __FUNCT__
696 #define __FUNCT__ "SNESLineSearchQuadratic_LS"
697 /*@C
698    SNESLineSearchQuadratic_LS - Performs a quadratic line search.
699 
700    Collective on SNES and Vec
701 
702    Input Parameters:
703 +  snes - the SNES context
704 .  lsctx - optional context for line search (not used here)
705 .  x - current iterate
706 .  f - residual evaluated at x
707 .  y - search direction
708 .  fnorm - 2-norm of f
709 -  xnorm - norm of x if known, otherwise 0
710 
711    Output Parameters:
712 +  g - residual evaluated at new iterate w
713 .  w - new iterate (x + lambda*y)
714 .  gnorm - 2-norm of g
715 .  ynorm - 2-norm of search length
716 -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
717 
718    Options Database Keys:
719 +  -snes_ls quadratic - Activates SNESLineSearchQuadratic()
720 .   -snes_ls_alpha <alpha> - Sets alpha
721 .   -snes_ls_maxstep <maxstep> - Sets the maximum stepsize the line search will use (if the 2-norm(y) > maxstep then scale y to be y = (maxstep/2-norm(y)) *y)
722 -   -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda/ max_i ( y[i]/x[i] )
723 
724    Notes:
725    Use SNESLineSearchSet() to set this routine within the SNESLS method.
726 
727    Level: advanced
728 
729 .keywords: SNES, nonlinear, quadratic, line search
730 
731 .seealso: SNESLineSearchCubic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms()
732 @*/
733 PetscErrorCode  SNESLineSearchQuadratic_LS(SNES snes,void *lsctx,Vec x,Vec f,Vec y,PetscReal fnorm,PetscReal xnorm,Vec g,Vec w,PetscReal *ynorm,PetscReal *gnorm,PetscBool  *flag)
734 {
735   /*
736      Note that for line search purposes we work with with the related
737      minimization problem:
738         min  z(x):  R^n -> R,
739      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
740    */
741   PetscReal      initslope,minlambda,lambda,lambdatemp,rellength;
742 #if defined(PETSC_USE_COMPLEX)
743   PetscScalar    cinitslope;
744 #endif
745   PetscErrorCode ierr;
746   PetscInt       count;
747   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
748 
749   PetscFunctionBegin;
750   ierr    = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
751   *flag   = PETSC_TRUE;
752 
753   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
754   if (*ynorm == 0.0) {
755     if (snes->ls_monitor) {
756       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
757       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"Line search: Direction and size is 0\n");CHKERRQ(ierr);
758       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
759     }
760     *gnorm = fnorm;
761     ierr   = VecCopy(x,w);CHKERRQ(ierr);
762     ierr   = VecCopy(f,g);CHKERRQ(ierr);
763     *flag  = PETSC_FALSE;
764     goto theend2;
765   }
766   if (*ynorm > snes->maxstep) {	/* Step too big, so scale back */
767     ierr   = VecScale(y,snes->maxstep/(*ynorm));CHKERRQ(ierr);
768     *ynorm = snes->maxstep;
769   }
770   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
771   minlambda = snes->steptol/rellength;
772   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
773 #if defined(PETSC_USE_COMPLEX)
774   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
775   initslope = PetscRealPart(cinitslope);
776 #else
777   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
778 #endif
779   if (initslope > 0.0)  initslope = -initslope;
780   if (initslope == 0.0) initslope = -1.0;
781   ierr = PetscInfo1(snes,"Initslope %14.12e \n",(double)initslope);CHKERRQ(ierr);
782 
783   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
784   if (snes->nfuncs >= snes->max_funcs) {
785     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
786     *flag = PETSC_FALSE;
787     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
788     goto theend2;
789   }
790   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
791   if (snes->domainerror) {
792     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
793     PetscFunctionReturn(0);
794   }
795   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
796   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
797   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + snes->ls_alpha*initslope) { /* Sufficient reduction */
798     if (snes->ls_monitor) {
799       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
800       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line Search: Using full step\n");CHKERRQ(ierr);
801       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
802     }
803     goto theend2;
804   }
805 
806   /* Fit points with quadratic */
807   lambda = 1.0;
808   count = 1;
809   while (PETSC_TRUE) {
810     if (lambda <= minlambda) { /* bad luck; use full step */
811       if (snes->ls_monitor) {
812         ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
813         ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"Line search: Unable to find good step length! %D \n",count);CHKERRQ(ierr);
814         ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"Line search: fnorm=%14.12e, gnorm=%14.12e, ynorm=%14.12e, lambda=%14.12e, initial slope=%14.12e\n",(double)fnorm,(double)*gnorm,(double)*ynorm,(double)lambda,(double)initslope);CHKERRQ(ierr);
815         ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
816       }
817       ierr = VecCopy(x,w);CHKERRQ(ierr);
818       *flag = PETSC_FALSE;
819       break;
820     }
821     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
822     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
823     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
824     else                         lambda     = lambdatemp;
825 
826     ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
827     if (snes->nfuncs >= snes->max_funcs) {
828       ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
829       ierr  = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",(double)fnorm,(double)*gnorm,(double)*ynorm,(double)lambda,(double)initslope);CHKERRQ(ierr);
830       *flag = PETSC_FALSE;
831       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
832       break;
833     }
834     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
835     if (snes->domainerror) {
836       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
837       PetscFunctionReturn(0);
838     }
839     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
840     if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
841     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*snes->ls_alpha*initslope) { /* sufficient reduction */
842       if (snes->ls_monitor) {
843         ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
844         ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line Search: Quadratically determined step, lambda=%14.12e\n",(double)lambda);CHKERRQ(ierr);
845         ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
846       }
847       break;
848     }
849     count++;
850   }
851   theend2:
852   /* Optional user-defined check for line search step validity */
853   if (snes->ops->postcheckstep) {
854     ierr = (*snes->ops->postcheckstep)(snes,x,y,w,snes->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
855     if (changed_y) {
856       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
857     }
858     if (changed_y || changed_w) { /* recompute the function if the step has changed */
859       ierr = SNESComputeFunction(snes,w,g);
860       if (snes->domainerror) {
861         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
862         PetscFunctionReturn(0);
863       }
864       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
865       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
866       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
867       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
868       if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
869     }
870   }
871   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
872   PetscFunctionReturn(0);
873 }
874 /* -------------------------------------------------------------------------- */
875 
876 /*
877    SNESView_LS - Prints info from the SNESLS data structure.
878 
879    Input Parameters:
880 .  SNES - the SNES context
881 .  viewer - visualization context
882 
883    Application Interface Routine: SNESView()
884 */
885 #undef __FUNCT__
886 #define __FUNCT__ "SNESView_LS"
887 static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer)
888 {
889   const char     *cstr;
890   PetscErrorCode ierr;
891   PetscBool      iascii;
892 
893   PetscFunctionBegin;
894   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
895   if (iascii) {
896     if (snes->ops->linesearch == SNESLineSearchNo_LS)             cstr = "SNESLineSearchNo";
897     if (snes->ops->linesearch == SNESLineSearchNoNorms_LS)        cstr = "SNESLineSearchNoNorms";
898     else if (snes->ops->linesearch == SNESLineSearchQuadratic_LS) cstr = "SNESLineSearchQuadratic";
899     else if (snes->ops->linesearch == SNESLineSearchCubic_LS)     cstr = "SNESLineSearchCubic";
900     else                                                   cstr = "unknown";
901     ierr = PetscViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
902     ierr = PetscViewerASCIIPrintf(viewer,"  alpha=%14.12e, maxstep=%14.12e, minlambda=%14.12e\n",(double)snes->ls_alpha,(double)snes->maxstep,(double)snes->steptol);CHKERRQ(ierr);
903     ierr = PetscViewerASCIIPrintf(viewer,"  damping factor=%14.12e\n",(double)snes->damping);CHKERRQ(ierr);
904   }
905   PetscFunctionReturn(0);
906 }
907 
908 /* -------------------------------------------------------------------------- */
909 /*
910    SNESSetFromOptions_LS - Sets various parameters for the SNESLS method.
911 
912    Input Parameter:
913 .  snes - the SNES context
914 
915    Application Interface Routine: SNESSetFromOptions()
916 */
917 #undef __FUNCT__
918 #define __FUNCT__ "SNESSetFromOptions_LS"
919 static PetscErrorCode SNESSetFromOptions_LS(SNES snes)
920 {
921   PetscErrorCode ierr;
922 
923   PetscFunctionBegin;
924   ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr);
925   ierr = PetscOptionsTail();CHKERRQ(ierr);
926   PetscFunctionReturn(0);
927 }
928 
929 EXTERN_C_BEGIN
930 extern PetscErrorCode  SNESLineSearchSetParams_LS(SNES,PetscReal,PetscReal,PetscReal);
931 EXTERN_C_END
932 
933 /* -------------------------------------------------------------------------- */
934 /*MC
935       SNESLS - Newton based nonlinear solver that uses a line search
936 
937    Options Database:
938 +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
939 .   -snes_ls_alpha <alpha> - Sets alpha used in determining if reduction in function norm is sufficient
940 .   -snes_ls_maxstep <maxstep> - Sets the maximum stepsize the line search will use (if the 2-norm(y) > maxstep then scale y to be y = (maxstep/2-norm(y)) *y)
941 .   -snes_ls_minlambda <minlambda>  - Sets the minimum lambda the line search will use  minlambda / max_i ( y[i]/x[i] )
942 .   -snes_ls_monitor - print information about progress of line searches
943 -   -snes_ls_damping - damping factor used if -snes_ls is basic or basicnonorms
944 
945 
946     Notes: This is the default nonlinear solver in SNES
947 
948    Level: beginner
949 
950 .seealso:  SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(),
951            SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(),
952           SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams()
953 
954 M*/
955 EXTERN_C_BEGIN
956 #undef __FUNCT__
957 #define __FUNCT__ "SNESCreate_LS"
958 PetscErrorCode  SNESCreate_LS(SNES snes)
959 {
960   PetscErrorCode ierr;
961   SNES_LS        *neP;
962 
963   PetscFunctionBegin;
964   snes->ops->setup           = SNESSetUp_LS;
965   snes->ops->solve           = SNESSolve_LS;
966   snes->ops->destroy         = SNESDestroy_LS;
967   snes->ops->setfromoptions  = SNESSetFromOptions_LS;
968   snes->ops->view            = SNESView_LS;
969   snes->ops->reset           = SNESReset_LS;
970 
971   snes->usesksp                      = PETSC_TRUE;
972   snes->usespc                       = PETSC_FALSE;
973   ierr                               = PetscNewLog(snes,SNES_LS,&neP);CHKERRQ(ierr);
974   snes->data                         = (void*)neP;
975   snes->ops->linesearch              = SNESLineSearchCubic_LS;
976   snes->ops->linesearchno            = SNESLineSearchNo_LS;
977   snes->ops->linesearchquadratic     = SNESLineSearchQuadratic_LS;
978   snes->ops->linesearchcubic         = SNESLineSearchCubic_LS;
979   snes->ops->linesearchexact         = PETSC_NULL;
980   snes->ops->linesearchtest          = PETSC_NULL;
981   snes->lsP                          = PETSC_NULL;
982   snes->ops->postcheckstep           = PETSC_NULL;
983   snes->postcheck                    = PETSC_NULL;
984   snes->ops->precheckstep            = PETSC_NULL;
985   snes->precheck                     = PETSC_NULL;
986 
987   PetscFunctionReturn(0);
988 }
989 EXTERN_C_END
990