Actual source code: ls.c
1: /*$Id: ls.c,v 1.170 2001/04/10 19:36:55 bsmith Exp $*/
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 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: int SNESLSCheckLocalMin_Private(Mat A,Vec F,Vec W,PetscReal fnorm,PetscTruth *ismin)
12: {
13: PetscReal a1;
14: int ierr;
15: PetscTruth hastranspose;
18: *ismin = PETSC_FALSE;
19: MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);
20: if (hastranspose) {
21: /* Compute || J^T F|| */
22: MatMultTranspose(A,F,W);
23: VecNorm(W,NORM_2,&a1);
24: PetscLogInfo(0,"SNESSolve_EQ_LS: || J^T F|| %g near zero implies found a local minimumn",a1/fnorm);
25: if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
26: } else {
27: Vec work;
28: Scalar result;
29: PetscReal wnorm;
31: VecSetRandom(PETSC_NULL,W);
32: VecNorm(W,NORM_2,&wnorm);
33: VecDuplicate(W,&work);
34: MatMult(A,W,work);
35: VecDot(F,work,&result);
36: VecDestroy(work);
37: a1 = PetscAbsScalar(result)/(fnorm*wnorm);
38: PetscLogInfo(0,"SNESSolve_EQ_LS: (F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimumn",a1);
39: if (a1 < 1.e-4) *ismin = PETSC_TRUE;
40: }
41: return(0);
42: }
44: /*
45: Checks if J^T(F - AX) = 0
46: */
47: int SNESLSCheckResidual_Private(Mat A,Vec F,Vec X,Vec W1,Vec W2)
48: {
49: PetscReal a1,a2;
50: int ierr;
51: PetscTruth hastranspose;
52: Scalar mone = -1.0;
55: MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);
56: if (hastranspose) {
57: MatMult(A,X,W1);
58: VecAXPY(&mone,F,W1);
60: /* Compute || J^T W|| */
61: MatMultTranspose(A,W1,W2);
62: VecNorm(W1,NORM_2,&a1);
63: VecNorm(W2,NORM_2,&a2);
64: if (a1 != 0) {
65: PetscLogInfo(0,"SNESSolve_EQ_LS: ||J^T(F-Ax)||/||F-AX|| %g near zero implies inconsistent rhsn",a2/a1);
66: }
67: }
68: return(0);
69: }
71: /* --------------------------------------------------------------------
73: This file implements a truncated Newton method with a line search,
74: for solving a system of nonlinear equations, using the SLES, Vec,
75: and Mat interfaces for linear solvers, vectors, and matrices,
76: respectively.
78: The following basic routines are required for each nonlinear solver:
79: SNESCreate_XXX() - Creates a nonlinear solver context
80: SNESSetFromOptions_XXX() - Sets runtime options
81: SNESSolve_XXX() - Solves the nonlinear system
82: SNESDestroy_XXX() - Destroys the nonlinear solver context
83: The suffix "_XXX" denotes a particular implementation, in this case
84: we use _EQ_LS (e.g., SNESCreate_EQ_LS, SNESSolve_EQ_LS) for solving
85: systems of nonlinear equations (EQ) with a line search (LS) method.
86: These routines are actually called via the common user interface
87: routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and
88: SNESDestroy(), so the application code interface remains identical
89: for all nonlinear solvers.
91: Another key routine is:
92: SNESSetUp_XXX() - Prepares for the use of a nonlinear solver
93: by setting data structures and options. The interface routine SNESSetUp()
94: is not usually called directly by the user, but instead is called by
95: SNESSolve() if necessary.
97: Additional basic routines are:
98: SNESView_XXX() - Prints details of runtime options that
99: have actually been used.
100: These are called by application codes via the interface routines
101: SNESView().
103: The various types of solvers (preconditioners, Krylov subspace methods,
104: nonlinear solvers, timesteppers) are all organized similarly, so the
105: above description applies to these categories also.
107: -------------------------------------------------------------------- */
108: /*
109: SNESSolve_EQ_LS - Solves a nonlinear system with a truncated Newton
110: method with a line search.
112: Input Parameters:
113: . snes - the SNES context
115: Output Parameter:
116: . outits - number of iterations until termination
118: Application Interface Routine: SNESSolve()
120: Notes:
121: This implements essentially a truncated Newton method with a
122: line search. By default a cubic backtracking line search
123: is employed, as described in the text "Numerical Methods for
124: Unconstrained Optimization and Nonlinear Equations" by Dennis
125: and Schnabel.
126: */
127: int SNESSolve_EQ_LS(SNES snes,int *outits)
128: {
129: SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data;
130: int maxits,i,ierr,lits,lsfail;
131: MatStructure flg = DIFFERENT_NONZERO_PATTERN;
132: PetscReal fnorm,gnorm,xnorm,ynorm;
133: Vec Y,X,F,G,W,TMP;
136: snes->reason = SNES_CONVERGED_ITERATING;
138: maxits = snes->max_its; /* maximum number of iterations */
139: X = snes->vec_sol; /* solution vector */
140: F = snes->vec_func; /* residual vector */
141: Y = snes->work[0]; /* work vectors */
142: G = snes->work[1];
143: W = snes->work[2];
145: PetscObjectTakeAccess(snes);
146: snes->iter = 0;
147: PetscObjectGrantAccess(snes);
148: SNESComputeFunction(snes,X,F); /* F(X) */
149: VecNorm(F,NORM_2,&fnorm); /* fnorm <- ||F|| */
150: PetscObjectTakeAccess(snes);
151: snes->norm = fnorm;
152: PetscObjectGrantAccess(snes);
153: SNESLogConvHistory(snes,fnorm,0);
154: SNESMonitor(snes,0,fnorm);
156: if (fnorm < snes->atol) {*outits = 0; snes->reason = SNES_CONVERGED_FNORM_ABS; return(0);}
158: /* set parameter for default relative tolerance convergence test */
159: snes->ttol = fnorm*snes->rtol;
161: for (i=0; i<maxits; i++) {
163: /* Solve J Y = F, where J is Jacobian matrix */
164: SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);
165: SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,flg);
166: SLESSolve(snes->sles,F,Y,&lits);
168: if (PetscLogPrintInfo){
169: SNESLSCheckResidual_Private(snes->jacobian,F,Y,G,W);
170: }
172: /* should check what happened to the linear solve? */
173: snes->linear_its += lits;
174: PetscLogInfo(snes,"SNESSolve_EQ_LS: iter=%d, linear solve iterations=%dn",snes->iter,lits);
176: /* Compute a (scaled) negative update in the line search routine:
177: Y <- X - lambda*Y
178: and evaluate G(Y) = function(Y))
179: */
180: VecCopy(Y,snes->vec_sol_update_always);
181: (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lsfail);
182: PetscLogInfo(snes,"SNESSolve_EQ_LS: fnorm=%g, gnorm=%g, ynorm=%g, lsfail=%dn",fnorm,gnorm,ynorm,lsfail);
184: TMP = F; F = G; snes->vec_func_always = F; G = TMP;
185: TMP = X; X = Y; snes->vec_sol_always = X; Y = TMP;
186: fnorm = gnorm;
188: if (lsfail) {
189: PetscTruth ismin;
190: snes->nfailures++;
191: snes->reason = SNES_DIVERGED_LS_FAILURE;
192: SNESLSCheckLocalMin_Private(snes->jacobian,F,W,fnorm,&ismin);
193: if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
194: break;
195: }
197: PetscObjectTakeAccess(snes);
198: snes->iter = i+1;
199: snes->norm = fnorm;
200: PetscObjectGrantAccess(snes);
201: SNESLogConvHistory(snes,fnorm,lits);
202: SNESMonitor(snes,i+1,fnorm);
204: /* Test for convergence */
205: if (snes->converged) {
206: VecNorm(X,NORM_2,&xnorm); /* xnorm = || X || */
207: (*snes->converged)(snes,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);
208: if (snes->reason) {
209: break;
210: }
211: }
212: }
213: if (X != snes->vec_sol) {
214: VecCopy(X,snes->vec_sol);
215: }
216: if (F != snes->vec_func) {
217: VecCopy(F,snes->vec_func);
218: }
219: snes->vec_sol_always = snes->vec_sol;
220: snes->vec_func_always = snes->vec_func;
221: if (i == maxits) {
222: PetscLogInfo(snes,"SNESSolve_EQ_LS: Maximum number of iterations has been reached: %dn",maxits);
223: i--;
224: snes->reason = SNES_DIVERGED_MAX_IT;
225: }
226: PetscObjectTakeAccess(snes);
227: PetscObjectGrantAccess(snes);
228: *outits = i+1;
229: return(0);
230: }
231: /* -------------------------------------------------------------------------- */
232: /*
233: SNESSetUp_EQ_LS - Sets up the internal data structures for the later use
234: of the SNESEQLS nonlinear solver.
236: Input Parameter:
237: . snes - the SNES context
238: . x - the solution vector
240: Application Interface Routine: SNESSetUp()
242: Notes:
243: For basic use of the SNES solvers the user need not explicitly call
244: SNESSetUp(), since these actions will automatically occur during
245: the call to SNESSolve().
246: */
247: int SNESSetUp_EQ_LS(SNES snes)
248: {
252: snes->nwork = 4;
253: VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);
254: PetscLogObjectParents(snes,snes->nwork,snes->work);
255: snes->vec_sol_update_always = snes->work[3];
256: return(0);
257: }
258: /* -------------------------------------------------------------------------- */
259: /*
260: SNESDestroy_EQ_LS - Destroys the private SNESEQLS context that was created
261: with SNESCreate_EQ_LS().
263: Input Parameter:
264: . snes - the SNES context
266: Application Interface Routine: SNESDestroy()
267: */
268: int SNESDestroy_EQ_LS(SNES snes)
269: {
270: int ierr;
273: if (snes->nwork) {
274: VecDestroyVecs(snes->work,snes->nwork);
275: }
276: PetscFree(snes->data);
277: return(0);
278: }
279: /* -------------------------------------------------------------------------- */
281: /*@C
282: SNESNoLineSearch - This routine is not a line search at all;
283: it simply uses the full Newton step. Thus, this routine is intended
284: to serve as a template and is not recommended for general use.
286: Collective on SNES and Vec
288: Input Parameters:
289: + snes - nonlinear context
290: . lsctx - optional context for line search (not used here)
291: . x - current iterate
292: . f - residual evaluated at x
293: . y - search direction (contains new iterate on output)
294: . w - work vector
295: - fnorm - 2-norm of f
297: Output Parameters:
298: + g - residual evaluated at new iterate y
299: . y - new iterate (contains search direction on input)
300: . gnorm - 2-norm of g
301: . ynorm - 2-norm of search length
302: - flag - set to 0, indicating a successful line search
304: Options Database Key:
305: . -snes_eq_ls basic - Activates SNESNoLineSearch()
307: Level: advanced
309: .keywords: SNES, nonlinear, line search, cubic
311: .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
312: SNESSetLineSearch(), SNESNoLineSearchNoNorms()
313: @*/
314: int SNESNoLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag)
315: {
316: int ierr;
317: Scalar mone = -1.0;
318: SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data;
319: PetscTruth change_y = PETSC_FALSE;
322: *flag = 0;
323: PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);
324: VecNorm(y,NORM_2,ynorm); /* ynorm = || y || */
325: VecAYPX(&mone,x,y); /* y <- y - x */
326: if (neP->CheckStep) {
327: (*neP->CheckStep)(snes,neP->checkP,y,&change_y);
328: }
329: SNESComputeFunction(snes,y,g); /* Compute F(y) */
330: VecNorm(g,NORM_2,gnorm); /* gnorm = || g || */
331: PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
332: return(0);
333: }
334: /* -------------------------------------------------------------------------- */
336: /*@C
337: SNESNoLineSearchNoNorms - This routine is not a line search at
338: all; it simply uses the full Newton step. This version does not
339: even compute the norm of the function or search direction; this
340: is intended only when you know the full step is fine and are
341: not checking for convergence of the nonlinear iteration (for
342: example, you are running always for a fixed number of Newton steps).
344: Collective on SNES and Vec
346: Input Parameters:
347: + snes - nonlinear context
348: . lsctx - optional context for line search (not used here)
349: . x - current iterate
350: . f - residual evaluated at x
351: . y - search direction (contains new iterate on output)
352: . w - work vector
353: - fnorm - 2-norm of f
355: Output Parameters:
356: + g - residual evaluated at new iterate y
357: . gnorm - not changed
358: . ynorm - not changed
359: - flag - set to 0, indicating a successful line search
361: Options Database Key:
362: . -snes_eq_ls basicnonorms - Activates SNESNoLineSearchNoNorms()
364: Notes:
365: SNESNoLineSearchNoNorms() must be used in conjunction with
366: either the options
367: $ -snes_no_convergence_test -snes_max_it <its>
368: or alternatively a user-defined custom test set via
369: SNESSetConvergenceTest(); or a -snes_max_it of 1,
370: otherwise, the SNES solver will generate an error.
372: During the final iteration this will not evaluate the function at
373: the solution point. This is to save a function evaluation while
374: using pseudo-timestepping.
376: The residual norms printed by monitoring routines such as
377: SNESDefaultMonitor() (as activated via -snes_monitor) will not be
378: correct, since they are not computed.
380: Level: advanced
382: .keywords: SNES, nonlinear, line search, cubic
384: .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
385: SNESSetLineSearch(), SNESNoLineSearch()
386: @*/
387: int SNESNoLineSearchNoNorms(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag)
388: {
389: int ierr;
390: Scalar mone = -1.0;
391: SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data;
392: PetscTruth change_y = PETSC_FALSE;
395: *flag = 0;
396: PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);
397: VecAYPX(&mone,x,y); /* y <- y - x */
398: if (neP->CheckStep) {
399: (*neP->CheckStep)(snes,neP->checkP,y,&change_y);
400: }
401:
402: /* don't evaluate function the last time through */
403: if (snes->iter < snes->max_its-1) {
404: SNESComputeFunction(snes,y,g); /* Compute F(y) */
405: }
406: PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
407: return(0);
408: }
409: /* -------------------------------------------------------------------------- */
410: /*@C
411: SNESCubicLineSearch - Performs a cubic line search (default line search method).
413: Collective on SNES
415: Input Parameters:
416: + snes - nonlinear context
417: . lsctx - optional context for line search (not used here)
418: . x - current iterate
419: . f - residual evaluated at x
420: . y - search direction (contains new iterate on output)
421: . w - work vector
422: - fnorm - 2-norm of f
424: Output Parameters:
425: + g - residual evaluated at new iterate y
426: . y - new iterate (contains search direction on input)
427: . gnorm - 2-norm of g
428: . ynorm - 2-norm of search length
429: - flag - 0 if line search succeeds; -1 on failure.
431: Options Database Key:
432: $ -snes_eq_ls cubic - Activates SNESCubicLineSearch()
434: Notes:
435: This line search is taken from "Numerical Methods for Unconstrained
436: Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
438: Level: advanced
440: .keywords: SNES, nonlinear, line search, cubic
442: .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
443: @*/
444: int SNESCubicLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag)
445: {
446: /*
447: Note that for line search purposes we work with with the related
448: minimization problem:
449: min z(x): R^n -> R,
450: where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
451: */
452:
453: PetscReal steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2;
454: PetscReal maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg;
455: #if defined(PETSC_USE_COMPLEX)
456: Scalar cinitslope,clambda;
457: #endif
458: int ierr,count;
459: SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data;
460: Scalar mone = -1.0,scale;
461: PetscTruth change_y = PETSC_FALSE;
464: PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);
465: *flag = 0;
466: alpha = neP->alpha;
467: maxstep = neP->maxstep;
468: steptol = neP->steptol;
470: VecNorm(y,NORM_2,ynorm);
471: if (*ynorm < snes->atol) {
472: PetscLogInfo(snes,"SNESCubicLineSearch: Search direction and size are nearly 0n");
473: *gnorm = fnorm;
474: VecCopy(x,y);
475: VecCopy(f,g);
476: goto theend1;
477: }
478: if (*ynorm > maxstep) { /* Step too big, so scale back */
479: scale = maxstep/(*ynorm);
480: #if defined(PETSC_USE_COMPLEX)
481: PetscLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %gn",PetscRealPart(scale),*ynorm);
482: #else
483: PetscLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %gn",scale,*ynorm);
484: #endif
485: VecScale(&scale,y);
486: *ynorm = maxstep;
487: }
488: minlambda = steptol/(*ynorm);
489: MatMult(snes->jacobian,y,w);
490: #if defined(PETSC_USE_COMPLEX)
491: VecDot(f,w,&cinitslope);
492: initslope = PetscRealPart(cinitslope);
493: #else
494: VecDot(f,w,&initslope);
495: #endif
496: if (initslope > 0.0) initslope = -initslope;
497: if (initslope == 0.0) initslope = -1.0;
499: VecCopy(y,w);
500: VecAYPX(&mone,x,w);
501: SNESComputeFunction(snes,w,g);
502: VecNorm(g,NORM_2,gnorm);
503: if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
504: VecCopy(w,y);
505: PetscLogInfo(snes,"SNESCubicLineSearch: Using full stepn");
506: goto theend1;
507: }
509: /* Fit points with quadratic */
510: lambda = 1.0;
511: lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
512: lambdaprev = lambda;
513: gnormprev = *gnorm;
514: if (lambdatemp > .5*lambda) lambdatemp = .5*lambda;
515: if (lambdatemp <= .1*lambda) lambda = .1*lambda;
516: else lambda = lambdatemp;
517: ierr = VecCopy(x,w);
518: lambdaneg = -lambda;
519: #if defined(PETSC_USE_COMPLEX)
520: clambda = lambdaneg; VecAXPY(&clambda,y,w);
521: #else
522: VecAXPY(&lambdaneg,y,w);
523: #endif
524: SNESComputeFunction(snes,w,g);
525: VecNorm(g,NORM_2,gnorm);
526: if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
527: VecCopy(w,y);
528: PetscLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%gn",lambda);
529: goto theend1;
530: }
532: /* Fit points with cubic */
533: count = 1;
534: while (1) {
535: if (lambda <= minlambda) { /* bad luck; use full step */
536: PetscLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d n",count);
537: PetscLogInfo(snes,"SNESCubicLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%gn",fnorm,*gnorm,*ynorm,lambda,initslope);
538: VecCopy(w,y);
539: *flag = -1; break;
540: }
541: t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
542: t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope;
543: a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
544: b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
545: d = b*b - 3*a*initslope;
546: if (d < 0.0) d = 0.0;
547: if (a == 0.0) {
548: lambdatemp = -initslope/(2.0*b);
549: } else {
550: lambdatemp = (-b + sqrt(d))/(3.0*a);
551: }
552: lambdaprev = lambda;
553: gnormprev = *gnorm;
554: if (lambdatemp > .5*lambda) lambdatemp = .5*lambda;
555: if (lambdatemp <= .1*lambda) lambda = .1*lambda;
556: else lambda = lambdatemp;
557: VecCopy(x,w);
558: lambdaneg = -lambda;
559: #if defined(PETSC_USE_COMPLEX)
560: clambda = lambdaneg;
561: VecAXPY(&clambda,y,w);
562: #else
563: VecAXPY(&lambdaneg,y,w);
564: #endif
565: SNESComputeFunction(snes,w,g);
566: VecNorm(g,NORM_2,gnorm);
567: if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */
568: VecCopy(w,y);
569: PetscLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%gn",lambda);
570: goto theend1;
571: } else {
572: PetscLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda, lambda=%gn",lambda);
573: }
574: count++;
575: }
576: theend1:
577: /* Optional user-defined check for line search step validity */
578: if (neP->CheckStep) {
579: (*neP->CheckStep)(snes,neP->checkP,y,&change_y);
580: if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */
581: SNESComputeFunction(snes,y,g);
582: VecNormBegin(y,NORM_2,ynorm);
583: VecNormBegin(g,NORM_2,gnorm);
584: VecNormEnd(y,NORM_2,ynorm);
585: VecNormEnd(g,NORM_2,gnorm);
586: }
587: }
588: PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
589: return(0);
590: }
591: /* -------------------------------------------------------------------------- */
592: /*@C
593: SNESQuadraticLineSearch - Performs a quadratic line search.
595: Collective on SNES and Vec
597: Input Parameters:
598: + snes - the SNES context
599: . lsctx - optional context for line search (not used here)
600: . x - current iterate
601: . f - residual evaluated at x
602: . y - search direction (contains new iterate on output)
603: . w - work vector
604: - fnorm - 2-norm of f
606: Output Parameters:
607: + g - residual evaluated at new iterate y
608: . y - new iterate (contains search direction on input)
609: . gnorm - 2-norm of g
610: . ynorm - 2-norm of search length
611: - flag - 0 if line search succeeds; -1 on failure.
613: Options Database Key:
614: . -snes_eq_ls quadratic - Activates SNESQuadraticLineSearch()
616: Notes:
617: Use SNESSetLineSearch() to set this routine within the SNESEQLS method.
619: Level: advanced
621: .keywords: SNES, nonlinear, quadratic, line search
623: .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
624: @*/
625: int SNESQuadraticLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag)
626: {
627: /*
628: Note that for line search purposes we work with with the related
629: minimization problem:
630: min z(x): R^n -> R,
631: where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
632: */
633: PetscReal steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg;
634: #if defined(PETSC_USE_COMPLEX)
635: Scalar cinitslope,clambda;
636: #endif
637: int ierr,count;
638: SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data;
639: Scalar mone = -1.0,scale;
640: PetscTruth change_y = PETSC_FALSE;
643: PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);
644: *flag = 0;
645: alpha = neP->alpha;
646: maxstep = neP->maxstep;
647: steptol = neP->steptol;
649: VecNorm(y,NORM_2,ynorm);
650: if (*ynorm < snes->atol) {
651: PetscLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0n");
652: *gnorm = fnorm;
653: VecCopy(x,y);
654: VecCopy(f,g);
655: goto theend2;
656: }
657: if (*ynorm > maxstep) { /* Step too big, so scale back */
658: scale = maxstep/(*ynorm);
659: VecScale(&scale,y);
660: *ynorm = maxstep;
661: }
662: minlambda = steptol/(*ynorm);
663: MatMult(snes->jacobian,y,w);
664: #if defined(PETSC_USE_COMPLEX)
665: VecDot(f,w,&cinitslope);
666: initslope = PetscRealPart(cinitslope);
667: #else
668: VecDot(f,w,&initslope);
669: #endif
670: if (initslope > 0.0) initslope = -initslope;
671: if (initslope == 0.0) initslope = -1.0;
673: VecCopy(y,w);
674: VecAYPX(&mone,x,w);
675: SNESComputeFunction(snes,w,g);
676: VecNorm(g,NORM_2,gnorm);
677: if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
678: VecCopy(w,y);
679: PetscLogInfo(snes,"SNESQuadraticLineSearch: Using full stepn");
680: goto theend2;
681: }
683: /* Fit points with quadratic */
684: lambda = 1.0;
685: count = 1;
686: while (1) {
687: if (lambda <= minlambda) { /* bad luck; use full step */
688: PetscLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d n",count);
689: PetscLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%gn",fnorm,*gnorm,*ynorm,lambda,initslope);
690: VecCopy(w,y);
691: *flag = -1; break;
692: }
693: lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
694: if (lambdatemp > .5*lambda) lambdatemp = .5*lambda;
695: if (lambdatemp <= .1*lambda) lambda = .1*lambda;
696: else lambda = lambdatemp;
697: VecCopy(x,w);
698: lambdaneg = -lambda;
699: #if defined(PETSC_USE_COMPLEX)
700: clambda = lambdaneg; VecAXPY(&clambda,y,w);
701: #else
702: VecAXPY(&lambdaneg,y,w);
703: #endif
704: SNESComputeFunction(snes,w,g);
705: VecNorm(g,NORM_2,gnorm);
706: if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
707: VecCopy(w,y);
708: PetscLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%gn",lambda);
709: break;
710: }
711: count++;
712: }
713: theend2:
714: /* Optional user-defined check for line search step validity */
715: if (neP->CheckStep) {
716: (*neP->CheckStep)(snes,neP->checkP,y,&change_y);
717: if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */
718: SNESComputeFunction(snes,y,g);
719: VecNormBegin(y,NORM_2,ynorm);
720: VecNormBegin(g,NORM_2,gnorm);
721: VecNormEnd(y,NORM_2,ynorm);
722: VecNormEnd(g,NORM_2,gnorm);
723: }
724: }
725: PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
726: return(0);
727: }
728: /* -------------------------------------------------------------------------- */
729: /*@C
730: SNESSetLineSearch - Sets the line search routine to be used
731: by the method SNESEQLS.
733: Input Parameters:
734: + snes - nonlinear context obtained from SNESCreate()
735: . lsctx - optional user-defined context for use by line search
736: - func - pointer to int function
738: Collective on SNES
740: Available Routines:
741: + SNESCubicLineSearch() - default line search
742: . SNESQuadraticLineSearch() - quadratic line search
743: . SNESNoLineSearch() - the full Newton step (actually not a line search)
744: - SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel)
746: Options Database Keys:
747: + -snes_eq_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
748: . -snes_eq_ls_alpha <alpha> - Sets alpha
749: . -snes_eq_ls_maxstep <max> - Sets maxstep
750: - -snes_eq_ls_steptol <steptol> - Sets steptol
752: Calling sequence of func:
753: .vb
754: func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
755: PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,*flag)
756: .ve
758: Input parameters for func:
759: + snes - nonlinear context
760: . lsctx - optional user-defined context for line search
761: . x - current iterate
762: . f - residual evaluated at x
763: . y - search direction (contains new iterate on output)
764: . w - work vector
765: - fnorm - 2-norm of f
767: Output parameters for func:
768: + g - residual evaluated at new iterate y
769: . y - new iterate (contains search direction on input)
770: . gnorm - 2-norm of g
771: . ynorm - 2-norm of search length
772: - flag - set to 0 if the line search succeeds; a nonzero integer
773: on failure.
775: Level: advanced
777: .keywords: SNES, nonlinear, set, line search, routine
779: .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(),
780: SNESSetLineSearchCheck(), SNESSetLineSearchParams(), SNESGetLineSearchParams()
781: @*/
782: int SNESSetLineSearch(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void *lsctx)
783: {
784: int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void*);
787: PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void (**)())&f);
788: if (f) {
789: (*f)(snes,func,lsctx);
790: }
791: return(0);
792: }
793: /* -------------------------------------------------------------------------- */
794: EXTERN_C_BEGIN
795: int SNESSetLineSearch_LS(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,
796: double,double*,double*,int*),void *lsctx)
797: {
799: ((SNES_EQ_LS *)(snes->data))->LineSearch = func;
800: ((SNES_EQ_LS *)(snes->data))->lsP = lsctx;
801: return(0);
802: }
803: EXTERN_C_END
804: /* -------------------------------------------------------------------------- */
805: /*@C
806: SNESSetLineSearchCheck - Sets a routine to check the validity of new iterate computed
807: by the line search routine in the Newton-based method SNESEQLS.
809: Input Parameters:
810: + snes - nonlinear context obtained from SNESCreate()
811: . func - pointer to int function
812: - checkctx - optional user-defined context for use by step checking routine
814: Collective on SNES
816: Calling sequence of func:
817: .vb
818: int func (SNES snes, void *checkctx, Vec x, PetscTruth *flag)
819: .ve
820: where func returns an error code of 0 on success and a nonzero
821: on failure.
823: Input parameters for func:
824: + snes - nonlinear context
825: . checkctx - optional user-defined context for use by step checking routine
826: - x - current candidate iterate
828: Output parameters for func:
829: + x - current iterate (possibly modified)
830: - flag - flag indicating whether x has been modified (either
831: PETSC_TRUE of PETSC_FALSE)
833: Level: advanced
835: Notes:
836: SNESNoLineSearch() and SNESNoLineSearchNoNorms() accept the new
837: iterate computed by the line search checking routine. In particular,
838: these routines (1) compute a candidate iterate u_{i+1}, (2) pass control
839: to the checking routine, and then (3) compute the corresponding nonlinear
840: function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}.
842: SNESQuadraticLineSearch() and SNESCubicLineSearch() also accept the
843: new iterate computed by the line search checking routine. In particular,
844: these routines (1) compute a candidate iterate u_{i+1} as well as a
845: candidate nonlinear function f(u_{i+1}), (2) pass control to the checking
846: routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes
847: were made to the candidate iterate in the checking routine (as indicated
848: by flag=PETSC_TRUE). The overhead of this function re-evaluation can be
849: very costly, so use this feature with caution!
851: .keywords: SNES, nonlinear, set, line search check, step check, routine
853: .seealso: SNESSetLineSearch()
854: @*/
855: int SNESSetLineSearchCheck(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx)
856: {
857: int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,PetscTruth*),void*);
860: PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearchCheck_C",(void (**)())&f);
861: if (f) {
862: (*f)(snes,func,checkctx);
863: }
864: return(0);
865: }
866: /* -------------------------------------------------------------------------- */
867: EXTERN_C_BEGIN
868: int SNESSetLineSearchCheck_LS(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx)
869: {
871: ((SNES_EQ_LS *)(snes->data))->CheckStep = func;
872: ((SNES_EQ_LS *)(snes->data))->checkP = checkctx;
873: return(0);
874: }
875: EXTERN_C_END
876: /* -------------------------------------------------------------------------- */
877: /*
878: SNESView_EQ_LS - Prints info from the SNESEQLS data structure.
880: Input Parameters:
881: . SNES - the SNES context
882: . viewer - visualization context
884: Application Interface Routine: SNESView()
885: */
886: static int SNESView_EQ_LS(SNES snes,PetscViewer viewer)
887: {
888: SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data;
889: char *cstr;
890: int ierr;
891: PetscTruth isascii;
894: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
895: if (isascii) {
896: if (ls->LineSearch == SNESNoLineSearch) cstr = "SNESNoLineSearch";
897: else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch";
898: else if (ls->LineSearch == SNESCubicLineSearch) cstr = "SNESCubicLineSearch";
899: else cstr = "unknown";
900: PetscViewerASCIIPrintf(viewer," line search variant: %sn",cstr);
901: PetscViewerASCIIPrintf(viewer," alpha=%g, maxstep=%g, steptol=%gn",ls->alpha,ls->maxstep,ls->steptol);
902: } else {
903: SETERRQ1(1,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name);
904: }
905: return(0);
906: }
907: /* -------------------------------------------------------------------------- */
908: /*
909: SNESSetFromOptions_EQ_LS - Sets various parameters for the SNESEQLS method.
911: Input Parameter:
912: . snes - the SNES context
914: Application Interface Routine: SNESSetFromOptions()
915: */
916: static int SNESSetFromOptions_EQ_LS(SNES snes)
917: {
918: SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data;
919: char ver[16],*lses[] = {"basic","basicnonorms","quadratic","cubic"};
920: int ierr;
921: PetscTruth flg;
924: PetscOptionsHead("SNES Line search options");
925: PetscOptionsDouble("-snes_eq_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);
926: PetscOptionsDouble("-snes_eq_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);
927: PetscOptionsDouble("-snes_eq_ls_steptol","Step must be greater than","None",ls->steptol,&ls->steptol,0);
929: PetscOptionsEList("-snes_eq_ls","Line search used","SNESSetLineSearch",lses,4,"cubic",ver,16,&flg);
930: if (flg) {
931: PetscTruth isbasic,isnonorms,isquad,iscubic;
933: PetscStrcmp(ver,lses[0],&isbasic);
934: PetscStrcmp(ver,lses[1],&isnonorms);
935: PetscStrcmp(ver,lses[2],&isquad);
936: PetscStrcmp(ver,lses[3],&iscubic);
938: if (isbasic) {
939: SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);
940: } else if (isnonorms) {
941: SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);
942: } else if (isquad) {
943: SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);
944: } else if (iscubic) {
945: SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);
946: }
947: else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Unknown line search");}
948: }
949: PetscOptionsTail();
950: return(0);
951: }
952: /* -------------------------------------------------------------------------- */
953: /*
954: SNESCreate_EQ_LS - Creates a nonlinear solver context for the SNESEQLS method,
955: SNES_EQ_LS, and sets this as the private data within the generic nonlinear solver
956: context, SNES, that was created within SNESCreate().
959: Input Parameter:
960: . snes - the SNES context
962: Application Interface Routine: SNESCreate()
963: */
964: EXTERN_C_BEGIN
965: int SNESCreate_EQ_LS(SNES snes)
966: {
967: int ierr;
968: SNES_EQ_LS *neP;
971: if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
972: SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only");
973: }
975: snes->setup = SNESSetUp_EQ_LS;
976: snes->solve = SNESSolve_EQ_LS;
977: snes->destroy = SNESDestroy_EQ_LS;
978: snes->converged = SNESConverged_EQ_LS;
979: snes->setfromoptions = SNESSetFromOptions_EQ_LS;
980: snes->view = SNESView_EQ_LS;
981: snes->nwork = 0;
983: ierr = PetscNew(SNES_EQ_LS,&neP);
984: PetscLogObjectMemory(snes,sizeof(SNES_EQ_LS));
985: snes->data = (void*)neP;
986: neP->alpha = 1.e-4;
987: neP->maxstep = 1.e8;
988: neP->steptol = 1.e-12;
989: neP->LineSearch = SNESCubicLineSearch;
990: neP->lsP = PETSC_NULL;
991: neP->CheckStep = PETSC_NULL;
992: neP->checkP = PETSC_NULL;
994: PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS",SNESSetLineSearch_LS);
995: PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS",SNESSetLineSearchCheck_LS);
997: return(0);
998: }
999: EXTERN_C_END