Actual source code: ls.c
1: /*$Id: ls.c,v 1.172 2001/08/07 03:04:11 balay Exp $*/
3: #include src/snes/impls/ls/ls.h
6: /*
7: max { p[i]/x[i] }
8: */
9: int VecMaxScale_SNES(Vec p,Vec x,PetscReal *m)
10: {
11: int ierr,i,n;
12: PetscScalar *pa,*xa;
13: PetscReal t;
14: MPI_Comm comm;
17: VecGetLocalSize(p,&n);
18: PetscObjectGetComm((PetscObject)p,&comm);
20: VecGetArray(p,&pa);
21: VecGetArray(x,&xa);
22: t = 0.0;
23: for ( i=0; i<n; i++) {
24: if (xa[i] != 0.0) {
25: t = PetscMax(PetscAbsScalar(pa[i]/xa[i]),t);
26: } else {
27: t = PetscMax(PetscAbsScalar(pa[i]),t);
28: }
29: }
30: MPI_Allreduce(&t,m,1,MPI_DOUBLE,MPI_MAX,comm);
31: VecRestoreArray(p,&pa);
32: VecRestoreArray(x,&xa);
33: return(0);
34: }
36: /*
37: Checks if J^T F = 0 which implies we've found a local minimum of the function,
38: but not a zero. In the case when one cannot compute J^T F we use the fact that
39: 0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
40: for this trick.
41: */
42: int SNESLSCheckLocalMin_Private(Mat A,Vec F,Vec W,PetscReal fnorm,PetscTruth *ismin)
43: {
44: PetscReal a1;
45: int ierr;
46: PetscTruth hastranspose;
49: *ismin = PETSC_FALSE;
50: MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);
51: if (hastranspose) {
52: /* Compute || J^T F|| */
53: MatMultTranspose(A,F,W);
54: VecNorm(W,NORM_2,&a1);
55: PetscLogInfo(0,"SNESSolve_LS: || J^T F|| %g near zero implies found a local minimumn",a1/fnorm);
56: if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
57: } else {
58: Vec work;
59: PetscScalar result;
60: PetscReal wnorm;
62: VecSetRandom(PETSC_NULL,W);
63: VecNorm(W,NORM_2,&wnorm);
64: VecDuplicate(W,&work);
65: MatMult(A,W,work);
66: VecDot(F,work,&result);
67: VecDestroy(work);
68: a1 = PetscAbsScalar(result)/(fnorm*wnorm);
69: PetscLogInfo(0,"SNESSolve_LS: (F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimumn",a1);
70: if (a1 < 1.e-4) *ismin = PETSC_TRUE;
71: }
72: return(0);
73: }
75: /*
76: Checks if J^T(F - J*X) = 0
77: */
78: int SNESLSCheckResidual_Private(Mat A,Vec F,Vec X,Vec W1,Vec W2)
79: {
80: PetscReal a1,a2;
81: int ierr;
82: PetscTruth hastranspose;
83: PetscScalar mone = -1.0;
86: MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);
87: if (hastranspose) {
88: MatMult(A,X,W1);
89: VecAXPY(&mone,F,W1);
91: /* Compute || J^T W|| */
92: MatMultTranspose(A,W1,W2);
93: VecNorm(W1,NORM_2,&a1);
94: VecNorm(W2,NORM_2,&a2);
95: if (a1 != 0) {
96: PetscLogInfo(0,"SNESSolve_LS: ||J^T(F-Ax)||/||F-AX|| %g near zero implies inconsistent rhsn",a2/a1);
97: }
98: }
99: return(0);
100: }
102: /* --------------------------------------------------------------------
104: This file implements a truncated Newton method with a line search,
105: for solving a system of nonlinear equations, using the SLES, Vec,
106: and Mat interfaces for linear solvers, vectors, and matrices,
107: respectively.
109: The following basic routines are required for each nonlinear solver:
110: SNESCreate_XXX() - Creates a nonlinear solver context
111: SNESSetFromOptions_XXX() - Sets runtime options
112: SNESSolve_XXX() - Solves the nonlinear system
113: SNESDestroy_XXX() - Destroys the nonlinear solver context
114: The suffix "_XXX" denotes a particular implementation, in this case
115: we use _LS (e.g., SNESCreate_LS, SNESSolve_LS) for solving
116: systems of nonlinear equations with a line search (LS) method.
117: These routines are actually called via the common user interface
118: routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and
119: SNESDestroy(), so the application code interface remains identical
120: for all nonlinear solvers.
122: Another key routine is:
123: SNESSetUp_XXX() - Prepares for the use of a nonlinear solver
124: by setting data structures and options. The interface routine SNESSetUp()
125: is not usually called directly by the user, but instead is called by
126: SNESSolve() if necessary.
128: Additional basic routines are:
129: SNESView_XXX() - Prints details of runtime options that
130: have actually been used.
131: These are called by application codes via the interface routines
132: SNESView().
134: The various types of solvers (preconditioners, Krylov subspace methods,
135: nonlinear solvers, timesteppers) are all organized similarly, so the
136: above description applies to these categories also.
138: -------------------------------------------------------------------- */
139: /*
140: SNESSolve_LS - Solves a nonlinear system with a truncated Newton
141: method with a line search.
143: Input Parameters:
144: . snes - the SNES context
146: Output Parameter:
147: . outits - number of iterations until termination
149: Application Interface Routine: SNESSolve()
151: Notes:
152: This implements essentially a truncated Newton method with a
153: line search. By default a cubic backtracking line search
154: is employed, as described in the text "Numerical Methods for
155: Unconstrained Optimization and Nonlinear Equations" by Dennis
156: and Schnabel.
157: */
158: int SNESSolve_LS(SNES snes,int *outits)
159: {
160: SNES_LS *neP = (SNES_LS*)snes->data;
161: int maxits,i,ierr,lits,lsfail;
162: MatStructure flg = DIFFERENT_NONZERO_PATTERN;
163: PetscReal fnorm,gnorm,xnorm,ynorm;
164: Vec Y,X,F,G,W,TMP;
167: snes->reason = SNES_CONVERGED_ITERATING;
169: maxits = snes->max_its; /* maximum number of iterations */
170: X = snes->vec_sol; /* solution vector */
171: F = snes->vec_func; /* residual vector */
172: Y = snes->work[0]; /* work vectors */
173: G = snes->work[1];
174: W = snes->work[2];
176: PetscObjectTakeAccess(snes);
177: snes->iter = 0;
178: PetscObjectGrantAccess(snes);
179: SNESComputeFunction(snes,X,F); /* F(X) */
180: VecNorm(F,NORM_2,&fnorm); /* fnorm <- ||F|| */
181: PetscObjectTakeAccess(snes);
182: snes->norm = fnorm;
183: PetscObjectGrantAccess(snes);
184: SNESLogConvHistory(snes,fnorm,0);
185: SNESMonitor(snes,0,fnorm);
187: if (fnorm < snes->atol) {*outits = 0; snes->reason = SNES_CONVERGED_FNORM_ABS; return(0);}
189: /* set parameter for default relative tolerance convergence test */
190: snes->ttol = fnorm*snes->rtol;
192: for (i=0; i<maxits; i++) {
194: /* Call general purpose update function */
195: if (snes->update != PETSC_NULL) {
196: (*snes->update)(snes, snes->iter);
197: }
199: /* Solve J Y = F, where J is Jacobian matrix */
200: SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);
201: SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,flg);
202: SLESSolve(snes->sles,F,Y,&lits);
204: if (PetscLogPrintInfo){
205: SNESLSCheckResidual_Private(snes->jacobian,F,Y,G,W);
206: }
208: /* should check what happened to the linear solve? */
209: snes->linear_its += lits;
210: PetscLogInfo(snes,"SNESSolve_LS: iter=%d, linear solve iterations=%dn",snes->iter,lits);
212: /* Compute a (scaled) negative update in the line search routine:
213: Y <- X - lambda*Y
214: and evaluate G(Y) = function(Y))
215: */
216: VecCopy(Y,snes->vec_sol_update_always);
217: (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lsfail);
218: PetscLogInfo(snes,"SNESSolve_LS: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lsfail=%dn",fnorm,gnorm,ynorm,lsfail);
220: TMP = F; F = G; snes->vec_func_always = F; G = TMP;
221: TMP = X; X = Y; snes->vec_sol_always = X; Y = TMP;
222: fnorm = gnorm;
224: PetscObjectTakeAccess(snes);
225: snes->iter = i+1;
226: snes->norm = fnorm;
227: PetscObjectGrantAccess(snes);
228: SNESLogConvHistory(snes,fnorm,lits);
229: SNESMonitor(snes,i+1,fnorm);
231: if (lsfail) {
232: PetscTruth ismin;
234: if (++snes->numFailures >= snes->maxFailures) {
235: snes->reason = SNES_DIVERGED_LS_FAILURE;
236: SNESLSCheckLocalMin_Private(snes->jacobian,F,W,fnorm,&ismin);
237: if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
238: break;
239: }
240: }
242: /* Test for convergence */
243: if (snes->converged) {
244: VecNorm(X,NORM_2,&xnorm); /* xnorm = || X || */
245: (*snes->converged)(snes,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);
246: if (snes->reason) {
247: break;
248: }
249: }
250: }
251: if (X != snes->vec_sol) {
252: VecCopy(X,snes->vec_sol);
253: }
254: if (F != snes->vec_func) {
255: VecCopy(F,snes->vec_func);
256: }
257: snes->vec_sol_always = snes->vec_sol;
258: snes->vec_func_always = snes->vec_func;
259: if (i == maxits) {
260: PetscLogInfo(snes,"SNESSolve_LS: Maximum number of iterations has been reached: %dn",maxits);
261: i--;
262: snes->reason = SNES_DIVERGED_MAX_IT;
263: }
264: PetscObjectTakeAccess(snes);
265: PetscObjectGrantAccess(snes);
266: *outits = i+1;
267: return(0);
268: }
269: /* -------------------------------------------------------------------------- */
270: /*
271: SNESSetUp_LS - Sets up the internal data structures for the later use
272: of the SNESLS nonlinear solver.
274: Input Parameter:
275: . snes - the SNES context
276: . x - the solution vector
278: Application Interface Routine: SNESSetUp()
280: Notes:
281: For basic use of the SNES solvers, the user need not explicitly call
282: SNESSetUp(), since these actions will automatically occur during
283: the call to SNESSolve().
284: */
285: int SNESSetUp_LS(SNES snes)
286: {
290: snes->nwork = 4;
291: VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);
292: PetscLogObjectParents(snes,snes->nwork,snes->work);
293: snes->vec_sol_update_always = snes->work[3];
294: return(0);
295: }
296: /* -------------------------------------------------------------------------- */
297: /*
298: SNESDestroy_LS - Destroys the private SNES_LS context that was created
299: with SNESCreate_LS().
301: Input Parameter:
302: . snes - the SNES context
304: Application Interface Routine: SNESDestroy()
305: */
306: int SNESDestroy_LS(SNES snes)
307: {
308: int ierr;
311: if (snes->nwork) {
312: VecDestroyVecs(snes->work,snes->nwork);
313: }
314: PetscFree(snes->data);
315: return(0);
316: }
317: /* -------------------------------------------------------------------------- */
319: /*@C
320: SNESNoLineSearch - This routine is not a line search at all;
321: it simply uses the full Newton step. Thus, this routine is intended
322: to serve as a template and is not recommended for general use.
324: Collective on SNES and Vec
326: Input Parameters:
327: + snes - nonlinear context
328: . lsctx - optional context for line search (not used here)
329: . x - current iterate
330: . f - residual evaluated at x
331: . y - search direction (contains new iterate on output)
332: . w - work vector
333: - fnorm - 2-norm of f
335: Output Parameters:
336: + g - residual evaluated at new iterate y
337: . y - new iterate (contains search direction on input)
338: . gnorm - 2-norm of g
339: . ynorm - 2-norm of search length
340: - flag - set to 0, indicating a successful line search
342: Options Database Key:
343: . -snes_ls basic - Activates SNESNoLineSearch()
345: Level: advanced
347: .keywords: SNES, nonlinear, line search, cubic
349: .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
350: SNESSetLineSearch(), SNESNoLineSearchNoNorms()
351: @*/
352: int SNESNoLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag)
353: {
354: int ierr;
355: PetscScalar mone = -1.0;
356: SNES_LS *neP = (SNES_LS*)snes->data;
357: PetscTruth change_y = PETSC_FALSE;
360: *flag = 0;
361: PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);
362: VecNorm(y,NORM_2,ynorm); /* ynorm = || y || */
363: VecAYPX(&mone,x,y); /* y <- y - x */
364: if (neP->CheckStep) {
365: (*neP->CheckStep)(snes,neP->checkP,y,&change_y);
366: }
367: SNESComputeFunction(snes,y,g); /* Compute F(y) */
368: VecNorm(g,NORM_2,gnorm); /* gnorm = || g || */
369: PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
370: return(0);
371: }
372: /* -------------------------------------------------------------------------- */
374: /*@C
375: SNESNoLineSearchNoNorms - This routine is not a line search at
376: all; it simply uses the full Newton step. This version does not
377: even compute the norm of the function or search direction; this
378: is intended only when you know the full step is fine and are
379: not checking for convergence of the nonlinear iteration (for
380: example, you are running always for a fixed number of Newton steps).
382: Collective on SNES and Vec
384: Input Parameters:
385: + snes - nonlinear context
386: . lsctx - optional context for line search (not used here)
387: . x - current iterate
388: . f - residual evaluated at x
389: . y - search direction (contains new iterate on output)
390: . w - work vector
391: - fnorm - 2-norm of f
393: Output Parameters:
394: + g - residual evaluated at new iterate y
395: . gnorm - not changed
396: . ynorm - not changed
397: - flag - set to 0, indicating a successful line search
399: Options Database Key:
400: . -snes_ls basicnonorms - Activates SNESNoLineSearchNoNorms()
402: Notes:
403: SNESNoLineSearchNoNorms() must be used in conjunction with
404: either the options
405: $ -snes_no_convergence_test -snes_max_it <its>
406: or alternatively a user-defined custom test set via
407: SNESSetConvergenceTest(); or a -snes_max_it of 1,
408: otherwise, the SNES solver will generate an error.
410: During the final iteration this will not evaluate the function at
411: the solution point. This is to save a function evaluation while
412: using pseudo-timestepping.
414: The residual norms printed by monitoring routines such as
415: SNESDefaultMonitor() (as activated via -snes_monitor) will not be
416: correct, since they are not computed.
418: Level: advanced
420: .keywords: SNES, nonlinear, line search, cubic
422: .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
423: SNESSetLineSearch(), SNESNoLineSearch()
424: @*/
425: int SNESNoLineSearchNoNorms(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag)
426: {
427: int ierr;
428: PetscScalar mone = -1.0;
429: SNES_LS *neP = (SNES_LS*)snes->data;
430: PetscTruth change_y = PETSC_FALSE;
433: *flag = 0;
434: PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);
435: VecAYPX(&mone,x,y); /* y <- y - x */
436: if (neP->CheckStep) {
437: (*neP->CheckStep)(snes,neP->checkP,y,&change_y);
438: }
439:
440: /* don't evaluate function the last time through */
441: if (snes->iter < snes->max_its-1) {
442: SNESComputeFunction(snes,y,g); /* Compute F(y) */
443: }
444: PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
445: return(0);
446: }
447: /* -------------------------------------------------------------------------- */
448: /*@C
449: SNESCubicLineSearch - Performs a cubic line search (default line search method).
451: Collective on SNES
453: Input Parameters:
454: + snes - nonlinear context
455: . lsctx - optional context for line search (not used here)
456: . x - current iterate
457: . f - residual evaluated at x
458: . y - search direction (contains new iterate on output)
459: . w - work vector
460: - fnorm - 2-norm of f
462: Output Parameters:
463: + g - residual evaluated at new iterate y
464: . y - new iterate (contains search direction on input)
465: . gnorm - 2-norm of g
466: . ynorm - 2-norm of search length
467: - flag - 0 if line search succeeds; -1 on failure.
469: Options Database Key:
470: $ -snes_ls cubic - Activates SNESCubicLineSearch()
472: Notes:
473: This line search is taken from "Numerical Methods for Unconstrained
474: Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
476: Level: advanced
478: .keywords: SNES, nonlinear, line search, cubic
480: .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
481: @*/
482: int SNESCubicLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag)
483: {
484: /*
485: Note that for line search purposes we work with with the related
486: minimization problem:
487: min z(x): R^n -> R,
488: where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
489: */
490:
491: PetscReal steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
492: PetscReal maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg;
493: #if defined(PETSC_USE_COMPLEX)
494: PetscScalar cinitslope,clambda;
495: #endif
496: int ierr,count;
497: SNES_LS *neP = (SNES_LS*)snes->data;
498: PetscScalar mone = -1.0,scale;
499: PetscTruth change_y = PETSC_FALSE;
502: PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);
503: *flag = 0;
504: alpha = neP->alpha;
505: maxstep = neP->maxstep;
506: steptol = neP->steptol;
508: VecNorm(y,NORM_2,ynorm);
509: if (*ynorm < snes->atol) {
510: PetscLogInfo(snes,"SNESCubicLineSearch: Search direction and size are nearly 0n");
511: *gnorm = fnorm;
512: VecCopy(x,y);
513: VecCopy(f,g);
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: PetscLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %gn",PetscRealPart(scale),*ynorm);
520: #else
521: PetscLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %gn",scale,*ynorm);
522: #endif
523: VecScale(&scale,y);
524: *ynorm = maxstep;
525: }
526: ierr = VecMaxScale_SNES(y,x,&rellength);
527: minlambda = steptol/rellength;
528: MatMult(snes->jacobian,y,w);
529: #if defined(PETSC_USE_COMPLEX)
530: VecDot(f,w,&cinitslope);
531: initslope = PetscRealPart(cinitslope);
532: #else
533: VecDot(f,w,&initslope);
534: #endif
535: if (initslope > 0.0) initslope = -initslope;
536: if (initslope == 0.0) initslope = -1.0;
538: VecCopy(y,w);
539: VecAYPX(&mone,x,w);
540: SNESComputeFunction(snes,w,g);
541: VecNorm(g,NORM_2,gnorm);
542: if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
543: VecCopy(w,y);
544: PetscLogInfo(snes,"SNESCubicLineSearch: Using full stepn");
545: goto theend1;
546: }
548: /* Fit points with quadratic */
549: lambda = 1.0;
550: lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
551: lambdaprev = lambda;
552: gnormprev = *gnorm;
553: if (lambdatemp > .5*lambda) lambdatemp = .5*lambda;
554: if (lambdatemp <= .1*lambda) lambda = .1*lambda;
555: else lambda = lambdatemp;
556: ierr = VecCopy(x,w);
557: lambdaneg = -lambda;
558: #if defined(PETSC_USE_COMPLEX)
559: clambda = lambdaneg; VecAXPY(&clambda,y,w);
560: #else
561: VecAXPY(&lambdaneg,y,w);
562: #endif
563: SNESComputeFunction(snes,w,g);
564: VecNorm(g,NORM_2,gnorm);
565: if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
566: VecCopy(w,y);
567: PetscLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%18.16en",lambda);
568: goto theend1;
569: }
571: /* Fit points with cubic */
572: count = 1;
573: while (PETSC_TRUE) {
574: if (lambda <= minlambda) { /* bad luck; use full step */
575: PetscLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d n",count);
576: PetscLogInfo(snes,"SNESCubicLineSearch:fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16en",fnorm,*gnorm,*ynorm,lambda,initslope);
577: VecCopy(x,y);
578: *flag = -1; break;
579: }
580: t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
581: t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope;
582: a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
583: b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
584: d = b*b - 3*a*initslope;
585: if (d < 0.0) d = 0.0;
586: if (a == 0.0) {
587: lambdatemp = -initslope/(2.0*b);
588: } else {
589: lambdatemp = (-b + sqrt(d))/(3.0*a);
590: }
591: lambdaprev = lambda;
592: gnormprev = *gnorm;
593: if (lambdatemp > .5*lambda) lambdatemp = .5*lambda;
594: if (lambdatemp <= .1*lambda) lambda = .1*lambda;
595: else lambda = lambdatemp;
596: VecCopy(x,w);
597: lambdaneg = -lambda;
598: #if defined(PETSC_USE_COMPLEX)
599: clambda = lambdaneg;
600: VecAXPY(&clambda,y,w);
601: #else
602: VecAXPY(&lambdaneg,y,w);
603: #endif
604: SNESComputeFunction(snes,w,g);
605: VecNorm(g,NORM_2,gnorm);
606: if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */
607: VecCopy(w,y);
608: PetscLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%18.16en",lambda);
609: goto theend1;
610: } else {
611: PetscLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda, lambda=%18.16en",lambda);
612: }
613: count++;
614: }
615: theend1:
616: /* Optional user-defined check for line search step validity */
617: if (neP->CheckStep) {
618: (*neP->CheckStep)(snes,neP->checkP,y,&change_y);
619: if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */
620: SNESComputeFunction(snes,y,g);
621: VecNormBegin(y,NORM_2,ynorm);
622: VecNormBegin(g,NORM_2,gnorm);
623: VecNormEnd(y,NORM_2,ynorm);
624: VecNormEnd(g,NORM_2,gnorm);
625: }
626: }
627: PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
628: return(0);
629: }
630: /* -------------------------------------------------------------------------- */
631: /*@C
632: SNESQuadraticLineSearch - Performs a quadratic line search.
634: Collective on SNES and Vec
636: Input Parameters:
637: + snes - the SNES context
638: . lsctx - optional context for line search (not used here)
639: . x - current iterate
640: . f - residual evaluated at x
641: . y - search direction (contains new iterate on output)
642: . w - work vector
643: - fnorm - 2-norm of f
645: Output Parameters:
646: + g - residual evaluated at new iterate y
647: . y - new iterate (contains search direction on input)
648: . gnorm - 2-norm of g
649: . ynorm - 2-norm of search length
650: - flag - 0 if line search succeeds; -1 on failure.
652: Options Database Key:
653: . -snes_ls quadratic - Activates SNESQuadraticLineSearch()
655: Notes:
656: Use SNESSetLineSearch() to set this routine within the SNESLS method.
658: Level: advanced
660: .keywords: SNES, nonlinear, quadratic, line search
662: .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
663: @*/
664: int SNESQuadraticLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag)
665: {
666: /*
667: Note that for line search purposes we work with with the related
668: minimization problem:
669: min z(x): R^n -> R,
670: where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
671: */
672: PetscReal steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg,rellength;
673: #if defined(PETSC_USE_COMPLEX)
674: PetscScalar cinitslope,clambda;
675: #endif
676: int ierr,count;
677: SNES_LS *neP = (SNES_LS*)snes->data;
678: PetscScalar mone = -1.0,scale;
679: PetscTruth change_y = PETSC_FALSE;
682: PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);
683: *flag = 0;
684: alpha = neP->alpha;
685: maxstep = neP->maxstep;
686: steptol = neP->steptol;
688: VecNorm(y,NORM_2,ynorm);
689: if (*ynorm < snes->atol) {
690: PetscLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0n");
691: *gnorm = fnorm;
692: VecCopy(x,y);
693: VecCopy(f,g);
694: goto theend2;
695: }
696: if (*ynorm > maxstep) { /* Step too big, so scale back */
697: scale = maxstep/(*ynorm);
698: VecScale(&scale,y);
699: *ynorm = maxstep;
700: }
701: ierr = VecMaxScale_SNES(y,x,&rellength);
702: minlambda = steptol/rellength;
703: MatMult(snes->jacobian,y,w);
704: #if defined(PETSC_USE_COMPLEX)
705: VecDot(f,w,&cinitslope);
706: initslope = PetscRealPart(cinitslope);
707: #else
708: VecDot(f,w,&initslope);
709: #endif
710: if (initslope > 0.0) initslope = -initslope;
711: if (initslope == 0.0) initslope = -1.0;
713: VecCopy(y,w);
714: VecAYPX(&mone,x,w);
715: SNESComputeFunction(snes,w,g);
716: VecNorm(g,NORM_2,gnorm);
717: if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
718: VecCopy(w,y);
719: PetscLogInfo(snes,"SNESQuadraticLineSearch: Using full stepn");
720: goto theend2;
721: }
723: /* Fit points with quadratic */
724: lambda = 1.0;
725: count = 1;
726: while (PETSC_TRUE) {
727: if (lambda <= minlambda) { /* bad luck; use full step */
728: PetscLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d n",count);
729: PetscLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%gn",fnorm,*gnorm,*ynorm,lambda,initslope);
730: VecCopy(x,y);
731: *flag = -1; break;
732: }
733: lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
734: if (lambdatemp > .5*lambda) lambdatemp = .5*lambda;
735: if (lambdatemp <= .1*lambda) lambda = .1*lambda;
736: else lambda = lambdatemp;
737: VecCopy(x,w);
738: lambdaneg = -lambda;
739: #if defined(PETSC_USE_COMPLEX)
740: clambda = lambdaneg; VecAXPY(&clambda,y,w);
741: #else
742: VecAXPY(&lambdaneg,y,w);
743: #endif
744: SNESComputeFunction(snes,w,g);
745: VecNorm(g,NORM_2,gnorm);
746: if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
747: VecCopy(w,y);
748: PetscLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%gn",lambda);
749: break;
750: }
751: count++;
752: }
753: theend2:
754: /* Optional user-defined check for line search step validity */
755: if (neP->CheckStep) {
756: (*neP->CheckStep)(snes,neP->checkP,y,&change_y);
757: if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */
758: SNESComputeFunction(snes,y,g);
759: VecNormBegin(y,NORM_2,ynorm);
760: VecNormBegin(g,NORM_2,gnorm);
761: VecNormEnd(y,NORM_2,ynorm);
762: VecNormEnd(g,NORM_2,gnorm);
763: }
764: }
765: PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
766: return(0);
767: }
768: /* -------------------------------------------------------------------------- */
769: /*@C
770: SNESSetLineSearch - Sets the line search routine to be used
771: by the method SNESLS.
773: Input Parameters:
774: + snes - nonlinear context obtained from SNESCreate()
775: . lsctx - optional user-defined context for use by line search
776: - func - pointer to int function
778: Collective on SNES
780: Available Routines:
781: + SNESCubicLineSearch() - default line search
782: . SNESQuadraticLineSearch() - quadratic line search
783: . SNESNoLineSearch() - the full Newton step (actually not a line search)
784: - SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel)
786: Options Database Keys:
787: + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
788: . -snes_ls_alpha <alpha> - Sets alpha
789: . -snes_ls_maxstep <max> - Sets maxstep
790: - -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code
791: will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length
792: the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x)
794: Calling sequence of func:
795: .vb
796: func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
797: PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,*flag)
798: .ve
800: Input parameters for func:
801: + snes - nonlinear context
802: . lsctx - optional user-defined context for line search
803: . x - current iterate
804: . f - residual evaluated at x
805: . y - search direction (contains new iterate on output)
806: . w - work vector
807: - fnorm - 2-norm of f
809: Output parameters for func:
810: + g - residual evaluated at new iterate y
811: . y - new iterate (contains search direction on input)
812: . gnorm - 2-norm of g
813: . ynorm - 2-norm of search length
814: - flag - set to 0 if the line search succeeds; a nonzero integer
815: on failure.
817: Level: advanced
819: .keywords: SNES, nonlinear, set, line search, routine
821: .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(),
822: SNESSetLineSearchCheck(), SNESSetLineSearchParams(), SNESGetLineSearchParams()
823: @*/
824: int SNESSetLineSearch(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void *lsctx)
825: {
826: int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void*);
829: PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void (**)(void))&f);
830: if (f) {
831: (*f)(snes,func,lsctx);
832: }
833: return(0);
834: }
835: /* -------------------------------------------------------------------------- */
836: EXTERN_C_BEGIN
837: int SNESSetLineSearch_LS(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,
838: PetscReal,PetscReal*,PetscReal*,int*),void *lsctx)
839: {
841: ((SNES_LS *)(snes->data))->LineSearch = func;
842: ((SNES_LS *)(snes->data))->lsP = lsctx;
843: return(0);
844: }
845: EXTERN_C_END
846: /* -------------------------------------------------------------------------- */
847: /*@C
848: SNESSetLineSearchCheck - Sets a routine to check the validity of new iterate computed
849: by the line search routine in the Newton-based method SNESLS.
851: Input Parameters:
852: + snes - nonlinear context obtained from SNESCreate()
853: . func - pointer to int function
854: - checkctx - optional user-defined context for use by step checking routine
856: Collective on SNES
858: Calling sequence of func:
859: .vb
860: int func (SNES snes, void *checkctx, Vec x, PetscTruth *flag)
861: .ve
862: where func returns an error code of 0 on success and a nonzero
863: on failure.
865: Input parameters for func:
866: + snes - nonlinear context
867: . checkctx - optional user-defined context for use by step checking routine
868: - x - current candidate iterate
870: Output parameters for func:
871: + x - current iterate (possibly modified)
872: - flag - flag indicating whether x has been modified (either
873: PETSC_TRUE of PETSC_FALSE)
875: Level: advanced
877: Notes:
878: SNESNoLineSearch() and SNESNoLineSearchNoNorms() accept the new
879: iterate computed by the line search checking routine. In particular,
880: these routines (1) compute a candidate iterate u_{i+1}, (2) pass control
881: to the checking routine, and then (3) compute the corresponding nonlinear
882: function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}.
884: SNESQuadraticLineSearch() and SNESCubicLineSearch() also accept the
885: new iterate computed by the line search checking routine. In particular,
886: these routines (1) compute a candidate iterate u_{i+1} as well as a
887: candidate nonlinear function f(u_{i+1}), (2) pass control to the checking
888: routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes
889: were made to the candidate iterate in the checking routine (as indicated
890: by flag=PETSC_TRUE). The overhead of this function re-evaluation can be
891: very costly, so use this feature with caution!
893: .keywords: SNES, nonlinear, set, line search check, step check, routine
895: .seealso: SNESSetLineSearch()
896: @*/
897: int SNESSetLineSearchCheck(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx)
898: {
899: int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,PetscTruth*),void*);
902: PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearchCheck_C",(void (**)(void))&f);
903: if (f) {
904: (*f)(snes,func,checkctx);
905: }
906: return(0);
907: }
908: /* -------------------------------------------------------------------------- */
909: EXTERN_C_BEGIN
910: int SNESSetLineSearchCheck_LS(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx)
911: {
913: ((SNES_LS *)(snes->data))->CheckStep = func;
914: ((SNES_LS *)(snes->data))->checkP = checkctx;
915: return(0);
916: }
917: EXTERN_C_END
918: /* -------------------------------------------------------------------------- */
919: /*
920: SNESPrintHelp_LS - Prints all options for the SNES_LS method.
922: Input Parameter:
923: . snes - the SNES context
925: Application Interface Routine: SNESPrintHelp()
926: */
927: static int SNESPrintHelp_LS(SNES snes,char *p)
928: {
929: SNES_LS *ls = (SNES_LS *)snes->data;
932: (*PetscHelpPrintf)(snes->comm," method SNES_LS (ls) for systems of nonlinear equations:n");
933: (*PetscHelpPrintf)(snes->comm," %ssnes_ls [cubic,quadratic,basic,basicnonorms,...]n",p);
934: (*PetscHelpPrintf)(snes->comm," %ssnes_ls_alpha <alpha> (default %g)n",p,ls->alpha);
935: (*PetscHelpPrintf)(snes->comm," %ssnes_ls_maxstep <max> (default %g)n",p,ls->maxstep);
936: (*PetscHelpPrintf)(snes->comm," %ssnes_ls_steptol <tol> (default %g)n",p,ls->steptol);
937: return(0);
938: }
940: /*
941: SNESView_LS - Prints info from the SNESLS data structure.
943: Input Parameters:
944: . SNES - the SNES context
945: . viewer - visualization context
947: Application Interface Routine: SNESView()
948: */
949: static int SNESView_LS(SNES snes,PetscViewer viewer)
950: {
951: SNES_LS *ls = (SNES_LS *)snes->data;
952: char *cstr;
953: int ierr;
954: PetscTruth isascii;
957: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
958: if (isascii) {
959: if (ls->LineSearch == SNESNoLineSearch) cstr = "SNESNoLineSearch";
960: else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch";
961: else if (ls->LineSearch == SNESCubicLineSearch) cstr = "SNESCubicLineSearch";
962: else cstr = "unknown";
963: PetscViewerASCIIPrintf(viewer," line search variant: %sn",cstr);
964: PetscViewerASCIIPrintf(viewer," alpha=%g, maxstep=%g, steptol=%gn",ls->alpha,ls->maxstep,ls->steptol);
965: } else {
966: SETERRQ1(1,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name);
967: }
968: return(0);
969: }
970: /* -------------------------------------------------------------------------- */
971: /*
972: SNESSetFromOptions_LS - Sets various parameters for the SNESLS method.
974: Input Parameter:
975: . snes - the SNES context
977: Application Interface Routine: SNESSetFromOptions()
978: */
979: static int SNESSetFromOptions_LS(SNES snes)
980: {
981: SNES_LS *ls = (SNES_LS *)snes->data;
982: char ver[16],*lses[] = {"basic","basicnonorms","quadratic","cubic"};
983: int ierr;
984: PetscTruth flg;
987: PetscOptionsHead("SNES Line search options");
988: PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);
989: PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);
990: PetscOptionsReal("-snes_ls_steptol","Step must be greater than","None",ls->steptol,&ls->steptol,0);
992: PetscOptionsEList("-snes_ls","Line search used","SNESSetLineSearch",lses,4,"cubic",ver,16,&flg);
993: if (flg) {
994: PetscTruth isbasic,isnonorms,isquad,iscubic;
996: PetscStrcmp(ver,lses[0],&isbasic);
997: PetscStrcmp(ver,lses[1],&isnonorms);
998: PetscStrcmp(ver,lses[2],&isquad);
999: PetscStrcmp(ver,lses[3],&iscubic);
1001: if (isbasic) {
1002: SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);
1003: } else if (isnonorms) {
1004: SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);
1005: } else if (isquad) {
1006: SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);
1007: } else if (iscubic) {
1008: SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);
1009: }
1010: else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Unknown line search");}
1011: }
1012: PetscOptionsTail();
1013: return(0);
1014: }
1015: /* -------------------------------------------------------------------------- */
1016: /*
1017: SNESCreate_LS - Creates a nonlinear solver context for the SNESLS method,
1018: SNES_LS, and sets this as the private data within the generic nonlinear solver
1019: context, SNES, that was created within SNESCreate().
1022: Input Parameter:
1023: . snes - the SNES context
1025: Application Interface Routine: SNESCreate()
1026: */
1027: EXTERN_C_BEGIN
1028: int SNESCreate_LS(SNES snes)
1029: {
1030: int ierr;
1031: SNES_LS *neP;
1034: snes->setup = SNESSetUp_LS;
1035: snes->solve = SNESSolve_LS;
1036: snes->destroy = SNESDestroy_LS;
1037: snes->converged = SNESConverged_LS;
1038: snes->printhelp = SNESPrintHelp_LS;
1039: snes->setfromoptions = SNESSetFromOptions_LS;
1040: snes->view = SNESView_LS;
1041: snes->nwork = 0;
1043: ierr = PetscNew(SNES_LS,&neP);
1044: PetscLogObjectMemory(snes,sizeof(SNES_LS));
1045: snes->data = (void*)neP;
1046: neP->alpha = 1.e-4;
1047: neP->maxstep = 1.e8;
1048: neP->steptol = 1.e-12;
1049: neP->LineSearch = SNESCubicLineSearch;
1050: neP->lsP = PETSC_NULL;
1051: neP->CheckStep = PETSC_NULL;
1052: neP->checkP = PETSC_NULL;
1054: PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS",SNESSetLineSearch_LS);
1055: PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS",SNESSetLineSearchCheck_LS);
1057: return(0);
1058: }
1059: EXTERN_C_END