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