Actual source code: snesncg.c
petsc-3.7.7 2017-09-25
1: #include <../src/snes/impls/ncg/snesncgimpl.h> /*I "petscsnes.h" I*/
2: const char *const SNESNCGTypes[] = {"FR","PRP","HS","DY","CD","SNESNCGType","SNES_NCG_",0};
6: PetscErrorCode SNESReset_NCG(SNES snes)
7: {
9: return(0);
10: }
12: #define SNESLINESEARCHNCGLINEAR "ncglinear"
14: /*
15: SNESDestroy_NCG - Destroys the private SNES_NCG context that was created with SNESCreate_NCG().
17: Input Parameter:
18: . snes - the SNES context
20: Application Interface Routine: SNESDestroy()
21: */
24: PetscErrorCode SNESDestroy_NCG(SNES snes)
25: {
29: PetscFree(snes->data);
30: return(0);
31: }
33: /*
34: SNESSetUp_NCG - Sets up the internal data structures for the later use
35: of the SNESNCG nonlinear solver.
37: Input Parameters:
38: + snes - the SNES context
39: - x - the solution vector
41: Application Interface Routine: SNESSetUp()
42: */
44: PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch);
48: PetscErrorCode SNESSetUp_NCG(SNES snes)
49: {
53: SNESSetWorkVecs(snes,2);
54: if (snes->pcside == PC_RIGHT) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNESNCG only supports left preconditioning");
55: if (snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED;
56: return(0);
57: }
58: /*
59: SNESSetFromOptions_NCG - Sets various parameters for the SNESNCG method.
61: Input Parameter:
62: . snes - the SNES context
64: Application Interface Routine: SNESSetFromOptions()
65: */
68: static PetscErrorCode SNESSetFromOptions_NCG(PetscOptionItems *PetscOptionsObject,SNES snes)
69: {
70: SNES_NCG *ncg = (SNES_NCG*)snes->data;
72: PetscBool debug = PETSC_FALSE;
73: SNESLineSearch linesearch;
74: SNESNCGType ncgtype=ncg->type;
77: PetscOptionsHead(PetscOptionsObject,"SNES NCG options");
78: PetscOptionsBool("-snes_ncg_monitor","Monitor NCG iterations","SNES",ncg->monitor ? PETSC_TRUE : PETSC_FALSE, &debug, NULL);
79: if (debug) {
80: ncg->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
81: }
82: PetscOptionsEnum("-snes_ncg_type","NCG Beta type used","SNESNCGSetType",SNESNCGTypes,(PetscEnum)ncg->type,(PetscEnum*)&ncgtype,NULL);
83: SNESNCGSetType(snes, ncgtype);
84: PetscOptionsTail();
85: if (!snes->linesearch) {
86: SNESGetLineSearch(snes, &linesearch);
87: if (!snes->pc) {
88: SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);
89: } else {
90: SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);
91: }
92: }
93: SNESLineSearchRegister(SNESLINESEARCHNCGLINEAR, SNESLineSearchCreate_NCGLinear);
94: return(0);
95: }
97: /*
98: SNESView_NCG - Prints info from the SNESNCG data structure.
100: Input Parameters:
101: + SNES - the SNES context
102: - viewer - visualization context
104: Application Interface Routine: SNESView()
105: */
108: static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer)
109: {
110: SNES_NCG *ncg = (SNES_NCG *) snes->data;
111: PetscBool iascii;
115: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);
116: if (iascii) {
117: PetscViewerASCIIPrintf(viewer, "ncg type: %s\n", SNESNCGTypes[ncg->type]);
118: }
119: return(0);
120: }
124: PetscErrorCode SNESLineSearchApply_NCGLinear(SNESLineSearch linesearch)
125: {
126: PetscScalar alpha, ptAp;
127: Vec X, Y, F, W;
128: SNES snes;
130: PetscReal *fnorm, *xnorm, *ynorm;
133: SNESLineSearchGetSNES(linesearch, &snes);
134: X = linesearch->vec_sol;
135: W = linesearch->vec_sol_new;
136: F = linesearch->vec_func;
137: Y = linesearch->vec_update;
138: fnorm = &linesearch->fnorm;
139: xnorm = &linesearch->xnorm;
140: ynorm = &linesearch->ynorm;
142: if (!snes->jacobian) {
143: SNESSetUpMatrices(snes);
144: }
146: /*
148: The exact step size for unpreconditioned linear CG is just:
149: alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy)
150: */
151: SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre);
152: VecDot(F, F, &alpha);
153: MatMult(snes->jacobian, Y, W);
154: VecDot(Y, W, &ptAp);
155: alpha = alpha / ptAp;
156: VecAXPY(X,-alpha,Y);
157: SNESComputeFunction(snes,X,F);
159: VecNorm(F,NORM_2,fnorm);
160: VecNorm(X,NORM_2,xnorm);
161: VecNorm(Y,NORM_2,ynorm);
162: return(0);
163: }
167: PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch linesearch)
168: {
170: linesearch->ops->apply = SNESLineSearchApply_NCGLinear;
171: linesearch->ops->destroy = NULL;
172: linesearch->ops->setfromoptions = NULL;
173: linesearch->ops->reset = NULL;
174: linesearch->ops->view = NULL;
175: linesearch->ops->setup = NULL;
176: return(0);
177: }
181: /*
183: Assuming F = SNESComputeFunction(X) compute Y^tJ^tF using a simple secant approximation of the jacobian.
185: */
186: PetscErrorCode SNESNCGComputeYtJtF_Private(SNES snes, Vec X, Vec F, Vec Y, Vec W, Vec G, PetscScalar * ytJtf)
187: {
189: PetscScalar ftf, ftg, fty, h;
192: VecDot(F, F, &ftf);
193: VecDot(F, Y, &fty);
194: h = 1e-5*fty / fty;
195: VecCopy(X, W);
196: VecAXPY(W, -h, Y); /* this is arbitrary */
197: SNESComputeFunction(snes, W, G);
198: VecDot(G, F, &ftg);
199: *ytJtf = (ftg - ftf) / h;
200: return(0);
201: }
205: /*@
206: SNESNCGSetType - Sets the conjugate update type for SNESNCG.
208: Logically Collective on SNES
210: Input Parameters:
211: + snes - the iterative context
212: - btype - update type
214: Options Database:
215: . -snes_ncg_type<prp,fr,hs,dy,cd>
217: Level: intermediate
219: SNESNCGSelectTypes:
220: + SNES_NCG_FR - Fletcher-Reeves update
221: . SNES_NCG_PRP - Polak-Ribiere-Polyak update
222: . SNES_NCG_HS - Hestenes-Steifel update
223: . SNES_NCG_DY - Dai-Yuan update
224: - SNES_NCG_CD - Conjugate Descent update
226: Notes:
227: PRP is the default, and the only one that tolerates generalized search directions.
229: .keywords: SNES, SNESNCG, selection, type, set
230: @*/
231: PetscErrorCode SNESNCGSetType(SNES snes, SNESNCGType btype)
232: {
237: PetscTryMethod(snes,"SNESNCGSetType_C",(SNES,SNESNCGType),(snes,btype));
238: return(0);
239: }
243: PetscErrorCode SNESNCGSetType_NCG(SNES snes, SNESNCGType btype)
244: {
245: SNES_NCG *ncg = (SNES_NCG*)snes->data;
248: ncg->type = btype;
249: return(0);
250: }
252: /*
253: SNESSolve_NCG - Solves a nonlinear system with the Nonlinear Conjugate Gradient method.
255: Input Parameters:
256: . snes - the SNES context
258: Output Parameter:
259: . outits - number of iterations until termination
261: Application Interface Routine: SNESSolve()
262: */
265: PetscErrorCode SNESSolve_NCG(SNES snes)
266: {
267: SNES_NCG *ncg = (SNES_NCG*)snes->data;
268: Vec X,dX,lX,F,dXold;
269: PetscReal fnorm, ynorm, xnorm, beta = 0.0;
270: PetscScalar dXdotdX, dXolddotdXold, dXdotdXold, lXdotdX, lXdotdXold;
271: PetscInt maxits, i;
272: PetscErrorCode ierr;
273: SNESLineSearchReason lsresult = SNES_LINESEARCH_SUCCEEDED;
274: SNESLineSearch linesearch;
275: SNESConvergedReason reason;
278: 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);
280: PetscCitationsRegister(SNESCitation,&SNEScite);
281: snes->reason = SNES_CONVERGED_ITERATING;
283: maxits = snes->max_its; /* maximum number of iterations */
284: X = snes->vec_sol; /* X^n */
285: dXold = snes->work[0]; /* The previous iterate of X */
286: dX = snes->work[1]; /* the preconditioned direction */
287: lX = snes->vec_sol_update; /* the conjugate direction */
288: F = snes->vec_func; /* residual vector */
290: SNESGetLineSearch(snes, &linesearch);
292: PetscObjectSAWsTakeAccess((PetscObject)snes);
293: snes->iter = 0;
294: snes->norm = 0.;
295: PetscObjectSAWsGrantAccess((PetscObject)snes);
297: /* compute the initial function and preconditioned update dX */
299: if (snes->pc && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
300: SNESApplyNPC(snes,X,NULL,dX);
301: SNESGetConvergedReason(snes->pc,&reason);
302: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
303: snes->reason = SNES_DIVERGED_INNER;
304: return(0);
305: }
306: VecCopy(dX,F);
307: VecNorm(F,NORM_2,&fnorm);
308: } else {
309: if (!snes->vec_func_init_set) {
310: SNESComputeFunction(snes,X,F);
311: } else snes->vec_func_init_set = PETSC_FALSE;
313: /* convergence test */
314: VecNorm(F,NORM_2,&fnorm);
315: SNESCheckFunctionNorm(snes,fnorm);
316: VecCopy(F,dX);
317: }
318: if (snes->pc) {
319: if (snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
320: SNESApplyNPC(snes,X,F,dX);
321: SNESGetConvergedReason(snes->pc,&reason);
322: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
323: snes->reason = SNES_DIVERGED_INNER;
324: return(0);
325: }
326: }
327: }
328: VecCopy(dX,lX);
329: VecDot(dX, dX, &dXdotdX);
332: PetscObjectSAWsTakeAccess((PetscObject)snes);
333: snes->norm = fnorm;
334: PetscObjectSAWsGrantAccess((PetscObject)snes);
335: SNESLogConvergenceHistory(snes,fnorm,0);
336: SNESMonitor(snes,0,fnorm);
338: /* test convergence */
339: (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
340: if (snes->reason) return(0);
342: /* Call general purpose update function */
343: if (snes->ops->update) {
344: (*snes->ops->update)(snes, snes->iter);
345: }
347: /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */
349: for (i = 1; i < maxits + 1; i++) {
350: /* some update types require the old update direction or conjugate direction */
351: if (ncg->type != SNES_NCG_FR) {
352: VecCopy(dX, dXold);
353: }
354: SNESLineSearchApply(linesearch,X,F,&fnorm,lX);
355: SNESLineSearchGetReason(linesearch, &lsresult);
356: SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);
357: if (lsresult) {
358: if (++snes->numFailures >= snes->maxFailures) {
359: snes->reason = SNES_DIVERGED_LINE_SEARCH;
360: return(0);
361: }
362: }
363: if (snes->nfuncs >= snes->max_funcs) {
364: snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
365: return(0);
366: }
367: /* Monitor convergence */
368: PetscObjectSAWsTakeAccess((PetscObject)snes);
369: snes->iter = i;
370: snes->norm = fnorm;
371: PetscObjectSAWsGrantAccess((PetscObject)snes);
372: SNESLogConvergenceHistory(snes,snes->norm,0);
373: SNESMonitor(snes,snes->iter,snes->norm);
375: /* Test for convergence */
376: (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);
377: if (snes->reason) return(0);
379: /* Call general purpose update function */
380: if (snes->ops->update) {
381: (*snes->ops->update)(snes, snes->iter);
382: }
383: if (snes->pc) {
384: if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
385: SNESApplyNPC(snes,X,NULL,dX);
386: SNESGetConvergedReason(snes->pc,&reason);
387: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
388: snes->reason = SNES_DIVERGED_INNER;
389: return(0);
390: }
391: VecCopy(dX,F);
392: } else {
393: SNESApplyNPC(snes,X,F,dX);
394: SNESGetConvergedReason(snes->pc,&reason);
395: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
396: snes->reason = SNES_DIVERGED_INNER;
397: return(0);
398: }
399: }
400: } else {
401: VecCopy(F,dX);
402: }
404: /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/
405: switch (ncg->type) {
406: case SNES_NCG_FR: /* Fletcher-Reeves */
407: dXolddotdXold= dXdotdX;
408: VecDot(dX, dX, &dXdotdX);
409: beta = PetscRealPart(dXdotdX / dXolddotdXold);
410: break;
411: case SNES_NCG_PRP: /* Polak-Ribiere-Poylak */
412: dXolddotdXold= dXdotdX;
413: VecDotBegin(dX, dX, &dXdotdXold);
414: VecDotBegin(dXold, dX, &dXdotdXold);
415: VecDotEnd(dX, dX, &dXdotdX);
416: VecDotEnd(dXold, dX, &dXdotdXold);
417: beta = PetscRealPart(((dXdotdX - dXdotdXold) / dXolddotdXold));
418: if (beta < 0.0) beta = 0.0; /* restart */
419: break;
420: case SNES_NCG_HS: /* Hestenes-Stiefel */
421: VecDotBegin(dX, dX, &dXdotdX);
422: VecDotBegin(dX, dXold, &dXdotdXold);
423: VecDotBegin(lX, dX, &lXdotdX);
424: VecDotBegin(lX, dXold, &lXdotdXold);
425: VecDotEnd(dX, dX, &dXdotdX);
426: VecDotEnd(dX, dXold, &dXdotdXold);
427: VecDotEnd(lX, dX, &lXdotdX);
428: VecDotEnd(lX, dXold, &lXdotdXold);
429: beta = PetscRealPart((dXdotdX - dXdotdXold) / (lXdotdX - lXdotdXold));
430: break;
431: case SNES_NCG_DY: /* Dai-Yuan */
432: VecDotBegin(dX, dX, &dXdotdX);
433: VecDotBegin(lX, dX, &lXdotdX);
434: VecDotBegin(lX, dXold, &lXdotdXold);
435: VecDotEnd(dX, dX, &dXdotdX);
436: VecDotEnd(lX, dX, &lXdotdX);
437: VecDotEnd(lX, dXold, &lXdotdXold);
438: beta = PetscRealPart(dXdotdX / (lXdotdXold - lXdotdX));
439: break;
440: case SNES_NCG_CD: /* Conjugate Descent */
441: VecDotBegin(dX, dX, &dXdotdX);
442: VecDotBegin(lX, dXold, &lXdotdXold);
443: VecDotEnd(dX, dX, &dXdotdX);
444: VecDotEnd(lX, dXold, &lXdotdXold);
445: beta = PetscRealPart(dXdotdX / lXdotdXold);
446: break;
447: }
448: if (ncg->monitor) {
449: PetscViewerASCIIPrintf(ncg->monitor, "beta = %e\n", (double)beta);
450: }
451: VecAYPX(lX, beta, dX);
452: }
453: PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);
454: if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
455: return(0);
456: }
460: /*MC
461: SNESNCG - Nonlinear Conjugate-Gradient method for the solution of general nonlinear systems.
463: Level: beginner
465: Options Database:
466: + -snes_ncg_type <fr, prp, dy, hs, cd> - Choice of conjugate-gradient update parameter, default is prp.
467: . -snes_linesearch_type <cp,l2,basic> - Line search type.
468: - -snes_ncg_monitor - Print relevant information about the ncg iteration.
470: Notes: This solves the nonlinear system of equations F(x) = 0 using the nonlinear generalization of the conjugate
471: gradient method. This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise
472: chooses the initial search direction as F(x) for the initial guess x.
474: Only supports left non-linear preconditioning.
476: References:
477: . 1. - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu,"Composing Scalable Nonlinear Algebraic Solvers",
478: SIAM Review, 57(4), 2015
481: .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR, SNESNGMRES, SNESNQN
482: M*/
485: PETSC_EXTERN PetscErrorCode SNESCreate_NCG(SNES snes)
486: {
488: SNES_NCG * neP;
491: snes->ops->destroy = SNESDestroy_NCG;
492: snes->ops->setup = SNESSetUp_NCG;
493: snes->ops->setfromoptions = SNESSetFromOptions_NCG;
494: snes->ops->view = SNESView_NCG;
495: snes->ops->solve = SNESSolve_NCG;
496: snes->ops->reset = SNESReset_NCG;
498: snes->usesksp = PETSC_FALSE;
499: snes->usespc = PETSC_TRUE;
500: snes->pcside = PC_LEFT;
502: if (!snes->tolerancesset) {
503: snes->max_funcs = 30000;
504: snes->max_its = 10000;
505: snes->stol = 1e-20;
506: }
508: PetscNewLog(snes,&neP);
509: snes->data = (void*) neP;
510: neP->monitor = NULL;
511: neP->type = SNES_NCG_PRP;
512: PetscObjectComposeFunction((PetscObject)snes,"SNESNCGSetType_C", SNESNCGSetType_NCG);
513: return(0);
514: }