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