xref: /petsc/src/snes/impls/richardson/snesrichardson.c (revision 420bcc1b905230dede3c88f397d1a4e60493adde)
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   /* set parameter for default relative tolerance convergence test */
107b7281c8aSPeter Brune   snes->ttol = fnorm * snes->rtol;
108b7281c8aSPeter Brune   /* test convergence */
1092d157150SStefano Zampini   PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm));
1103ba16761SJacob Faibussowitsch   if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
111b7281c8aSPeter Brune 
112b7281c8aSPeter Brune   for (i = 1; i < maxits + 1; i++) {
1139566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y));
1149566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetReason(snes->linesearch, &lsresult));
1159566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm));
116422a814eSBarry Smith     if (lsresult) {
117028b6e76SBarry Smith       if (++snes->numFailures >= snes->maxFailures) {
118028b6e76SBarry Smith         snes->reason = SNES_DIVERGED_LINE_SEARCH;
119028b6e76SBarry Smith         break;
120028b6e76SBarry Smith       }
121028b6e76SBarry Smith     }
122e71169deSBarry Smith     if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
123028b6e76SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
124028b6e76SBarry Smith       break;
125028b6e76SBarry Smith     }
126b7281c8aSPeter Brune 
127028b6e76SBarry Smith     /* Monitor convergence */
1289566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
129b7281c8aSPeter Brune     snes->iter  = i;
130028b6e76SBarry Smith     snes->norm  = fnorm;
131c1e67a49SFande Kong     snes->xnorm = xnorm;
132c1e67a49SFande Kong     snes->ynorm = ynorm;
1339566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
1349566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
135028b6e76SBarry Smith     /* Test for convergence */
1362d157150SStefano Zampini     PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm));
1372d157150SStefano Zampini     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
138028b6e76SBarry Smith     if (snes->reason) break;
139b7281c8aSPeter Brune 
140b7281c8aSPeter Brune     /* Call general purpose update function */
141dbbe0bcdSBarry Smith     PetscTryTypeMethod(snes, update, snes->iter);
142b7281c8aSPeter Brune 
143efd4aadfSBarry Smith     if (snes->npc) {
144b7281c8aSPeter Brune       if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
1459566063dSJacob Faibussowitsch         PetscCall(SNESApplyNPC(snes, X, NULL, Y));
1469566063dSJacob Faibussowitsch         PetscCall(VecNorm(F, NORM_2, &fnorm));
1479566063dSJacob Faibussowitsch         PetscCall(VecCopy(Y, F));
148b7281c8aSPeter Brune       } else {
1499566063dSJacob Faibussowitsch         PetscCall(SNESApplyNPC(snes, X, F, Y));
150b7281c8aSPeter Brune       }
1519566063dSJacob Faibussowitsch       PetscCall(SNESGetConvergedReason(snes->npc, &reason));
152b7281c8aSPeter Brune       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
153b7281c8aSPeter Brune         snes->reason = SNES_DIVERGED_INNER;
1543ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
155b7281c8aSPeter Brune       }
156b7281c8aSPeter Brune     } else {
1579566063dSJacob Faibussowitsch       PetscCall(VecCopy(F, Y));
158b7281c8aSPeter Brune     }
159b7281c8aSPeter Brune   }
160b7281c8aSPeter Brune   if (i == maxits + 1) {
16163a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits));
162028b6e76SBarry Smith     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
163028b6e76SBarry Smith   }
1643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
165028b6e76SBarry Smith }
166028b6e76SBarry Smith 
167028b6e76SBarry Smith /*MC
168d5c3842bSBarry Smith    SNESNRICHARDSON - Richardson nonlinear solver that uses successive substitutions, also sometimes known as Picard iteration.
169028b6e76SBarry Smith 
170f6dfbefdSBarry Smith    Options Database Keys:
17167b8a455SSatish Balay +  -snes_linesearch_type <l2,cp,basic> - Line search type.
17267b8a455SSatish Balay -  -snes_linesearch_damping<1.0> - Damping for the line search.
173028b6e76SBarry Smith 
174f6dfbefdSBarry Smith    Level: beginner
175f6dfbefdSBarry Smith 
17695452b02SPatrick Sanan    Notes:
17795452b02SPatrick Sanan    If no inner nonlinear preconditioner is provided then solves F(x) - b = 0 using x^{n+1} = x^{n} - lambda
178f6dfbefdSBarry Smith    (F(x^n) - b) where lambda is obtained either `SNESLineSearchSetDamping()`, -snes_damping or a line search.  If
1790462cc06SPierre Jolivet    an inner nonlinear preconditioner is provided (either with -npc_snes_type or `SNESSetNPC()`) then the inner
180af60355fSPeter Brune    solver is called an initial solution x^n and the nonlinear Richardson uses x^{n+1} = x^{n} + lambda d^{n}
181af60355fSPeter Brune    where d^{n} = \hat{x}^{n} - x^{n} where \hat{x}^{n} is the solution returned from the inner solver.
182028b6e76SBarry Smith 
18372835e02SPeter Brune    The update, especially without inner nonlinear preconditioner, may be ill-scaled.  If using the basic
18472835e02SPeter Brune    linesearch, one may have to scale the update with -snes_linesearch_damping
18572835e02SPeter Brune 
186028b6e76SBarry Smith    This uses no derivative information thus will be much slower then Newton's method obtained with -snes_type ls
187028b6e76SBarry Smith 
1882d547940SBarry Smith    Only supports left non-linear preconditioning.
1892d547940SBarry Smith 
190*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNEWTONTR`, `SNESNGMRES`, `SNESQN`, `SNESNCG`,
191*420bcc1bSBarry Smith           `SNESLineSearchSetDamping()`
192028b6e76SBarry Smith M*/
193d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NRichardson(SNES snes)
194d71ae5a4SJacob Faibussowitsch {
195bf7f4e0aSPeter Brune   SNES_NRichardson *neP;
196d8d34be6SBarry Smith   SNESLineSearch    linesearch;
1975fd66863SKarl Rupp 
198028b6e76SBarry Smith   PetscFunctionBegin;
199d5c3842bSBarry Smith   snes->ops->destroy        = SNESDestroy_NRichardson;
200d5c3842bSBarry Smith   snes->ops->setup          = SNESSetUp_NRichardson;
201d5c3842bSBarry Smith   snes->ops->setfromoptions = SNESSetFromOptions_NRichardson;
202d5c3842bSBarry Smith   snes->ops->view           = SNESView_NRichardson;
203d5c3842bSBarry Smith   snes->ops->solve          = SNESSolve_NRichardson;
204d5c3842bSBarry Smith   snes->ops->reset          = SNESReset_NRichardson;
205028b6e76SBarry Smith 
206028b6e76SBarry Smith   snes->usesksp = PETSC_FALSE;
207efd4aadfSBarry Smith   snes->usesnpc = PETSC_TRUE;
208028b6e76SBarry Smith 
209efd4aadfSBarry Smith   snes->npcside = PC_LEFT;
210c6b63b32SPeter Brune 
2119566063dSJacob Faibussowitsch   PetscCall(SNESGetLineSearch(snes, &linesearch));
21248a46eb9SPierre Jolivet   if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHL2));
213d8d34be6SBarry Smith 
2144fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
2154fc747eaSLawrence Mitchell 
2164dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&neP));
217bf7f4e0aSPeter Brune   snes->data = (void *)neP;
218bf7f4e0aSPeter Brune 
21988976e71SPeter Brune   if (!snes->tolerancesset) {
2200e444f03SPeter Brune     snes->max_funcs = 30000;
2210e444f03SPeter Brune     snes->max_its   = 10000;
222c522fa08SPeter Brune     snes->stol      = 1e-20;
22388976e71SPeter Brune   }
2243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
225028b6e76SBarry Smith }
226