142f4f86dSBarry Smith #include <../src/snes/impls/richardson/snesrichardsonimpl.h>
2028b6e76SBarry Smith
SNESDestroy_NRichardson(SNES snes)366976f2fSJacob Faibussowitsch static PetscErrorCode SNESDestroy_NRichardson(SNES snes)
4d71ae5a4SJacob Faibussowitsch {
5028b6e76SBarry Smith PetscFunctionBegin;
69566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data));
73ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
8028b6e76SBarry Smith }
9028b6e76SBarry Smith
SNESSetUp_NRichardson(SNES snes)1066976f2fSJacob Faibussowitsch static PetscErrorCode SNESSetUp_NRichardson(SNES snes)
11d71ae5a4SJacob Faibussowitsch {
12028b6e76SBarry Smith PetscFunctionBegin;
1308401ef6SPierre Jolivet PetscCheck(snes->npcside != PC_RIGHT, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "NRichardson only supports left preconditioning");
146c67d002SPeter Brune if (snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED;
153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
16028b6e76SBarry Smith }
17028b6e76SBarry Smith
SNESSetFromOptions_NRichardson(SNES snes,PetscOptionItems PetscOptionsObject)18ce78bad3SBarry Smith static PetscErrorCode SNESSetFromOptions_NRichardson(SNES snes, PetscOptionItems PetscOptionsObject)
19d71ae5a4SJacob Faibussowitsch {
20028b6e76SBarry Smith PetscFunctionBegin;
21d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SNES Richardson options");
22d0609cedSBarry Smith PetscOptionsHeadEnd();
233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
24028b6e76SBarry Smith }
25028b6e76SBarry Smith
SNESSolve_NRichardson(SNES snes)2666976f2fSJacob Faibussowitsch static PetscErrorCode SNESSolve_NRichardson(SNES snes)
27d71ae5a4SJacob Faibussowitsch {
28bf7f4e0aSPeter Brune Vec X, Y, F;
29e7058c64SPeter Brune PetscReal xnorm, fnorm, ynorm;
30028b6e76SBarry Smith PetscInt maxits, i;
31028b6e76SBarry Smith SNESConvergedReason reason;
32028b6e76SBarry Smith
33028b6e76SBarry Smith PetscFunctionBegin;
340b121fc5SBarry Smith 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);
35c579b300SPatrick Farrell
36028b6e76SBarry Smith snes->reason = SNES_CONVERGED_ITERATING;
37028b6e76SBarry Smith
38028b6e76SBarry Smith maxits = snes->max_its; /* maximum number of iterations */
39028b6e76SBarry Smith X = snes->vec_sol; /* X^n */
40028b6e76SBarry Smith Y = snes->vec_sol_update; /* \tilde X */
41028b6e76SBarry Smith F = snes->vec_func; /* residual vector */
42028b6e76SBarry Smith
439566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
44028b6e76SBarry Smith snes->iter = 0;
45028b6e76SBarry Smith snes->norm = 0.;
469566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
47b7281c8aSPeter Brune
48efd4aadfSBarry Smith if (snes->npc && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
499566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes, X, NULL, F));
509566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason));
51b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
52b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER;
533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
54b7281c8aSPeter Brune }
559566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm));
56b7281c8aSPeter Brune } else {
57ac530a7eSPierre Jolivet if (!snes->vec_func_init_set) PetscCall(SNESComputeFunction(snes, X, F));
58ac530a7eSPierre Jolivet else snes->vec_func_init_set = PETSC_FALSE;
59c1c75074SPeter Brune
609566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm));
61*76c63389SBarry Smith SNESCheckFunctionDomainError(snes, fnorm);
62b7281c8aSPeter Brune }
63efd4aadfSBarry Smith if (snes->npc && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
649566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes, X, F, Y));
659566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason));
66b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
67b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER;
683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
69b7281c8aSPeter Brune }
70b7281c8aSPeter Brune } else {
719566063dSJacob Faibussowitsch PetscCall(VecCopy(F, Y));
72b7281c8aSPeter Brune }
73e4ed7901SPeter Brune
749566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
75028b6e76SBarry Smith snes->norm = fnorm;
769566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
779566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
78028b6e76SBarry Smith
79028b6e76SBarry Smith /* test convergence */
802d157150SStefano Zampini PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm));
812d157150SStefano Zampini PetscCall(SNESMonitor(snes, 0, fnorm));
823ba16761SJacob Faibussowitsch if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
83028b6e76SBarry Smith
84028b6e76SBarry Smith /* Call general purpose update function */
85dbbe0bcdSBarry Smith PetscTryTypeMethod(snes, update, snes->iter);
86b7281c8aSPeter Brune
87b7281c8aSPeter Brune for (i = 1; i < maxits + 1; i++) {
889566063dSJacob Faibussowitsch PetscCall(SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y));
89*76c63389SBarry Smith if (snes->reason) break;
90*76c63389SBarry Smith SNESCheckLineSearchFailure(snes);
919566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm));
92e71169deSBarry Smith if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
93028b6e76SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
94028b6e76SBarry Smith break;
95028b6e76SBarry Smith }
96b7281c8aSPeter Brune
97028b6e76SBarry Smith /* Monitor convergence */
989566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
99b7281c8aSPeter Brune snes->iter = i;
100028b6e76SBarry Smith snes->norm = fnorm;
101c1e67a49SFande Kong snes->xnorm = xnorm;
102c1e67a49SFande Kong snes->ynorm = ynorm;
1039566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
1049566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
105028b6e76SBarry Smith /* Test for convergence */
1062d157150SStefano Zampini PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm));
1072d157150SStefano Zampini PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
108028b6e76SBarry Smith if (snes->reason) break;
109b7281c8aSPeter Brune
110b7281c8aSPeter Brune /* Call general purpose update function */
111dbbe0bcdSBarry Smith PetscTryTypeMethod(snes, update, snes->iter);
112b7281c8aSPeter Brune
113efd4aadfSBarry Smith if (snes->npc) {
114b7281c8aSPeter Brune if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
1159566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes, X, NULL, Y));
1169566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm));
1179566063dSJacob Faibussowitsch PetscCall(VecCopy(Y, F));
118b7281c8aSPeter Brune } else {
1199566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes, X, F, Y));
120b7281c8aSPeter Brune }
1219566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason));
122b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
123b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER;
1243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
125b7281c8aSPeter Brune }
126b7281c8aSPeter Brune } else {
1279566063dSJacob Faibussowitsch PetscCall(VecCopy(F, Y));
128b7281c8aSPeter Brune }
129b7281c8aSPeter Brune }
130b7281c8aSPeter Brune if (i == maxits + 1) {
13163a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits));
132028b6e76SBarry Smith if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
133028b6e76SBarry Smith }
1343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
135028b6e76SBarry Smith }
136028b6e76SBarry Smith
137028b6e76SBarry Smith /*MC
138d5c3842bSBarry Smith SNESNRICHARDSON - Richardson nonlinear solver that uses successive substitutions, also sometimes known as Picard iteration.
139028b6e76SBarry Smith
140f6dfbefdSBarry Smith Options Database Keys:
14167b8a455SSatish Balay + -snes_linesearch_type <l2,cp,basic> - Line search type.
14267b8a455SSatish Balay - -snes_linesearch_damping <1.0> - Damping for the line search.
143028b6e76SBarry Smith
144f6dfbefdSBarry Smith Level: beginner
145f6dfbefdSBarry Smith
14695452b02SPatrick Sanan Notes:
1470b4b7b1cSBarry Smith If no inner nonlinear preconditioner is provided then solves $F(x) - b = 0 $ using $x^{n+1} = x^{n} - \lambda
1480b4b7b1cSBarry Smith (F(x^n) - b) $ where $ \lambda$ is obtained with either `SNESLineSearchSetDamping()`, `-snes_damping` or a line search. If
1490b4b7b1cSBarry Smith an inner nonlinear preconditioner is provided (either with `-npc_snes_typ`e or `SNESSetNPC()`) then the inner
1500b4b7b1cSBarry Smith solver is called on the initial solution $x^n$ and the nonlinear Richardson uses $ x^{n+1} = x^{n} + \lambda d^{n}$
1510b4b7b1cSBarry Smith where $d^{n} = \hat{x}^{n} - x^{n} $ where $\hat{x}^{n} $ is the solution returned from the inner solver.
152028b6e76SBarry Smith
15372835e02SPeter Brune The update, especially without inner nonlinear preconditioner, may be ill-scaled. If using the basic
1540b4b7b1cSBarry Smith linesearch, one may have to scale the update with `-snes_linesearch_damping`
15572835e02SPeter Brune
156e87b5d96SPierre Jolivet This uses no derivative information provided with `SNESSetJacobian()` thus it will be much slower than Newton's method obtained with `-snes_type ls`
157028b6e76SBarry Smith
1582d547940SBarry Smith Only supports left non-linear preconditioning.
1592d547940SBarry Smith
160420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNEWTONTR`, `SNESNGMRES`, `SNESQN`, `SNESNCG`,
161420bcc1bSBarry Smith `SNESLineSearchSetDamping()`
162028b6e76SBarry Smith M*/
SNESCreate_NRichardson(SNES snes)163d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NRichardson(SNES snes)
164d71ae5a4SJacob Faibussowitsch {
165bf7f4e0aSPeter Brune SNES_NRichardson *neP;
166d8d34be6SBarry Smith SNESLineSearch linesearch;
1675fd66863SKarl Rupp
168028b6e76SBarry Smith PetscFunctionBegin;
169d5c3842bSBarry Smith snes->ops->destroy = SNESDestroy_NRichardson;
170d5c3842bSBarry Smith snes->ops->setup = SNESSetUp_NRichardson;
171d5c3842bSBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_NRichardson;
172d5c3842bSBarry Smith snes->ops->solve = SNESSolve_NRichardson;
173028b6e76SBarry Smith
174028b6e76SBarry Smith snes->usesksp = PETSC_FALSE;
175efd4aadfSBarry Smith snes->usesnpc = PETSC_TRUE;
176028b6e76SBarry Smith
177efd4aadfSBarry Smith snes->npcside = PC_LEFT;
178c6b63b32SPeter Brune
1799566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(snes, &linesearch));
180a99ef635SJonas Heinzmann if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHSECANT));
181d8d34be6SBarry Smith
1824fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_TRUE;
1834fc747eaSLawrence Mitchell
18477e5a1f9SBarry Smith PetscCall(SNESParametersInitialize(snes));
18577e5a1f9SBarry Smith PetscObjectParameterSetDefault(snes, max_funcs, 30000);
18677e5a1f9SBarry Smith PetscObjectParameterSetDefault(snes, max_its, 10000);
18777e5a1f9SBarry Smith PetscObjectParameterSetDefault(snes, stol, 1e-20);
18877e5a1f9SBarry Smith
1894dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&neP));
190bf7f4e0aSPeter Brune snes->data = (void *)neP;
1913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
192028b6e76SBarry Smith }
193