xref: /petsc/src/snes/impls/richardson/snesrichardson.c (revision b2ccae6bdc8edea944f1c160ca3b2eb32c69ecb2)
1 #include <../src/snes/impls/richardson/snesrichardsonimpl.h>
2 
3 static PetscErrorCode SNESDestroy_NRichardson(SNES snes)
4 {
5   PetscFunctionBegin;
6   PetscCall(PetscFree(snes->data));
7   PetscFunctionReturn(PETSC_SUCCESS);
8 }
9 
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 }
17 
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 }
25 
26 static PetscErrorCode SNESView_NRichardson(SNES snes, PetscViewer viewer)
27 {
28   PetscBool isascii;
29 
30   PetscFunctionBegin;
31   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
32   if (isascii) { }
33   PetscFunctionReturn(PETSC_SUCCESS);
34 }
35 
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;
43 
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);
46 
47   snes->reason = SNES_CONVERGED_ITERATING;
48 
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 */
53 
54   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
55   snes->iter = 0;
56   snes->norm = 0.;
57   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
58 
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) PetscCall(SNESComputeFunction(snes, X, F));
69     else snes->vec_func_init_set = PETSC_FALSE;
70 
71     PetscCall(VecNorm(F, NORM_2, &fnorm));
72     SNESCheckFunctionNorm(snes, fnorm);
73   }
74   if (snes->npc && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
75     PetscCall(SNESApplyNPC(snes, X, F, Y));
76     PetscCall(SNESGetConvergedReason(snes->npc, &reason));
77     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
78       snes->reason = SNES_DIVERGED_INNER;
79       PetscFunctionReturn(PETSC_SUCCESS);
80     }
81   } else {
82     PetscCall(VecCopy(F, Y));
83   }
84 
85   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
86   snes->norm = fnorm;
87   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
88   PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
89 
90   /* test convergence */
91   PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm));
92   PetscCall(SNESMonitor(snes, 0, fnorm));
93   if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
94 
95   /* Call general purpose update function */
96   PetscTryTypeMethod(snes, update, snes->iter);
97 
98   for (i = 1; i < maxits + 1; i++) {
99     PetscCall(SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y));
100     PetscCall(SNESLineSearchGetReason(snes->linesearch, &lsresult));
101     PetscCall(SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm));
102     if (lsresult) {
103       if (++snes->numFailures >= snes->maxFailures) {
104         snes->reason = SNES_DIVERGED_LINE_SEARCH;
105         break;
106       }
107     }
108     if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
109       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
110       break;
111     }
112 
113     /* Monitor convergence */
114     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
115     snes->iter  = i;
116     snes->norm  = fnorm;
117     snes->xnorm = xnorm;
118     snes->ynorm = ynorm;
119     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
120     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
121     /* Test for convergence */
122     PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm));
123     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
124     if (snes->reason) break;
125 
126     /* Call general purpose update function */
127     PetscTryTypeMethod(snes, update, snes->iter);
128 
129     if (snes->npc) {
130       if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
131         PetscCall(SNESApplyNPC(snes, X, NULL, Y));
132         PetscCall(VecNorm(F, NORM_2, &fnorm));
133         PetscCall(VecCopy(Y, F));
134       } else {
135         PetscCall(SNESApplyNPC(snes, X, F, Y));
136       }
137       PetscCall(SNESGetConvergedReason(snes->npc, &reason));
138       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
139         snes->reason = SNES_DIVERGED_INNER;
140         PetscFunctionReturn(PETSC_SUCCESS);
141       }
142     } else {
143       PetscCall(VecCopy(F, Y));
144     }
145   }
146   if (i == maxits + 1) {
147     PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits));
148     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
149   }
150   PetscFunctionReturn(PETSC_SUCCESS);
151 }
152 
153 /*MC
154    SNESNRICHARDSON - Richardson nonlinear solver that uses successive substitutions, also sometimes known as Picard iteration.
155 
156    Options Database Keys:
157 +  -snes_linesearch_type <l2,cp,basic> - Line search type.
158 -  -snes_linesearch_damping <1.0>      - Damping for the line search.
159 
160    Level: beginner
161 
162    Notes:
163    If no inner nonlinear preconditioner is provided then solves $F(x) - b = 0 $ using $x^{n+1} = x^{n} - \lambda
164    (F(x^n) - b) $ where $ \lambda$ is obtained with either `SNESLineSearchSetDamping()`, `-snes_damping` or a line search.  If
165    an inner nonlinear preconditioner is provided (either with `-npc_snes_typ`e or `SNESSetNPC()`) then the inner
166    solver is called on the initial solution $x^n$ and the nonlinear Richardson uses $ x^{n+1} = x^{n} + \lambda d^{n}$
167    where $d^{n} = \hat{x}^{n} - x^{n} $ where $\hat{x}^{n} $ is the solution returned from the inner solver.
168 
169    The update, especially without inner nonlinear preconditioner, may be ill-scaled.  If using the basic
170    linesearch, one may have to scale the update with `-snes_linesearch_damping`
171 
172    This uses no derivative information provided with `SNESSetJacobian()` thus it will be much slower than Newton's method obtained with `-snes_type ls`
173 
174    Only supports left non-linear preconditioning.
175 
176 .seealso: [](ch_snes), `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNEWTONTR`, `SNESNGMRES`, `SNESQN`, `SNESNCG`,
177           `SNESLineSearchSetDamping()`
178 M*/
179 PETSC_EXTERN PetscErrorCode SNESCreate_NRichardson(SNES snes)
180 {
181   SNES_NRichardson *neP;
182   SNESLineSearch    linesearch;
183 
184   PetscFunctionBegin;
185   snes->ops->destroy        = SNESDestroy_NRichardson;
186   snes->ops->setup          = SNESSetUp_NRichardson;
187   snes->ops->setfromoptions = SNESSetFromOptions_NRichardson;
188   snes->ops->view           = SNESView_NRichardson;
189   snes->ops->solve          = SNESSolve_NRichardson;
190 
191   snes->usesksp = PETSC_FALSE;
192   snes->usesnpc = PETSC_TRUE;
193 
194   snes->npcside = PC_LEFT;
195 
196   PetscCall(SNESGetLineSearch(snes, &linesearch));
197   if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHSECANT));
198 
199   snes->alwayscomputesfinalresidual = PETSC_TRUE;
200 
201   PetscCall(SNESParametersInitialize(snes));
202   PetscObjectParameterSetDefault(snes, max_funcs, 30000);
203   PetscObjectParameterSetDefault(snes, max_its, 10000);
204   PetscObjectParameterSetDefault(snes, stol, 1e-20);
205 
206   PetscCall(PetscNew(&neP));
207   snes->data = (void *)neP;
208   PetscFunctionReturn(PETSC_SUCCESS);
209 }
210