xref: /petsc/src/snes/impls/richardson/snesrichardson.c (revision 2b8d69ca7ea5fe9190df62c1dce3bbd66fce84dd)
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(PetscOptionItems *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   SNESLineSearchReason lsresult;
120   SNESConvergedReason  reason;
121 
122   PetscFunctionBegin;
123   if (snes->xl || snes->xu || snes->ops->computevariablebounds) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
124 
125   snes->reason = SNES_CONVERGED_ITERATING;
126 
127   maxits = snes->max_its;        /* maximum number of iterations */
128   X      = snes->vec_sol;        /* X^n */
129   Y      = snes->vec_sol_update; /* \tilde X */
130   F      = snes->vec_func;       /* residual vector */
131 
132   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
133   snes->iter = 0;
134   snes->norm = 0.;
135   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
136 
137   if (snes->pc && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
138     ierr = SNESApplyNPC(snes,X,NULL,F);CHKERRQ(ierr);
139     ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
140     if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
141       snes->reason = SNES_DIVERGED_INNER;
142       PetscFunctionReturn(0);
143     }
144     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
145   } else {
146     if (!snes->vec_func_init_set) {
147       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
148     } else snes->vec_func_init_set = PETSC_FALSE;
149 
150     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
151     SNESCheckFunctionNorm(snes,fnorm);
152   }
153   if (snes->pc && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
154       ierr = SNESApplyNPC(snes,X,F,Y);CHKERRQ(ierr);
155       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
156       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
157         snes->reason = SNES_DIVERGED_INNER;
158         PetscFunctionReturn(0);
159       }
160   } else {
161     ierr = VecCopy(F,Y);CHKERRQ(ierr);
162   }
163 
164   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
165   snes->norm = fnorm;
166   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
167   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
168   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
169 
170   /* test convergence */
171   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
172   if (snes->reason) PetscFunctionReturn(0);
173 
174   /* Call general purpose update function */
175   if (snes->ops->update) {
176     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
177   }
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   for (i = 1; i < maxits+1; i++) {
186     ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);CHKERRQ(ierr);
187     ierr = SNESLineSearchGetReason(snes->linesearch, &lsresult);CHKERRQ(ierr);
188     ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
189     if (lsresult) {
190       if (++snes->numFailures >= snes->maxFailures) {
191         snes->reason = SNES_DIVERGED_LINE_SEARCH;
192         break;
193       }
194     }
195     if (snes->nfuncs >= snes->max_funcs) {
196       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
197       break;
198     }
199 
200     /* Monitor convergence */
201     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
202     snes->iter = i;
203     snes->norm = fnorm;
204     ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
205     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
206     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
207     /* Test for convergence */
208     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
209     if (snes->reason) break;
210 
211     /* Call general purpose update function */
212     if (snes->ops->update) {
213       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
214     }
215 
216     if (snes->pc) {
217       if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
218         ierr = SNESApplyNPC(snes,X,NULL,Y);CHKERRQ(ierr);
219         ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
220         ierr = VecCopy(Y,F);CHKERRQ(ierr);
221       } else {
222         ierr = SNESApplyNPC(snes,X,F,Y);CHKERRQ(ierr);
223       }
224       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
225       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
226         snes->reason = SNES_DIVERGED_INNER;
227         PetscFunctionReturn(0);
228       }
229     } else {
230       ierr = VecCopy(F,Y);CHKERRQ(ierr);
231     }
232   }
233   if (i == maxits+1) {
234     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
235     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
236   }
237   PetscFunctionReturn(0);
238 }
239 
240 /*MC
241   SNESNRICHARDSON - Richardson nonlinear solver that uses successive substitutions, also sometimes known as Picard iteration.
242 
243   Level: beginner
244 
245   Options Database:
246 +   -snes_linesearch_type <l2,cp,basic> Line search type.
247 -   -snes_linesearch_damping<1.0> Damping for the line search.
248 
249   Notes: If no inner nonlinear preconditioner is provided then solves F(x) - b = 0 using x^{n+1} = x^{n} - lambda
250             (F(x^n) - b) where lambda is obtained either SNESLineSearchSetDamping(), -snes_damping or a line search.  If
251             an inner nonlinear preconditioner is provided (either with -npc_snes_type or SNESSetNPC()) then the inner
252             solver is called an initial solution x^n and the nonlinear Richardson uses x^{n+1} = x^{n} + lambda d^{n}
253             where d^{n} = \hat{x}^{n} - x^{n} where \hat{x}^{n} is the solution returned from the inner solver.
254 
255             The update, especially without inner nonlinear preconditioner, may be ill-scaled.  If using the basic
256             linesearch, one may have to scale the update with -snes_linesearch_damping
257 
258      This uses no derivative information thus will be much slower then Newton's method obtained with -snes_type ls
259 
260 .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR, SNESNGMRES, SNESQN, SNESNCG
261 M*/
262 #undef __FUNCT__
263 #define __FUNCT__ "SNESCreate_NRichardson"
264 PETSC_EXTERN PetscErrorCode SNESCreate_NRichardson(SNES snes)
265 {
266   PetscErrorCode   ierr;
267   SNES_NRichardson *neP;
268 
269   PetscFunctionBegin;
270   snes->ops->destroy        = SNESDestroy_NRichardson;
271   snes->ops->setup          = SNESSetUp_NRichardson;
272   snes->ops->setfromoptions = SNESSetFromOptions_NRichardson;
273   snes->ops->view           = SNESView_NRichardson;
274   snes->ops->solve          = SNESSolve_NRichardson;
275   snes->ops->reset          = SNESReset_NRichardson;
276 
277   snes->usesksp = PETSC_FALSE;
278   snes->usespc  = PETSC_TRUE;
279 
280   snes->pcside = PC_LEFT;
281 
282   ierr       = PetscNewLog(snes,&neP);CHKERRQ(ierr);
283   snes->data = (void*) neP;
284 
285   if (!snes->tolerancesset) {
286     snes->max_funcs = 30000;
287     snes->max_its   = 10000;
288     snes->stol      = 1e-20;
289   }
290   PetscFunctionReturn(0);
291 }
292