xref: /petsc/src/snes/impls/fas/fas.c (revision 84df9cb40eca90ea9b18a456fab7a4ecc7f6c1a4)
1 /* Defines the basic SNES object */
2 #include <../src/snes/impls/fas/fasimpls.h>
3 
4 /*MC
5 Full Approximation Scheme nonlinear multigrid solver.
6 
7 The nonlinear problem is solved via the repeated application of nonlinear preconditioners and coarse-grid corrections
8 
9 .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
10 M*/
11 
12 extern PetscErrorCode SNESDestroy_FAS(SNES snes);
13 extern PetscErrorCode SNESSetUp_FAS(SNES snes);
14 extern PetscErrorCode SNESSetFromOptions_FAS(SNES snes);
15 extern PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer);
16 extern PetscErrorCode SNESSolve_FAS(SNES snes);
17 extern PetscErrorCode SNESReset_FAS(SNES snes);
18 
19 EXTERN_C_BEGIN
20 
21 #undef __FUNCT__
22 #define __FUNCT__ "SNESCreate_FAS"
23 PetscErrorCode SNESCreate_FAS(SNES snes)
24 {
25   SNES_FAS * fas;
26   PetscErrorCode ierr;
27 
28   PetscFunctionBegin;
29   snes->ops->destroy        = SNESDestroy_FAS;
30   snes->ops->setup          = SNESSetUp_FAS;
31   snes->ops->setfromoptions = SNESSetFromOptions_FAS;
32   snes->ops->view           = SNESView_FAS;
33   snes->ops->solve          = SNESSolve_FAS;
34   snes->ops->reset          = SNESReset_FAS;
35 
36   ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr);
37   snes->data = (void*) fas;
38   fas->level = 0;
39   fas->presmooth  = PETSC_NULL;
40   fas->postsmooth = PETSC_NULL;
41   fas->next = PETSC_NULL;
42   fas->interpolate = PETSC_NULL;
43   fas->restrct = PETSC_NULL;
44 
45   PetscFunctionReturn(0);
46 }
47 EXTERN_C_END
48 
49 #undef __FUNCT__
50 #define __FUNCT__ "SNESFASGetLevels"
51 PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt * levels) {
52   SNES_FAS * fas = (SNES_FAS *)snes->data;
53   PetscFunctionBegin;
54   *levels = fas->level;
55   PetscFunctionReturn(0);
56 }
57 
58 #undef __FUNCT__
59 #define __FUNCT__ "SNESFASGetSNES"
60 PetscErrorCode SNESFASGetSNES(SNES snes, PetscInt level, SNES * lsnes) {
61   SNES_FAS * fas = (SNES_FAS *)snes->data;
62   PetscInt levels = fas->level;
63   PetscInt i;
64   PetscFunctionBegin;
65   *lsnes = snes;
66   if (fas->level < level) {
67     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES should only be called on a finer SNESFAS instance than the level.");
68   }
69   if (level > levels - 1) {
70     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Level %d doesn't exist in the SNESFAS.");
71   }
72   for (i = fas->level; i > level; i--) {
73     *lsnes = fas->next;
74     fas = (SNES_FAS *)(*lsnes)->data;
75   }
76   if (fas->level != level) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES didn't return the right level!");
77   PetscFunctionReturn(0);
78 }
79 
80 #undef __FUNCT__
81 #define __FUNCT__ "SNESFASSetLevels"
82 PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) {
83   PetscErrorCode ierr;
84   PetscInt i;
85   SNES_FAS * fas = (SNES_FAS *)snes->data;
86   MPI_Comm comm;
87 
88   PetscFunctionBegin;
89   comm = PETSC_COMM_WORLD;
90   /* user has changed the number of levels; reset */
91   ierr = SNESReset(snes);CHKERRQ(ierr);
92   /* destroy any coarser levels if necessary */
93   if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr);
94   /* setup the finest level */
95   for (i = levels - 1; i >= 0; i--) {
96     if (comms) comm = comms[i];
97     fas->level = i;
98     fas->levels = levels;
99     if (i > 0) {
100       ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr);
101       ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr);
102       ierr = SNESSetOptionsPrefix(fas->next, "fas_");CHKERRQ(ierr);
103       ierr = SNESSetFromOptions(fas->next);CHKERRQ(ierr);
104       fas = (SNES_FAS *)fas->next->data;
105     }
106   }
107   PetscFunctionReturn(0);
108 }
109 
110 #undef __FUNCT__
111 #define __FUNCT__ "SNESFASSetInterpolation"
112 PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) {
113   SNES_FAS * fas =  (SNES_FAS *)snes->data;
114   PetscInt top_level = fas->level,i;
115 
116   PetscFunctionBegin;
117   if (level > top_level)
118     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInterpolation", level);
119   /* get to the correct level */
120   for (i = fas->level; i > level; i--) {
121     fas = (SNES_FAS *)fas->next->data;
122   }
123   if (fas->level != level)
124     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInterpolation");
125   fas->interpolate = mat;
126   PetscFunctionReturn(0);
127 }
128 
129 #undef __FUNCT__
130 #define __FUNCT__ "SNESFASSetRestriction"
131 PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) {
132   SNES_FAS * fas =  (SNES_FAS *)snes->data;
133   PetscInt top_level = fas->level,i;
134 
135   PetscFunctionBegin;
136   if (level > top_level)
137     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level);
138   /* get to the correct level */
139   for (i = fas->level; i > level; i--) {
140     fas = (SNES_FAS *)fas->next->data;
141   }
142   if (fas->level != level)
143     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction");
144   fas->restrct = mat;
145   PetscFunctionReturn(0);
146 }
147 
148 #undef __FUNCT__
149 #define __FUNCT__ "SNESFASSetRScale"
150 PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) {
151   SNES_FAS * fas =  (SNES_FAS *)snes->data;
152   PetscInt top_level = fas->level,i;
153 
154   PetscFunctionBegin;
155   if (level > top_level)
156     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level);
157   /* get to the correct level */
158   for (i = fas->level; i > level; i--) {
159     fas = (SNES_FAS *)fas->next->data;
160   }
161   if (fas->level != level)
162     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction");
163   fas->rscale = rscale;
164   PetscFunctionReturn(0);
165 }
166 
167 #undef __FUNCT__
168 #define __FUNCT__ "SNESReset_FAS"
169 PetscErrorCode SNESReset_FAS(SNES snes)
170 {
171   PetscErrorCode ierr = 0;
172   SNES_FAS * fas = (SNES_FAS *)snes->data;
173 
174   PetscFunctionBegin;
175   /* destroy local data created in SNESSetup_FAS */
176   if (fas->presmooth)    ierr = SNESDestroy(&fas->presmooth);CHKERRQ(ierr);
177   if (fas->postsmooth)   ierr = SNESDestroy(&fas->postsmooth);CHKERRQ(ierr);
178   if (fas->interpolate)  ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
179   if (fas->restrct)      ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
180   if (fas->rscale)       ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
181 
182   /* recurse -- reset should NOT destroy the structures -- destroy should destroy the structures recursively */
183   if (fas->next) ierr = SNESReset(fas->next);CHKERRQ(ierr);
184   if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);}
185   PetscFunctionReturn(0);
186 }
187 
188 #undef __FUNCT__
189 #define __FUNCT__ "SNESDestroy_FAS"
190 PetscErrorCode SNESDestroy_FAS(SNES snes)
191 {
192   SNES_FAS * fas = (SNES_FAS *)snes->data;
193   PetscErrorCode ierr;
194 
195   PetscFunctionBegin;
196   /* recursively resets and then destroys */
197   ierr = SNESReset_FAS(snes);CHKERRQ(ierr);
198   if (fas->next) ierr = SNESDestroy(&fas->next);CHKERRQ(ierr);
199   ierr = PetscFree(fas);CHKERRQ(ierr);
200   PetscFunctionReturn(0);
201 }
202 
203 #undef __FUNCT__
204 #define __FUNCT__ "SNESSetUp_FAS"
205 PetscErrorCode SNESSetUp_FAS(SNES snes)
206 {
207   SNES_FAS   *fas = (SNES_FAS *) snes->data;
208   PetscErrorCode ierr;
209 
210   PetscFunctionBegin;
211   ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* the four work vectors are used to transfer stuff BACK */
212   /* gets the solver ready for solution */
213   if (snes->dm) {
214     /* construct EVERYTHING from the DM -- including the progressive set of smoothers */
215     if (fas->next) {
216       /* for now -- assume the DM and the evaluation functions have been set externally */
217       if (!fas->interpolate) {
218         if (!fas->next->dm) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONGSTATE, "All levels must be assigned a DM");
219         /* set the interpolation and restriction from the DM */
220         ierr = DMGetInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr);
221         fas->restrct = fas->interpolate;
222       }
223       /* TODO LATER: Preconditioner setup goes here */
224     }
225   }
226   if (fas->next) {
227     /* gotta set up the solution vector for this to work */
228     ierr = VecDuplicate(fas->rscale, &fas->next->vec_sol);CHKERRQ(ierr);
229     ierr = SNESSetUp(fas->next);CHKERRQ(ierr);
230   }
231   /* got to set them all up at once */
232   PetscFunctionReturn(0);
233 }
234 
235 #undef __FUNCT__
236 #define __FUNCT__ "SNESSetFromOptions_FAS"
237 PetscErrorCode SNESSetFromOptions_FAS(SNES snes)
238 {
239 
240   PetscErrorCode ierr;
241 
242   PetscFunctionBegin;
243   /* types for the pre and postsmoothers */
244   ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr);
245   ierr = PetscOptionsTail();CHKERRQ(ierr);
246   PetscFunctionReturn(0);
247 }
248 
249 #undef __FUNCT__
250 #define __FUNCT__ "SNESView_FAS"
251 PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
252 {
253   SNES_FAS   *fas = (SNES_FAS *) snes->data;
254   PetscBool      iascii;
255   PetscErrorCode ierr;
256   PetscInt levels = fas->levels;
257   PetscInt i;
258 
259   PetscFunctionBegin;
260   ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
261   if (iascii) {
262     ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %d\n",  fas->levels);CHKERRQ(ierr);
263     ierr = PetscViewerASCIIPushTab(viewer);
264     for (i = levels - 1; i >= 0; i--) {
265       ierr = PetscViewerASCIIPrintf(viewer, "level: %d\n",  fas->level);CHKERRQ(ierr);
266       if (fas->presmooth) {
267         ierr = PetscViewerASCIIPrintf(viewer, "pre-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
268         ierr = PetscViewerASCIIPushTab(viewer);
269         ierr = SNESView(fas->presmooth, viewer);CHKERRQ(ierr);
270         ierr = PetscViewerASCIIPopTab(viewer);
271       } else {
272         ierr = PetscViewerASCIIPrintf(viewer, "no pre-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
273       }
274       if (fas->postsmooth) {
275         ierr = PetscViewerASCIIPrintf(viewer, "post-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
276         ierr = PetscViewerASCIIPushTab(viewer);
277         ierr = SNESView(fas->postsmooth, viewer);CHKERRQ(ierr);
278         ierr = PetscViewerASCIIPopTab(viewer);
279       } else {
280         ierr = PetscViewerASCIIPrintf(viewer, "no post-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
281       }
282       if (fas->next) fas = (SNES_FAS *)fas->next->data;
283     }
284     ierr = PetscViewerASCIIPopTab(viewer);
285   } else {
286     SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name);
287   }
288   PetscFunctionReturn(0);
289 }
290 
291 /*
292 
293 Defines the FAS cycle as:
294 
295 fine problem: F(x) = 0
296 coarse problem: F^c(x) = b^c
297 
298 b^c = F^c(I^c_fx^f - I^c_fF(x))
299 
300 correction:
301 
302 x = x + I(x^c - Rx)
303 
304  */
305 
306 #undef __FUNCT__
307 #define __FUNCT__ "FASCycle_Private"
308 PetscErrorCode FASCycle_Private(SNES snes, Vec B, Vec X, Vec F) {
309 
310   PetscErrorCode ierr;
311   Vec X_c, Xo_c, F_c, B_c;
312   SNES_FAS * fas = (SNES_FAS *)snes->data;
313 
314   PetscFunctionBegin;
315   /* pre-smooth -- just update using the pre-smoother */
316   if (fas->presmooth) {
317     ierr = SNESSolve(fas->presmooth, B, X);CHKERRQ(ierr);
318   } else if (snes->pc) {
319     ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
320   }
321   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
322   if (fas->next) {
323     X_c  = fas->next->vec_sol;
324     Xo_c = fas->next->work[1];
325     F_c  = fas->next->vec_func;
326     B_c  = fas->next->work[2];
327     /* project to coarse */
328     ierr = MatRestrict(fas->restrct, X, Xo_c);CHKERRQ(ierr);
329     ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr);
330     ierr = MatRestrict(fas->restrct, F, F_c);CHKERRQ(ierr);
331     ierr = VecPointwiseMult(F_c,  fas->rscale, F_c);CHKERRQ(ierr);
332     if (B) {
333       ierr = MatRestrict(fas->restrct, B, B_c);CHKERRQ(ierr);
334       ierr = VecPointwiseMult(B_c,  fas->rscale, B_c);CHKERRQ(ierr);
335     } else {
336       ierr = VecSet(B_c, 0.0);CHKERRQ(ierr);
337     }
338     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
339     fas->next->vec_rhs = PETSC_NULL; /*unset the RHS for the next problem so we may evaluate the function rather than the residual */
340     ierr = VecAXPY(B_c, -1.0, F_c);CHKERRQ(ierr);                     /* B_c = RB + F(X_c) - F_C */
341     ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr);
342     ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr);
343     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
344 
345     /* test -- initial residuals with and without the RHS set */
346     fas->next->vec_rhs = B_c;
347     ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr);
348     fas->next->vec_rhs = PETSC_NULL;
349     ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr);
350 
351     /* recurse to the next level */
352     ierr = FASCycle_Private(fas->next, B_c, X_c, F_c);CHKERRQ(ierr);
353     /* correct as x <- x + I(x^c - Rx)*/
354     ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
355     ierr = MatInterpolate(fas->interpolate, X_c, F);CHKERRQ(ierr);
356     ierr = VecAXPY(X, 1.0, F);CHKERRQ(ierr);
357   }
358   /* post-smooth -- just update using the post-smoother */
359   if (fas->level != 0) {
360     if (fas->postsmooth) {
361       ierr = SNESSolve(fas->postsmooth, B, X);CHKERRQ(ierr);
362     } else if (snes->pc) {
363       ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
364     }
365   }
366   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
367 
368   PetscFunctionReturn(0);
369 }
370 
371 #undef __FUNCT__
372 #define __FUNCT__ "FASInitialGuess_Private"
373 
374 
375 #undef __FUNCT__
376 #define __FUNCT__ "SNESSolve_FAS"
377 
378 PetscErrorCode SNESSolve_FAS(SNES snes)
379 {
380   PetscErrorCode ierr;
381   PetscInt i, maxits;
382   Vec X, F, B;
383   PetscReal fnorm;
384   PetscFunctionBegin;
385   maxits = snes->max_its;	     /* maximum number of iterations */
386   snes->reason = SNES_CONVERGED_ITERATING;
387   X = snes->vec_sol;
388   F = snes->vec_func;
389   B = snes->vec_rhs;
390   /* initial iteration */
391   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
392   snes->iter = 0;
393   snes->norm = 0.;
394   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
395   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
396   if (snes->domainerror) {
397     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
398     PetscFunctionReturn(0);
399   }
400   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
401   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
402   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
403   snes->norm = fnorm;
404   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
405   SNESLogConvHistory(snes,fnorm,0);
406   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
407 
408   /* set parameter for default relative tolerance convergence test */
409   snes->ttol = fnorm*snes->rtol;
410   /* test convergence */
411   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
412   if (snes->reason) PetscFunctionReturn(0);
413   for (i = 0; i < maxits; i++) {
414     /* Call general purpose update function */
415     if (snes->ops->update) {
416       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
417     }
418     ierr = FASCycle_Private(snes, B, X, F);CHKERRQ(ierr);
419     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
420     /* Monitor convergence */
421     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
422     snes->iter = i+1;
423     snes->norm = fnorm;
424     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
425     SNESLogConvHistory(snes,snes->norm,0);
426     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
427     /* Test for convergence */
428     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
429     if (snes->reason) break;
430   }
431   if (i == maxits) {
432     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
433     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
434   }
435   PetscFunctionReturn(0);
436 }
437