Actual source code: snesncg.c
1: #include <../src/snes/impls/ncg/snesncgimpl.h>
2: const char *const SNESNCGTypes[] = {"FR","PRP","HS","DY","CD","SNESNCGType","SNES_NCG_",NULL};
4: static PetscErrorCode SNESReset_NCG(SNES snes)
5: {
7: return(0);
8: }
10: /*
11: SNESDestroy_NCG - Destroys the private SNES_NCG context that was created with SNESCreate_NCG().
13: Input Parameter:
14: . snes - the SNES context
16: Application Interface Routine: SNESDestroy()
17: */
18: static PetscErrorCode SNESDestroy_NCG(SNES snes)
19: {
23: PetscFree(snes->data);
24: return(0);
25: }
27: /*
28: SNESSetUp_NCG - Sets up the internal data structures for the later use
29: of the SNESNCG nonlinear solver.
31: Input Parameters:
32: + snes - the SNES context
33: - x - the solution vector
35: Application Interface Routine: SNESSetUp()
36: */
38: static PetscErrorCode SNESSetUp_NCG(SNES snes)
39: {
43: SNESSetWorkVecs(snes,2);
44: if (snes->npcside== PC_RIGHT) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNESNCG only supports left preconditioning");
45: if (snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED;
46: return(0);
47: }
49: static PetscErrorCode SNESLineSearchApply_NCGLinear(SNESLineSearch linesearch)
50: {
51: PetscScalar alpha, ptAp;
52: Vec X, Y, F, W;
53: SNES snes;
55: PetscReal *fnorm, *xnorm, *ynorm;
58: SNESLineSearchGetSNES(linesearch, &snes);
59: X = linesearch->vec_sol;
60: W = linesearch->vec_sol_new;
61: F = linesearch->vec_func;
62: Y = linesearch->vec_update;
63: fnorm = &linesearch->fnorm;
64: xnorm = &linesearch->xnorm;
65: ynorm = &linesearch->ynorm;
67: if (!snes->jacobian) {
68: SNESSetUpMatrices(snes);
69: }
71: /*
73: The exact step size for unpreconditioned linear CG is just:
74: alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy)
75: */
76: SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre);
77: VecDot(F, F, &alpha);
78: MatMult(snes->jacobian, Y, W);
79: VecDot(Y, W, &ptAp);
80: alpha = alpha / ptAp;
81: VecAXPY(X,-alpha,Y);
82: SNESComputeFunction(snes,X,F);
84: VecNorm(F,NORM_2,fnorm);
85: VecNorm(X,NORM_2,xnorm);
86: VecNorm(Y,NORM_2,ynorm);
87: return(0);
88: }
90: /*MC
91: SNESLINESEARCHNCGLINEAR - Special line search only for SNESNCG
93: This line search uses the length "as if" the problem is linear (that is what is computed by the linear CG method) using the Jacobian of the function.
94: alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy) where r (f) is the current residual (function value), p (y) is the current search direction.
96: Notes: This requires a Jacobian-vector product but does not require the solution of a linear system with the Jacobian
98: This is a "odd-ball" line search, we don't know if it is in the literature or used in practice by anyone.
100: Level: advanced
102: .seealso: SNESLineSearchCreate(), SNESLineSearchSetType()
103: M*/
105: PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch linesearch)
106: {
108: linesearch->ops->apply = SNESLineSearchApply_NCGLinear;
109: linesearch->ops->destroy = NULL;
110: linesearch->ops->setfromoptions = NULL;
111: linesearch->ops->reset = NULL;
112: linesearch->ops->view = NULL;
113: linesearch->ops->setup = NULL;
114: return(0);
115: }
117: /*
118: SNESSetFromOptions_NCG - Sets various parameters for the SNESNCG method.
120: Input Parameter:
121: . snes - the SNES context
123: Application Interface Routine: SNESSetFromOptions()
124: */
125: static PetscErrorCode SNESSetFromOptions_NCG(PetscOptionItems *PetscOptionsObject,SNES snes)
126: {
127: SNES_NCG *ncg = (SNES_NCG*)snes->data;
129: PetscBool debug = PETSC_FALSE;
130: SNESNCGType ncgtype=ncg->type;
131: SNESLineSearch linesearch;
134: PetscOptionsHead(PetscOptionsObject,"SNES NCG options");
135: PetscOptionsBool("-snes_ncg_monitor","Monitor the beta values used in the NCG iterations","SNES",ncg->monitor ? PETSC_TRUE : PETSC_FALSE, &debug, NULL);
136: if (debug) {
137: ncg->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
138: }
139: PetscOptionsEnum("-snes_ncg_type","NCG Beta type used","SNESNCGSetType",SNESNCGTypes,(PetscEnum)ncg->type,(PetscEnum*)&ncgtype,NULL);
140: SNESNCGSetType(snes, ncgtype);
141: PetscOptionsTail();
142: if (!snes->linesearch) {
143: SNESGetLineSearch(snes, &linesearch);
144: if (!((PetscObject)linesearch)->type_name) {
145: if (!snes->npc) {
146: SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);
147: } else {
148: SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);
149: }
150: }
151: }
152: return(0);
153: }
155: /*
156: SNESView_NCG - Prints info from the SNESNCG data structure.
158: Input Parameters:
159: + SNES - the SNES context
160: - viewer - visualization context
162: Application Interface Routine: SNESView()
163: */
164: static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer)
165: {
166: SNES_NCG *ncg = (SNES_NCG *) snes->data;
167: PetscBool iascii;
171: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);
172: if (iascii) {
173: PetscViewerASCIIPrintf(viewer, " type: %s\n", SNESNCGTypes[ncg->type]);
174: }
175: return(0);
176: }
178: /*
180: Assuming F = SNESComputeFunction(X) compute Y^tJ^tF using a simple secant approximation of the jacobian.
182: This routine is not currently used. I don't know what its intended purpose is.
184: Note it has a hardwired differencing parameter of 1e-5
186: */
187: PetscErrorCode SNESNCGComputeYtJtF_Private(SNES snes, Vec X, Vec F, Vec Y, Vec W, Vec G, PetscScalar * ytJtf)
188: {
190: PetscScalar ftf, ftg, fty, h;
193: VecDot(F, F, &ftf);
194: VecDot(F, Y, &fty);
195: h = 1e-5*fty / fty;
196: VecCopy(X, W);
197: VecAXPY(W, -h, Y); /* this is arbitrary */
198: SNESComputeFunction(snes, W, G);
199: VecDot(G, F, &ftg);
200: *ytJtf = (ftg - ftf) / h;
201: return(0);
202: }
204: /*@
205: SNESNCGSetType - Sets the conjugate update type for SNESNCG.
207: Logically Collective on SNES
209: Input Parameters:
210: + snes - the iterative context
211: - btype - update type
213: Options Database:
214: . -snes_ncg_type <prp,fr,hs,dy,cd> -strategy for selecting algorithm for computing beta
216: Level: intermediate
218: SNESNCGSelectTypes:
219: + SNES_NCG_FR - Fletcher-Reeves update
220: . SNES_NCG_PRP - Polak-Ribiere-Polyak update
221: . SNES_NCG_HS - Hestenes-Steifel update
222: . SNES_NCG_DY - Dai-Yuan update
223: - SNES_NCG_CD - Conjugate Descent update
225: Notes:
226: SNES_NCG_PRP is the default, and the only one that tolerates generalized search directions.
228: It is not clear what "generalized search directions" means, does it mean use with a nonlinear preconditioner,
229: that is using -npc_snes_type <type>, SNESSetNPC(), or SNESGetNPC()?
231: @*/
232: PetscErrorCode SNESNCGSetType(SNES snes, SNESNCGType btype)
233: {
238: PetscTryMethod(snes,"SNESNCGSetType_C",(SNES,SNESNCGType),(snes,btype));
239: return(0);
240: }
242: static PetscErrorCode SNESNCGSetType_NCG(SNES snes, SNESNCGType btype)
243: {
244: SNES_NCG *ncg = (SNES_NCG*)snes->data;
247: ncg->type = btype;
248: return(0);
249: }
251: /*
252: SNESSolve_NCG - Solves a nonlinear system with the Nonlinear Conjugate Gradient method.
254: Input Parameters:
255: . snes - the SNES context
257: Output Parameter:
258: . outits - number of iterations until termination
260: Application Interface Routine: SNESSolve()
261: */
262: static PetscErrorCode SNESSolve_NCG(SNES snes)
263: {
264: SNES_NCG *ncg = (SNES_NCG*)snes->data;
265: Vec X,dX,lX,F,dXold;
266: PetscReal fnorm, ynorm, xnorm, beta = 0.0;
267: PetscScalar dXdotdX, dXolddotdXold, dXdotdXold, lXdotdX, lXdotdXold;
268: PetscInt maxits, i;
269: PetscErrorCode ierr;
270: SNESLineSearchReason lsresult = SNES_LINESEARCH_SUCCEEDED;
271: SNESLineSearch linesearch;
272: SNESConvergedReason reason;
275: if (snes->xl || snes->xu || snes->ops->computevariablebounds) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
277: PetscCitationsRegister(SNESCitation,&SNEScite);
278: snes->reason = SNES_CONVERGED_ITERATING;
280: maxits = snes->max_its; /* maximum number of iterations */
281: X = snes->vec_sol; /* X^n */
282: dXold = snes->work[0]; /* The previous iterate of X */
283: dX = snes->work[1]; /* the preconditioned direction */
284: lX = snes->vec_sol_update; /* the conjugate direction */
285: F = snes->vec_func; /* residual vector */
287: SNESGetLineSearch(snes, &linesearch);
289: PetscObjectSAWsTakeAccess((PetscObject)snes);
290: snes->iter = 0;
291: snes->norm = 0.;
292: PetscObjectSAWsGrantAccess((PetscObject)snes);
294: /* compute the initial function and preconditioned update dX */
296: if (snes->npc && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
297: SNESApplyNPC(snes,X,NULL,dX);
298: SNESGetConvergedReason(snes->npc,&reason);
299: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
300: snes->reason = SNES_DIVERGED_INNER;
301: return(0);
302: }
303: VecCopy(dX,F);
304: VecNorm(F,NORM_2,&fnorm);
305: } else {
306: if (!snes->vec_func_init_set) {
307: SNESComputeFunction(snes,X,F);
308: } else snes->vec_func_init_set = PETSC_FALSE;
310: /* convergence test */
311: VecNorm(F,NORM_2,&fnorm);
312: SNESCheckFunctionNorm(snes,fnorm);
313: VecCopy(F,dX);
314: }
315: if (snes->npc) {
316: if (snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
317: SNESApplyNPC(snes,X,F,dX);
318: SNESGetConvergedReason(snes->npc,&reason);
319: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
320: snes->reason = SNES_DIVERGED_INNER;
321: return(0);
322: }
323: }
324: }
325: VecCopy(dX,lX);
326: VecDot(dX, dX, &dXdotdX);
329: PetscObjectSAWsTakeAccess((PetscObject)snes);
330: snes->norm = fnorm;
331: PetscObjectSAWsGrantAccess((PetscObject)snes);
332: SNESLogConvergenceHistory(snes,fnorm,0);
333: SNESMonitor(snes,0,fnorm);
335: /* test convergence */
336: (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
337: if (snes->reason) return(0);
339: /* Call general purpose update function */
340: if (snes->ops->update) {
341: (*snes->ops->update)(snes, snes->iter);
342: }
344: /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */
346: for (i = 1; i < maxits + 1; i++) {
347: /* some update types require the old update direction or conjugate direction */
348: if (ncg->type != SNES_NCG_FR) {
349: VecCopy(dX, dXold);
350: }
351: SNESLineSearchApply(linesearch,X,F,&fnorm,lX);
352: SNESLineSearchGetReason(linesearch, &lsresult);
353: SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);
354: if (lsresult) {
355: if (++snes->numFailures >= snes->maxFailures) {
356: snes->reason = SNES_DIVERGED_LINE_SEARCH;
357: return(0);
358: }
359: }
360: if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
361: snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
362: return(0);
363: }
364: /* Monitor convergence */
365: PetscObjectSAWsTakeAccess((PetscObject)snes);
366: snes->iter = i;
367: snes->norm = fnorm;
368: snes->xnorm = xnorm;
369: snes->ynorm = ynorm;
370: PetscObjectSAWsGrantAccess((PetscObject)snes);
371: SNESLogConvergenceHistory(snes,snes->norm,0);
372: SNESMonitor(snes,snes->iter,snes->norm);
374: /* Test for convergence */
375: (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);
376: if (snes->reason) return(0);
378: /* Call general purpose update function */
379: if (snes->ops->update) {
380: (*snes->ops->update)(snes, snes->iter);
381: }
382: if (snes->npc) {
383: if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
384: SNESApplyNPC(snes,X,NULL,dX);
385: SNESGetConvergedReason(snes->npc,&reason);
386: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
387: snes->reason = SNES_DIVERGED_INNER;
388: return(0);
389: }
390: VecCopy(dX,F);
391: } else {
392: SNESApplyNPC(snes,X,F,dX);
393: SNESGetConvergedReason(snes->npc,&reason);
394: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
395: snes->reason = SNES_DIVERGED_INNER;
396: return(0);
397: }
398: }
399: } else {
400: VecCopy(F,dX);
401: }
403: /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/
404: switch (ncg->type) {
405: case SNES_NCG_FR: /* Fletcher-Reeves */
406: dXolddotdXold= dXdotdX;
407: VecDot(dX, dX, &dXdotdX);
408: beta = PetscRealPart(dXdotdX / dXolddotdXold);
409: break;
410: case SNES_NCG_PRP: /* Polak-Ribiere-Poylak */
411: dXolddotdXold= dXdotdX;
412: VecDotBegin(dX, dX, &dXdotdXold);
413: VecDotBegin(dXold, dX, &dXdotdXold);
414: VecDotEnd(dX, dX, &dXdotdX);
415: VecDotEnd(dXold, dX, &dXdotdXold);
416: beta = PetscRealPart(((dXdotdX - dXdotdXold) / dXolddotdXold));
417: if (beta < 0.0) beta = 0.0; /* restart */
418: break;
419: case SNES_NCG_HS: /* Hestenes-Stiefel */
420: VecDotBegin(dX, dX, &dXdotdX);
421: VecDotBegin(dX, dXold, &dXdotdXold);
422: VecDotBegin(lX, dX, &lXdotdX);
423: VecDotBegin(lX, dXold, &lXdotdXold);
424: VecDotEnd(dX, dX, &dXdotdX);
425: VecDotEnd(dX, dXold, &dXdotdXold);
426: VecDotEnd(lX, dX, &lXdotdX);
427: VecDotEnd(lX, dXold, &lXdotdXold);
428: beta = PetscRealPart((dXdotdX - dXdotdXold) / (lXdotdX - lXdotdXold));
429: break;
430: case SNES_NCG_DY: /* Dai-Yuan */
431: VecDotBegin(dX, dX, &dXdotdX);
432: VecDotBegin(lX, dX, &lXdotdX);
433: VecDotBegin(lX, dXold, &lXdotdXold);
434: VecDotEnd(dX, dX, &dXdotdX);
435: VecDotEnd(lX, dX, &lXdotdX);
436: VecDotEnd(lX, dXold, &lXdotdXold);
437: beta = PetscRealPart(dXdotdX / (lXdotdXold - lXdotdX));
438: break;
439: case SNES_NCG_CD: /* Conjugate Descent */
440: VecDotBegin(dX, dX, &dXdotdX);
441: VecDotBegin(lX, dXold, &lXdotdXold);
442: VecDotEnd(dX, dX, &dXdotdX);
443: VecDotEnd(lX, dXold, &lXdotdXold);
444: beta = PetscRealPart(dXdotdX / lXdotdXold);
445: break;
446: }
447: if (ncg->monitor) {
448: PetscViewerASCIIPrintf(ncg->monitor, "beta = %e\n", (double)beta);
449: }
450: VecAYPX(lX, beta, dX);
451: }
452: PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);
453: if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
454: return(0);
455: }
457: /*MC
458: SNESNCG - Nonlinear Conjugate-Gradient method for the solution of general nonlinear systems.
460: Level: beginner
462: Options Database:
463: + -snes_ncg_type <fr, prp, dy, hs, cd> - Choice of conjugate-gradient update parameter, default is prp.
464: . -snes_linesearch_type <cp,l2,basic> - Line search type.
465: - -snes_ncg_monitor - Print the beta values used in the ncg iteration, .
467: Notes:
468: This solves the nonlinear system of equations F(x) = 0 using the nonlinear generalization of the conjugate
469: gradient method. This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise
470: chooses the initial search direction as F(x) for the initial guess x.
472: Only supports left non-linear preconditioning.
474: Default line search is SNESLINESEARCHCP, unless a nonlinear preconditioner is used with -npc_snes_type <type>, SNESSetNPC(), or SNESGetNPC() then
475: SNESLINESEARCHL2 is used. Also supports the special purpose line search SNESLINESEARCHNCGLINEAR
477: References:
478: . 1. - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu,"Composing Scalable Nonlinear Algebraic Solvers",
479: SIAM Review, 57(4), 2015
482: .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR, SNESNGMRES, SNESQN, SNESLINESEARCHNCGLINEAR, SNESNCGSetType(), SNESNCGGetType(), SNESLineSearchSetType()
483: M*/
484: PETSC_EXTERN PetscErrorCode SNESCreate_NCG(SNES snes)
485: {
487: SNES_NCG * neP;
490: snes->ops->destroy = SNESDestroy_NCG;
491: snes->ops->setup = SNESSetUp_NCG;
492: snes->ops->setfromoptions = SNESSetFromOptions_NCG;
493: snes->ops->view = SNESView_NCG;
494: snes->ops->solve = SNESSolve_NCG;
495: snes->ops->reset = SNESReset_NCG;
497: snes->usesksp = PETSC_FALSE;
498: snes->usesnpc = PETSC_TRUE;
499: snes->npcside = PC_LEFT;
501: snes->alwayscomputesfinalresidual = PETSC_TRUE;
503: if (!snes->tolerancesset) {
504: snes->max_funcs = 30000;
505: snes->max_its = 10000;
506: snes->stol = 1e-20;
507: }
509: PetscNewLog(snes,&neP);
510: snes->data = (void*) neP;
511: neP->monitor = NULL;
512: neP->type = SNES_NCG_PRP;
513: PetscObjectComposeFunction((PetscObject)snes,"SNESNCGSetType_C", SNESNCGSetType_NCG);
514: return(0);
515: }