xref: /petsc/src/snes/impls/richardson/snesrichardson.c (revision 8895d3f2aaf765dc4e8d62d73f37241d00bd6bbd)
1 #include <../src/snes/impls/richardson/snesrichardsonimpl.h>
2 
3 
4 #undef __FUNCT__
5 #define __FUNCT__ "SNESReset_NRichardson"
6 PetscErrorCode SNESReset_NRichardson(SNES snes)
7 {
8   PetscFunctionBegin;
9   PetscFunctionReturn(0);
10 }
11 
12 /*
13   SNESDestroy_NRichardson - Destroys the private SNES_NRichardson context that was created with SNESCreate_NRichardson().
14 
15   Input Parameter:
16 . snes - the SNES context
17 
18   Application Interface Routine: SNESDestroy()
19 */
20 #undef __FUNCT__
21 #define __FUNCT__ "SNESDestroy_NRichardson"
22 PetscErrorCode SNESDestroy_NRichardson(SNES snes)
23 {
24   PetscErrorCode ierr;
25 
26   PetscFunctionBegin;
27   ierr = SNESReset_NRichardson(snes);CHKERRQ(ierr);
28   ierr = PetscFree(snes->data);CHKERRQ(ierr);
29   PetscFunctionReturn(0);
30 }
31 
32 /*
33    SNESSetUp_NRichardson - Sets up the internal data structures for the later use
34    of the SNESNRICHARDSON nonlinear solver.
35 
36    Input Parameters:
37 +  snes - the SNES context
38 -  x - the solution vector
39 
40    Application Interface Routine: SNESSetUp()
41  */
42 #undef __FUNCT__
43 #define __FUNCT__ "SNESSetUp_NRichardson"
44 PetscErrorCode SNESSetUp_NRichardson(SNES snes)
45 {
46   PetscFunctionBegin;
47   if (snes->pcside == PC_RIGHT) {SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"NRichardson only supports left preconditioning");}
48   if (snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED;
49   PetscFunctionReturn(0);
50 }
51 
52 /*
53   SNESSetFromOptions_NRichardson - Sets various parameters for the SNESNEWTONLS method.
54 
55   Input Parameter:
56 . snes - the SNES context
57 
58   Application Interface Routine: SNESSetFromOptions()
59 */
60 #undef __FUNCT__
61 #define __FUNCT__ "SNESSetFromOptions_NRichardson"
62 static PetscErrorCode SNESSetFromOptions_NRichardson(PetscOptions *PetscOptionsObject,SNES snes)
63 {
64   PetscErrorCode ierr;
65   SNESLineSearch linesearch;
66 
67   PetscFunctionBegin;
68   ierr = PetscOptionsHead(PetscOptionsObject,"SNES Richardson options");CHKERRQ(ierr);
69   ierr = PetscOptionsTail();CHKERRQ(ierr);
70   if (!snes->linesearch) {
71     ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr);
72     ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr);
73   }
74   PetscFunctionReturn(0);
75 }
76 
77 /*
78   SNESView_NRichardson - Prints info from the SNESRichardson data structure.
79 
80   Input Parameters:
81 + SNES - the SNES context
82 - viewer - visualization context
83 
84   Application Interface Routine: SNESView()
85 */
86 #undef __FUNCT__
87 #define __FUNCT__ "SNESView_NRichardson"
88 static PetscErrorCode SNESView_NRichardson(SNES snes, PetscViewer viewer)
89 {
90   PetscBool      iascii;
91   PetscErrorCode ierr;
92 
93   PetscFunctionBegin;
94   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
95   if (iascii) {
96   }
97   PetscFunctionReturn(0);
98 }
99 
100 /*
101   SNESSolve_NRichardson - Solves a nonlinear system with the Richardson method.
102 
103   Input Parameters:
104 . snes - the SNES context
105 
106   Output Parameter:
107 . outits - number of iterations until termination
108 
109   Application Interface Routine: SNESSolve()
110 */
111 #undef __FUNCT__
112 #define __FUNCT__ "SNESSolve_NRichardson"
113 PetscErrorCode SNESSolve_NRichardson(SNES snes)
114 {
115   Vec                 X, Y, F;
116   PetscReal           xnorm, fnorm, ynorm;
117   PetscInt            maxits, i;
118   PetscErrorCode      ierr;
119   PetscBool           lsSuccess;
120   SNESConvergedReason reason;
121 
122   PetscFunctionBegin;
123 
124   if (snes->xl || snes->xu || snes->ops->computevariablebounds) {
125     SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
126   }
127 
128   snes->reason = SNES_CONVERGED_ITERATING;
129 
130   maxits = snes->max_its;        /* maximum number of iterations */
131   X      = snes->vec_sol;        /* X^n */
132   Y      = snes->vec_sol_update; /* \tilde X */
133   F      = snes->vec_func;       /* residual vector */
134 
135   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
136   snes->iter = 0;
137   snes->norm = 0.;
138   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
139 
140   if (snes->pc && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
141     ierr = SNESApplyNPC(snes,X,NULL,F);CHKERRQ(ierr);
142     ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
143     if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
144       snes->reason = SNES_DIVERGED_INNER;
145       PetscFunctionReturn(0);
146     }
147     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
148   } else {
149     if (!snes->vec_func_init_set) {
150       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
151       if (snes->domainerror) {
152         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
153         PetscFunctionReturn(0);
154       }
155     } else snes->vec_func_init_set = PETSC_FALSE;
156 
157     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
158     if (PetscIsInfOrNanReal(fnorm)) {
159       snes->reason = SNES_DIVERGED_FNORM_NAN;
160       PetscFunctionReturn(0);
161     }
162   }
163   if (snes->pc && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
164       ierr = SNESApplyNPC(snes,X,F,Y);CHKERRQ(ierr);
165       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
166       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
167         snes->reason = SNES_DIVERGED_INNER;
168         PetscFunctionReturn(0);
169       }
170   } else {
171     ierr = VecCopy(F,Y);CHKERRQ(ierr);
172   }
173 
174   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
175   snes->norm = fnorm;
176   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
177   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
178   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
179 
180   /* test convergence */
181   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
182   if (snes->reason) PetscFunctionReturn(0);
183 
184   /* Call general purpose update function */
185   if (snes->ops->update) {
186     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
187   }
188 
189   /* set parameter for default relative tolerance convergence test */
190   snes->ttol = fnorm*snes->rtol;
191   /* test convergence */
192   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
193   if (snes->reason) PetscFunctionReturn(0);
194 
195   for (i = 1; i < maxits+1; i++) {
196     lsSuccess = PETSC_TRUE;
197 
198     ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);CHKERRQ(ierr);
199     ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
200     ierr = SNESLineSearchGetSuccess(snes->linesearch, &lsSuccess);CHKERRQ(ierr);
201     if (!lsSuccess) {
202       if (++snes->numFailures >= snes->maxFailures) {
203         snes->reason = SNES_DIVERGED_LINE_SEARCH;
204         break;
205       }
206     }
207     if (snes->nfuncs >= snes->max_funcs) {
208       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
209       break;
210     }
211     if (snes->domainerror) {
212       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
213       PetscFunctionReturn(0);
214     }
215 
216     /* Monitor convergence */
217     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
218     snes->iter = i;
219     snes->norm = fnorm;
220     ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
221     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
222     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
223     /* Test for convergence */
224     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
225     if (snes->reason) break;
226 
227     /* Call general purpose update function */
228     if (snes->ops->update) {
229       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
230     }
231 
232     if (snes->pc) {
233       if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
234         ierr = SNESApplyNPC(snes,X,NULL,Y);CHKERRQ(ierr);
235         ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
236         ierr = VecCopy(Y,F);CHKERRQ(ierr);
237       } else {
238         ierr = SNESApplyNPC(snes,X,F,Y);CHKERRQ(ierr);
239       }
240       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
241       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
242         snes->reason = SNES_DIVERGED_INNER;
243         PetscFunctionReturn(0);
244       }
245     } else {
246       ierr = VecCopy(F,Y);CHKERRQ(ierr);
247     }
248   }
249   if (i == maxits+1) {
250     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
251     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
252   }
253   PetscFunctionReturn(0);
254 }
255 
256 /*MC
257   SNESNRICHARDSON - Richardson nonlinear solver that uses successive substitutions, also sometimes known as Picard iteration.
258 
259   Level: beginner
260 
261   Options Database:
262 +   -snes_linesearch_type <l2,cp,basic> Line search type.
263 -   -snes_linesearch_damping<1.0> Damping for the line search.
264 
265   Notes: If no inner nonlinear preconditioner is provided then solves F(x) - b = 0 using x^{n+1} = x^{n} - lambda
266             (F(x^n) - b) where lambda is obtained either SNESLineSearchSetDamping(), -snes_damping or a line search.  If
267             an inner nonlinear preconditioner is provided (either with -npc_snes_type or SNESSetNPC()) then the inner
268             solver is called an initial solution x^n and the nonlinear Richardson uses x^{n+1} = x^{n} + lambda d^{n}
269             where d^{n} = \hat{x}^{n} - x^{n} where \hat{x}^{n} is the solution returned from the inner solver.
270 
271             The update, especially without inner nonlinear preconditioner, may be ill-scaled.  If using the basic
272             linesearch, one may have to scale the update with -snes_linesearch_damping
273 
274      This uses no derivative information thus will be much slower then Newton's method obtained with -snes_type ls
275 
276 .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR, SNESNGMRES, SNESQN, SNESNCG
277 M*/
278 #undef __FUNCT__
279 #define __FUNCT__ "SNESCreate_NRichardson"
280 PETSC_EXTERN PetscErrorCode SNESCreate_NRichardson(SNES snes)
281 {
282   PetscErrorCode   ierr;
283   SNES_NRichardson *neP;
284 
285   PetscFunctionBegin;
286   snes->ops->destroy        = SNESDestroy_NRichardson;
287   snes->ops->setup          = SNESSetUp_NRichardson;
288   snes->ops->setfromoptions = SNESSetFromOptions_NRichardson;
289   snes->ops->view           = SNESView_NRichardson;
290   snes->ops->solve          = SNESSolve_NRichardson;
291   snes->ops->reset          = SNESReset_NRichardson;
292 
293   snes->usesksp = PETSC_FALSE;
294   snes->usespc  = PETSC_TRUE;
295 
296   snes->pcside = PC_LEFT;
297 
298   ierr       = PetscNewLog(snes,&neP);CHKERRQ(ierr);
299   snes->data = (void*) neP;
300 
301   if (!snes->tolerancesset) {
302     snes->max_funcs = 30000;
303     snes->max_its   = 10000;
304     snes->stol      = 1e-20;
305   }
306   PetscFunctionReturn(0);
307 }
308