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