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: {
6: PetscFunctionBegin;
7: PetscFunctionReturn(PETSC_SUCCESS);
8: }
10: static PetscErrorCode SNESDestroy_NCG(SNES snes)
11: {
12: PetscFunctionBegin;
13: PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNCGSetType_C", NULL));
14: PetscCall(PetscFree(snes->data));
15: PetscFunctionReturn(PETSC_SUCCESS);
16: }
18: static PetscErrorCode SNESSetUp_NCG(SNES snes)
19: {
20: PetscFunctionBegin;
21: PetscCall(SNESSetWorkVecs(snes, 2));
22: PetscCheck(snes->npcside != PC_RIGHT, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNESNCG only supports left preconditioning");
23: if (snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED;
24: PetscFunctionReturn(PETSC_SUCCESS);
25: }
27: static PetscErrorCode SNESLineSearchApply_NCGLinear(SNESLineSearch linesearch)
28: {
29: PetscScalar alpha, ptAp;
30: Vec X, Y, F, W;
31: SNES snes;
32: PetscReal *fnorm, *xnorm, *ynorm;
34: PetscFunctionBegin;
35: PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
36: X = linesearch->vec_sol;
37: W = linesearch->vec_sol_new;
38: F = linesearch->vec_func;
39: Y = linesearch->vec_update;
40: fnorm = &linesearch->fnorm;
41: xnorm = &linesearch->xnorm;
42: ynorm = &linesearch->ynorm;
44: if (!snes->jacobian) PetscCall(SNESSetUpMatrices(snes));
46: /*
47: The exact step size for unpreconditioned linear CG is just:
48: alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy)
49: */
50: PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
51: PetscCall(VecDot(F, F, &alpha));
52: PetscCall(MatMult(snes->jacobian, Y, W));
53: PetscCall(VecDot(Y, W, &ptAp));
54: alpha = alpha / ptAp;
55: PetscCall(VecAXPY(X, -alpha, Y));
56: PetscCall(SNESComputeFunction(snes, X, F));
58: PetscCall(VecNorm(F, NORM_2, fnorm));
59: PetscCall(VecNorm(X, NORM_2, xnorm));
60: PetscCall(VecNorm(Y, NORM_2, ynorm));
61: PetscFunctionReturn(PETSC_SUCCESS);
62: }
64: /*MC
65: SNESLINESEARCHNCGLINEAR - Special line search only for the nonlinear CG solver `SNESNCG`
67: 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.
68: 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.
70: Level: advanced
72: Notes:
73: This requires a Jacobian-vector product but does not require the solution of a linear system with the Jacobian
75: This is a "odd-ball" line search, we don't know if it is in the literature or used in practice by anyone.
77: .seealso: [](ch_snes), `SNES`, `SNESNCG`, `SNESLineSearchCreate()`, `SNESLineSearchSetType()`
78: M*/
80: PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch linesearch)
81: {
82: PetscFunctionBegin;
83: linesearch->ops->apply = SNESLineSearchApply_NCGLinear;
84: linesearch->ops->destroy = NULL;
85: linesearch->ops->setfromoptions = NULL;
86: linesearch->ops->reset = NULL;
87: linesearch->ops->view = NULL;
88: linesearch->ops->setup = NULL;
89: PetscFunctionReturn(PETSC_SUCCESS);
90: }
92: static PetscErrorCode SNESSetFromOptions_NCG(SNES snes, PetscOptionItems *PetscOptionsObject)
93: {
94: SNES_NCG *ncg = (SNES_NCG *)snes->data;
95: PetscBool debug = PETSC_FALSE;
96: SNESNCGType ncgtype = ncg->type;
97: SNESLineSearch linesearch;
99: PetscFunctionBegin;
100: PetscOptionsHeadBegin(PetscOptionsObject, "SNES NCG options");
101: PetscCall(PetscOptionsBool("-snes_ncg_monitor", "Monitor the beta values used in the NCG iterations", "SNES", ncg->monitor ? PETSC_TRUE : PETSC_FALSE, &debug, NULL));
102: if (debug) ncg->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
103: PetscCall(PetscOptionsEnum("-snes_ncg_type", "NCG Beta type used", "SNESNCGSetType", SNESNCGTypes, (PetscEnum)ncg->type, (PetscEnum *)&ncgtype, NULL));
104: PetscCall(SNESNCGSetType(snes, ncgtype));
105: PetscOptionsHeadEnd();
106: if (!snes->linesearch) {
107: PetscCall(SNESGetLineSearch(snes, &linesearch));
108: if (!((PetscObject)linesearch)->type_name) {
109: if (!snes->npc) {
110: PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHCP));
111: } else {
112: PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHL2));
113: }
114: }
115: }
116: PetscFunctionReturn(PETSC_SUCCESS);
117: }
119: static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer)
120: {
121: SNES_NCG *ncg = (SNES_NCG *)snes->data;
122: PetscBool iascii;
124: PetscFunctionBegin;
125: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
126: if (iascii) PetscCall(PetscViewerASCIIPrintf(viewer, " type: %s\n", SNESNCGTypes[ncg->type]));
127: PetscFunctionReturn(PETSC_SUCCESS);
128: }
130: /*@
131: SNESNCGSetType - Sets the conjugate update type for nonlinear CG `SNESNCG`.
133: Logically Collective
135: Input Parameters:
136: + snes - the iterative context
137: - btype - update type, see `SNESNCGType`
139: Options Database Key:
140: . -snes_ncg_type <prp,fr,hs,dy,cd> - strategy for selecting algorithm for computing beta
142: Level: intermediate
144: Notes:
145: `SNES_NCG_PRP` is the default, and the only one that tolerates generalized search directions.
147: It is not clear what "generalized search directions" means, does it mean use with a nonlinear preconditioner,
148: that is using -npc_snes_type <type>, `SNESSetNPC()`, or `SNESGetNPC()`?
150: .seealso: [](ch_snes), `SNES`, `SNESNCG`, `SNESNCGType`, `SNES_NCG_FR`, `SNES_NCG_PRP`, `SNES_NCG_HS`, `SNES_NCG_DY`, `SNES_NCG_CD`
151: @*/
152: PetscErrorCode SNESNCGSetType(SNES snes, SNESNCGType btype)
153: {
154: PetscFunctionBegin;
156: PetscTryMethod(snes, "SNESNCGSetType_C", (SNES, SNESNCGType), (snes, btype));
157: PetscFunctionReturn(PETSC_SUCCESS);
158: }
160: static PetscErrorCode SNESNCGSetType_NCG(SNES snes, SNESNCGType btype)
161: {
162: SNES_NCG *ncg = (SNES_NCG *)snes->data;
164: PetscFunctionBegin;
165: ncg->type = btype;
166: PetscFunctionReturn(PETSC_SUCCESS);
167: }
169: /*
170: SNESSolve_NCG - Solves a nonlinear system with the Nonlinear Conjugate Gradient method.
172: Input Parameter:
173: . snes - the `SNES` context
175: Output Parameter:
176: . outits - number of iterations until termination
178: Application Interface Routine: SNESSolve()
179: */
180: static PetscErrorCode SNESSolve_NCG(SNES snes)
181: {
182: SNES_NCG *ncg = (SNES_NCG *)snes->data;
183: Vec X, dX, lX, F, dXold;
184: PetscReal fnorm, ynorm, xnorm, beta = 0.0;
185: PetscScalar dXdotdX, dXolddotdXold, dXdotdXold, lXdotdX, lXdotdXold;
186: PetscInt maxits, i;
187: SNESLineSearchReason lsresult = SNES_LINESEARCH_SUCCEEDED;
188: SNESLineSearch linesearch;
189: SNESConvergedReason reason;
191: PetscFunctionBegin;
192: PetscCheck(!snes->xl && !snes->xu && !snes->ops->computevariablebounds, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
194: PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite));
195: snes->reason = SNES_CONVERGED_ITERATING;
197: maxits = snes->max_its; /* maximum number of iterations */
198: X = snes->vec_sol; /* X^n */
199: dXold = snes->work[0]; /* The previous iterate of X */
200: dX = snes->work[1]; /* the preconditioned direction */
201: lX = snes->vec_sol_update; /* the conjugate direction */
202: F = snes->vec_func; /* residual vector */
204: PetscCall(SNESGetLineSearch(snes, &linesearch));
206: PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
207: snes->iter = 0;
208: snes->norm = 0.;
209: PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
211: /* compute the initial function and preconditioned update dX */
213: if (snes->npc && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
214: PetscCall(SNESApplyNPC(snes, X, NULL, dX));
215: PetscCall(SNESGetConvergedReason(snes->npc, &reason));
216: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
217: snes->reason = SNES_DIVERGED_INNER;
218: PetscFunctionReturn(PETSC_SUCCESS);
219: }
220: PetscCall(VecCopy(dX, F));
221: PetscCall(VecNorm(F, NORM_2, &fnorm));
222: } else {
223: if (!snes->vec_func_init_set) {
224: PetscCall(SNESComputeFunction(snes, X, F));
225: } else snes->vec_func_init_set = PETSC_FALSE;
227: /* convergence test */
228: PetscCall(VecNorm(F, NORM_2, &fnorm));
229: SNESCheckFunctionNorm(snes, fnorm);
230: PetscCall(VecCopy(F, dX));
231: }
232: if (snes->npc) {
233: if (snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
234: PetscCall(SNESApplyNPC(snes, X, F, dX));
235: PetscCall(SNESGetConvergedReason(snes->npc, &reason));
236: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
237: snes->reason = SNES_DIVERGED_INNER;
238: PetscFunctionReturn(PETSC_SUCCESS);
239: }
240: }
241: }
242: PetscCall(VecCopy(dX, lX));
243: PetscCall(VecDot(dX, dX, &dXdotdX));
245: PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
246: snes->norm = fnorm;
247: PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
248: PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
250: /* test convergence */
251: PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm));
252: PetscCall(SNESMonitor(snes, 0, fnorm));
253: if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
255: /* Call general purpose update function */
256: PetscTryTypeMethod(snes, update, snes->iter);
258: /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */
260: for (i = 1; i < maxits + 1; i++) {
261: /* some update types require the old update direction or conjugate direction */
262: if (ncg->type != SNES_NCG_FR) PetscCall(VecCopy(dX, dXold));
263: PetscCall(SNESLineSearchApply(linesearch, X, F, &fnorm, lX));
264: PetscCall(SNESLineSearchGetReason(linesearch, &lsresult));
265: PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm));
266: if (lsresult) {
267: if (++snes->numFailures >= snes->maxFailures) {
268: snes->reason = SNES_DIVERGED_LINE_SEARCH;
269: PetscFunctionReturn(PETSC_SUCCESS);
270: }
271: }
272: if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
273: snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
274: PetscFunctionReturn(PETSC_SUCCESS);
275: }
276: /* Monitor convergence */
277: PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
278: snes->iter = i;
279: snes->norm = fnorm;
280: snes->xnorm = xnorm;
281: snes->ynorm = ynorm;
282: PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
283: PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
285: /* Test for convergence */
286: PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm));
287: PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
288: if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
290: /* Call general purpose update function */
291: PetscTryTypeMethod(snes, update, snes->iter);
292: if (snes->npc) {
293: if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
294: PetscCall(SNESApplyNPC(snes, X, NULL, dX));
295: PetscCall(SNESGetConvergedReason(snes->npc, &reason));
296: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
297: snes->reason = SNES_DIVERGED_INNER;
298: PetscFunctionReturn(PETSC_SUCCESS);
299: }
300: PetscCall(VecCopy(dX, F));
301: } else {
302: PetscCall(SNESApplyNPC(snes, X, F, dX));
303: PetscCall(SNESGetConvergedReason(snes->npc, &reason));
304: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
305: snes->reason = SNES_DIVERGED_INNER;
306: PetscFunctionReturn(PETSC_SUCCESS);
307: }
308: }
309: } else {
310: PetscCall(VecCopy(F, dX));
311: }
313: /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/
314: switch (ncg->type) {
315: case SNES_NCG_FR: /* Fletcher-Reeves */
316: dXolddotdXold = dXdotdX;
317: PetscCall(VecDot(dX, dX, &dXdotdX));
318: beta = PetscRealPart(dXdotdX / dXolddotdXold);
319: break;
320: case SNES_NCG_PRP: /* Polak-Ribiere-Poylak */
321: dXolddotdXold = dXdotdX;
322: PetscCall(VecDotBegin(dX, dX, &dXdotdX));
323: PetscCall(VecDotBegin(dXold, dX, &dXdotdXold));
324: PetscCall(VecDotEnd(dX, dX, &dXdotdX));
325: PetscCall(VecDotEnd(dXold, dX, &dXdotdXold));
326: beta = PetscRealPart((dXdotdX - dXdotdXold) / dXolddotdXold);
327: if (beta < 0.0) beta = 0.0; /* restart */
328: break;
329: case SNES_NCG_HS: /* Hestenes-Stiefel */
330: PetscCall(VecDotBegin(dX, dX, &dXdotdX));
331: PetscCall(VecDotBegin(dX, dXold, &dXdotdXold));
332: PetscCall(VecDotBegin(lX, dX, &lXdotdX));
333: PetscCall(VecDotBegin(lX, dXold, &lXdotdXold));
334: PetscCall(VecDotEnd(dX, dX, &dXdotdX));
335: PetscCall(VecDotEnd(dX, dXold, &dXdotdXold));
336: PetscCall(VecDotEnd(lX, dX, &lXdotdX));
337: PetscCall(VecDotEnd(lX, dXold, &lXdotdXold));
338: beta = PetscRealPart((dXdotdX - dXdotdXold) / (lXdotdX - lXdotdXold));
339: break;
340: case SNES_NCG_DY: /* Dai-Yuan */
341: PetscCall(VecDotBegin(dX, dX, &dXdotdX));
342: PetscCall(VecDotBegin(lX, dX, &lXdotdX));
343: PetscCall(VecDotBegin(lX, dXold, &lXdotdXold));
344: PetscCall(VecDotEnd(dX, dX, &dXdotdX));
345: PetscCall(VecDotEnd(lX, dX, &lXdotdX));
346: PetscCall(VecDotEnd(lX, dXold, &lXdotdXold));
347: beta = PetscRealPart(dXdotdX / (lXdotdXold - lXdotdX));
348: break;
349: case SNES_NCG_CD: /* Conjugate Descent */
350: PetscCall(VecDotBegin(dX, dX, &dXdotdX));
351: PetscCall(VecDotBegin(lX, dXold, &lXdotdXold));
352: PetscCall(VecDotEnd(dX, dX, &dXdotdX));
353: PetscCall(VecDotEnd(lX, dXold, &lXdotdXold));
354: beta = PetscRealPart(dXdotdX / lXdotdXold);
355: break;
356: }
357: if (ncg->monitor) PetscCall(PetscViewerASCIIPrintf(ncg->monitor, "beta = %e\n", (double)beta));
358: PetscCall(VecAYPX(lX, beta, dX));
359: }
360: PetscFunctionReturn(PETSC_SUCCESS);
361: }
363: /*MC
364: SNESNCG - Nonlinear Conjugate-Gradient method for the solution of nonlinear systems {cite}`bruneknepleysmithtu15`.
366: Level: beginner
368: Options Database Keys:
369: + -snes_ncg_type <fr, prp, dy, hs, cd> - Choice of conjugate-gradient update parameter, default is prp.
370: . -snes_linesearch_type <cp,l2,basic> - Line search type.
371: - -snes_ncg_monitor - Print the beta values nonlinear Conjugate-Gradient used in the iteration, .
373: Notes:
374: This solves the nonlinear system of equations $ F(x) = 0 $ using the nonlinear generalization of the conjugate
375: gradient method. This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise
376: chooses the initial search direction as $ F(x) $ for the initial guess $x$.
378: Only supports left non-linear preconditioning.
380: Default line search is `SNESLINESEARCHCP`, unless a nonlinear preconditioner is used with -npc_snes_type <type>, `SNESSetNPC()`, or `SNESGetNPC()` then
381: `SNESLINESEARCHL2` is used. Also supports the special purpose line search `SNESLINESEARCHNCGLINEAR`
383: .seealso: [](ch_snes), `SNES`, `SNESNCG`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNEWTONTR`, `SNESNGMRES`, `SNESQN`, `SNESLINESEARCHNCGLINEAR`, `SNESNCGSetType()`, `SNESLineSearchSetType()`
384: M*/
385: PETSC_EXTERN PetscErrorCode SNESCreate_NCG(SNES snes)
386: {
387: SNES_NCG *neP;
389: PetscFunctionBegin;
390: snes->ops->destroy = SNESDestroy_NCG;
391: snes->ops->setup = SNESSetUp_NCG;
392: snes->ops->setfromoptions = SNESSetFromOptions_NCG;
393: snes->ops->view = SNESView_NCG;
394: snes->ops->solve = SNESSolve_NCG;
395: snes->ops->reset = SNESReset_NCG;
397: snes->usesksp = PETSC_FALSE;
398: snes->usesnpc = PETSC_TRUE;
399: snes->npcside = PC_LEFT;
401: snes->alwayscomputesfinalresidual = PETSC_TRUE;
403: PetscCall(SNESParametersInitialize(snes));
404: PetscObjectParameterSetDefault(snes, max_funcs, 30000);
405: PetscObjectParameterSetDefault(snes, max_its, 10000);
406: PetscObjectParameterSetDefault(snes, stol, 1e-20);
408: PetscCall(PetscNew(&neP));
409: snes->data = (void *)neP;
410: neP->monitor = NULL;
411: neP->type = SNES_NCG_PRP;
412: PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNCGSetType_C", SNESNCGSetType_NCG));
413: PetscFunctionReturn(PETSC_SUCCESS);
414: }