xref: /petsc/src/snes/impls/fas/fas.c (revision 421d9b3201840ccac4284b98ab49dc267dcc958d)
1*421d9b32SPeter Brune /* Defines the basic SNES object */
2*421d9b32SPeter Brune #include <../src/snes/impls/fas/fasimpls.h>
3*421d9b32SPeter Brune 
4*421d9b32SPeter Brune /*MC
5*421d9b32SPeter Brune Full Approximation Scheme nonlinear multigrid solver.
6*421d9b32SPeter Brune 
7*421d9b32SPeter Brune The nonlinear problem is solved via the repeated application of nonlinear preconditioners and coarse-grid corrections
8*421d9b32SPeter Brune 
9*421d9b32SPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
10*421d9b32SPeter Brune M*/
11*421d9b32SPeter Brune 
12*421d9b32SPeter Brune extern PetscErrorCode SNESDestroy_FAS(SNES snes);
13*421d9b32SPeter Brune extern PetscErrorCode SNESSetUp_FAS(SNES snes);
14*421d9b32SPeter Brune extern PetscErrorCode SNESSetFromOptions_FAS(SNES snes);
15*421d9b32SPeter Brune extern PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer);
16*421d9b32SPeter Brune extern PetscErrorCode SNESSolve_FAS(SNES snes);
17*421d9b32SPeter Brune extern PetscErrorCode SNESReset_FAS(SNES snes);
18*421d9b32SPeter Brune 
19*421d9b32SPeter Brune EXTERN_C_BEGIN
20*421d9b32SPeter Brune 
21*421d9b32SPeter Brune #undef __FUNCT__
22*421d9b32SPeter Brune #define __FUNCT__ "SNESCreate_FAS"
23*421d9b32SPeter Brune PetscErrorCode SNESCreate_FAS(SNES snes)
24*421d9b32SPeter Brune {
25*421d9b32SPeter Brune   SNES_FAS * fas;
26*421d9b32SPeter Brune   PetscErrorCode ierr;
27*421d9b32SPeter Brune 
28*421d9b32SPeter Brune   PetscFunctionBegin;
29*421d9b32SPeter Brune   snes->ops->destroy        = SNESDestroy_FAS;
30*421d9b32SPeter Brune   snes->ops->setup          = SNESSetUp_FAS;
31*421d9b32SPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_FAS;
32*421d9b32SPeter Brune   snes->ops->view           = SNESView_FAS;
33*421d9b32SPeter Brune   snes->ops->solve          = SNESSolve_FAS;
34*421d9b32SPeter Brune   snes->ops->reset          = SNESReset_FAS;
35*421d9b32SPeter Brune 
36*421d9b32SPeter Brune   ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr);
37*421d9b32SPeter Brune   snes->data = (void*) fas;
38*421d9b32SPeter Brune   fas->level = 0;
39*421d9b32SPeter Brune   fas->presmooth  = PETSC_NULL;
40*421d9b32SPeter Brune   fas->postsmooth = PETSC_NULL;
41*421d9b32SPeter Brune   fas->next = PETSC_NULL;
42*421d9b32SPeter Brune   fas->interpolate = PETSC_NULL;
43*421d9b32SPeter Brune   fas->restrct = PETSC_NULL;
44*421d9b32SPeter Brune 
45*421d9b32SPeter Brune   PetscFunctionReturn(0);
46*421d9b32SPeter Brune }
47*421d9b32SPeter Brune EXTERN_C_END
48*421d9b32SPeter Brune 
49*421d9b32SPeter Brune #undef __FUNCT__
50*421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetLevels"
51*421d9b32SPeter Brune PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt * levels) {
52*421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
53*421d9b32SPeter Brune   PetscFunctionBegin;
54*421d9b32SPeter Brune   *levels = fas->level;
55*421d9b32SPeter Brune   PetscFunctionReturn(0);
56*421d9b32SPeter Brune }
57*421d9b32SPeter Brune 
58*421d9b32SPeter Brune #undef __FUNCT__
59*421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetSNES"
60*421d9b32SPeter Brune PetscErrorCode SNESFASGetSNES(SNES snes, PetscInt level, SNES * lsnes) {
61*421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
62*421d9b32SPeter Brune   PetscInt levels = fas->level;
63*421d9b32SPeter Brune   PetscInt i;
64*421d9b32SPeter Brune   PetscFunctionBegin;
65*421d9b32SPeter Brune   *lsnes = snes;
66*421d9b32SPeter Brune   if (fas->level < level) {
67*421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES should only be called on a finer SNESFAS instance than the level.");
68*421d9b32SPeter Brune   }
69*421d9b32SPeter Brune   if (level > levels - 1) {
70*421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Level %d doesn't exist in the SNESFAS.");
71*421d9b32SPeter Brune   }
72*421d9b32SPeter Brune   for (i = fas->level; i > level; i--) {
73*421d9b32SPeter Brune     *lsnes = fas->next;
74*421d9b32SPeter Brune     fas = (SNES_FAS *)(*lsnes)->data;
75*421d9b32SPeter Brune   }
76*421d9b32SPeter Brune   if (fas->level != level) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES didn't return the right level!");
77*421d9b32SPeter Brune   PetscFunctionReturn(0);
78*421d9b32SPeter Brune }
79*421d9b32SPeter Brune 
80*421d9b32SPeter Brune #undef __FUNCT__
81*421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetLevels"
82*421d9b32SPeter Brune PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) {
83*421d9b32SPeter Brune   PetscErrorCode ierr;
84*421d9b32SPeter Brune   PetscInt i;
85*421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
86*421d9b32SPeter Brune   MPI_Comm comm;
87*421d9b32SPeter Brune 
88*421d9b32SPeter Brune   PetscFunctionBegin;
89*421d9b32SPeter Brune   comm = PETSC_COMM_WORLD;
90*421d9b32SPeter Brune   /* user has changed the number of levels; reset */
91*421d9b32SPeter Brune   ierr = SNESReset(snes);CHKERRQ(ierr);
92*421d9b32SPeter Brune   /* destroy any coarser levels if necessary */
93*421d9b32SPeter Brune   if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr);
94*421d9b32SPeter Brune   /* setup the finest level */
95*421d9b32SPeter Brune   for (i = levels - 1; i >= 0; i--) {
96*421d9b32SPeter Brune     if (comms) comm = comms[i];
97*421d9b32SPeter Brune     fas->level = i;
98*421d9b32SPeter Brune     fas->levels = levels;
99*421d9b32SPeter Brune     if (i > 0) {
100*421d9b32SPeter Brune       ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr);
101*421d9b32SPeter Brune       ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr);
102*421d9b32SPeter Brune       fas = (SNES_FAS *)fas->next->data;
103*421d9b32SPeter Brune     }
104*421d9b32SPeter Brune   }
105*421d9b32SPeter Brune   PetscFunctionReturn(0);
106*421d9b32SPeter Brune }
107*421d9b32SPeter Brune 
108*421d9b32SPeter Brune #undef __FUNCT__
109*421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetInterpolation"
110*421d9b32SPeter Brune PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) {
111*421d9b32SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
112*421d9b32SPeter Brune   PetscInt top_level = fas->level;
113*421d9b32SPeter Brune 
114*421d9b32SPeter Brune   PetscFunctionBegin;
115*421d9b32SPeter Brune   if (level > top_level)
116*421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInterpolation", level);
117*421d9b32SPeter Brune   /* get to the correct level */
118*421d9b32SPeter Brune   for (int i = fas->level; i > level; i--) {
119*421d9b32SPeter Brune     fas = (SNES_FAS *)fas->next->data;
120*421d9b32SPeter Brune   }
121*421d9b32SPeter Brune   if (fas->level != level)
122*421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInterpolation");
123*421d9b32SPeter Brune   fas->interpolate = mat;
124*421d9b32SPeter Brune   PetscFunctionReturn(0);
125*421d9b32SPeter Brune }
126*421d9b32SPeter Brune 
127*421d9b32SPeter Brune #undef __FUNCT__
128*421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRestriction"
129*421d9b32SPeter Brune PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) {
130*421d9b32SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
131*421d9b32SPeter Brune   PetscInt top_level = fas->level;
132*421d9b32SPeter Brune 
133*421d9b32SPeter Brune   PetscFunctionBegin;
134*421d9b32SPeter Brune   if (level > top_level)
135*421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level);
136*421d9b32SPeter Brune   /* get to the correct level */
137*421d9b32SPeter Brune   for (int i = fas->level; i > level; i--) {
138*421d9b32SPeter Brune     fas = (SNES_FAS *)fas->next->data;
139*421d9b32SPeter Brune   }
140*421d9b32SPeter Brune   if (fas->level != level)
141*421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction");
142*421d9b32SPeter Brune   fas->restrct = mat;
143*421d9b32SPeter Brune   PetscFunctionReturn(0);
144*421d9b32SPeter Brune }
145*421d9b32SPeter Brune 
146*421d9b32SPeter Brune #undef __FUNCT__
147*421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRScale"
148*421d9b32SPeter Brune PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) {
149*421d9b32SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
150*421d9b32SPeter Brune   PetscInt top_level = fas->level;
151*421d9b32SPeter Brune 
152*421d9b32SPeter Brune   PetscFunctionBegin;
153*421d9b32SPeter Brune   if (level > top_level)
154*421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level);
155*421d9b32SPeter Brune   /* get to the correct level */
156*421d9b32SPeter Brune   for (int i = fas->level; i > level; i--) {
157*421d9b32SPeter Brune     fas = (SNES_FAS *)fas->next->data;
158*421d9b32SPeter Brune   }
159*421d9b32SPeter Brune   if (fas->level != level)
160*421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction");
161*421d9b32SPeter Brune   fas->rscale = rscale;
162*421d9b32SPeter Brune   PetscFunctionReturn(0);
163*421d9b32SPeter Brune }
164*421d9b32SPeter Brune 
165*421d9b32SPeter Brune #undef __FUNCT__
166*421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS"
167*421d9b32SPeter Brune PetscErrorCode SNESReset_FAS(SNES snes)
168*421d9b32SPeter Brune {
169*421d9b32SPeter Brune   PetscErrorCode ierr;
170*421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
171*421d9b32SPeter Brune 
172*421d9b32SPeter Brune   PetscFunctionBegin;
173*421d9b32SPeter Brune   /* destroy local data created in SNESSetup_FAS */
174*421d9b32SPeter Brune   if (fas->presmooth)    ierr = SNESDestroy(&fas->presmooth);CHKERRQ(ierr);
175*421d9b32SPeter Brune   if (fas->postsmooth)   ierr = SNESDestroy(&fas->postsmooth);CHKERRQ(ierr);
176*421d9b32SPeter Brune   if (fas->interpolate)  ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
177*421d9b32SPeter Brune   if (fas->restrct)      ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
178*421d9b32SPeter Brune   if (fas->rscale)       ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
179*421d9b32SPeter Brune 
180*421d9b32SPeter Brune     /* recurse -- reset should NOT destroy the structures -- destroy should destroy the structures recursively */
181*421d9b32SPeter Brune   if (fas->next) ierr = SNESReset(fas->next);CHKERRQ(ierr);
182*421d9b32SPeter Brune 
183*421d9b32SPeter Brune   PetscFunctionReturn(0);
184*421d9b32SPeter Brune }
185*421d9b32SPeter Brune 
186*421d9b32SPeter Brune #undef __FUNCT__
187*421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS"
188*421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes)
189*421d9b32SPeter Brune {
190*421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
191*421d9b32SPeter Brune   PetscErrorCode ierr;
192*421d9b32SPeter Brune 
193*421d9b32SPeter Brune   PetscFunctionBegin;
194*421d9b32SPeter Brune   /* recursively resets and then destroys */
195*421d9b32SPeter Brune   ierr = SNESReset_FAS(snes);CHKERRQ(ierr);
196*421d9b32SPeter Brune   if (fas->next) ierr = SNESDestroy(&fas->next);CHKERRQ(ierr);
197*421d9b32SPeter Brune   ierr = PetscFree(fas);CHKERRQ(ierr);
198*421d9b32SPeter Brune   PetscFunctionReturn(0);
199*421d9b32SPeter Brune }
200*421d9b32SPeter Brune 
201*421d9b32SPeter Brune #undef __FUNCT__
202*421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS"
203*421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes)
204*421d9b32SPeter Brune {
205*421d9b32SPeter Brune   SNES_FAS   *fas = (SNES_FAS *) snes->data;
206*421d9b32SPeter Brune   PetscErrorCode ierr;
207*421d9b32SPeter Brune 
208*421d9b32SPeter Brune   PetscFunctionBegin;
209*421d9b32SPeter Brune 
210*421d9b32SPeter Brune   /* gets the solver ready for solution */
211*421d9b32SPeter Brune   if (snes->dm) {
212*421d9b32SPeter Brune     /* construct EVERYTHING from the DM -- including the progressive set of smoothers */
213*421d9b32SPeter Brune     if (fas->next) {
214*421d9b32SPeter Brune       /* for now -- assume the DM and the evaluation functions have been set externally */
215*421d9b32SPeter Brune       if (!fas->interpolate) {
216*421d9b32SPeter Brune         if (!fas->next->dm) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONGSTATE, "All levels must be assigned a DM");
217*421d9b32SPeter Brune         /* set the interpolation and restriction from the DM */
218*421d9b32SPeter Brune         ierr = DMGetInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr);
219*421d9b32SPeter Brune         fas->restrct = fas->interpolate;
220*421d9b32SPeter Brune       }
221*421d9b32SPeter Brune       /* TODO LATER: Preconditioner setup goes here */
222*421d9b32SPeter Brune     }
223*421d9b32SPeter Brune   } else {
224*421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_SUP, "SNESSetup_FAS presently only works with DM Coarsening.  This will be fixed. ");
225*421d9b32SPeter Brune   }
226*421d9b32SPeter Brune 
227*421d9b32SPeter Brune   PetscFunctionReturn(0);
228*421d9b32SPeter Brune }
229*421d9b32SPeter Brune 
230*421d9b32SPeter Brune #undef __FUNCT__
231*421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS"
232*421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(SNES snes)
233*421d9b32SPeter Brune {
234*421d9b32SPeter Brune 
235*421d9b32SPeter Brune   PetscErrorCode ierr;
236*421d9b32SPeter Brune 
237*421d9b32SPeter Brune   PetscFunctionBegin;
238*421d9b32SPeter Brune   /* types for the pre and postsmoothers */
239*421d9b32SPeter Brune   ierr = PetscOptionsHead("SNESFAS Options");CHKERRQ(ierr);
240*421d9b32SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
241*421d9b32SPeter Brune   PetscFunctionReturn(0);
242*421d9b32SPeter Brune }
243*421d9b32SPeter Brune 
244*421d9b32SPeter Brune #undef __FUNCT__
245*421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS"
246*421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
247*421d9b32SPeter Brune {
248*421d9b32SPeter Brune   SNES_FAS   *fas = (SNES_FAS *) snes->data;
249*421d9b32SPeter Brune   PetscBool      iascii;
250*421d9b32SPeter Brune   PetscErrorCode ierr;
251*421d9b32SPeter Brune   PetscInt levels = fas->levels;
252*421d9b32SPeter Brune   PetscInt i;
253*421d9b32SPeter Brune 
254*421d9b32SPeter Brune   PetscFunctionBegin;
255*421d9b32SPeter Brune   ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
256*421d9b32SPeter Brune   if (iascii) {
257*421d9b32SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %d\n",  fas->levels);CHKERRQ(ierr);
258*421d9b32SPeter Brune     ierr = PetscViewerASCIIPushTab(viewer);
259*421d9b32SPeter Brune     for (i = levels - 1; i >= 0; i--) {
260*421d9b32SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "level: %d\n",  fas->level);CHKERRQ(ierr);
261*421d9b32SPeter Brune       if (fas->presmooth) {
262*421d9b32SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "pre-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
263*421d9b32SPeter Brune         ierr = PetscViewerASCIIPushTab(viewer);
264*421d9b32SPeter Brune         ierr = SNESView(fas->presmooth, viewer);CHKERRQ(ierr);
265*421d9b32SPeter Brune         ierr = PetscViewerASCIIPopTab(viewer);
266*421d9b32SPeter Brune       } else {
267*421d9b32SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "no pre-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
268*421d9b32SPeter Brune       }
269*421d9b32SPeter Brune       if (fas->postsmooth) {
270*421d9b32SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "post-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
271*421d9b32SPeter Brune         ierr = PetscViewerASCIIPushTab(viewer);
272*421d9b32SPeter Brune         ierr = SNESView(fas->postsmooth, viewer);CHKERRQ(ierr);
273*421d9b32SPeter Brune         ierr = PetscViewerASCIIPopTab(viewer);
274*421d9b32SPeter Brune       } else {
275*421d9b32SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "no post-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
276*421d9b32SPeter Brune       }
277*421d9b32SPeter Brune       if (fas->next) fas = (SNES_FAS *)fas->next->data;
278*421d9b32SPeter Brune     }
279*421d9b32SPeter Brune     ierr = PetscViewerASCIIPopTab(viewer);
280*421d9b32SPeter Brune   } else {
281*421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name);
282*421d9b32SPeter Brune   }
283*421d9b32SPeter Brune   PetscFunctionReturn(0);
284*421d9b32SPeter Brune }
285*421d9b32SPeter Brune 
286*421d9b32SPeter Brune /*
287*421d9b32SPeter Brune  */
288*421d9b32SPeter Brune #undef __FUNCT__
289*421d9b32SPeter Brune #define __FUNCT__ "FASCycle_Private"
290*421d9b32SPeter Brune PetscErrorCode FASCycle_Private(SNES snes) {
291*421d9b32SPeter Brune 
292*421d9b32SPeter Brune   PetscFunctionBegin;
293*421d9b32SPeter Brune   PetscFunctionReturn(0);
294*421d9b32SPeter Brune 
295*421d9b32SPeter Brune }
296*421d9b32SPeter Brune 
297*421d9b32SPeter Brune #undef __FUNCT__
298*421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS"
299*421d9b32SPeter Brune 
300*421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes)
301*421d9b32SPeter Brune {
302*421d9b32SPeter Brune   PetscFunctionBegin;
303*421d9b32SPeter Brune 
304*421d9b32SPeter Brune   PetscFunctionReturn(0);
305*421d9b32SPeter Brune }
306