xref: /petsc/src/snes/impls/fas/fas.c (revision 6fc7ef2bf116ea7e508976bc51f847a86697e13b)
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       fas = (SNES_FAS *)fas->next->data;
103     }
104   }
105   PetscFunctionReturn(0);
106 }
107 
108 #undef __FUNCT__
109 #define __FUNCT__ "SNESFASSetInterpolation"
110 PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) {
111   SNES_FAS * fas =  (SNES_FAS *)snes->data;
112   PetscInt top_level = fas->level,i;
113 
114   PetscFunctionBegin;
115   if (level > top_level)
116     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInterpolation", level);
117   /* get to the correct level */
118   for (i = fas->level; i > level; i--) {
119     fas = (SNES_FAS *)fas->next->data;
120   }
121   if (fas->level != level)
122     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInterpolation");
123   fas->interpolate = mat;
124   PetscFunctionReturn(0);
125 }
126 
127 #undef __FUNCT__
128 #define __FUNCT__ "SNESFASSetRestriction"
129 PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) {
130   SNES_FAS * fas =  (SNES_FAS *)snes->data;
131   PetscInt top_level = fas->level,i;
132 
133   PetscFunctionBegin;
134   if (level > top_level)
135     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level);
136   /* get to the correct level */
137   for (i = fas->level; i > level; i--) {
138     fas = (SNES_FAS *)fas->next->data;
139   }
140   if (fas->level != level)
141     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction");
142   fas->restrct = mat;
143   PetscFunctionReturn(0);
144 }
145 
146 #undef __FUNCT__
147 #define __FUNCT__ "SNESFASSetRScale"
148 PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) {
149   SNES_FAS * fas =  (SNES_FAS *)snes->data;
150   PetscInt top_level = fas->level,i;
151 
152   PetscFunctionBegin;
153   if (level > top_level)
154     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level);
155   /* get to the correct level */
156   for (i = fas->level; i > level; i--) {
157     fas = (SNES_FAS *)fas->next->data;
158   }
159   if (fas->level != level)
160     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction");
161   fas->rscale = rscale;
162   PetscFunctionReturn(0);
163 }
164 
165 #undef __FUNCT__
166 #define __FUNCT__ "SNESReset_FAS"
167 PetscErrorCode SNESReset_FAS(SNES snes)
168 {
169   PetscErrorCode ierr = 0;
170   SNES_FAS * fas = (SNES_FAS *)snes->data;
171 
172   PetscFunctionBegin;
173   /* destroy local data created in SNESSetup_FAS */
174   if (fas->presmooth)    ierr = SNESDestroy(&fas->presmooth);CHKERRQ(ierr);
175   if (fas->postsmooth)   ierr = SNESDestroy(&fas->postsmooth);CHKERRQ(ierr);
176   if (fas->interpolate)  ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
177   if (fas->restrct)      ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
178   if (fas->rscale)       ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
179 
180     /* recurse -- reset should NOT destroy the structures -- destroy should destroy the structures recursively */
181   if (fas->next) ierr = SNESReset(fas->next);CHKERRQ(ierr);
182 
183   PetscFunctionReturn(0);
184 }
185 
186 #undef __FUNCT__
187 #define __FUNCT__ "SNESDestroy_FAS"
188 PetscErrorCode SNESDestroy_FAS(SNES snes)
189 {
190   SNES_FAS * fas = (SNES_FAS *)snes->data;
191   PetscErrorCode ierr;
192 
193   PetscFunctionBegin;
194   /* recursively resets and then destroys */
195   ierr = SNESReset_FAS(snes);CHKERRQ(ierr);
196   if (fas->next) ierr = SNESDestroy(&fas->next);CHKERRQ(ierr);
197   ierr = PetscFree(fas);CHKERRQ(ierr);
198   PetscFunctionReturn(0);
199 }
200 
201 #undef __FUNCT__
202 #define __FUNCT__ "SNESSetUp_FAS"
203 PetscErrorCode SNESSetUp_FAS(SNES snes)
204 {
205   SNES_FAS   *fas = (SNES_FAS *) snes->data;
206   PetscErrorCode ierr;
207 
208   PetscFunctionBegin;
209 
210   /* gets the solver ready for solution */
211   if (snes->dm) {
212     /* construct EVERYTHING from the DM -- including the progressive set of smoothers */
213     if (fas->next) {
214       /* for now -- assume the DM and the evaluation functions have been set externally */
215       if (!fas->interpolate) {
216         if (!fas->next->dm) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONGSTATE, "All levels must be assigned a DM");
217         /* set the interpolation and restriction from the DM */
218         ierr = DMGetInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr);
219         fas->restrct = fas->interpolate;
220       }
221       /* TODO LATER: Preconditioner setup goes here */
222     }
223   } else {
224     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_SUP, "SNESSetup_FAS presently only works with DM Coarsening.  This will be fixed. ");
225   }
226 
227   PetscFunctionReturn(0);
228 }
229 
230 #undef __FUNCT__
231 #define __FUNCT__ "SNESSetFromOptions_FAS"
232 PetscErrorCode SNESSetFromOptions_FAS(SNES snes)
233 {
234 
235   PetscErrorCode ierr;
236 
237   PetscFunctionBegin;
238   /* types for the pre and postsmoothers */
239   ierr = PetscOptionsHead("SNESFAS Options");CHKERRQ(ierr);
240   ierr = PetscOptionsTail();CHKERRQ(ierr);
241   PetscFunctionReturn(0);
242 }
243 
244 #undef __FUNCT__
245 #define __FUNCT__ "SNESView_FAS"
246 PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
247 {
248   SNES_FAS   *fas = (SNES_FAS *) snes->data;
249   PetscBool      iascii;
250   PetscErrorCode ierr;
251   PetscInt levels = fas->levels;
252   PetscInt i;
253 
254   PetscFunctionBegin;
255   ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
256   if (iascii) {
257     ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %d\n",  fas->levels);CHKERRQ(ierr);
258     ierr = PetscViewerASCIIPushTab(viewer);
259     for (i = levels - 1; i >= 0; i--) {
260       ierr = PetscViewerASCIIPrintf(viewer, "level: %d\n",  fas->level);CHKERRQ(ierr);
261       if (fas->presmooth) {
262         ierr = PetscViewerASCIIPrintf(viewer, "pre-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
263         ierr = PetscViewerASCIIPushTab(viewer);
264         ierr = SNESView(fas->presmooth, viewer);CHKERRQ(ierr);
265         ierr = PetscViewerASCIIPopTab(viewer);
266       } else {
267         ierr = PetscViewerASCIIPrintf(viewer, "no pre-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
268       }
269       if (fas->postsmooth) {
270         ierr = PetscViewerASCIIPrintf(viewer, "post-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
271         ierr = PetscViewerASCIIPushTab(viewer);
272         ierr = SNESView(fas->postsmooth, viewer);CHKERRQ(ierr);
273         ierr = PetscViewerASCIIPopTab(viewer);
274       } else {
275         ierr = PetscViewerASCIIPrintf(viewer, "no post-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
276       }
277       if (fas->next) fas = (SNES_FAS *)fas->next->data;
278     }
279     ierr = PetscViewerASCIIPopTab(viewer);
280   } else {
281     SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name);
282   }
283   PetscFunctionReturn(0);
284 }
285 
286 /*
287  */
288 #undef __FUNCT__
289 #define __FUNCT__ "FASCycle_Private"
290 PetscErrorCode FASCycle_Private(SNES snes) {
291 
292   PetscFunctionBegin;
293   PetscFunctionReturn(0);
294 
295 }
296 
297 #undef __FUNCT__
298 #define __FUNCT__ "SNESSolve_FAS"
299 
300 PetscErrorCode SNESSolve_FAS(SNES snes)
301 {
302   PetscFunctionBegin;
303 
304   PetscFunctionReturn(0);
305 }
306