Actual source code: snesncg.c
petsc-3.4.5 2014-06-29
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: return(0);
55: }
56: /*
57: SNESSetFromOptions_NCG - Sets various parameters for the SNESNCG method.
59: Input Parameter:
60: . snes - the SNES context
62: Application Interface Routine: SNESSetFromOptions()
63: */
66: static PetscErrorCode SNESSetFromOptions_NCG(SNES snes)
67: {
68: SNES_NCG *ncg = (SNES_NCG*)snes->data;
70: PetscBool debug;
71: SNESLineSearch linesearch;
72: SNESNCGType ncgtype=ncg->type;
75: PetscOptionsHead("SNES NCG options");
76: PetscOptionsBool("-snes_ncg_monitor","Monitor NCG iterations","SNES",ncg->monitor ? PETSC_TRUE : PETSC_FALSE, &debug, NULL);
77: PetscOptionsEnum("-snes_ncg_type","NCG Beta type used","SNESNCGSetType",SNESNCGTypes,(PetscEnum)ncg->type,(PetscEnum*)&ncgtype,NULL);
78: SNESNCGSetType(snes, ncgtype);
79: if (debug) {
80: ncg->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
81: }
82: PetscOptionsTail();
83: if (!snes->linesearch) {
84: SNESGetLineSearch(snes, &linesearch);
85: if (!snes->pc) {
86: SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);
87: } else {
88: SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);
89: }
90: }
91: SNESLineSearchRegister(SNESLINESEARCHNCGLINEAR, SNESLineSearchCreate_NCGLinear);
92: return(0);
93: }
95: /*
96: SNESView_NCG - Prints info from the SNESNCG data structure.
98: Input Parameters:
99: + SNES - the SNES context
100: - viewer - visualization context
102: Application Interface Routine: SNESView()
103: */
106: static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer)
107: {
108: PetscBool iascii;
112: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);
113: if (iascii) {
114: }
115: return(0);
116: }
120: PetscErrorCode SNESLineSearchApply_NCGLinear(SNESLineSearch linesearch)
121: {
122: PetscScalar alpha, ptAp;
123: Vec X, Y, F, W;
124: SNES snes;
126: PetscReal *fnorm, *xnorm, *ynorm;
127: MatStructure flg = DIFFERENT_NONZERO_PATTERN;
130: SNESLineSearchGetSNES(linesearch, &snes);
131: X = linesearch->vec_sol;
132: W = linesearch->vec_sol_new;
133: F = linesearch->vec_func;
134: Y = linesearch->vec_update;
135: fnorm = &linesearch->fnorm;
136: xnorm = &linesearch->xnorm;
137: ynorm = &linesearch->ynorm;
139: if (!snes->jacobian) {
140: SNESSetUpMatrices(snes);
141: }
143: /*
145: The exact step size for unpreconditioned linear CG is just:
146: alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy)
147: */
148: SNESComputeJacobian(snes, X, &snes->jacobian, &snes->jacobian_pre, &flg);
149: VecDot(F, F, &alpha);
150: MatMult(snes->jacobian, Y, W);
151: VecDot(Y, W, &ptAp);
152: alpha = alpha / ptAp;
153: VecAXPY(X,-alpha,Y);
154: SNESComputeFunction(snes,X,F);
156: VecNorm(F,NORM_2,fnorm);
157: VecNorm(X,NORM_2,xnorm);
158: VecNorm(Y,NORM_2,ynorm);
159: return(0);
160: }
164: PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch linesearch)
165: {
167: linesearch->ops->apply = SNESLineSearchApply_NCGLinear;
168: linesearch->ops->destroy = NULL;
169: linesearch->ops->setfromoptions = NULL;
170: linesearch->ops->reset = NULL;
171: linesearch->ops->view = NULL;
172: linesearch->ops->setup = NULL;
173: return(0);
174: }
178: /*
180: Assuming F = SNESComputeFunction(X) compute Y^tJ^tF using a simple secant approximation of the jacobian.
182: */
183: PetscErrorCode SNESNCGComputeYtJtF_Private(SNES snes, Vec X, Vec F, Vec Y, Vec W, Vec G, PetscScalar * ytJtf)
184: {
186: PetscScalar ftf, ftg, fty, h;
189: VecDot(F, F, &ftf);
190: VecDot(F, Y, &fty);
191: h = 1e-5*fty / fty;
192: VecCopy(X, W);
193: VecAXPY(W, -h, Y); /* this is arbitrary */
194: SNESComputeFunction(snes, W, G);
195: VecDot(G, F, &ftg);
196: *ytJtf = (ftg - ftf) / h;
197: return(0);
198: }
202: /*@
203: SNESNCGSetType - Sets the conjugate update type for SNESNCG.
205: Logically Collective on SNES
207: Input Parameters:
208: + snes - the iterative context
209: - btype - update type
211: Options Database:
212: . -snes_ncg_type<prp,fr,hs,dy,cd>
214: Level: intermediate
216: SNESNCGSelectTypes:
217: + SNES_NCG_FR - Fletcher-Reeves update
218: . SNES_NCG_PRP - Polak-Ribiere-Polyak update
219: . SNES_NCG_HS - Hestenes-Steifel update
220: . SNES_NCG_DY - Dai-Yuan update
221: - SNES_NCG_CD - Conjugate Descent update
223: Notes:
224: PRP is the default, and the only one that tolerates generalized search directions.
226: .keywords: SNES, SNESNCG, selection, type, set
227: @*/
228: PetscErrorCode SNESNCGSetType(SNES snes, SNESNCGType btype)
229: {
234: PetscTryMethod(snes,"SNESNCGSetType_C",(SNES,SNESNCGType),(snes,btype));
235: return(0);
236: }
240: PetscErrorCode SNESNCGSetType_NCG(SNES snes, SNESNCGType btype)
241: {
242: SNES_NCG *ncg = (SNES_NCG*)snes->data;
245: ncg->type = btype;
246: return(0);
247: }
249: /*
250: SNESSolve_NCG - Solves a nonlinear system with the Nonlinear Conjugate Gradient method.
252: Input Parameters:
253: . snes - the SNES context
255: Output Parameter:
256: . outits - number of iterations until termination
258: Application Interface Routine: SNESSolve()
259: */
262: PetscErrorCode SNESSolve_NCG(SNES snes)
263: {
264: SNES_NCG *ncg = (SNES_NCG*)snes->data;
265: Vec X, dX, lX, F, B, Fold;
266: PetscReal fnorm, ynorm, xnorm, beta = 0.0;
267: PetscScalar dXdotF, dXolddotFold, dXdotFold, lXdotF, lXdotFold;
268: PetscInt maxits, i;
269: PetscErrorCode ierr;
270: SNESConvergedReason reason;
271: PetscBool lsSuccess = PETSC_TRUE;
272: SNESLineSearch linesearch;
275: snes->reason = SNES_CONVERGED_ITERATING;
277: maxits = snes->max_its; /* maximum number of iterations */
278: X = snes->vec_sol; /* X^n */
279: Fold = snes->work[0]; /* The previous iterate of X */
280: dX = snes->work[1]; /* the preconditioned direction */
281: lX = snes->vec_sol_update; /* the conjugate direction */
282: F = snes->vec_func; /* residual vector */
283: B = snes->vec_rhs; /* the right hand side */
285: SNESGetLineSearch(snes, &linesearch);
287: PetscObjectAMSTakeAccess((PetscObject)snes);
288: snes->iter = 0;
289: snes->norm = 0.;
290: PetscObjectAMSGrantAccess((PetscObject)snes);
292: /* compute the initial function and preconditioned update dX */
293: if (!snes->vec_func_init_set) {
294: SNESComputeFunction(snes,X,F);
295: if (snes->domainerror) {
296: snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
297: return(0);
298: }
299: } else snes->vec_func_init_set = PETSC_FALSE;
301: if (!snes->norm_init_set) {
302: /* convergence test */
303: VecNorm(F, NORM_2, &fnorm); /* fnorm <- ||F|| */
304: if (PetscIsInfOrNanReal(fnorm)) {
305: snes->reason = SNES_DIVERGED_FNORM_NAN;
306: return(0);
307: }
308: } else {
309: fnorm = snes->norm_init;
310: snes->norm_init_set = PETSC_FALSE;
311: }
312: PetscObjectAMSTakeAccess((PetscObject)snes);
313: snes->norm = fnorm;
314: PetscObjectAMSGrantAccess((PetscObject)snes);
315: SNESLogConvergenceHistory(snes,fnorm,0);
316: SNESMonitor(snes,0,fnorm);
318: /* set parameter for default relative tolerance convergence test */
319: snes->ttol = fnorm*snes->rtol;
320: /* test convergence */
321: (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
322: if (snes->reason) return(0);
324: /* Call general purpose update function */
325: if (snes->ops->update) {
326: (*snes->ops->update)(snes, snes->iter);
327: }
329: /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */
331: if (snes->pc && snes->pcside == PC_RIGHT) {
332: VecCopy(X, dX);
333: SNESSetInitialFunction(snes->pc, F);
334: SNESSetInitialFunctionNorm(snes->pc, fnorm);
336: PetscLogEventBegin(SNES_NPCSolve,snes->pc,dX,B,0);
337: SNESSolve(snes->pc, B, dX);
338: PetscLogEventEnd(SNES_NPCSolve,snes->pc,dX,B,0);
340: SNESGetConvergedReason(snes->pc,&reason);
341: if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
342: snes->reason = SNES_DIVERGED_INNER;
343: return(0);
344: }
345: VecAYPX(dX,-1.0,X);
346: } else {
347: VecCopy(F, dX);
348: }
349: VecCopy(dX, lX);
350: VecDot(F, dX, &dXdotF);
351: /*
352: } else {
353: SNESNCGComputeYtJtF_Private(snes, X, F, dX, W, G, &dXdotF);
354: }
355: */
356: for (i = 1; i < maxits + 1; i++) {
357: lsSuccess = PETSC_TRUE;
358: /* some update types require the old update direction or conjugate direction */
359: if (ncg->type != SNES_NCG_FR) {
360: VecCopy(F, Fold);
361: }
362: SNESLineSearchApply(linesearch, X, F, &fnorm, lX);
363: SNESLineSearchGetSuccess(linesearch, &lsSuccess);
364: if (!lsSuccess) {
365: if (++snes->numFailures >= snes->maxFailures) {
366: snes->reason = SNES_DIVERGED_LINE_SEARCH;
367: return(0);
368: }
369: }
370: if (snes->nfuncs >= snes->max_funcs) {
371: snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
372: return(0);
373: }
374: if (snes->domainerror) {
375: snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
376: return(0);
377: }
378: SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);
379: /* Monitor convergence */
380: PetscObjectAMSTakeAccess((PetscObject)snes);
381: snes->iter = i;
382: snes->norm = fnorm;
383: PetscObjectAMSGrantAccess((PetscObject)snes);
384: SNESLogConvergenceHistory(snes,snes->norm,0);
385: SNESMonitor(snes,snes->iter,snes->norm);
387: /* Test for convergence */
388: (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);
389: if (snes->reason) return(0);
391: /* Call general purpose update function */
392: if (snes->ops->update) {
393: (*snes->ops->update)(snes, snes->iter);
394: }
395: if (snes->pc && snes->pcside == PC_RIGHT) {
396: VecCopy(X,dX);
397: SNESSetInitialFunction(snes->pc, F);
398: SNESSetInitialFunctionNorm(snes->pc, fnorm);
399: PetscLogEventBegin(SNES_NPCSolve,snes->pc,dX,B,0);
400: SNESSolve(snes->pc, B, dX);
401: PetscLogEventEnd(SNES_NPCSolve,snes->pc,dX,B,0);
402: SNESGetConvergedReason(snes->pc,&reason);
403: if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
404: snes->reason = SNES_DIVERGED_INNER;
405: return(0);
406: }
407: VecAYPX(dX,-1.0,X);
408: } else {
409: VecCopy(F, dX);
410: }
412: /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/
413: switch (ncg->type) {
414: case SNES_NCG_FR: /* Fletcher-Reeves */
415: dXolddotFold = dXdotF;
416: VecDot(dX, dX, &dXdotF);
417: beta = PetscRealPart(dXdotF / dXolddotFold);
418: break;
419: case SNES_NCG_PRP: /* Polak-Ribiere-Poylak */
420: dXolddotFold = dXdotF;
421: VecDotBegin(F, dX, &dXdotF);
422: VecDotBegin(Fold, dX, &dXdotFold);
423: VecDotEnd(F, dX, &dXdotF);
424: VecDotEnd(Fold, dX, &dXdotFold);
425: beta = PetscRealPart(((dXdotF - dXdotFold) / dXolddotFold));
426: if (beta < 0.0) beta = 0.0; /* restart */
427: break;
428: case SNES_NCG_HS: /* Hestenes-Stiefel */
429: VecDotBegin(dX, F, &dXdotF);
430: VecDotBegin(dX, Fold, &dXdotFold);
431: VecDotBegin(lX, F, &lXdotF);
432: VecDotBegin(lX, Fold, &lXdotFold);
433: VecDotEnd(dX, F, &dXdotF);
434: VecDotEnd(dX, Fold, &dXdotFold);
435: VecDotEnd(lX, F, &lXdotF);
436: VecDotEnd(lX, Fold, &lXdotFold);
437: beta = PetscRealPart((dXdotF - dXdotFold) / (lXdotF - lXdotFold));
438: break;
439: case SNES_NCG_DY: /* Dai-Yuan */
440: VecDotBegin(dX, F, &dXdotF);
441: VecDotBegin(lX, F, &lXdotF);
442: VecDotBegin(lX, Fold, &lXdotFold);
443: VecDotEnd(dX, F, &dXdotF);
444: VecDotEnd(lX, F, &lXdotF);
445: VecDotEnd(lX, Fold, &lXdotFold);
446: beta = PetscRealPart(dXdotF / (lXdotFold - lXdotF));
447: break;
448: case SNES_NCG_CD: /* Conjugate Descent */
449: VecDotBegin(dX, F, &dXdotF);
450: VecDotBegin(lX, Fold, &lXdotFold);
451: VecDotEnd(dX, F, &dXdotF);
452: VecDotEnd(lX, Fold, &lXdotFold);
453: beta = PetscRealPart(dXdotF / lXdotFold);
454: break;
455: }
456: if (ncg->monitor) {
457: PetscViewerASCIIPrintf(ncg->monitor, "beta = %e\n", beta);
458: }
459: VecAYPX(lX, beta, dX);
460: }
461: PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);
462: if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
463: return(0);
464: }
468: /*MC
469: SNESNCG - Nonlinear Conjugate-Gradient method for the solution of general nonlinear systems.
471: Level: beginner
473: Options Database:
474: + -snes_ncg_type <fr, prp, dy, hs, cd> - Choice of conjugate-gradient update parameter.
475: . -snes_linesearch_type <cp,l2,basic> - Line search type.
476: - -snes_ncg_monitor - Print relevant information about the ncg iteration.
478: Notes: This solves the nonlinear system of equations F(x) = 0 using the nonlinear generalization of the conjugate
479: gradient method. This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise
480: chooses the initial search direction as F(x) for the initial guess x.
483: .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR, SNESNGMRES, SNESNQN
484: M*/
487: PETSC_EXTERN PetscErrorCode SNESCreate_NCG(SNES snes)
488: {
490: SNES_NCG * neP;
493: snes->ops->destroy = SNESDestroy_NCG;
494: snes->ops->setup = SNESSetUp_NCG;
495: snes->ops->setfromoptions = SNESSetFromOptions_NCG;
496: snes->ops->view = SNESView_NCG;
497: snes->ops->solve = SNESSolve_NCG;
498: snes->ops->reset = SNESReset_NCG;
500: snes->usesksp = PETSC_FALSE;
501: snes->usespc = 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, SNES_NCG, &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: }