Actual source code: ls.c
2: #include src/snes/impls/ls/ls.h
4: /*
5: Checks if J^T F = 0 which implies we've found a local minimum of the function,
6: but not a zero. In the case when one cannot compute J^T F we use the fact that
7: 0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
8: for this trick.
9: */
12: PetscErrorCode SNESLSCheckLocalMin_Private(Mat A,Vec F,Vec W,PetscReal fnorm,PetscTruth *ismin)
13: {
14: PetscReal a1;
16: PetscTruth hastranspose;
19: *ismin = PETSC_FALSE;
20: MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);
21: if (hastranspose) {
22: /* Compute || J^T F|| */
23: MatMultTranspose(A,F,W);
24: VecNorm(W,NORM_2,&a1);
25: PetscLogInfo(0,"SNESSolve_LS: || J^T F|| %g near zero implies found a local minimum\n",a1/fnorm);
26: if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
27: } else {
28: Vec work;
29: PetscScalar result;
30: PetscReal wnorm;
32: VecSetRandom(PETSC_NULL,W);
33: VecNorm(W,NORM_2,&wnorm);
34: VecDuplicate(W,&work);
35: MatMult(A,W,work);
36: VecDot(F,work,&result);
37: VecDestroy(work);
38: a1 = PetscAbsScalar(result)/(fnorm*wnorm);
39: PetscLogInfo(0,"SNESSolve_LS: (F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n",a1);
40: if (a1 < 1.e-4) *ismin = PETSC_TRUE;
41: }
42: return(0);
43: }
45: /*
46: Checks if J^T(F - J*X) = 0
47: */
50: PetscErrorCode SNESLSCheckResidual_Private(Mat A,Vec F,Vec X,Vec W1,Vec W2)
51: {
52: PetscReal a1,a2;
54: PetscTruth hastranspose;
55: PetscScalar mone = -1.0;
58: MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);
59: if (hastranspose) {
60: MatMult(A,X,W1);
61: VecAXPY(&mone,F,W1);
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 != 0) {
68: PetscLogInfo(0,"SNESSolve_LS: ||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;
136: PetscInt maxits,i,lits;
137: PetscTruth lssucceed;
138: MatStructure flg = DIFFERENT_NONZERO_PATTERN;
139: PetscReal fnorm,gnorm,xnorm,ynorm;
140: Vec Y,X,F,G,W,TMP;
141: KSP ksp;
144: SNESGetKSP(snes,&ksp);
145: snes->reason = SNES_CONVERGED_ITERATING;
147: maxits = snes->max_its; /* maximum number of iterations */
148: X = snes->vec_sol; /* solution vector */
149: F = snes->vec_func; /* residual vector */
150: Y = snes->work[0]; /* work vectors */
151: G = snes->work[1];
152: W = snes->work[2];
154: PetscObjectTakeAccess(snes);
155: snes->iter = 0;
156: PetscObjectGrantAccess(snes);
157: SNESComputeFunction(snes,X,F); /* F(X) */
158: VecNorm(F,NORM_2,&fnorm); /* fnorm <- ||F|| */
159: PetscObjectTakeAccess(snes);
160: snes->norm = fnorm;
161: PetscObjectGrantAccess(snes);
162: SNESLogConvHistory(snes,fnorm,0);
163: SNESMonitor(snes,0,fnorm);
165: if (fnorm < snes->abstol) {snes->reason = SNES_CONVERGED_FNORM_ABS; return(0);}
167: /* set parameter for default relative tolerance convergence test */
168: snes->ttol = fnorm*snes->rtol;
170: for (i=0; i<maxits; i++) {
172: /* Call general purpose update function */
173: if (snes->update != PETSC_NULL) {
174: (*snes->update)(snes, snes->iter);
175: }
177: /* Solve J Y = F, where J is Jacobian matrix */
178: SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);
179: KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);
180: KSPSolve(snes->ksp,F,Y);
181: KSPGetIterationNumber(ksp,&lits);
183: if (PetscLogPrintInfo){
184: SNESLSCheckResidual_Private(snes->jacobian,F,Y,G,W);
185: }
187: /* should check what happened to the linear solve? */
188: snes->linear_its += lits;
189: PetscLogInfo(snes,"SNESSolve_LS: iter=%D, linear solve iterations=%D\n",snes->iter,lits);
191: /* Compute a (scaled) negative update in the line search routine:
192: Y <- X - lambda*Y
193: and evaluate G(Y) = function(Y))
194: */
195: VecCopy(Y,snes->vec_sol_update_always);
196: (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lssucceed);
197: PetscLogInfo(snes,"SNESSolve_LS: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",fnorm,gnorm,ynorm,(int)lssucceed);
199: TMP = F; F = G; snes->vec_func_always = F; G = TMP;
200: TMP = X; X = Y; snes->vec_sol_always = X; Y = TMP;
201: fnorm = gnorm;
203: PetscObjectTakeAccess(snes);
204: snes->iter = i+1;
205: snes->norm = fnorm;
206: PetscObjectGrantAccess(snes);
207: SNESLogConvHistory(snes,fnorm,lits);
208: SNESMonitor(snes,i+1,fnorm);
210: if (!lssucceed) {
211: PetscTruth ismin;
213: if (++snes->numFailures >= snes->maxFailures) {
214: snes->reason = SNES_DIVERGED_LS_FAILURE;
215: SNESLSCheckLocalMin_Private(snes->jacobian,F,W,fnorm,&ismin);
216: if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
217: break;
218: }
219: }
221: /* Test for convergence */
222: if (snes->converged) {
223: VecNorm(X,NORM_2,&xnorm); /* xnorm = || X || */
224: (*snes->converged)(snes,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);
225: if (snes->reason) {
226: break;
227: }
228: }
229: }
230: if (X != snes->vec_sol) {
231: VecCopy(X,snes->vec_sol);
232: }
233: if (F != snes->vec_func) {
234: VecCopy(F,snes->vec_func);
235: }
236: snes->vec_sol_always = snes->vec_sol;
237: snes->vec_func_always = snes->vec_func;
238: if (i == maxits) {
239: PetscLogInfo(snes,"SNESSolve_LS: Maximum number of iterations has been reached: %D\n",maxits);
240: snes->reason = SNES_DIVERGED_MAX_IT;
241: }
242: return(0);
243: }
244: /* -------------------------------------------------------------------------- */
245: /*
246: SNESSetUp_LS - Sets up the internal data structures for the later use
247: of the SNESLS nonlinear solver.
249: Input Parameter:
250: . snes - the SNES context
251: . x - the solution vector
253: Application Interface Routine: SNESSetUp()
255: Notes:
256: For basic use of the SNES solvers, the user need not explicitly call
257: SNESSetUp(), since these actions will automatically occur during
258: the call to SNESSolve().
259: */
262: PetscErrorCode SNESSetUp_LS(SNES snes)
263: {
267: snes->nwork = 4;
268: VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);
269: PetscLogObjectParents(snes,snes->nwork,snes->work);
270: snes->vec_sol_update_always = snes->work[3];
271: return(0);
272: }
273: /* -------------------------------------------------------------------------- */
274: /*
275: SNESDestroy_LS - Destroys the private SNES_LS context that was created
276: with SNESCreate_LS().
278: Input Parameter:
279: . snes - the SNES context
281: Application Interface Routine: SNESDestroy()
282: */
285: PetscErrorCode SNESDestroy_LS(SNES snes)
286: {
290: if (snes->nwork) {
291: VecDestroyVecs(snes->work,snes->nwork);
292: }
293: PetscFree(snes->data);
294: return(0);
295: }
296: /* -------------------------------------------------------------------------- */
300: /*@C
301: SNESNoLineSearch - This routine is not a line search at all;
302: it simply uses the full Newton step. Thus, this routine is intended
303: to serve as a template and is not recommended for general use.
305: Collective on SNES and Vec
307: Input Parameters:
308: + snes - nonlinear context
309: . lsctx - optional context for line search (not used here)
310: . x - current iterate
311: . f - residual evaluated at x
312: . y - search direction (contains new iterate on output)
313: . w - work vector
314: - fnorm - 2-norm of f
316: Output Parameters:
317: + g - residual evaluated at new iterate y
318: . y - new iterate (contains search direction on input)
319: . gnorm - 2-norm of g
320: . ynorm - 2-norm of search length
321: - flag - PETSC_TRUE on success, PETSC_FALSE on failure
323: Options Database Key:
324: . -snes_ls basic - Activates SNESNoLineSearch()
326: Level: advanced
328: .keywords: SNES, nonlinear, line search, cubic
330: .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
331: SNESSetLineSearch(), SNESNoLineSearchNoNorms()
332: @*/
333: PetscErrorCode SNESNoLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
334: {
336: PetscScalar mone = -1.0;
337: SNES_LS *neP = (SNES_LS*)snes->data;
338: PetscTruth change_y = PETSC_FALSE;
341: *flag = PETSC_TRUE;
342: PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);
343: VecNorm(y,NORM_2,ynorm); /* ynorm = || y || */
344: VecAYPX(&mone,x,y); /* y <- y - x */
345: if (neP->CheckStep) {
346: (*neP->CheckStep)(snes,neP->checkP,y,&change_y);
347: }
348: SNESComputeFunction(snes,y,g); /* Compute F(y) */
349: VecNorm(g,NORM_2,gnorm); /* gnorm = || g || */
350: PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
351: return(0);
352: }
353: /* -------------------------------------------------------------------------- */
357: /*@C
358: SNESNoLineSearchNoNorms - This routine is not a line search at
359: all; it simply uses the full Newton step. This version does not
360: even compute the norm of the function or search direction; this
361: is intended only when you know the full step is fine and are
362: not checking for convergence of the nonlinear iteration (for
363: example, you are running always for a fixed number of Newton steps).
365: Collective on SNES and Vec
367: Input Parameters:
368: + snes - nonlinear context
369: . lsctx - optional context for line search (not used here)
370: . x - current iterate
371: . f - residual evaluated at x
372: . y - search direction (contains new iterate on output)
373: . w - work vector
374: - fnorm - 2-norm of f
376: Output Parameters:
377: + g - residual evaluated at new iterate y
378: . gnorm - not changed
379: . ynorm - not changed
380: - flag - set to PETSC_TRUE indicating a successful line search
382: Options Database Key:
383: . -snes_ls basicnonorms - Activates SNESNoLineSearchNoNorms()
385: Notes:
386: SNESNoLineSearchNoNorms() must be used in conjunction with
387: either the options
388: $ -snes_no_convergence_test -snes_max_it <its>
389: or alternatively a user-defined custom test set via
390: SNESSetConvergenceTest(); or a -snes_max_it of 1,
391: otherwise, the SNES solver will generate an error.
393: During the final iteration this will not evaluate the function at
394: the solution point. This is to save a function evaluation while
395: using pseudo-timestepping.
397: The residual norms printed by monitoring routines such as
398: SNESDefaultMonitor() (as activated via -snes_monitor) will not be
399: correct, since they are not computed.
401: Level: advanced
403: .keywords: SNES, nonlinear, line search, cubic
405: .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
406: SNESSetLineSearch(), SNESNoLineSearch()
407: @*/
408: PetscErrorCode SNESNoLineSearchNoNorms(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
409: {
411: PetscScalar mone = -1.0;
412: SNES_LS *neP = (SNES_LS*)snes->data;
413: PetscTruth change_y = PETSC_FALSE;
416: *flag = PETSC_TRUE;
417: PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);
418: VecAYPX(&mone,x,y); /* y <- y - x */
419: if (neP->CheckStep) {
420: (*neP->CheckStep)(snes,neP->checkP,y,&change_y);
421: }
422:
423: /* don't evaluate function the last time through */
424: if (snes->iter < snes->max_its-1) {
425: SNESComputeFunction(snes,y,g); /* Compute F(y) */
426: }
427: PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
428: return(0);
429: }
430: /* -------------------------------------------------------------------------- */
433: /*@C
434: SNESCubicLineSearch - Performs a cubic line search (default line search method).
436: Collective on SNES
438: Input Parameters:
439: + snes - nonlinear context
440: . lsctx - optional context for line search (not used here)
441: . x - current iterate
442: . f - residual evaluated at x
443: . y - search direction (contains new iterate on output)
444: . w - work vector
445: - fnorm - 2-norm of f
447: Output Parameters:
448: + g - residual evaluated at new iterate y
449: . y - new iterate (contains search direction on input)
450: . gnorm - 2-norm of g
451: . ynorm - 2-norm of search length
452: - flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
454: Options Database Key:
455: $ -snes_ls cubic - Activates SNESCubicLineSearch()
457: Notes:
458: This line search is taken from "Numerical Methods for Unconstrained
459: Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
461: Level: advanced
463: .keywords: SNES, nonlinear, line search, cubic
465: .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
466: @*/
467: PetscErrorCode SNESCubicLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
468: {
469: /*
470: Note that for line search purposes we work with with the related
471: minimization problem:
472: min z(x): R^n -> R,
473: where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
474: */
475:
476: PetscReal steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
477: PetscReal maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg;
478: #if defined(PETSC_USE_COMPLEX)
479: PetscScalar cinitslope,clambda;
480: #endif
482: PetscInt count;
483: SNES_LS *neP = (SNES_LS*)snes->data;
484: PetscScalar mone = -1.0,scale;
485: PetscTruth change_y = PETSC_FALSE;
488: PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);
489: *flag = PETSC_TRUE;
490: alpha = neP->alpha;
491: maxstep = neP->maxstep;
492: steptol = neP->steptol;
494: VecNorm(y,NORM_2,ynorm);
495: if (*ynorm == 0.0) {
496: PetscLogInfo(snes,"SNESCubicLineSearch: Search direction and size is 0\n");
497: *gnorm = fnorm;
498: VecCopy(x,y);
499: VecCopy(f,g);
500: *flag = PETSC_FALSE;
501: goto theend1;
502: }
503: if (*ynorm > maxstep) { /* Step too big, so scale back */
504: scale = maxstep/(*ynorm);
505: #if defined(PETSC_USE_COMPLEX)
506: PetscLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",PetscRealPart(scale),*ynorm);
507: #else
508: PetscLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",scale,*ynorm);
509: #endif
510: VecScale(&scale,y);
511: *ynorm = maxstep;
512: }
513: VecMaxPointwiseDivide(y,x,&rellength);
514: minlambda = steptol/rellength;
515: MatMult(snes->jacobian,y,w);
516: #if defined(PETSC_USE_COMPLEX)
517: VecDot(f,w,&cinitslope);
518: initslope = PetscRealPart(cinitslope);
519: #else
520: VecDot(f,w,&initslope);
521: #endif
522: if (initslope > 0.0) initslope = -initslope;
523: if (initslope == 0.0) initslope = -1.0;
525: VecCopy(y,w);
526: VecAYPX(&mone,x,w);
527: SNESComputeFunction(snes,w,g);
528: VecNorm(g,NORM_2,gnorm);
529: if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
530: VecCopy(w,y);
531: PetscLogInfo(snes,"SNESCubicLineSearch: Using full step\n");
532: goto theend1;
533: }
535: /* Fit points with quadratic */
536: lambda = 1.0;
537: lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
538: lambdaprev = lambda;
539: gnormprev = *gnorm;
540: if (lambdatemp > .5*lambda) lambdatemp = .5*lambda;
541: if (lambdatemp <= .1*lambda) lambda = .1*lambda;
542: else lambda = lambdatemp;
543: VecCopy(x,w);
544: lambdaneg = -lambda;
545: #if defined(PETSC_USE_COMPLEX)
546: clambda = lambdaneg; VecAXPY(&clambda,y,w);
547: #else
548: VecAXPY(&lambdaneg,y,w);
549: #endif
550: SNESComputeFunction(snes,w,g);
551: VecNorm(g,NORM_2,gnorm);
552: if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
553: VecCopy(w,y);
554: PetscLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%18.16e\n",lambda);
555: goto theend1;
556: }
558: /* Fit points with cubic */
559: count = 1;
560: while (count < 10000) {
561: if (lambda <= minlambda) { /* bad luck; use full step */
562: PetscLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %D \n",count);
563: PetscLogInfo(snes,"SNESCubicLineSearch:fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope);
564: VecCopy(x,y);
565: *flag = PETSC_FALSE; break;
566: }
567: t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
568: t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope;
569: a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
570: b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
571: d = b*b - 3*a*initslope;
572: if (d < 0.0) d = 0.0;
573: if (a == 0.0) {
574: lambdatemp = -initslope/(2.0*b);
575: } else {
576: lambdatemp = (-b + sqrt(d))/(3.0*a);
577: }
578: lambdaprev = lambda;
579: gnormprev = *gnorm;
580: if (lambdatemp > .5*lambda) lambdatemp = .5*lambda;
581: if (lambdatemp <= .1*lambda) lambda = .1*lambda;
582: else lambda = lambdatemp;
583: VecCopy(x,w);
584: lambdaneg = -lambda;
585: #if defined(PETSC_USE_COMPLEX)
586: clambda = lambdaneg;
587: VecAXPY(&clambda,y,w);
588: #else
589: VecAXPY(&lambdaneg,y,w);
590: #endif
591: SNESComputeFunction(snes,w,g);
592: VecNorm(g,NORM_2,gnorm);
593: if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */
594: VecCopy(w,y);
595: PetscLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%18.16e\n",lambda);
596: goto theend1;
597: } else {
598: PetscLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda, lambda=%18.16e\n",lambda);
599: }
600: count++;
601: }
602: if (count >= 10000) {
603: SETERRQ(PETSC_ERR_LIB, "Lambda was decreased more than 10,000 times, so something is probably wrong with the function evaluation");
604: }
605: theend1:
606: /* Optional user-defined check for line search step validity */
607: if (neP->CheckStep) {
608: (*neP->CheckStep)(snes,neP->checkP,y,&change_y);
609: if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */
610: SNESComputeFunction(snes,y,g);
611: VecNormBegin(y,NORM_2,ynorm);
612: VecNormBegin(g,NORM_2,gnorm);
613: VecNormEnd(y,NORM_2,ynorm);
614: VecNormEnd(g,NORM_2,gnorm);
615: }
616: }
617: PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
618: return(0);
619: }
620: /* -------------------------------------------------------------------------- */
623: /*@C
624: SNESQuadraticLineSearch - Performs a quadratic line search.
626: Collective on SNES and Vec
628: Input Parameters:
629: + snes - the SNES context
630: . lsctx - optional context for line search (not used here)
631: . x - current iterate
632: . f - residual evaluated at x
633: . y - search direction (contains new iterate on output)
634: . w - work vector
635: - fnorm - 2-norm of f
637: Output Parameters:
638: + g - residual evaluated at new iterate y
639: . y - new iterate (contains search direction on input)
640: . gnorm - 2-norm of g
641: . ynorm - 2-norm of search length
642: - flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
644: Options Database Key:
645: . -snes_ls quadratic - Activates SNESQuadraticLineSearch()
647: Notes:
648: Use SNESSetLineSearch() to set this routine within the SNESLS method.
650: Level: advanced
652: .keywords: SNES, nonlinear, quadratic, line search
654: .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
655: @*/
656: PetscErrorCode SNESQuadraticLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
657: {
658: /*
659: Note that for line search purposes we work with with the related
660: minimization problem:
661: min z(x): R^n -> R,
662: where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
663: */
664: PetscReal steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg,rellength;
665: #if defined(PETSC_USE_COMPLEX)
666: PetscScalar cinitslope,clambda;
667: #endif
669: PetscInt count;
670: SNES_LS *neP = (SNES_LS*)snes->data;
671: PetscScalar mone = -1.0,scale;
672: PetscTruth change_y = PETSC_FALSE;
675: PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);
676: *flag = PETSC_TRUE;
677: alpha = neP->alpha;
678: maxstep = neP->maxstep;
679: steptol = neP->steptol;
681: VecNorm(y,NORM_2,ynorm);
682: if (*ynorm == 0.0) {
683: PetscLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n");
684: *gnorm = fnorm;
685: VecCopy(x,y);
686: VecCopy(f,g);
687: *flag = PETSC_FALSE;
688: goto theend2;
689: }
690: if (*ynorm > maxstep) { /* Step too big, so scale back */
691: scale = maxstep/(*ynorm);
692: VecScale(&scale,y);
693: *ynorm = maxstep;
694: }
695: VecMaxPointwiseDivide(y,x,&rellength);
696: minlambda = steptol/rellength;
697: MatMult(snes->jacobian,y,w);
698: #if defined(PETSC_USE_COMPLEX)
699: VecDot(f,w,&cinitslope);
700: initslope = PetscRealPart(cinitslope);
701: #else
702: VecDot(f,w,&initslope);
703: #endif
704: if (initslope > 0.0) initslope = -initslope;
705: if (initslope == 0.0) initslope = -1.0;
707: VecCopy(y,w);
708: VecAYPX(&mone,x,w);
709: SNESComputeFunction(snes,w,g);
710: VecNorm(g,NORM_2,gnorm);
711: if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
712: VecCopy(w,y);
713: PetscLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n");
714: goto theend2;
715: }
717: /* Fit points with quadratic */
718: lambda = 1.0;
719: count = 1;
720: while (PETSC_TRUE) {
721: if (lambda <= minlambda) { /* bad luck; use full step */
722: PetscLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %D \n",count);
723: PetscLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",fnorm,*gnorm,*ynorm,lambda,initslope);
724: VecCopy(x,y);
725: *flag = PETSC_FALSE; break;
726: }
727: lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
728: if (lambdatemp > .5*lambda) lambdatemp = .5*lambda;
729: if (lambdatemp <= .1*lambda) lambda = .1*lambda;
730: else lambda = lambdatemp;
731: VecCopy(x,w);
732: lambdaneg = -lambda;
733: #if defined(PETSC_USE_COMPLEX)
734: clambda = lambdaneg; VecAXPY(&clambda,y,w);
735: #else
736: VecAXPY(&lambdaneg,y,w);
737: #endif
738: SNESComputeFunction(snes,w,g);
739: VecNorm(g,NORM_2,gnorm);
740: if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
741: VecCopy(w,y);
742: PetscLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda);
743: break;
744: }
745: count++;
746: }
747: theend2:
748: /* Optional user-defined check for line search step validity */
749: if (neP->CheckStep) {
750: (*neP->CheckStep)(snes,neP->checkP,y,&change_y);
751: if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */
752: SNESComputeFunction(snes,y,g);
753: VecNormBegin(y,NORM_2,ynorm);
754: VecNormBegin(g,NORM_2,gnorm);
755: VecNormEnd(y,NORM_2,ynorm);
756: VecNormEnd(g,NORM_2,gnorm);
757: }
758: }
759: PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
760: return(0);
761: }
763: /* -------------------------------------------------------------------------- */
766: /*@C
767: SNESSetLineSearch - Sets the line search routine to be used
768: by the method SNESLS.
770: Input Parameters:
771: + snes - nonlinear context obtained from SNESCreate()
772: . lsctx - optional user-defined context for use by line search
773: - func - pointer to int function
775: Collective on SNES
777: Available Routines:
778: + SNESCubicLineSearch() - default line search
779: . SNESQuadraticLineSearch() - quadratic line search
780: . SNESNoLineSearch() - the full Newton step (actually not a line search)
781: - SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel)
783: Options Database Keys:
784: + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
785: . -snes_ls_alpha <alpha> - Sets alpha
786: . -snes_ls_maxstep <max> - Sets maxstep
787: - -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code
788: will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length
789: the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x)
791: Calling sequence of func:
792: .vb
793: func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
794: PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
795: .ve
797: Input parameters for func:
798: + snes - nonlinear context
799: . lsctx - optional user-defined context for line search
800: . x - current iterate
801: . f - residual evaluated at x
802: . y - search direction (contains new iterate on output)
803: . w - work vector
804: - fnorm - 2-norm of f
806: Output parameters for func:
807: + g - residual evaluated at new iterate y
808: . y - new iterate (contains search direction on input)
809: . gnorm - 2-norm of g
810: . ynorm - 2-norm of search length
811: - flag - set to PETSC_TRUE if the line search succeeds; PETSC_FALSE on failure.
813: Level: advanced
815: .keywords: SNES, nonlinear, set, line search, routine
817: .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(),
818: SNESSetLineSearchCheck(), SNESSetLineSearchParams(), SNESGetLineSearchParams()
819: @*/
820: PetscErrorCode SNESSetLineSearch(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void *lsctx)
821: {
822: PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void*);
825: PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void (**)(void))&f);
826: if (f) {
827: (*f)(snes,func,lsctx);
828: }
829: return(0);
830: }
833: /* -------------------------------------------------------------------------- */
837: PetscErrorCode SNESSetLineSearch_LS(SNES snes,FCN2 func,void *lsctx)
838: {
840: ((SNES_LS *)(snes->data))->LineSearch = func;
841: ((SNES_LS *)(snes->data))->lsP = lsctx;
842: return(0);
843: }
845: /* -------------------------------------------------------------------------- */
848: /*@C
849: SNESSetLineSearchCheck - Sets a routine to check the validity of new iterate computed
850: by the line search routine in the Newton-based method SNESLS.
852: Input Parameters:
853: + snes - nonlinear context obtained from SNESCreate()
854: . func - pointer to int function
855: - checkctx - optional user-defined context for use by step checking routine
857: Collective on SNES
859: Calling sequence of func:
860: .vb
861: int func (SNES snes, void *checkctx, Vec x, PetscTruth *flag)
862: .ve
863: where func returns an error code of 0 on success and a nonzero
864: on failure.
866: Input parameters for func:
867: + snes - nonlinear context
868: . checkctx - optional user-defined context for use by step checking routine
869: - x - current candidate iterate
871: Output parameters for func:
872: + x - current iterate (possibly modified)
873: - flag - flag indicating whether x has been modified (either
874: PETSC_TRUE of PETSC_FALSE)
876: Level: advanced
878: Notes:
879: SNESNoLineSearch() and SNESNoLineSearchNoNorms() accept the new
880: iterate computed by the line search checking routine. In particular,
881: these routines (1) compute a candidate iterate u_{i+1}, (2) pass control
882: to the checking routine, and then (3) compute the corresponding nonlinear
883: function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}.
885: SNESQuadraticLineSearch() and SNESCubicLineSearch() also accept the
886: new iterate computed by the line search checking routine. In particular,
887: these routines (1) compute a candidate iterate u_{i+1} as well as a
888: candidate nonlinear function f(u_{i+1}), (2) pass control to the checking
889: routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes
890: were made to the candidate iterate in the checking routine (as indicated
891: by flag=PETSC_TRUE). The overhead of this function re-evaluation can be
892: very costly, so use this feature with caution!
894: .keywords: SNES, nonlinear, set, line search check, step check, routine
896: .seealso: SNESSetLineSearch()
897: @*/
898: PetscErrorCode SNESSetLineSearchCheck(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx)
899: {
900: PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,void*,Vec,PetscTruth*),void*);
903: PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearchCheck_C",(void (**)(void))&f);
904: if (f) {
905: (*f)(snes,func,checkctx);
906: }
907: return(0);
908: }
909: /* -------------------------------------------------------------------------- */
914: PetscErrorCode SNESSetLineSearchCheck_LS(SNES snes,FCN func,void *checkctx)
915: {
917: ((SNES_LS *)(snes->data))->CheckStep = func;
918: ((SNES_LS *)(snes->data))->checkP = checkctx;
919: return(0);
920: }
922: /* -------------------------------------------------------------------------- */
923: /*
924: SNESPrintHelp_LS - Prints all options for the SNES_LS method.
926: Input Parameter:
927: . snes - the SNES context
929: Application Interface Routine: SNESPrintHelp()
930: */
933: static PetscErrorCode SNESPrintHelp_LS(SNES snes,char *p)
934: {
935: SNES_LS *ls = (SNES_LS *)snes->data;
938: (*PetscHelpPrintf)(snes->comm," method SNES_LS (ls) for systems of nonlinear equations:\n");
939: (*PetscHelpPrintf)(snes->comm," %ssnes_ls [cubic,quadratic,basic,basicnonorms,...]\n",p);
940: (*PetscHelpPrintf)(snes->comm," %ssnes_ls_alpha <alpha> (default %g)\n",p,ls->alpha);
941: (*PetscHelpPrintf)(snes->comm," %ssnes_ls_maxstep <max> (default %g)\n",p,ls->maxstep);
942: (*PetscHelpPrintf)(snes->comm," %ssnes_ls_steptol <tol> (default %g)\n",p,ls->steptol);
943: return(0);
944: }
946: /*
947: SNESView_LS - Prints info from the SNESLS data structure.
949: Input Parameters:
950: . SNES - the SNES context
951: . viewer - visualization context
953: Application Interface Routine: SNESView()
954: */
957: static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer)
958: {
959: SNES_LS *ls = (SNES_LS *)snes->data;
960: const char *cstr;
962: PetscTruth iascii;
965: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
966: if (iascii) {
967: if (ls->LineSearch == SNESNoLineSearch) cstr = "SNESNoLineSearch";
968: else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch";
969: else if (ls->LineSearch == SNESCubicLineSearch) cstr = "SNESCubicLineSearch";
970: else cstr = "unknown";
971: PetscViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);
972: PetscViewerASCIIPrintf(viewer," alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);
973: } else {
974: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name);
975: }
976: return(0);
977: }
978: /* -------------------------------------------------------------------------- */
979: /*
980: SNESSetFromOptions_LS - Sets various parameters for the SNESLS method.
982: Input Parameter:
983: . snes - the SNES context
985: Application Interface Routine: SNESSetFromOptions()
986: */
989: static PetscErrorCode SNESSetFromOptions_LS(SNES snes)
990: {
991: SNES_LS *ls = (SNES_LS *)snes->data;
992: const char *lses[] = {"basic","basicnonorms","quadratic","cubic"};
994: PetscInt indx;
995: PetscTruth flg;
998: PetscOptionsHead("SNES Line search options");
999: PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);
1000: PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);
1001: PetscOptionsReal("-snes_ls_steptol","Step must be greater than","None",ls->steptol,&ls->steptol,0);
1003: PetscOptionsEList("-snes_ls","Line search used","SNESSetLineSearch",lses,4,"cubic",&indx,&flg);
1004: if (flg) {
1005: switch (indx) {
1006: case 0:
1007: SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);
1008: break;
1009: case 1:
1010: SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);
1011: break;
1012: case 2:
1013: SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);
1014: break;
1015: case 3:
1016: SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);
1017: break;
1018: }
1019: }
1020: PetscOptionsTail();
1021: return(0);
1022: }
1023: /* -------------------------------------------------------------------------- */
1024: /*MC
1025: SNESLS - Newton based nonlinear solver that uses a line search
1027: Options Database:
1028: + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
1029: . -snes_ls_alpha <alpha> - Sets alpha
1030: . -snes_ls_maxstep <max> - Sets maxstep
1031: - -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code
1032: will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length
1033: the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x)
1035: Notes: This is the default nonlinear solver in SNES
1037: Level: beginner
1039: .seealso: SNESCreate(), SNES, SNESSetType(), SNESTR, SNESSetLineSearch(),
1040: SNESSetLineSearchCheck(), SNESNoLineSearch(), SNESCubicLineSearch(), SNESQuadraticLineSearch(),
1041: SNESSetLineSearch(), SNESNoLineSearchNoNorms()
1043: M*/
1047: PetscErrorCode SNESCreate_LS(SNES snes)
1048: {
1050: SNES_LS *neP;
1053: snes->setup = SNESSetUp_LS;
1054: snes->solve = SNESSolve_LS;
1055: snes->destroy = SNESDestroy_LS;
1056: snes->converged = SNESConverged_LS;
1057: snes->printhelp = SNESPrintHelp_LS;
1058: snes->setfromoptions = SNESSetFromOptions_LS;
1059: snes->view = SNESView_LS;
1060: snes->nwork = 0;
1062: PetscNew(SNES_LS,&neP);
1063: PetscLogObjectMemory(snes,sizeof(SNES_LS));
1064: snes->data = (void*)neP;
1065: neP->alpha = 1.e-4;
1066: neP->maxstep = 1.e8;
1067: neP->steptol = 1.e-12;
1068: neP->LineSearch = SNESCubicLineSearch;
1069: neP->lsP = PETSC_NULL;
1070: neP->CheckStep = PETSC_NULL;
1071: neP->checkP = PETSC_NULL;
1073: PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS",SNESSetLineSearch_LS);
1074: PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS",SNESSetLineSearchCheck_LS);
1076: return(0);
1077: }