Actual source code: snesrichardson.c
1: #include <../src/snes/impls/richardson/snesrichardsonimpl.h>
3: static PetscErrorCode SNESDestroy_NRichardson(SNES snes)
4: {
5: PetscFunctionBegin;
6: PetscCall(PetscFree(snes->data));
7: PetscFunctionReturn(PETSC_SUCCESS);
8: }
10: static PetscErrorCode SNESSetUp_NRichardson(SNES snes)
11: {
12: PetscFunctionBegin;
13: PetscCheck(snes->npcside != PC_RIGHT, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "NRichardson only supports left preconditioning");
14: if (snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED;
15: PetscFunctionReturn(PETSC_SUCCESS);
16: }
18: static PetscErrorCode SNESSetFromOptions_NRichardson(SNES snes, PetscOptionItems PetscOptionsObject)
19: {
20: PetscFunctionBegin;
21: PetscOptionsHeadBegin(PetscOptionsObject, "SNES Richardson options");
22: PetscOptionsHeadEnd();
23: PetscFunctionReturn(PETSC_SUCCESS);
24: }
26: static PetscErrorCode SNESView_NRichardson(SNES snes, PetscViewer viewer)
27: {
28: PetscBool iascii;
30: PetscFunctionBegin;
31: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
32: if (iascii) { }
33: PetscFunctionReturn(PETSC_SUCCESS);
34: }
36: static PetscErrorCode SNESSolve_NRichardson(SNES snes)
37: {
38: Vec X, Y, F;
39: PetscReal xnorm, fnorm, ynorm;
40: PetscInt maxits, i;
41: SNESLineSearchReason lsresult;
42: SNESConvergedReason reason;
44: PetscFunctionBegin;
45: 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);
47: snes->reason = SNES_CONVERGED_ITERATING;
49: maxits = snes->max_its; /* maximum number of iterations */
50: X = snes->vec_sol; /* X^n */
51: Y = snes->vec_sol_update; /* \tilde X */
52: F = snes->vec_func; /* residual vector */
54: PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
55: snes->iter = 0;
56: snes->norm = 0.;
57: PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
59: if (snes->npc && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
60: PetscCall(SNESApplyNPC(snes, X, NULL, F));
61: PetscCall(SNESGetConvergedReason(snes->npc, &reason));
62: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
63: snes->reason = SNES_DIVERGED_INNER;
64: PetscFunctionReturn(PETSC_SUCCESS);
65: }
66: PetscCall(VecNorm(F, NORM_2, &fnorm));
67: } else {
68: if (!snes->vec_func_init_set) {
69: PetscCall(SNESComputeFunction(snes, X, F));
70: } else snes->vec_func_init_set = PETSC_FALSE;
72: PetscCall(VecNorm(F, NORM_2, &fnorm));
73: SNESCheckFunctionNorm(snes, fnorm);
74: }
75: if (snes->npc && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
76: PetscCall(SNESApplyNPC(snes, X, F, Y));
77: PetscCall(SNESGetConvergedReason(snes->npc, &reason));
78: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
79: snes->reason = SNES_DIVERGED_INNER;
80: PetscFunctionReturn(PETSC_SUCCESS);
81: }
82: } else {
83: PetscCall(VecCopy(F, Y));
84: }
86: PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
87: snes->norm = fnorm;
88: PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
89: PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
91: /* test convergence */
92: PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm));
93: PetscCall(SNESMonitor(snes, 0, fnorm));
94: if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
96: /* Call general purpose update function */
97: PetscTryTypeMethod(snes, update, snes->iter);
99: for (i = 1; i < maxits + 1; i++) {
100: PetscCall(SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y));
101: PetscCall(SNESLineSearchGetReason(snes->linesearch, &lsresult));
102: PetscCall(SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm));
103: if (lsresult) {
104: if (++snes->numFailures >= snes->maxFailures) {
105: snes->reason = SNES_DIVERGED_LINE_SEARCH;
106: break;
107: }
108: }
109: if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
110: snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
111: break;
112: }
114: /* Monitor convergence */
115: PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
116: snes->iter = i;
117: snes->norm = fnorm;
118: snes->xnorm = xnorm;
119: snes->ynorm = ynorm;
120: PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
121: PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
122: /* Test for convergence */
123: PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm));
124: PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
125: if (snes->reason) break;
127: /* Call general purpose update function */
128: PetscTryTypeMethod(snes, update, snes->iter);
130: if (snes->npc) {
131: if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
132: PetscCall(SNESApplyNPC(snes, X, NULL, Y));
133: PetscCall(VecNorm(F, NORM_2, &fnorm));
134: PetscCall(VecCopy(Y, F));
135: } else {
136: PetscCall(SNESApplyNPC(snes, X, F, Y));
137: }
138: PetscCall(SNESGetConvergedReason(snes->npc, &reason));
139: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
140: snes->reason = SNES_DIVERGED_INNER;
141: PetscFunctionReturn(PETSC_SUCCESS);
142: }
143: } else {
144: PetscCall(VecCopy(F, Y));
145: }
146: }
147: if (i == maxits + 1) {
148: PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits));
149: if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
150: }
151: PetscFunctionReturn(PETSC_SUCCESS);
152: }
154: /*MC
155: SNESNRICHARDSON - Richardson nonlinear solver that uses successive substitutions, also sometimes known as Picard iteration.
157: Options Database Keys:
158: + -snes_linesearch_type <l2,cp,basic> - Line search type.
159: - -snes_linesearch_damping <1.0> - Damping for the line search.
161: Level: beginner
163: Notes:
164: If no inner nonlinear preconditioner is provided then solves $F(x) - b = 0 $ using $x^{n+1} = x^{n} - \lambda
165: (F(x^n) - b) $ where $ \lambda$ is obtained with either `SNESLineSearchSetDamping()`, `-snes_damping` or a line search. If
166: an inner nonlinear preconditioner is provided (either with `-npc_snes_typ`e or `SNESSetNPC()`) then the inner
167: solver is called on the initial solution $x^n$ and the nonlinear Richardson uses $ x^{n+1} = x^{n} + \lambda d^{n}$
168: where $d^{n} = \hat{x}^{n} - x^{n} $ where $\hat{x}^{n} $ is the solution returned from the inner solver.
170: The update, especially without inner nonlinear preconditioner, may be ill-scaled. If using the basic
171: linesearch, one may have to scale the update with `-snes_linesearch_damping`
173: This uses no derivative information provided with `SNESSetJacobian()` thus it will be much slower then Newton's method obtained with `-snes_type ls`
175: Only supports left non-linear preconditioning.
177: .seealso: [](ch_snes), `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNEWTONTR`, `SNESNGMRES`, `SNESQN`, `SNESNCG`,
178: `SNESLineSearchSetDamping()`
179: M*/
180: PETSC_EXTERN PetscErrorCode SNESCreate_NRichardson(SNES snes)
181: {
182: SNES_NRichardson *neP;
183: SNESLineSearch linesearch;
185: PetscFunctionBegin;
186: snes->ops->destroy = SNESDestroy_NRichardson;
187: snes->ops->setup = SNESSetUp_NRichardson;
188: snes->ops->setfromoptions = SNESSetFromOptions_NRichardson;
189: snes->ops->view = SNESView_NRichardson;
190: snes->ops->solve = SNESSolve_NRichardson;
192: snes->usesksp = PETSC_FALSE;
193: snes->usesnpc = PETSC_TRUE;
195: snes->npcside = PC_LEFT;
197: PetscCall(SNESGetLineSearch(snes, &linesearch));
198: if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHL2));
200: snes->alwayscomputesfinalresidual = PETSC_TRUE;
202: PetscCall(SNESParametersInitialize(snes));
203: PetscObjectParameterSetDefault(snes, max_funcs, 30000);
204: PetscObjectParameterSetDefault(snes, max_its, 10000);
205: PetscObjectParameterSetDefault(snes, stol, 1e-20);
207: PetscCall(PetscNew(&neP));
208: snes->data = (void *)neP;
209: PetscFunctionReturn(PETSC_SUCCESS);
210: }