xref: /petsc/src/snes/impls/richardson/snesrichardson.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
142f4f86dSBarry Smith #include <../src/snes/impls/richardson/snesrichardsonimpl.h>
2028b6e76SBarry Smith 
39371c9d4SSatish Balay PetscErrorCode SNESReset_NRichardson(SNES snes) {
4028b6e76SBarry Smith   PetscFunctionBegin;
5028b6e76SBarry Smith   PetscFunctionReturn(0);
6028b6e76SBarry Smith }
7028b6e76SBarry Smith 
8028b6e76SBarry Smith /*
9d5c3842bSBarry Smith   SNESDestroy_NRichardson - Destroys the private SNES_NRichardson context that was created with SNESCreate_NRichardson().
10028b6e76SBarry Smith 
11028b6e76SBarry Smith   Input Parameter:
12028b6e76SBarry Smith . snes - the SNES context
13028b6e76SBarry Smith 
14028b6e76SBarry Smith   Application Interface Routine: SNESDestroy()
15028b6e76SBarry Smith */
169371c9d4SSatish Balay PetscErrorCode SNESDestroy_NRichardson(SNES snes) {
17028b6e76SBarry Smith   PetscFunctionBegin;
189566063dSJacob Faibussowitsch   PetscCall(SNESReset_NRichardson(snes));
199566063dSJacob Faibussowitsch   PetscCall(PetscFree(snes->data));
20028b6e76SBarry Smith   PetscFunctionReturn(0);
21028b6e76SBarry Smith }
22028b6e76SBarry Smith 
23028b6e76SBarry Smith /*
24d5c3842bSBarry Smith    SNESSetUp_NRichardson - Sets up the internal data structures for the later use
25d5c3842bSBarry Smith    of the SNESNRICHARDSON nonlinear solver.
26028b6e76SBarry Smith 
27028b6e76SBarry Smith    Input Parameters:
28028b6e76SBarry Smith +  snes - the SNES context
29028b6e76SBarry Smith -  x - the solution vector
30028b6e76SBarry Smith 
31028b6e76SBarry Smith    Application Interface Routine: SNESSetUp()
32028b6e76SBarry Smith  */
339371c9d4SSatish Balay PetscErrorCode SNESSetUp_NRichardson(SNES snes) {
34028b6e76SBarry Smith   PetscFunctionBegin;
3508401ef6SPierre Jolivet   PetscCheck(snes->npcside != PC_RIGHT, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "NRichardson only supports left preconditioning");
366c67d002SPeter Brune   if (snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED;
37028b6e76SBarry Smith   PetscFunctionReturn(0);
38028b6e76SBarry Smith }
39028b6e76SBarry Smith 
40028b6e76SBarry Smith /*
4104d7464bSBarry Smith   SNESSetFromOptions_NRichardson - Sets various parameters for the SNESNEWTONLS method.
42028b6e76SBarry Smith 
43028b6e76SBarry Smith   Input Parameter:
44028b6e76SBarry Smith . snes - the SNES context
45028b6e76SBarry Smith 
46028b6e76SBarry Smith   Application Interface Routine: SNESSetFromOptions()
47028b6e76SBarry Smith */
489371c9d4SSatish Balay static PetscErrorCode SNESSetFromOptions_NRichardson(SNES snes, PetscOptionItems *PetscOptionsObject) {
49028b6e76SBarry Smith   PetscFunctionBegin;
50d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SNES Richardson options");
51d0609cedSBarry Smith   PetscOptionsHeadEnd();
52028b6e76SBarry Smith   PetscFunctionReturn(0);
53028b6e76SBarry Smith }
54028b6e76SBarry Smith 
55028b6e76SBarry Smith /*
56d5c3842bSBarry Smith   SNESView_NRichardson - Prints info from the SNESRichardson data structure.
57028b6e76SBarry Smith 
58028b6e76SBarry Smith   Input Parameters:
59028b6e76SBarry Smith + SNES - the SNES context
60028b6e76SBarry Smith - viewer - visualization context
61028b6e76SBarry Smith 
62028b6e76SBarry Smith   Application Interface Routine: SNESView()
63028b6e76SBarry Smith */
649371c9d4SSatish Balay static PetscErrorCode SNESView_NRichardson(SNES snes, PetscViewer viewer) {
65028b6e76SBarry Smith   PetscBool iascii;
66028b6e76SBarry Smith 
67028b6e76SBarry Smith   PetscFunctionBegin;
689566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
699371c9d4SSatish Balay   if (iascii) { }
70028b6e76SBarry Smith   PetscFunctionReturn(0);
71028b6e76SBarry Smith }
72028b6e76SBarry Smith 
73028b6e76SBarry Smith /*
74d5c3842bSBarry Smith   SNESSolve_NRichardson - Solves a nonlinear system with the Richardson method.
75028b6e76SBarry Smith 
76028b6e76SBarry Smith   Input Parameters:
77028b6e76SBarry Smith . snes - the SNES context
78028b6e76SBarry Smith 
79028b6e76SBarry Smith   Output Parameter:
80028b6e76SBarry Smith . outits - number of iterations until termination
81028b6e76SBarry Smith 
82028b6e76SBarry Smith   Application Interface Routine: SNESSolve()
83028b6e76SBarry Smith */
849371c9d4SSatish Balay PetscErrorCode SNESSolve_NRichardson(SNES snes) {
85bf7f4e0aSPeter Brune   Vec                  X, Y, F;
86e7058c64SPeter Brune   PetscReal            xnorm, fnorm, ynorm;
87028b6e76SBarry Smith   PetscInt             maxits, i;
88422a814eSBarry Smith   SNESLineSearchReason lsresult;
89028b6e76SBarry Smith   SNESConvergedReason  reason;
90028b6e76SBarry Smith 
91028b6e76SBarry Smith   PetscFunctionBegin;
920b121fc5SBarry 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);
93c579b300SPatrick Farrell 
94028b6e76SBarry Smith   snes->reason = SNES_CONVERGED_ITERATING;
95028b6e76SBarry Smith 
96028b6e76SBarry Smith   maxits = snes->max_its;        /* maximum number of iterations */
97028b6e76SBarry Smith   X      = snes->vec_sol;        /* X^n */
98028b6e76SBarry Smith   Y      = snes->vec_sol_update; /* \tilde X */
99028b6e76SBarry Smith   F      = snes->vec_func;       /* residual vector */
100028b6e76SBarry Smith 
1019566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
102028b6e76SBarry Smith   snes->iter = 0;
103028b6e76SBarry Smith   snes->norm = 0.;
1049566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
105b7281c8aSPeter Brune 
106efd4aadfSBarry Smith   if (snes->npc && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
1079566063dSJacob Faibussowitsch     PetscCall(SNESApplyNPC(snes, X, NULL, F));
1089566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(snes->npc, &reason));
109b7281c8aSPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
110b7281c8aSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
111b7281c8aSPeter Brune       PetscFunctionReturn(0);
112b7281c8aSPeter Brune     }
1139566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, &fnorm));
114b7281c8aSPeter Brune   } else {
115e4ed7901SPeter Brune     if (!snes->vec_func_init_set) {
1169566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, X, F));
1171aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
118c1c75074SPeter Brune 
1199566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, &fnorm));
120422a814eSBarry Smith     SNESCheckFunctionNorm(snes, fnorm);
121b7281c8aSPeter Brune   }
122efd4aadfSBarry Smith   if (snes->npc && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
1239566063dSJacob Faibussowitsch     PetscCall(SNESApplyNPC(snes, X, F, Y));
1249566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(snes->npc, &reason));
125b7281c8aSPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
126b7281c8aSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
127b7281c8aSPeter Brune       PetscFunctionReturn(0);
128b7281c8aSPeter Brune     }
129b7281c8aSPeter Brune   } else {
1309566063dSJacob Faibussowitsch     PetscCall(VecCopy(F, Y));
131b7281c8aSPeter Brune   }
132e4ed7901SPeter Brune 
1339566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
134028b6e76SBarry Smith   snes->norm = fnorm;
1359566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
1369566063dSJacob Faibussowitsch   PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
1379566063dSJacob Faibussowitsch   PetscCall(SNESMonitor(snes, 0, fnorm));
138028b6e76SBarry Smith 
139028b6e76SBarry Smith   /* test convergence */
140dbbe0bcdSBarry Smith   PetscUseTypeMethod(snes, converged, 0, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP);
141028b6e76SBarry Smith   if (snes->reason) PetscFunctionReturn(0);
142028b6e76SBarry Smith 
143028b6e76SBarry Smith   /* Call general purpose update function */
144dbbe0bcdSBarry Smith   PetscTryTypeMethod(snes, update, snes->iter);
145b7281c8aSPeter Brune 
146b7281c8aSPeter Brune   /* set parameter for default relative tolerance convergence test */
147b7281c8aSPeter Brune   snes->ttol = fnorm * snes->rtol;
148b7281c8aSPeter Brune   /* test convergence */
149dbbe0bcdSBarry Smith   PetscUseTypeMethod(snes, converged, 0, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP);
150b7281c8aSPeter Brune   if (snes->reason) PetscFunctionReturn(0);
151b7281c8aSPeter Brune 
152b7281c8aSPeter Brune   for (i = 1; i < maxits + 1; i++) {
1539566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y));
1549566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetReason(snes->linesearch, &lsresult));
1559566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm));
156422a814eSBarry Smith     if (lsresult) {
157028b6e76SBarry Smith       if (++snes->numFailures >= snes->maxFailures) {
158028b6e76SBarry Smith         snes->reason = SNES_DIVERGED_LINE_SEARCH;
159028b6e76SBarry Smith         break;
160028b6e76SBarry Smith       }
161028b6e76SBarry Smith     }
162e71169deSBarry Smith     if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
163028b6e76SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
164028b6e76SBarry Smith       break;
165028b6e76SBarry Smith     }
166b7281c8aSPeter Brune 
167028b6e76SBarry Smith     /* Monitor convergence */
1689566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
169b7281c8aSPeter Brune     snes->iter  = i;
170028b6e76SBarry Smith     snes->norm  = fnorm;
171c1e67a49SFande Kong     snes->xnorm = xnorm;
172c1e67a49SFande Kong     snes->ynorm = ynorm;
1739566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
1749566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
1759566063dSJacob Faibussowitsch     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
176028b6e76SBarry Smith     /* Test for convergence */
177dbbe0bcdSBarry Smith     PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &snes->reason, snes->cnvP);
178028b6e76SBarry Smith     if (snes->reason) break;
179b7281c8aSPeter Brune 
180b7281c8aSPeter Brune     /* Call general purpose update function */
181dbbe0bcdSBarry Smith     PetscTryTypeMethod(snes, update, snes->iter);
182b7281c8aSPeter Brune 
183efd4aadfSBarry Smith     if (snes->npc) {
184b7281c8aSPeter Brune       if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
1859566063dSJacob Faibussowitsch         PetscCall(SNESApplyNPC(snes, X, NULL, Y));
1869566063dSJacob Faibussowitsch         PetscCall(VecNorm(F, NORM_2, &fnorm));
1879566063dSJacob Faibussowitsch         PetscCall(VecCopy(Y, F));
188b7281c8aSPeter Brune       } else {
1899566063dSJacob Faibussowitsch         PetscCall(SNESApplyNPC(snes, X, F, Y));
190b7281c8aSPeter Brune       }
1919566063dSJacob Faibussowitsch       PetscCall(SNESGetConvergedReason(snes->npc, &reason));
192b7281c8aSPeter Brune       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
193b7281c8aSPeter Brune         snes->reason = SNES_DIVERGED_INNER;
194b7281c8aSPeter Brune         PetscFunctionReturn(0);
195b7281c8aSPeter Brune       }
196b7281c8aSPeter Brune     } else {
1979566063dSJacob Faibussowitsch       PetscCall(VecCopy(F, Y));
198b7281c8aSPeter Brune     }
199b7281c8aSPeter Brune   }
200b7281c8aSPeter Brune   if (i == maxits + 1) {
20163a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits));
202028b6e76SBarry Smith     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
203028b6e76SBarry Smith   }
204028b6e76SBarry Smith   PetscFunctionReturn(0);
205028b6e76SBarry Smith }
206028b6e76SBarry Smith 
207028b6e76SBarry Smith /*MC
208d5c3842bSBarry Smith    SNESNRICHARDSON - Richardson nonlinear solver that uses successive substitutions, also sometimes known as Picard iteration.
209028b6e76SBarry Smith 
210f6dfbefdSBarry Smith    Options Database Keys:
21167b8a455SSatish Balay +  -snes_linesearch_type <l2,cp,basic> - Line search type.
21267b8a455SSatish Balay -  -snes_linesearch_damping<1.0> - Damping for the line search.
213028b6e76SBarry Smith 
214f6dfbefdSBarry Smith    Level: beginner
215f6dfbefdSBarry Smith 
21695452b02SPatrick Sanan    Notes:
21795452b02SPatrick Sanan    If no inner nonlinear preconditioner is provided then solves F(x) - b = 0 using x^{n+1} = x^{n} - lambda
218f6dfbefdSBarry Smith    (F(x^n) - b) where lambda is obtained either `SNESLineSearchSetDamping()`, -snes_damping or a line search.  If
219f6dfbefdSBarry Smith    an inner nonlinear preconditioner is provided (either with -npc_snes_type or `SNESSetNPC())` then the inner
220af60355fSPeter Brune    solver is called an initial solution x^n and the nonlinear Richardson uses x^{n+1} = x^{n} + lambda d^{n}
221af60355fSPeter Brune    where d^{n} = \hat{x}^{n} - x^{n} where \hat{x}^{n} is the solution returned from the inner solver.
222028b6e76SBarry Smith 
22372835e02SPeter Brune    The update, especially without inner nonlinear preconditioner, may be ill-scaled.  If using the basic
22472835e02SPeter Brune    linesearch, one may have to scale the update with -snes_linesearch_damping
22572835e02SPeter Brune 
226028b6e76SBarry Smith    This uses no derivative information thus will be much slower then Newton's method obtained with -snes_type ls
227028b6e76SBarry Smith 
2282d547940SBarry Smith    Only supports left non-linear preconditioning.
2292d547940SBarry Smith 
230db781477SPatrick Sanan .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNEWTONTR`, `SNESNGMRES`, `SNESQN`, `SNESNCG`
231028b6e76SBarry Smith M*/
2329371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode SNESCreate_NRichardson(SNES snes) {
233bf7f4e0aSPeter Brune   SNES_NRichardson *neP;
234d8d34be6SBarry Smith   SNESLineSearch    linesearch;
2355fd66863SKarl Rupp 
236028b6e76SBarry Smith   PetscFunctionBegin;
237d5c3842bSBarry Smith   snes->ops->destroy        = SNESDestroy_NRichardson;
238d5c3842bSBarry Smith   snes->ops->setup          = SNESSetUp_NRichardson;
239d5c3842bSBarry Smith   snes->ops->setfromoptions = SNESSetFromOptions_NRichardson;
240d5c3842bSBarry Smith   snes->ops->view           = SNESView_NRichardson;
241d5c3842bSBarry Smith   snes->ops->solve          = SNESSolve_NRichardson;
242d5c3842bSBarry Smith   snes->ops->reset          = SNESReset_NRichardson;
243028b6e76SBarry Smith 
244028b6e76SBarry Smith   snes->usesksp = PETSC_FALSE;
245efd4aadfSBarry Smith   snes->usesnpc = PETSC_TRUE;
246028b6e76SBarry Smith 
247efd4aadfSBarry Smith   snes->npcside = PC_LEFT;
248c6b63b32SPeter Brune 
2499566063dSJacob Faibussowitsch   PetscCall(SNESGetLineSearch(snes, &linesearch));
25048a46eb9SPierre Jolivet   if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHL2));
251d8d34be6SBarry Smith 
2524fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
2534fc747eaSLawrence Mitchell 
254*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&neP));
255bf7f4e0aSPeter Brune   snes->data = (void *)neP;
256bf7f4e0aSPeter Brune 
25788976e71SPeter Brune   if (!snes->tolerancesset) {
2580e444f03SPeter Brune     snes->max_funcs = 30000;
2590e444f03SPeter Brune     snes->max_its   = 10000;
260c522fa08SPeter Brune     snes->stol      = 1e-20;
26188976e71SPeter Brune   }
262028b6e76SBarry Smith   PetscFunctionReturn(0);
263028b6e76SBarry Smith }
264