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