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