Actual source code: ls.c
1: #define PETSCSNES_DLL
3: #include ../src/snes/impls/ls/ls.h
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: */
13: PetscErrorCode SNESLSCheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscTruth *ismin)
14: {
15: PetscReal a1;
17: PetscTruth hastranspose;
20: *ismin = PETSC_FALSE;
21: MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);
22: if (hastranspose) {
23: /* Compute || J^T F|| */
24: MatMultTranspose(A,F,W);
25: VecNorm(W,NORM_2,&a1);
26: PetscInfo1(snes,"|| J^T F|| %G near zero implies found a local minimum\n",a1/fnorm);
27: if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
28: } else {
29: Vec work;
30: PetscScalar result;
31: PetscReal wnorm;
33: VecSetRandom(W,PETSC_NULL);
34: VecNorm(W,NORM_2,&wnorm);
35: VecDuplicate(W,&work);
36: MatMult(A,W,work);
37: VecDot(F,work,&result);
38: VecDestroy(work);
39: a1 = PetscAbsScalar(result)/(fnorm*wnorm);
40: PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %G near zero implies found a local minimum\n",a1);
41: if (a1 < 1.e-4) *ismin = PETSC_TRUE;
42: }
43: return(0);
44: }
46: /*
47: Checks if J^T(F - J*X) = 0
48: */
51: PetscErrorCode SNESLSCheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2)
52: {
53: PetscReal a1,a2;
55: PetscTruth hastranspose;
58: MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);
59: if (hastranspose) {
60: MatMult(A,X,W1);
61: VecAXPY(W1,-1.0,F);
63: /* Compute || J^T W|| */
64: MatMultTranspose(A,W1,W2);
65: VecNorm(W1,NORM_2,&a1);
66: VecNorm(W2,NORM_2,&a2);
67: if (a1) {
68: PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %G near zero implies inconsistent rhs\n",a2/a1);
69: }
70: }
71: return(0);
72: }
74: /* --------------------------------------------------------------------
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.
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.
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.
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().
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.
110: -------------------------------------------------------------------- */
111: /*
112: SNESSolve_LS - Solves a nonlinear system with a truncated Newton
113: method with a line search.
115: Input Parameters:
116: . snes - the SNES context
118: Output Parameter:
119: . outits - number of iterations until termination
121: Application Interface Routine: SNESSolve()
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: */
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;
144: snes->numFailures = 0;
145: snes->numLinearSolveFailures = 0;
146: snes->reason = SNES_CONVERGED_ITERATING;
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];
155: PetscObjectTakeAccess(snes);
156: snes->iter = 0;
157: snes->norm = 0;
158: PetscObjectGrantAccess(snes);
159: SNESComputeFunction(snes,X,F);
160: if (snes->domainerror) {
161: snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
162: return(0);
163: }
164: VecNormBegin(F,NORM_2,&fnorm); /* fnorm <- ||F|| */
165: VecNormBegin(X,NORM_2,&xnorm); /* xnorm <- ||x|| */
166: VecNormEnd(F,NORM_2,&fnorm);
167: VecNormEnd(X,NORM_2,&xnorm);
168: if PetscIsInfOrNanReal(fnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
169: PetscObjectTakeAccess(snes);
170: snes->norm = fnorm;
171: PetscObjectGrantAccess(snes);
172: SNESLogConvHistory(snes,fnorm,0);
173: SNESMonitor(snes,0,fnorm);
175: /* set parameter for default relative tolerance convergence test */
176: snes->ttol = fnorm*snes->rtol;
177: /* test convergence */
178: (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
179: if (snes->reason) return(0);
181: for (i=0; i<maxits; i++) {
183: /* Call general purpose update function */
184: if (snes->ops->update) {
185: (*snes->ops->update)(snes, snes->iter);
186: }
188: /* Solve J Y = F, where J is Jacobian matrix */
189: SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);
190: KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);
191: SNES_KSPSolve(snes,snes->ksp,F,Y);
192: KSPGetConvergedReason(snes->ksp,&kspreason);
193: if (kspreason < 0) {
194: if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
195: PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);
196: snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
197: break;
198: }
199: }
200: KSPGetIterationNumber(snes->ksp,&lits);
201: snes->linear_its += lits;
202: PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);
204: if (neP->precheckstep) {
205: PetscTruth changed_y = PETSC_FALSE;
206: (*neP->precheckstep)(snes,X,Y,neP->precheck,&changed_y);
207: }
209: if (PetscLogPrintInfo){
210: SNESLSCheckResidual_Private(snes,snes->jacobian,F,Y,G,W);
211: }
213: /* Compute a (scaled) negative update in the line search routine:
214: Y <- X - lambda*Y
215: and evaluate G = function(Y) (depends on the line search).
216: */
217: VecCopy(Y,snes->vec_sol_update);
218: ynorm = 1; gnorm = fnorm;
219: (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,xnorm,&ynorm,&gnorm,&lssucceed);
220: PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",fnorm,gnorm,ynorm,(int)lssucceed);
221: if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
222: if (snes->domainerror) {
223: snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
224: return(0);
225: }
226: if (!lssucceed) {
227: if (++snes->numFailures >= snes->maxFailures) {
228: PetscTruth ismin;
229: snes->reason = SNES_DIVERGED_LS_FAILURE;
230: SNESLSCheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);
231: if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
232: break;
233: }
234: }
235: /* Update function and solution vectors */
236: fnorm = gnorm;
237: VecCopy(G,F);
238: VecCopy(W,X);
239: /* Monitor convergence */
240: PetscObjectTakeAccess(snes);
241: snes->iter = i+1;
242: snes->norm = fnorm;
243: PetscObjectGrantAccess(snes);
244: SNESLogConvHistory(snes,snes->norm,lits);
245: SNESMonitor(snes,snes->iter,snes->norm);
246: /* Test for convergence, xnorm = || X || */
247: if (snes->ops->converged != SNESSkipConverged) { VecNorm(X,NORM_2,&xnorm); }
248: (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);
249: if (snes->reason) break;
250: }
251: if (i == maxits) {
252: PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);
253: if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
254: }
255: return(0);
256: }
257: /* -------------------------------------------------------------------------- */
258: /*
259: SNESSetUp_LS - Sets up the internal data structures for the later use
260: of the SNESLS nonlinear solver.
262: Input Parameter:
263: . snes - the SNES context
264: . x - the solution vector
266: Application Interface Routine: SNESSetUp()
268: Notes:
269: For basic use of the SNES solvers, the user need not explicitly call
270: SNESSetUp(), since these actions will automatically occur during
271: the call to SNESSolve().
272: */
275: PetscErrorCode SNESSetUp_LS(SNES snes)
276: {
280: if (!snes->vec_sol_update) {
281: VecDuplicate(snes->vec_sol,&snes->vec_sol_update);
282: PetscLogObjectParent(snes,snes->vec_sol_update);
283: }
284: if (!snes->work) {
285: snes->nwork = 3;
286: VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);
287: PetscLogObjectParents(snes,snes->nwork,snes->work);
288: }
289: return(0);
290: }
291: /* -------------------------------------------------------------------------- */
292: /*
293: SNESDestroy_LS - Destroys the private SNES_LS context that was created
294: with SNESCreate_LS().
296: Input Parameter:
297: . snes - the SNES context
299: Application Interface Routine: SNESDestroy()
300: */
303: PetscErrorCode SNESDestroy_LS(SNES snes)
304: {
308: if (snes->vec_sol_update) {
309: VecDestroy(snes->vec_sol_update);
310: snes->vec_sol_update = PETSC_NULL;
311: }
312: if (snes->nwork) {
313: VecDestroyVecs(snes->work,snes->nwork);
314: snes->nwork = 0;
315: snes->work = PETSC_NULL;
316: }
317: PetscFree(snes->data);
319: /* clear composed functions */
320: PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);
321: PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","",PETSC_NULL);
322: PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C","",PETSC_NULL);
324: return(0);
325: }
326: /* -------------------------------------------------------------------------- */
330: /*@C
331: SNESLineSearchNo - This routine is not a line search at all;
332: it simply uses the full Newton step. Thus, this routine is intended
333: to serve as a template and is not recommended for general use.
335: Collective on SNES and Vec
337: Input Parameters:
338: + snes - nonlinear context
339: . lsctx - optional context for line search (not used here)
340: . x - current iterate
341: . f - residual evaluated at x
342: . y - search direction
343: . w - work vector
344: . fnorm - 2-norm of f
345: - xnorm - norm of x if known, otherwise 0
347: Output Parameters:
348: + g - residual evaluated at new iterate y
349: . w - new iterate
350: . gnorm - 2-norm of g
351: . ynorm - 2-norm of search length
352: - flag - PETSC_TRUE on success, PETSC_FALSE on failure
354: Options Database Key:
355: . -snes_ls basic - Activates SNESLineSearchNo()
357: Level: advanced
359: .keywords: SNES, nonlinear, line search, cubic
361: .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(),
362: SNESLineSearchSet(), SNESLineSearchNoNorms()
363: @*/
364: PetscErrorCode SNESLineSearchNo(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
365: {
367: SNES_LS *neP = (SNES_LS*)snes->data;
368: PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
371: *flag = PETSC_TRUE;
372: PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);
373: VecNorm(y,NORM_2,ynorm); /* ynorm = || y || */
374: VecWAXPY(w,-1.0,y,x); /* w <- x - y */
375: if (neP->postcheckstep) {
376: (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);
377: }
378: if (changed_y) {
379: VecWAXPY(w,-1.0,y,x); /* w <- x - y */
380: }
381: SNESComputeFunction(snes,w,g);
382: if (!snes->domainerror) {
383: VecNorm(g,NORM_2,gnorm); /* gnorm = || g || */
384: if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
385: }
386: PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
387: return(0);
388: }
389: /* -------------------------------------------------------------------------- */
393: /*@C
394: SNESLineSearchNoNorms - This routine is not a line search at
395: all; it simply uses the full Newton step. This version does not
396: even compute the norm of the function or search direction; this
397: is intended only when you know the full step is fine and are
398: not checking for convergence of the nonlinear iteration (for
399: example, you are running always for a fixed number of Newton steps).
401: Collective on SNES and Vec
403: Input Parameters:
404: + snes - nonlinear context
405: . lsctx - optional context for line search (not used here)
406: . x - current iterate
407: . f - residual evaluated at x
408: . y - search direction
409: . w - work vector
410: . fnorm - 2-norm of f
411: - xnorm - norm of x if known, otherwise 0
413: Output Parameters:
414: + g - residual evaluated at new iterate y
415: . w - new iterate
416: . gnorm - not changed
417: . ynorm - not changed
418: - flag - set to PETSC_TRUE indicating a successful line search
420: Options Database Key:
421: . -snes_ls basicnonorms - Activates SNESLineSearchNoNorms()
423: Notes:
424: SNESLineSearchNoNorms() must be used in conjunction with
425: either the options
426: $ -snes_no_convergence_test -snes_max_it <its>
427: or alternatively a user-defined custom test set via
428: SNESSetConvergenceTest(); or a -snes_max_it of 1,
429: otherwise, the SNES solver will generate an error.
431: During the final iteration this will not evaluate the function at
432: the solution point. This is to save a function evaluation while
433: using pseudo-timestepping.
435: The residual norms printed by monitoring routines such as
436: SNESMonitorDefault() (as activated via -snes_monitor) will not be
437: correct, since they are not computed.
439: Level: advanced
441: .keywords: SNES, nonlinear, line search, cubic
443: .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(),
444: SNESLineSearchSet(), SNESLineSearchNo()
445: @*/
446: PetscErrorCode SNESLineSearchNoNorms(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
447: {
449: SNES_LS *neP = (SNES_LS*)snes->data;
450: PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
453: *flag = PETSC_TRUE;
454: PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);
455: VecWAXPY(w,-1.0,y,x); /* w <- x - y */
456: if (neP->postcheckstep) {
457: (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);
458: }
459: if (changed_y) {
460: VecWAXPY(w,-1.0,y,x); /* w <- x - y */
461: }
462:
463: /* don't evaluate function the last time through */
464: if (snes->iter < snes->max_its-1) {
465: SNESComputeFunction(snes,w,g);
466: }
467: PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
468: return(0);
469: }
470: /* -------------------------------------------------------------------------- */
473: /*@C
474: SNESLineSearchCubic - Performs a cubic line search (default line search method).
476: Collective on SNES
478: Input Parameters:
479: + snes - nonlinear context
480: . lsctx - optional context for line search (not used here)
481: . x - current iterate
482: . f - residual evaluated at x
483: . y - search direction
484: . w - work vector
485: . fnorm - 2-norm of f
486: - xnorm - norm of x if known, otherwise 0
488: Output Parameters:
489: + g - residual evaluated at new iterate y
490: . w - new iterate
491: . gnorm - 2-norm of g
492: . ynorm - 2-norm of search length
493: - flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
495: Options Database Key:
496: + -snes_ls cubic - Activates SNESLineSearchCubic()
497: . -snes_ls_alpha <alpha> - Sets alpha
498: . -snes_ls_maxstep <maxstep> - Sets the maximum stepsize the line search will use (if the 2-norm(y) > maxstep then scale y to be y = (maxstep/2-norm(y)) *y)
499: - -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda/ max_i ( y[i]/x[i] )
501:
502: Notes:
503: This line search is taken from "Numerical Methods for Unconstrained
504: Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
506: Level: advanced
508: .keywords: SNES, nonlinear, line search, cubic
510: .seealso: SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms()
511: @*/
512: PetscErrorCode SNESLineSearchCubic(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
513: {
514: /*
515: Note that for line search purposes we work with with the related
516: minimization problem:
517: min z(x): R^n -> R,
518: where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
519: */
520:
521: PetscReal initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
522: PetscReal minlambda,lambda,lambdatemp;
523: #if defined(PETSC_USE_COMPLEX)
524: PetscScalar cinitslope;
525: #endif
527: PetscInt count;
528: SNES_LS *neP = (SNES_LS*)snes->data;
529: PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
532: PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);
533: *flag = PETSC_TRUE;
535: VecNorm(y,NORM_2,ynorm);
536: if (!*ynorm) {
537: PetscInfo(snes,"Search direction and size is 0\n");
538: *gnorm = fnorm;
539: VecCopy(x,w);
540: VecCopy(f,g);
541: *flag = PETSC_FALSE;
542: goto theend1;
543: }
544: if (*ynorm > neP->maxstep) { /* Step too big, so scale back */
545: PetscInfo2(snes,"Scaling step by %G old ynorm %G\n",neP->maxstep/(*ynorm),*ynorm);
546: VecScale(y,neP->maxstep/(*ynorm));
547: *ynorm = neP->maxstep;
548: }
549: VecMaxPointwiseDivide(y,x,&rellength);
550: minlambda = neP->minlambda/rellength;
551: MatMult(snes->jacobian,y,w);
552: #if defined(PETSC_USE_COMPLEX)
553: VecDot(f,w,&cinitslope);
554: initslope = PetscRealPart(cinitslope);
555: #else
556: VecDot(f,w,&initslope);
557: #endif
558: if (initslope > 0.0) initslope = -initslope;
559: if (initslope == 0.0) initslope = -1.0;
561: VecWAXPY(w,-1.0,y,x);
562: if (snes->nfuncs >= snes->max_funcs) {
563: PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");
564: *flag = PETSC_FALSE;
565: snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
566: goto theend1;
567: }
568: SNESComputeFunction(snes,w,g);
569: if (snes->domainerror) {
570: PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
571: return(0);
572: }
573: VecNorm(g,NORM_2,gnorm);
574: if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
575: PetscInfo2(snes,"Initial fnorm %G gnorm %G\n",fnorm,*gnorm);
576: if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + neP->alpha*initslope) { /* Sufficient reduction */
577: PetscInfo2(snes,"Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);
578: goto theend1;
579: }
581: /* Fit points with quadratic */
582: lambda = 1.0;
583: lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
584: lambdaprev = lambda;
585: gnormprev = *gnorm;
586: if (lambdatemp > .5*lambda) lambdatemp = .5*lambda;
587: if (lambdatemp <= .1*lambda) lambda = .1*lambda;
588: else lambda = lambdatemp;
590: VecWAXPY(w,-lambda,y,x);
591: if (snes->nfuncs >= snes->max_funcs) {
592: PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);
593: *flag = PETSC_FALSE;
594: snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
595: goto theend1;
596: }
597: SNESComputeFunction(snes,w,g);
598: if (snes->domainerror) {
599: PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
600: return(0);
601: }
602: VecNorm(g,NORM_2,gnorm);
603: if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
604: PetscInfo1(snes,"gnorm after quadratic fit %G\n",*gnorm);
605: if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*neP->alpha*initslope) { /* sufficient reduction */
606: PetscInfo1(snes,"Quadratically determined step, lambda=%18.16e\n",lambda);
607: goto theend1;
608: }
610: /* Fit points with cubic */
611: count = 1;
612: while (PETSC_TRUE) {
613: if (lambda <= minlambda) {
614: PetscInfo1(snes,"Unable to find good step length! After %D tries \n",count);
615: PetscInfo6(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,minlambda,lambda,initslope);
616: *flag = PETSC_FALSE;
617: break;
618: }
619: t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
620: t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope;
621: a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
622: b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
623: d = b*b - 3*a*initslope;
624: if (d < 0.0) d = 0.0;
625: if (a == 0.0) {
626: lambdatemp = -initslope/(2.0*b);
627: } else {
628: lambdatemp = (-b + sqrt(d))/(3.0*a);
629: }
630: lambdaprev = lambda;
631: gnormprev = *gnorm;
632: if (lambdatemp > .5*lambda) lambdatemp = .5*lambda;
633: if (lambdatemp <= .1*lambda) lambda = .1*lambda;
634: else lambda = lambdatemp;
635: VecWAXPY(w,-lambda,y,x);
636: if (snes->nfuncs >= snes->max_funcs) {
637: PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);
638: PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope);
639: *flag = PETSC_FALSE;
640: snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
641: break;
642: }
643: SNESComputeFunction(snes,w,g);
644: if (snes->domainerror) {
645: PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
646: return(0);
647: }
648: VecNorm(g,NORM_2,gnorm);
649: if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
650: if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*neP->alpha*initslope) { /* is reduction enough? */
651: PetscInfo2(snes,"Cubically determined step, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);
652: break;
653: } else {
654: PetscInfo2(snes,"Cubic step no good, shrinking lambda, current gnorem %G lambda=%18.16e\n",*gnorm,lambda);
655: }
656: count++;
657: }
658: theend1:
659: /* Optional user-defined check for line search step validity */
660: if (neP->postcheckstep && *flag) {
661: (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);
662: if (changed_y) {
663: VecWAXPY(w,-1.0,y,x);
664: }
665: if (changed_y || changed_w) { /* recompute the function if the step has changed */
666: SNESComputeFunction(snes,w,g);
667: if (snes->domainerror) {
668: PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
669: return(0);
670: }
671: VecNormBegin(g,NORM_2,gnorm);
672: if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
673: VecNormBegin(y,NORM_2,ynorm);
674: VecNormEnd(g,NORM_2,gnorm);
675: VecNormEnd(y,NORM_2,ynorm);
676: }
677: }
678: PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
679: return(0);
680: }
681: /* -------------------------------------------------------------------------- */
684: /*@C
685: SNESLineSearchQuadratic - Performs a quadratic line search.
687: Collective on SNES and Vec
689: Input Parameters:
690: + snes - the SNES context
691: . lsctx - optional context for line search (not used here)
692: . x - current iterate
693: . f - residual evaluated at x
694: . y - search direction
695: . w - work vector
696: . fnorm - 2-norm of f
697: - xnorm - norm of x if known, otherwise 0
699: Output Parameters:
700: + g - residual evaluated at new iterate w
701: . w - new iterate (x + lambda*y)
702: . gnorm - 2-norm of g
703: . ynorm - 2-norm of search length
704: - flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
706: Options Database Keys:
707: + -snes_ls quadratic - Activates SNESLineSearchQuadratic()
708: . -snes_ls_alpha <alpha> - Sets alpha
709: . -snes_ls_maxstep <maxstep> - Sets the maximum stepsize the line search will use (if the 2-norm(y) > maxstep then scale y to be y = (maxstep/2-norm(y)) *y)
710: - -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda/ max_i ( y[i]/x[i] )
712: Notes:
713: Use SNESLineSearchSet() to set this routine within the SNESLS method.
715: Level: advanced
717: .keywords: SNES, nonlinear, quadratic, line search
719: .seealso: SNESLineSearchCubic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms()
720: @*/
721: PetscErrorCode SNESLineSearchQuadratic(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
722: {
723: /*
724: Note that for line search purposes we work with with the related
725: minimization problem:
726: min z(x): R^n -> R,
727: where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
728: */
729: PetscReal initslope,minlambda,lambda,lambdatemp,rellength;
730: #if defined(PETSC_USE_COMPLEX)
731: PetscScalar cinitslope;
732: #endif
734: PetscInt count;
735: SNES_LS *neP = (SNES_LS*)snes->data;
736: PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
739: PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);
740: *flag = PETSC_TRUE;
742: VecNorm(y,NORM_2,ynorm);
743: if (*ynorm == 0.0) {
744: PetscInfo(snes,"Search direction and size is 0\n");
745: *gnorm = fnorm;
746: VecCopy(x,w);
747: VecCopy(f,g);
748: *flag = PETSC_FALSE;
749: goto theend2;
750: }
751: if (*ynorm > neP->maxstep) { /* Step too big, so scale back */
752: VecScale(y,neP->maxstep/(*ynorm));
753: *ynorm = neP->maxstep;
754: }
755: VecMaxPointwiseDivide(y,x,&rellength);
756: minlambda = neP->minlambda/rellength;
757: MatMult(snes->jacobian,y,w);
758: #if defined(PETSC_USE_COMPLEX)
759: VecDot(f,w,&cinitslope);
760: initslope = PetscRealPart(cinitslope);
761: #else
762: VecDot(f,w,&initslope);
763: #endif
764: if (initslope > 0.0) initslope = -initslope;
765: if (initslope == 0.0) initslope = -1.0;
766: PetscInfo1(snes,"Initslope %G \n",initslope);
768: VecWAXPY(w,-1.0,y,x);
769: if (snes->nfuncs >= snes->max_funcs) {
770: PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");
771: *flag = PETSC_FALSE;
772: snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
773: goto theend2;
774: }
775: SNESComputeFunction(snes,w,g);
776: if (snes->domainerror) {
777: PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
778: return(0);
779: }
780: VecNorm(g,NORM_2,gnorm);
781: if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
782: if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + neP->alpha*initslope) { /* Sufficient reduction */
783: PetscInfo(snes,"Using full step\n");
784: goto theend2;
785: }
787: /* Fit points with quadratic */
788: lambda = 1.0;
789: count = 1;
790: while (PETSC_TRUE) {
791: if (lambda <= minlambda) { /* bad luck; use full step */
792: PetscInfo1(snes,"Unable to find good step length! %D \n",count);
793: PetscInfo5(snes,"fnorm=%G, gnorm=%G, ynorm=%G, lambda=%G, initial slope=%G\n",fnorm,*gnorm,*ynorm,lambda,initslope);
794: VecCopy(x,w);
795: *flag = PETSC_FALSE;
796: break;
797: }
798: lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
799: if (lambdatemp > .5*lambda) lambdatemp = .5*lambda;
800: if (lambdatemp <= .1*lambda) lambda = .1*lambda;
801: else lambda = lambdatemp;
802:
803: VecWAXPY(w,-lambda,y,x);
804: if (snes->nfuncs >= snes->max_funcs) {
805: PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);
806: PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope);
807: *flag = PETSC_FALSE;
808: snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
809: break;
810: }
811: SNESComputeFunction(snes,w,g);
812: if (snes->domainerror) {
813: PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
814: return(0);
815: }
816: VecNorm(g,NORM_2,gnorm);
817: if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
818: if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*neP->alpha*initslope) { /* sufficient reduction */
819: PetscInfo1(snes,"Quadratically determined step, lambda=%G\n",lambda);
820: break;
821: }
822: count++;
823: }
824: theend2:
825: /* Optional user-defined check for line search step validity */
826: if (neP->postcheckstep) {
827: (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);
828: if (changed_y) {
829: VecWAXPY(w,-1.0,y,x);
830: }
831: if (changed_y || changed_w) { /* recompute the function if the step has changed */
832: SNESComputeFunction(snes,w,g);
833: if (snes->domainerror) {
834: PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
835: return(0);
836: }
837: VecNormBegin(g,NORM_2,gnorm);
838: VecNormBegin(y,NORM_2,ynorm);
839: VecNormEnd(g,NORM_2,gnorm);
840: VecNormEnd(y,NORM_2,ynorm);
841: if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
842: }
843: }
844: PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
845: return(0);
846: }
848: /* -------------------------------------------------------------------------- */
851: /*@C
852: SNESLineSearchSet - Sets the line search routine to be used
853: by the method SNESLS.
855: Input Parameters:
856: + snes - nonlinear context obtained from SNESCreate()
857: . lsctx - optional user-defined context for use by line search
858: - func - pointer to int function
860: Collective on SNES
862: Available Routines:
863: + SNESLineSearchCubic() - default line search
864: . SNESLineSearchQuadratic() - quadratic line search
865: . SNESLineSearchNo() - the full Newton step (actually not a line search)
866: - SNESLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel)
868: Options Database Keys:
869: + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
870: . -snes_ls_alpha <alpha> - Sets alpha
871: . -snes_ls_maxstep <maxstep> - Sets maximum step the line search will use (if the 2-norm(y) > maxstep then scale y to be y = (maxstep/2-norm(y)) *y)
872: - -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda / max_i ( y[i]/x[i] )
874: Calling sequence of func:
875: .vb
876: func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
877: .ve
879: Input parameters for func:
880: + snes - nonlinear context
881: . lsctx - optional user-defined context for line search
882: . x - current iterate
883: . f - residual evaluated at x
884: . y - search direction
885: . w - work vector
886: - fnorm - 2-norm of f
888: Output parameters for func:
889: + g - residual evaluated at new iterate y
890: . w - new iterate
891: . gnorm - 2-norm of g
892: . ynorm - 2-norm of search length
893: - flag - set to PETSC_TRUE if the line search succeeds; PETSC_FALSE on failure.
895: Level: advanced
897: .keywords: SNES, nonlinear, set, line search, routine
899: .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchNoNorms(),
900: SNESLineSearchSetPostCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams(), SNESLineSearchSetPreCheck()
901: @*/
902: PetscErrorCode SNESLineSearchSet(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void *lsctx)
903: {
904: PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void*);
907: PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSet_C",(void (**)(void))&f);
908: if (f) {
909: (*f)(snes,func,lsctx);
910: }
911: return(0);
912: }
915: /* -------------------------------------------------------------------------- */
919: PetscErrorCode SNESLineSearchSet_LS(SNES snes,FCN2 func,void *lsctx)
920: {
922: ((SNES_LS *)(snes->data))->LineSearch = func;
923: ((SNES_LS *)(snes->data))->lsP = lsctx;
924: return(0);
925: }
927: /* -------------------------------------------------------------------------- */
930: /*@C
931: SNESLineSearchSetPostCheck - Sets a routine to check the validity of new iterate computed
932: by the line search routine in the Newton-based method SNESLS.
934: Input Parameters:
935: + snes - nonlinear context obtained from SNESCreate()
936: . func - pointer to function
937: - checkctx - optional user-defined context for use by step checking routine
939: Collective on SNES
941: Calling sequence of func:
942: .vb
943: int func (SNES snes, Vec x,Vec y,Vec w,void *checkctx, PetscTruth *changed_y,PetscTruth *changed_w)
944: .ve
945: where func returns an error code of 0 on success and a nonzero
946: on failure.
948: Input parameters for func:
949: + snes - nonlinear context
950: . checkctx - optional user-defined context for use by step checking routine
951: . x - previous iterate
952: . y - new search direction and length
953: - w - current candidate iterate
955: Output parameters for func:
956: + y - search direction (possibly changed)
957: . w - current iterate (possibly modified)
958: . changed_y - indicates search direction was changed by this routine
959: - changed_w - indicates current iterate was changed by this routine
961: Level: advanced
963: Notes: All line searches accept the new iterate computed by the line search checking routine.
965: Only one of changed_y and changed_w can be PETSC_TRUE
967: On input w = x + y
969: SNESLineSearchNo() and SNESLineSearchNoNorms() (1) compute a candidate iterate u_{i+1}, (2) pass control
970: to the checking routine, and then (3) compute the corresponding nonlinear
971: function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}.
973: SNESLineSearchQuadratic() and SNESLineSearchCubic() (1) compute a candidate iterate u_{i+1} as well as a
974: candidate nonlinear function f(u_{i+1}), (2) pass control to the checking
975: routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes
976: were made to the candidate iterate in the checking routine (as indicated
977: by flag=PETSC_TRUE). The overhead of this extra function re-evaluation can be
978: very costly, so use this feature with caution!
980: .keywords: SNES, nonlinear, set, line search check, step check, routine
982: .seealso: SNESLineSearchSet(), SNESLineSearchSetPreCheck()
983: @*/
984: PetscErrorCode SNESLineSearchSetPostCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void *checkctx)
985: {
986: PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void*);
989: PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPostCheck_C",(void (**)(void))&f);
990: if (f) {
991: (*f)(snes,func,checkctx);
992: }
993: return(0);
994: }
997: /*@C
998: SNESLineSearchSetPreCheck - Sets a routine to check the validity of a new direction given by the linear solve
999: before the line search is called.
1001: Input Parameters:
1002: + snes - nonlinear context obtained from SNESCreate()
1003: . func - pointer to function
1004: - checkctx - optional user-defined context for use by step checking routine
1006: Collective on SNES
1008: Calling sequence of func:
1009: .vb
1010: int func (SNES snes, Vec x,Vec y,void *checkctx, PetscTruth *changed_y)
1011: .ve
1012: where func returns an error code of 0 on success and a nonzero
1013: on failure.
1015: Input parameters for func:
1016: + snes - nonlinear context
1017: . checkctx - optional user-defined context for use by step checking routine
1018: . x - previous iterate
1019: - y - new search direction and length
1021: Output parameters for func:
1022: + y - search direction (possibly changed)
1023: - changed_y - indicates search direction was changed by this routine
1025: Level: advanced
1027: Notes: All line searches accept the new iterate computed by the line search checking routine.
1029: .keywords: SNES, nonlinear, set, line search check, step check, routine
1031: .seealso: SNESLineSearchSet(), SNESLineSearchSetPostCheck(), SNESSetUpdate()
1032: @*/
1033: PetscErrorCode SNESLineSearchSetPreCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,void*,PetscTruth*),void *checkctx)
1034: {
1035: PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,void*,PetscTruth*),void*);
1038: PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPreCheck_C",(void (**)(void))&f);
1039: if (f) {
1040: (*f)(snes,func,checkctx);
1041: }
1042: return(0);
1043: }
1045: /* -------------------------------------------------------------------------- */
1051: PetscErrorCode SNESLineSearchSetPostCheck_LS(SNES snes,FCN1 func,void *checkctx)
1052: {
1054: ((SNES_LS *)(snes->data))->postcheckstep = func;
1055: ((SNES_LS *)(snes->data))->postcheck = checkctx;
1056: return(0);
1057: }
1063: PetscErrorCode SNESLineSearchSetPreCheck_LS(SNES snes,FCN3 func,void *checkctx)
1064: {
1066: ((SNES_LS *)(snes->data))->precheckstep = func;
1067: ((SNES_LS *)(snes->data))->precheck = checkctx;
1068: return(0);
1069: }
1072: /*
1073: SNESView_LS - Prints info from the SNESLS data structure.
1075: Input Parameters:
1076: . SNES - the SNES context
1077: . viewer - visualization context
1079: Application Interface Routine: SNESView()
1080: */
1083: static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer)
1084: {
1085: SNES_LS *ls = (SNES_LS *)snes->data;
1086: const char *cstr;
1088: PetscTruth iascii;
1091: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
1092: if (iascii) {
1093: if (ls->LineSearch == SNESLineSearchNo) cstr = "SNESLineSearchNo";
1094: else if (ls->LineSearch == SNESLineSearchQuadratic) cstr = "SNESLineSearchQuadratic";
1095: else if (ls->LineSearch == SNESLineSearchCubic) cstr = "SNESLineSearchCubic";
1096: else cstr = "unknown";
1097: PetscViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);
1098: PetscViewerASCIIPrintf(viewer," alpha=%G, maxstep=%G, minlambda=%G\n",ls->alpha,ls->maxstep,ls->minlambda);
1099: } else {
1100: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name);
1101: }
1102: return(0);
1103: }
1104: /* -------------------------------------------------------------------------- */
1105: /*
1106: SNESSetFromOptions_LS - Sets various parameters for the SNESLS method.
1108: Input Parameter:
1109: . snes - the SNES context
1111: Application Interface Routine: SNESSetFromOptions()
1112: */
1115: static PetscErrorCode SNESSetFromOptions_LS(SNES snes)
1116: {
1117: SNES_LS *ls = (SNES_LS *)snes->data;
1118: const char *lses[] = {"basic","basicnonorms","quadratic","cubic"};
1120: PetscInt indx;
1121: PetscTruth flg;
1124: PetscOptionsHead("SNES Line search options");
1125: PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);
1126: PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);
1127: PetscOptionsReal("-snes_ls_minlambda","Minimum lambda allowed","None",ls->minlambda,&ls->minlambda,0);
1129: PetscOptionsEList("-snes_ls","Line search used","SNESLineSearchSet",lses,4,"cubic",&indx,&flg);
1130: if (flg) {
1131: switch (indx) {
1132: case 0:
1133: SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);
1134: break;
1135: case 1:
1136: SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);
1137: break;
1138: case 2:
1139: SNESLineSearchSet(snes,SNESLineSearchQuadratic,PETSC_NULL);
1140: break;
1141: case 3:
1142: SNESLineSearchSet(snes,SNESLineSearchCubic,PETSC_NULL);
1143: break;
1144: }
1145: }
1146: PetscOptionsTail();
1147: return(0);
1148: }
1149: /* -------------------------------------------------------------------------- */
1150: /*MC
1151: SNESLS - Newton based nonlinear solver that uses a line search
1153: Options Database:
1154: + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
1155: . -snes_ls_alpha <alpha> - Sets alpha
1156: . -snes_ls_maxstep <maxstep> - Sets the maximum stepsize the line search will use (if the 2-norm(y) > maxstep then scale y to be y = (maxstep/2-norm(y)) *y)
1157: - -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda / max_i ( y[i]/x[i] )
1159: Notes: This is the default nonlinear solver in SNES
1161: Level: beginner
1163: .seealso: SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(),
1164: SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(),
1165: SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck()
1167: M*/
1171: PetscErrorCode SNESCreate_LS(SNES snes)
1172: {
1174: SNES_LS *neP;
1177: snes->ops->setup = SNESSetUp_LS;
1178: snes->ops->solve = SNESSolve_LS;
1179: snes->ops->destroy = SNESDestroy_LS;
1180: snes->ops->setfromoptions = SNESSetFromOptions_LS;
1181: snes->ops->view = SNESView_LS;
1183: PetscNewLog(snes,SNES_LS,&neP);
1184: snes->data = (void*)neP;
1185: neP->alpha = 1.e-4;
1186: neP->maxstep = 1.e8;
1187: neP->minlambda = 1.e-12;
1188: neP->LineSearch = SNESLineSearchCubic;
1189: neP->lsP = PETSC_NULL;
1190: neP->postcheckstep = PETSC_NULL;
1191: neP->postcheck = PETSC_NULL;
1192: neP->precheckstep = PETSC_NULL;
1193: neP->precheck = PETSC_NULL;
1195: PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C",
1196: "SNESLineSearchSet_LS",
1197: SNESLineSearchSet_LS);
1198: PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C",
1199: "SNESLineSearchSetPostCheck_LS",
1200: SNESLineSearchSetPostCheck_LS);
1201: PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C",
1202: "SNESLineSearchSetPreCheck_LS",
1203: SNESLineSearchSetPreCheck_LS);
1205: return(0);
1206: }