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