xref: /petsc/src/snes/impls/richardson/snesrichardson.c (revision 0b4b7b1c20c2ed4ade67e3d50a7710fe0ffbfca5)
142f4f86dSBarry Smith #include <../src/snes/impls/richardson/snesrichardsonimpl.h>
2028b6e76SBarry Smith 
366976f2fSJacob Faibussowitsch static PetscErrorCode SNESReset_NRichardson(SNES snes)
4d71ae5a4SJacob Faibussowitsch {
5028b6e76SBarry Smith   PetscFunctionBegin;
63ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7028b6e76SBarry Smith }
8028b6e76SBarry Smith 
966976f2fSJacob Faibussowitsch static PetscErrorCode SNESDestroy_NRichardson(SNES snes)
10d71ae5a4SJacob Faibussowitsch {
11028b6e76SBarry Smith   PetscFunctionBegin;
129566063dSJacob Faibussowitsch   PetscCall(SNESReset_NRichardson(snes));
139566063dSJacob Faibussowitsch   PetscCall(PetscFree(snes->data));
143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15028b6e76SBarry Smith }
16028b6e76SBarry Smith 
1766976f2fSJacob Faibussowitsch static PetscErrorCode SNESSetUp_NRichardson(SNES snes)
18d71ae5a4SJacob Faibussowitsch {
19028b6e76SBarry Smith   PetscFunctionBegin;
2008401ef6SPierre Jolivet   PetscCheck(snes->npcside != PC_RIGHT, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "NRichardson only supports left preconditioning");
216c67d002SPeter Brune   if (snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED;
223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
23028b6e76SBarry Smith }
24028b6e76SBarry Smith 
25d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetFromOptions_NRichardson(SNES snes, PetscOptionItems *PetscOptionsObject)
26d71ae5a4SJacob Faibussowitsch {
27028b6e76SBarry Smith   PetscFunctionBegin;
28d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SNES Richardson options");
29d0609cedSBarry Smith   PetscOptionsHeadEnd();
303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
31028b6e76SBarry Smith }
32028b6e76SBarry Smith 
33d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_NRichardson(SNES snes, PetscViewer viewer)
34d71ae5a4SJacob Faibussowitsch {
35028b6e76SBarry Smith   PetscBool iascii;
36028b6e76SBarry Smith 
37028b6e76SBarry Smith   PetscFunctionBegin;
389566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
399371c9d4SSatish Balay   if (iascii) { }
403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
41028b6e76SBarry Smith }
42028b6e76SBarry Smith 
4366976f2fSJacob Faibussowitsch static PetscErrorCode SNESSolve_NRichardson(SNES snes)
44d71ae5a4SJacob Faibussowitsch {
45bf7f4e0aSPeter Brune   Vec                  X, Y, F;
46e7058c64SPeter Brune   PetscReal            xnorm, fnorm, ynorm;
47028b6e76SBarry Smith   PetscInt             maxits, i;
48422a814eSBarry Smith   SNESLineSearchReason lsresult;
49028b6e76SBarry Smith   SNESConvergedReason  reason;
50028b6e76SBarry Smith 
51028b6e76SBarry Smith   PetscFunctionBegin;
520b121fc5SBarry 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);
53c579b300SPatrick Farrell 
54028b6e76SBarry Smith   snes->reason = SNES_CONVERGED_ITERATING;
55028b6e76SBarry Smith 
56028b6e76SBarry Smith   maxits = snes->max_its;        /* maximum number of iterations */
57028b6e76SBarry Smith   X      = snes->vec_sol;        /* X^n */
58028b6e76SBarry Smith   Y      = snes->vec_sol_update; /* \tilde X */
59028b6e76SBarry Smith   F      = snes->vec_func;       /* residual vector */
60028b6e76SBarry Smith 
619566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
62028b6e76SBarry Smith   snes->iter = 0;
63028b6e76SBarry Smith   snes->norm = 0.;
649566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
65b7281c8aSPeter Brune 
66efd4aadfSBarry Smith   if (snes->npc && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
679566063dSJacob Faibussowitsch     PetscCall(SNESApplyNPC(snes, X, NULL, F));
689566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(snes->npc, &reason));
69b7281c8aSPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
70b7281c8aSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
713ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
72b7281c8aSPeter Brune     }
739566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, &fnorm));
74b7281c8aSPeter Brune   } else {
75e4ed7901SPeter Brune     if (!snes->vec_func_init_set) {
769566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, X, F));
771aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
78c1c75074SPeter Brune 
799566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, &fnorm));
80422a814eSBarry Smith     SNESCheckFunctionNorm(snes, fnorm);
81b7281c8aSPeter Brune   }
82efd4aadfSBarry Smith   if (snes->npc && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
839566063dSJacob Faibussowitsch     PetscCall(SNESApplyNPC(snes, X, F, Y));
849566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(snes->npc, &reason));
85b7281c8aSPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
86b7281c8aSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
873ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
88b7281c8aSPeter Brune     }
89b7281c8aSPeter Brune   } else {
909566063dSJacob Faibussowitsch     PetscCall(VecCopy(F, Y));
91b7281c8aSPeter Brune   }
92e4ed7901SPeter Brune 
939566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
94028b6e76SBarry Smith   snes->norm = fnorm;
959566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
969566063dSJacob Faibussowitsch   PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
97028b6e76SBarry Smith 
98028b6e76SBarry Smith   /* test convergence */
992d157150SStefano Zampini   PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm));
1002d157150SStefano Zampini   PetscCall(SNESMonitor(snes, 0, fnorm));
1013ba16761SJacob Faibussowitsch   if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
102028b6e76SBarry Smith 
103028b6e76SBarry Smith   /* Call general purpose update function */
104dbbe0bcdSBarry Smith   PetscTryTypeMethod(snes, update, snes->iter);
105b7281c8aSPeter Brune 
106b7281c8aSPeter Brune   for (i = 1; i < maxits + 1; i++) {
1079566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y));
1089566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetReason(snes->linesearch, &lsresult));
1099566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm));
110422a814eSBarry Smith     if (lsresult) {
111028b6e76SBarry Smith       if (++snes->numFailures >= snes->maxFailures) {
112028b6e76SBarry Smith         snes->reason = SNES_DIVERGED_LINE_SEARCH;
113028b6e76SBarry Smith         break;
114028b6e76SBarry Smith       }
115028b6e76SBarry Smith     }
116e71169deSBarry Smith     if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
117028b6e76SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
118028b6e76SBarry Smith       break;
119028b6e76SBarry Smith     }
120b7281c8aSPeter Brune 
121028b6e76SBarry Smith     /* Monitor convergence */
1229566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
123b7281c8aSPeter Brune     snes->iter  = i;
124028b6e76SBarry Smith     snes->norm  = fnorm;
125c1e67a49SFande Kong     snes->xnorm = xnorm;
126c1e67a49SFande Kong     snes->ynorm = ynorm;
1279566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
1289566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
129028b6e76SBarry Smith     /* Test for convergence */
1302d157150SStefano Zampini     PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm));
1312d157150SStefano Zampini     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
132028b6e76SBarry Smith     if (snes->reason) break;
133b7281c8aSPeter Brune 
134b7281c8aSPeter Brune     /* Call general purpose update function */
135dbbe0bcdSBarry Smith     PetscTryTypeMethod(snes, update, snes->iter);
136b7281c8aSPeter Brune 
137efd4aadfSBarry Smith     if (snes->npc) {
138b7281c8aSPeter Brune       if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
1399566063dSJacob Faibussowitsch         PetscCall(SNESApplyNPC(snes, X, NULL, Y));
1409566063dSJacob Faibussowitsch         PetscCall(VecNorm(F, NORM_2, &fnorm));
1419566063dSJacob Faibussowitsch         PetscCall(VecCopy(Y, F));
142b7281c8aSPeter Brune       } else {
1439566063dSJacob Faibussowitsch         PetscCall(SNESApplyNPC(snes, X, F, Y));
144b7281c8aSPeter Brune       }
1459566063dSJacob Faibussowitsch       PetscCall(SNESGetConvergedReason(snes->npc, &reason));
146b7281c8aSPeter Brune       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
147b7281c8aSPeter Brune         snes->reason = SNES_DIVERGED_INNER;
1483ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
149b7281c8aSPeter Brune       }
150b7281c8aSPeter Brune     } else {
1519566063dSJacob Faibussowitsch       PetscCall(VecCopy(F, Y));
152b7281c8aSPeter Brune     }
153b7281c8aSPeter Brune   }
154b7281c8aSPeter Brune   if (i == maxits + 1) {
15563a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits));
156028b6e76SBarry Smith     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
157028b6e76SBarry Smith   }
1583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
159028b6e76SBarry Smith }
160028b6e76SBarry Smith 
161028b6e76SBarry Smith /*MC
162d5c3842bSBarry Smith    SNESNRICHARDSON - Richardson nonlinear solver that uses successive substitutions, also sometimes known as Picard iteration.
163028b6e76SBarry Smith 
164f6dfbefdSBarry Smith    Options Database Keys:
16567b8a455SSatish Balay +  -snes_linesearch_type <l2,cp,basic> - Line search type.
16667b8a455SSatish Balay -  -snes_linesearch_damping <1.0>      - Damping for the line search.
167028b6e76SBarry Smith 
168f6dfbefdSBarry Smith    Level: beginner
169f6dfbefdSBarry Smith 
17095452b02SPatrick Sanan    Notes:
171*0b4b7b1cSBarry Smith    If no inner nonlinear preconditioner is provided then solves $F(x) - b = 0 $ using $x^{n+1} = x^{n} - \lambda
172*0b4b7b1cSBarry Smith    (F(x^n) - b) $ where $ \lambda$ is obtained with either `SNESLineSearchSetDamping()`, `-snes_damping` or a line search.  If
173*0b4b7b1cSBarry Smith    an inner nonlinear preconditioner is provided (either with `-npc_snes_typ`e or `SNESSetNPC()`) then the inner
174*0b4b7b1cSBarry Smith    solver is called on the initial solution $x^n$ and the nonlinear Richardson uses $ x^{n+1} = x^{n} + \lambda d^{n}$
175*0b4b7b1cSBarry Smith    where $d^{n} = \hat{x}^{n} - x^{n} $ where $\hat{x}^{n} $ is the solution returned from the inner solver.
176028b6e76SBarry Smith 
17772835e02SPeter Brune    The update, especially without inner nonlinear preconditioner, may be ill-scaled.  If using the basic
178*0b4b7b1cSBarry Smith    linesearch, one may have to scale the update with `-snes_linesearch_damping`
17972835e02SPeter Brune 
180*0b4b7b1cSBarry Smith    This uses no derivative information provided with `SNESSetJacobian()` thus it will be much slower then Newton's method obtained with `-snes_type ls`
181028b6e76SBarry Smith 
1822d547940SBarry Smith    Only supports left non-linear preconditioning.
1832d547940SBarry Smith 
184420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNEWTONTR`, `SNESNGMRES`, `SNESQN`, `SNESNCG`,
185420bcc1bSBarry Smith           `SNESLineSearchSetDamping()`
186028b6e76SBarry Smith M*/
187d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NRichardson(SNES snes)
188d71ae5a4SJacob Faibussowitsch {
189bf7f4e0aSPeter Brune   SNES_NRichardson *neP;
190d8d34be6SBarry Smith   SNESLineSearch    linesearch;
1915fd66863SKarl Rupp 
192028b6e76SBarry Smith   PetscFunctionBegin;
193d5c3842bSBarry Smith   snes->ops->destroy        = SNESDestroy_NRichardson;
194d5c3842bSBarry Smith   snes->ops->setup          = SNESSetUp_NRichardson;
195d5c3842bSBarry Smith   snes->ops->setfromoptions = SNESSetFromOptions_NRichardson;
196d5c3842bSBarry Smith   snes->ops->view           = SNESView_NRichardson;
197d5c3842bSBarry Smith   snes->ops->solve          = SNESSolve_NRichardson;
198d5c3842bSBarry Smith   snes->ops->reset          = SNESReset_NRichardson;
199028b6e76SBarry Smith 
200028b6e76SBarry Smith   snes->usesksp = PETSC_FALSE;
201efd4aadfSBarry Smith   snes->usesnpc = PETSC_TRUE;
202028b6e76SBarry Smith 
203efd4aadfSBarry Smith   snes->npcside = PC_LEFT;
204c6b63b32SPeter Brune 
2059566063dSJacob Faibussowitsch   PetscCall(SNESGetLineSearch(snes, &linesearch));
20648a46eb9SPierre Jolivet   if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHL2));
207d8d34be6SBarry Smith 
2084fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
2094fc747eaSLawrence Mitchell 
21077e5a1f9SBarry Smith   PetscCall(SNESParametersInitialize(snes));
21177e5a1f9SBarry Smith   PetscObjectParameterSetDefault(snes, max_funcs, 30000);
21277e5a1f9SBarry Smith   PetscObjectParameterSetDefault(snes, max_its, 10000);
21377e5a1f9SBarry Smith   PetscObjectParameterSetDefault(snes, stol, 1e-20);
21477e5a1f9SBarry Smith 
2154dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&neP));
216bf7f4e0aSPeter Brune   snes->data = (void *)neP;
2173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
218028b6e76SBarry Smith }
219