xref: /petsc/src/snes/impls/richardson/snesrichardson.c (revision bccf9bb3a792d9a0cb91bcafafa02bff20973da7)
1 
2 #include <../src/snes/impls/richardson/snesrichardsonimpl.h>
3 
4 #undef __FUNCT__
5 #define __FUNCT__ "SNESReset_NRichardson"
6 PetscErrorCode SNESReset_NRichardson(SNES snes)
7 {
8 
9   PetscFunctionBegin;
10   PetscFunctionReturn(0);
11 }
12 
13 /*
14   SNESDestroy_NRichardson - Destroys the private SNES_NRichardson context that was created with SNESCreate_NRichardson().
15 
16   Input Parameter:
17 . snes - the SNES context
18 
19   Application Interface Routine: SNESDestroy()
20 */
21 #undef __FUNCT__
22 #define __FUNCT__ "SNESDestroy_NRichardson"
23 PetscErrorCode SNESDestroy_NRichardson(SNES snes)
24 {
25   PetscErrorCode   ierr;
26 
27   PetscFunctionBegin;
28   ierr = SNESReset_NRichardson(snes);CHKERRQ(ierr);
29   ierr = PetscFree(snes->data);CHKERRQ(ierr);
30   PetscFunctionReturn(0);
31 }
32 
33 /*
34    SNESSetUp_NRichardson - Sets up the internal data structures for the later use
35    of the SNESNRICHARDSON nonlinear solver.
36 
37    Input Parameters:
38 +  snes - the SNES context
39 -  x - the solution vector
40 
41    Application Interface Routine: SNESSetUp()
42  */
43 #undef __FUNCT__
44 #define __FUNCT__ "SNESSetUp_NRichardson"
45 PetscErrorCode SNESSetUp_NRichardson(SNES snes)
46 {
47   PetscErrorCode ierr;
48 
49   PetscFunctionBegin;
50   ierr = SNESDefaultGetWork(snes,2);CHKERRQ(ierr);
51   PetscFunctionReturn(0);
52 }
53 
54 /*
55   SNESSetFromOptions_NRichardson - Sets various parameters for the SNESLS method.
56 
57   Input Parameter:
58 . snes - the SNES context
59 
60   Application Interface Routine: SNESSetFromOptions()
61 */
62 #undef __FUNCT__
63 #define __FUNCT__ "SNESSetFromOptions_NRichardson"
64 static PetscErrorCode SNESSetFromOptions_NRichardson(SNES snes)
65 {
66   PetscErrorCode ierr;
67   PetscFunctionBegin;
68     ierr = PetscOptionsHead("SNES Richardson options");CHKERRQ(ierr);
69     ierr = PetscOptionsTail();CHKERRQ(ierr);
70   PetscFunctionReturn(0);
71 }
72 
73 /*
74   SNESView_NRichardson - Prints info from the SNESRichardson data structure.
75 
76   Input Parameters:
77 + SNES - the SNES context
78 - viewer - visualization context
79 
80   Application Interface Routine: SNESView()
81 */
82 #undef __FUNCT__
83 #define __FUNCT__ "SNESView_NRichardson"
84 static PetscErrorCode SNESView_NRichardson(SNES snes, PetscViewer viewer)
85 {
86   PetscBool        iascii;
87   PetscErrorCode   ierr;
88 
89   PetscFunctionBegin;
90   ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
91   if (iascii) {
92     ierr = PetscViewerASCIIPrintf(viewer,"  richardson variant: %s\n", SNESLineSearchTypeName(snes->ls_type));CHKERRQ(ierr);
93   }
94   PetscFunctionReturn(0);
95 }
96 
97 #undef __FUNCT__
98 #define __FUNCT__ "SNESLineSearchQuadratic_NRichardson"
99 PetscErrorCode SNESLineSearchQuadratic_NRichardson(SNES snes,void *lsctx,Vec X,Vec F,Vec Y,PetscReal fnorm,PetscReal dummyXnorm,Vec G,Vec W,PetscReal *dummyYnorm,PetscReal *gnorm,PetscBool *flag)
100 {
101   PetscInt         i;
102   PetscReal        alphas[3] = {0.0, 0.5, 1.0};
103   PetscReal        norms[3];
104   PetscReal        alpha,a,b;
105   PetscErrorCode   ierr;
106 
107   PetscFunctionBegin;
108   norms[0]  = fnorm;
109   for(i=1; i < 3; ++i) {
110     ierr = VecWAXPY(W, -alphas[i], Y, X);CHKERRQ(ierr);     /* W =  X^n - \alpha Y */
111     ierr = SNESComputeFunction(snes, W, G);CHKERRQ(ierr);
112     ierr = VecNorm(G, NORM_2, &norms[i]);CHKERRQ(ierr);
113   }
114   for(i = 0; i < 3; ++i) {
115     norms[i] = PetscSqr(norms[i]);
116   }
117   /* Fit a quadratic:
118        If we have x_{0,1,2} = 0, x_1, x_2 which generate norms y_{0,1,2}
119        a = (x_1 y_2 - x_2 y_1 + (x_2 - x_1) y_0)/(x^2_2 x_1 - x_2 x^2_1)
120        b = (x^2_1 y_2 - x^2_2 y_1 + (x^2_2 - x^2_1) y_0)/(x_2 x^2_1 - x^2_2 x_1)
121        c = y_0
122        x_min = -b/2a
123 
124        If we let x_{0,1,2} = 0, 0.5, 1.0
125        a = 2 y_2 - 4 y_1 + 2 y_0
126        b =  -y_2 + 4 y_1 - 3 y_0
127        c =   y_0
128   */
129   a = (alphas[1]*norms[2] - alphas[2]*norms[1] + (alphas[2] - alphas[1])*norms[0])/(PetscSqr(alphas[2])*alphas[1] - alphas[2]*PetscSqr(alphas[1]));
130   b = (PetscSqr(alphas[1])*norms[2] - PetscSqr(alphas[2])*norms[1] + (PetscSqr(alphas[2]) - PetscSqr(alphas[1]))*norms[0])/(alphas[2]*PetscSqr(alphas[1]) - PetscSqr(alphas[2])*alphas[1]);
131   /* Check for positive a (concave up) */
132   if (a >= 0.0) {
133     alpha = -b/(2.0*a);
134     alpha = PetscMin(alpha, alphas[2]);
135     alpha = PetscMax(alpha, alphas[0]);
136   } else {
137     alpha = 1.0;
138   }
139   if (snes->ls_monitor) {
140     ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
141     ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: norms[0] = %g, norms[1] = %g, norms[2] = %g alpha %g\n",
142                                   PetscSqrtReal(norms[0]),PetscSqrtReal(norms[1]),PetscSqrtReal(norms[2]),alpha);CHKERRQ(ierr);
143     ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
144   }
145   ierr = VecCopy(X, W);CHKERRQ(ierr);
146   ierr = VecAXPY(W, -alpha, Y);CHKERRQ(ierr);
147   if (alpha != 1.0) {
148     ierr = SNESComputeFunction(snes, W, G);CHKERRQ(ierr);
149     ierr = VecNorm(G, NORM_2, gnorm);CHKERRQ(ierr);
150   } else {
151     *gnorm = PetscSqrtReal(norms[2]);
152   }
153   if (alpha == 0.0) *flag = PETSC_FALSE;
154   else              *flag = PETSC_TRUE;
155   PetscFunctionReturn(0);
156 }
157 
158 /*
159   SNESSolve_NRichardson - Solves a nonlinear system with the Richardson method.
160 
161   Input Parameters:
162 . snes - the SNES context
163 
164   Output Parameter:
165 . outits - number of iterations until termination
166 
167   Application Interface Routine: SNESSolve()
168 */
169 #undef __FUNCT__
170 #define __FUNCT__ "SNESSolve_NRichardson"
171 PetscErrorCode SNESSolve_NRichardson(SNES snes)
172 {
173   Vec                 X, Y, F, W, G;
174   PetscReal           fnorm;
175   PetscInt            maxits, i;
176   PetscErrorCode      ierr;
177   SNESConvergedReason reason;
178 
179   PetscFunctionBegin;
180   snes->reason = SNES_CONVERGED_ITERATING;
181 
182   maxits = snes->max_its;        /* maximum number of iterations */
183   X      = snes->vec_sol;        /* X^n */
184   Y      = snes->vec_sol_update; /* \tilde X */
185   F      = snes->vec_func;       /* residual vector */
186   W      = snes->work[0];        /* work vector */
187   G      = snes->work[1];        /* line search function vector */
188 
189   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
190   snes->iter = 0;
191   snes->norm = 0.;
192   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
193   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
194   if (snes->domainerror) {
195     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
196     PetscFunctionReturn(0);
197   }
198   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
199   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
200   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
201   snes->norm = fnorm;
202   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
203   SNESLogConvHistory(snes,fnorm,0);
204   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
205 
206   /* set parameter for default relative tolerance convergence test */
207   snes->ttol = fnorm*snes->rtol;
208   /* test convergence */
209   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
210   if (snes->reason) PetscFunctionReturn(0);
211 
212   for(i = 0; i < maxits; i++) {
213     PetscBool  lsSuccess = PETSC_TRUE;
214     PetscReal  dummyNorm;
215 
216     /* Call general purpose update function */
217     if (snes->ops->update) {
218       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
219     }
220     else if (snes->pc) {
221       ierr = VecCopy(X,Y);CHKERRQ(ierr);
222       ierr = SNESSolve(snes->pc, snes->vec_rhs, Y);CHKERRQ(ierr);
223       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
224       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
225         snes->reason = SNES_DIVERGED_INNER;
226         PetscFunctionReturn(0);
227       }
228       ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr);
229     } else {
230       ierr = VecCopy(F,Y);CHKERRQ(ierr);
231     }
232     ierr = (*snes->ops->linesearch)(snes, snes->lsP, X, F, Y, fnorm, 0.0, G, W, &dummyNorm, &fnorm, &lsSuccess);CHKERRQ(ierr);
233     if (!lsSuccess) {
234       if (++snes->numFailures >= snes->maxFailures) {
235         snes->reason = SNES_DIVERGED_LINE_SEARCH;
236         break;
237       }
238     }
239     if (snes->nfuncs >= snes->max_funcs) {
240       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
241       break;
242     }
243     if (snes->domainerror) {
244       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
245       PetscFunctionReturn(0);
246     }
247     ierr = VecCopy(G, F);CHKERRQ(ierr);
248     ierr = VecCopy(W, X);CHKERRQ(ierr);
249     /* Monitor convergence */
250     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
251     snes->iter = i+1;
252     snes->norm = fnorm;
253     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
254     SNESLogConvHistory(snes,snes->norm,0);
255     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
256     /* Test for convergence */
257     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
258     if (snes->reason) break;
259   }
260   if (i == maxits) {
261     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
262     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
263   }
264   PetscFunctionReturn(0);
265 }
266 
267 
268 EXTERN_C_BEGIN
269 #undef __FUNCT__
270 #define __FUNCT__ "SNESLineSearchSetType_NRichardson"
271 PetscErrorCode  SNESLineSearchSetType_NRichardson(SNES snes, SNESLineSearchType type)
272 {
273   PetscErrorCode ierr;
274   PetscFunctionBegin;
275 
276   switch (type) {
277   case SNES_LS_BASIC:
278     ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr);
279     break;
280   case SNES_LS_BASIC_NONORMS:
281     ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
282     break;
283   case SNES_LS_QUADRATIC:
284     ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic_NRichardson,PETSC_NULL);CHKERRQ(ierr);
285     break;
286   case SNES_LS_SECANT:
287     ierr = SNESLineSearchSet(snes,SNESLineSearchQuadraticSecant,PETSC_NULL);CHKERRQ(ierr);
288     break;
289   default:
290     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type");
291     break;
292   }
293   snes->ls_type = type;
294   PetscFunctionReturn(0);
295 }
296 EXTERN_C_END
297 
298 /*MC
299   SNESNRICHARDSON - Richardson nonlinear solver that uses successive substitutions, also sometimes known as Picard iteration.
300 
301   Level: beginner
302 
303   Options Database:
304 +   -snes_ls_damping - damping factor to apply to F(x) (used only if -snes_ls is basic or basicnonorms)
305 -   -snes_ls <basic,basicnormnorms,quadratic>
306 
307   Notes: If no inner nonlinear preconditioner is provided then solves F(x) - b = 0 using x^{n+1} = x^{n} - lambda
308             (F(x^n) - b) where lambda is obtained either SNESLineSearchSetDamping(), -snes_damping or a line search.  If
309             an inner nonlinear preconditioner is provided (either with -npc_snes_type or SNESSetPC()) then the inner
310             solver is called an initial solution x^n and the nonlinear Richardson uses x^{n+1} = x^{n} + lambda d^{n}
311             where d^{n} = \hat{x}^{n} - x^{n} where \hat{x}^{n} is the solution returned from the inner solver.
312 
313      This uses no derivative information thus will be much slower then Newton's method obtained with -snes_type ls
314 
315 .seealso:  SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR, SNESNGMRES, SNESNQN
316 M*/
317 EXTERN_C_BEGIN
318 #undef __FUNCT__
319 #define __FUNCT__ "SNESCreate_NRichardson"
320 PetscErrorCode  SNESCreate_NRichardson(SNES snes)
321 {
322   PetscErrorCode ierr;
323   PetscFunctionBegin;
324   snes->ops->destroy         = SNESDestroy_NRichardson;
325   snes->ops->setup           = SNESSetUp_NRichardson;
326   snes->ops->setfromoptions  = SNESSetFromOptions_NRichardson;
327   snes->ops->view            = SNESView_NRichardson;
328   snes->ops->solve           = SNESSolve_NRichardson;
329   snes->ops->reset           = SNESReset_NRichardson;
330 
331   snes->usesksp              = PETSC_FALSE;
332   snes->usespc               = PETSC_TRUE;
333 
334   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","SNESLineSearchSetType_NRichardson",SNESLineSearchSetType_NRichardson);CHKERRQ(ierr);
335   ierr = SNESLineSearchSetType(snes, SNES_LS_SECANT);CHKERRQ(ierr);
336 
337   PetscFunctionReturn(0);
338 }
339 EXTERN_C_END
340