xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision 8bc8193efbc389280f83b3d41dffa9e2d23e2ace)
1 /*
2     Defines the multigrid preconditioner interface.
3 */
4 #include "src/ksp/pc/impls/mg/mgimpl.h"                    /*I "petscmg.h" I*/
5 
6 
7 /*
8        MGMCycle_Private - Given an MG structure created with MGCreate() runs
9                   one multiplicative cycle down through the levels and
10                   back up.
11 
12     Input Parameter:
13 .   mg - structure created with  MGCreate().
14 */
15 #undef __FUNCT__
16 #define __FUNCT__ "MGMCycle_Private"
17 PetscErrorCode MGMCycle_Private(MG *mglevels,PetscTruth *converged)
18 {
19   MG             mg = *mglevels,mgc;
20   PetscErrorCode ierr;
21   PetscInt       cycles = mg->cycles;
22   PetscScalar    zero = 0.0;
23 
24   PetscFunctionBegin;
25   if (converged) *converged = PETSC_FALSE;
26 
27   if (mg->eventsolve) {ierr = PetscLogEventBegin(mg->eventsolve,0,0,0,0);CHKERRQ(ierr);}
28   ierr = KSPSolve(mg->smoothd,mg->b,mg->x);CHKERRQ(ierr);
29   if (mg->eventsolve) {ierr = PetscLogEventEnd(mg->eventsolve,0,0,0,0);CHKERRQ(ierr);}
30   if (mg->level) {  /* not the coarsest grid */
31     ierr = (*mg->residual)(mg->A,mg->b,mg->x,mg->r);CHKERRQ(ierr);
32 
33     /* if on finest level and have convergence criteria set */
34     if (mg->level == mg->levels-1 && mg->ttol) {
35       PetscReal rnorm;
36       ierr = VecNorm(mg->r,NORM_2,&rnorm);CHKERRQ(ierr);
37       if (rnorm <= mg->ttol) {
38         *converged = PETSC_TRUE;
39         if (rnorm < mg->abstol) {
40           PetscLogInfo(0,"Linear solver has converged. Residual norm %g is less than absolute tolerance %g\n",rnorm,mg->abstol);
41         } else {
42           PetscLogInfo(0,"Linear solver has converged. Residual norm %g is less than relative tolerance times initial residual norm %g\n",rnorm,mg->ttol);
43         }
44         PetscFunctionReturn(0);
45       }
46     }
47 
48     mgc = *(mglevels - 1);
49     ierr = MatRestrict(mg->restrct,mg->r,mgc->b);CHKERRQ(ierr);
50     ierr = VecSet(&zero,mgc->x);CHKERRQ(ierr);
51     while (cycles--) {
52       ierr = MGMCycle_Private(mglevels-1,converged);CHKERRQ(ierr);
53     }
54     ierr = MatInterpolateAdd(mg->interpolate,mgc->x,mg->x,mg->x);CHKERRQ(ierr);
55     if (mg->eventsolve) {ierr = PetscLogEventBegin(mg->eventsolve,0,0,0,0);CHKERRQ(ierr);}
56     ierr = KSPSolve(mg->smoothu,mg->b,mg->x);CHKERRQ(ierr);
57     if (mg->eventsolve) {ierr = PetscLogEventEnd(mg->eventsolve,0,0,0,0);CHKERRQ(ierr);}
58   }
59   PetscFunctionReturn(0);
60 }
61 
62 /*
63        MGCreate_Private - Creates a MG structure for use with the
64                multigrid code. Level 0 is the coarsest. (But the
65                finest level is stored first in the array).
66 
67 */
68 #undef __FUNCT__
69 #define __FUNCT__ "MGCreate_Private"
70 static PetscErrorCode MGCreate_Private(MPI_Comm comm,PetscInt levels,PC pc,MPI_Comm *comms,MG **result)
71 {
72   MG             *mg;
73   PetscErrorCode ierr;
74   PetscInt       i;
75   PetscMPIInt    size;
76   char           *prefix;
77   PC             ipc;
78 
79   PetscFunctionBegin;
80   ierr = PetscMalloc(levels*sizeof(MG),&mg);CHKERRQ(ierr);
81   PetscLogObjectMemory(pc,levels*(sizeof(MG)+sizeof(struct _MG)));
82 
83   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
84 
85   for (i=0; i<levels; i++) {
86     ierr = PetscNew(struct _MG,&mg[i]);CHKERRQ(ierr);
87     ierr = PetscMemzero(mg[i],sizeof(struct _MG));CHKERRQ(ierr);
88     mg[i]->level  = i;
89     mg[i]->levels = levels;
90     mg[i]->cycles = 1;
91     mg[i]->default_smoothu = 1;
92     mg[i]->default_smoothd = 1;
93 
94     if (comms) comm = comms[i];
95     ierr = KSPCreate(comm,&mg[i]->smoothd);CHKERRQ(ierr);
96     ierr = KSPSetTolerances(mg[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg[i]->default_smoothd);CHKERRQ(ierr);
97     ierr = KSPSetOptionsPrefix(mg[i]->smoothd,prefix);CHKERRQ(ierr);
98 
99     /* do special stuff for coarse grid */
100     if (!i && levels > 1) {
101       ierr = KSPAppendOptionsPrefix(mg[0]->smoothd,"mg_coarse_");CHKERRQ(ierr);
102 
103       /* coarse solve is (redundant) LU by default */
104       ierr = KSPSetType(mg[0]->smoothd,KSPPREONLY);CHKERRQ(ierr);
105       ierr = KSPGetPC(mg[0]->smoothd,&ipc);CHKERRQ(ierr);
106       ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
107       if (size > 1) {
108         ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr);
109         ierr = PCRedundantGetPC(ipc,&ipc);CHKERRQ(ierr);
110       }
111       ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr);
112 
113     } else {
114       char tprefix[128];
115       sprintf(tprefix,"mg_levels_%d_",(int)i);
116       ierr = KSPAppendOptionsPrefix(mg[i]->smoothd,tprefix);CHKERRQ(ierr);
117     }
118     PetscLogObjectParent(pc,mg[i]->smoothd);
119     mg[i]->smoothu         = mg[i]->smoothd;
120     mg[i]->rtol = 0.0;
121     mg[i]->abstol = 0.0;
122     mg[i]->dtol = 0.0;
123     mg[i]->ttol = 0.0;
124     mg[i]->eventsetup = 0;
125     mg[i]->eventsolve = 0;
126   }
127   *result = mg;
128   PetscFunctionReturn(0);
129 }
130 
131 #undef __FUNCT__
132 #define __FUNCT__ "PCDestroy_MG"
133 static PetscErrorCode PCDestroy_MG(PC pc)
134 {
135   MG             *mg = (MG*)pc->data;
136   PetscErrorCode ierr;
137   PetscInt       i,n = mg[0]->levels;
138 
139   PetscFunctionBegin;
140   for (i=0; i<n; i++) {
141     if (mg[i]->smoothd != mg[i]->smoothu) {
142       ierr = KSPDestroy(mg[i]->smoothd);CHKERRQ(ierr);
143     }
144     ierr = KSPDestroy(mg[i]->smoothu);CHKERRQ(ierr);
145     ierr = PetscFree(mg[i]);CHKERRQ(ierr);
146   }
147   ierr = PetscFree(mg);CHKERRQ(ierr);
148   PetscFunctionReturn(0);
149 }
150 
151 
152 
153 EXTERN PetscErrorCode MGACycle_Private(MG*);
154 EXTERN PetscErrorCode MGFCycle_Private(MG*);
155 EXTERN PetscErrorCode MGKCycle_Private(MG*);
156 
157 /*
158    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
159              or full cycle of multigrid.
160 
161   Note:
162   A simple wrapper which calls MGMCycle(),MGACycle(), or MGFCycle().
163 */
164 #undef __FUNCT__
165 #define __FUNCT__ "PCApply_MG"
166 static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
167 {
168   MG             *mg = (MG*)pc->data;
169   PetscScalar    zero = 0.0;
170   PetscErrorCode ierr;
171   PetscInt       levels = mg[0]->levels;
172 
173   PetscFunctionBegin;
174   mg[levels-1]->b = b;
175   mg[levels-1]->x = x;
176   if (mg[0]->am == MGMULTIPLICATIVE) {
177     ierr = VecSet(&zero,x);CHKERRQ(ierr);
178     ierr = MGMCycle_Private(mg+levels-1,PETSC_NULL);CHKERRQ(ierr);
179   }
180   else if (mg[0]->am == MGADDITIVE) {
181     ierr = MGACycle_Private(mg);CHKERRQ(ierr);
182   }
183   else if (mg[0]->am == MGKASKADE) {
184     ierr = MGKCycle_Private(mg);CHKERRQ(ierr);
185   }
186   else {
187     ierr = MGFCycle_Private(mg);CHKERRQ(ierr);
188   }
189   PetscFunctionReturn(0);
190 }
191 
192 #undef __FUNCT__
193 #define __FUNCT__ "PCApplyRichardson_MG"
194 static PetscErrorCode PCApplyRichardson_MG(PC pc,Vec b,Vec x,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its)
195 {
196   MG             *mg = (MG*)pc->data;
197   PetscErrorCode ierr;
198   PetscInt       levels = mg[0]->levels;
199   PetscTruth     converged = PETSC_FALSE;
200 
201   PetscFunctionBegin;
202   mg[levels-1]->b    = b;
203   mg[levels-1]->x    = x;
204 
205   mg[levels-1]->rtol = rtol;
206   mg[levels-1]->abstol = abstol;
207   mg[levels-1]->dtol = dtol;
208   if (rtol) {
209     /* compute initial residual norm for relative convergence test */
210     PetscReal rnorm;
211     ierr               = (*mg[levels-1]->residual)(mg[levels-1]->A,b,x,w);CHKERRQ(ierr);
212     ierr               = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr);
213     mg[levels-1]->ttol = PetscMax(rtol*rnorm,abstol);
214   } else if (abstol) {
215     mg[levels-1]->ttol = abstol;
216   } else {
217     mg[levels-1]->ttol = 0.0;
218   }
219 
220   while (its-- && !converged) {
221     ierr = MGMCycle_Private(mg+levels-1,&converged);CHKERRQ(ierr);
222   }
223   PetscFunctionReturn(0);
224 }
225 
226 #undef __FUNCT__
227 #define __FUNCT__ "PCSetFromOptions_MG"
228 static PetscErrorCode PCSetFromOptions_MG(PC pc)
229 {
230   PetscErrorCode ierr;
231   PetscInt       indx,m,levels = 1;
232   PetscTruth     flg;
233   const char     *type[] = {"additive","multiplicative","full","cascade","kascade"};
234 
235   PetscFunctionBegin;
236 
237   ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr);
238     if (!pc->data) {
239       ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","MGSetLevels",levels,&levels,&flg);CHKERRQ(ierr);
240       ierr = MGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr);
241     }
242     ierr = PetscOptionsInt("-pc_mg_cycles","1 for V cycle, 2 for W-cycle","MGSetCycles",1,&m,&flg);CHKERRQ(ierr);
243     if (flg) {
244       ierr = MGSetCycles(pc,m);CHKERRQ(ierr);
245     }
246     ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","MGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr);
247     if (flg) {
248       ierr = MGSetNumberSmoothUp(pc,m);CHKERRQ(ierr);
249     }
250     ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","MGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr);
251     if (flg) {
252       ierr = MGSetNumberSmoothDown(pc,m);CHKERRQ(ierr);
253     }
254     ierr = PetscOptionsEList("-pc_mg_type","Multigrid type","MGSetType",type,5,type[1],&indx,&flg);CHKERRQ(ierr);
255     if (flg) {
256       MGType mg = (MGType) 0;
257       switch (indx) {
258       case 0:
259         mg = MGADDITIVE;
260         break;
261       case 1:
262         mg = MGMULTIPLICATIVE;
263         break;
264       case 2:
265         mg = MGFULL;
266         break;
267       case 3:
268         mg = MGKASKADE;
269         break;
270       case 4:
271         mg = MGKASKADE;
272         break;
273       }
274       ierr = MGSetType(pc,mg);CHKERRQ(ierr);
275     }
276     ierr = PetscOptionsName("-pc_mg_log","Log times for each multigrid level","None",&flg);CHKERRQ(ierr);
277     if (flg) {
278       MG   *mg = (MG*)pc->data;
279       int  i;
280       char eventname[128];
281       if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
282       levels = mg[0]->levels;
283       for (i=0; i<levels; i++) {
284         sprintf(eventname,"MSetup Level %d",(int)i);
285         ierr = PetscLogEventRegister(&mg[i]->eventsetup,eventname,pc->cookie);CHKERRQ(ierr);
286         sprintf(eventname,"MGSolve Level %d",(int)i);
287         ierr = PetscLogEventRegister(&mg[i]->eventsolve,eventname,pc->cookie);CHKERRQ(ierr);
288       }
289     }
290   ierr = PetscOptionsTail();CHKERRQ(ierr);
291   PetscFunctionReturn(0);
292 }
293 
294 #undef __FUNCT__
295 #define __FUNCT__ "PCView_MG"
296 static PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
297 {
298   MG             *mg = (MG*)pc->data;
299   PetscErrorCode ierr;
300   PetscInt       levels = mg[0]->levels,i;
301   const char     *cstring;
302   PetscTruth     iascii;
303 
304   PetscFunctionBegin;
305   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
306   if (iascii) {
307     if (mg[0]->am == MGMULTIPLICATIVE) cstring = "multiplicative";
308     else if (mg[0]->am == MGADDITIVE)  cstring = "additive";
309     else if (mg[0]->am == MGFULL)      cstring = "full";
310     else if (mg[0]->am == MGKASKADE)   cstring = "Kaskade";
311     else cstring = "unknown";
312     ierr = PetscViewerASCIIPrintf(viewer,"  MG: type is %s, levels=%D cycles=%D, pre-smooths=%D, post-smooths=%D\n",
313                       cstring,levels,mg[0]->cycles,mg[0]->default_smoothd,mg[0]->default_smoothu);CHKERRQ(ierr);
314     for (i=0; i<levels; i++) {
315       ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
316       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
317       ierr = KSPView(mg[i]->smoothd,viewer);CHKERRQ(ierr);
318       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
319       if (mg[i]->smoothd == mg[i]->smoothu) {
320         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
321       } else {
322         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
323         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
324         ierr = KSPView(mg[i]->smoothu,viewer);CHKERRQ(ierr);
325         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
326       }
327     }
328   } else {
329     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCMG",((PetscObject)viewer)->type_name);
330   }
331   PetscFunctionReturn(0);
332 }
333 
334 /*
335     Calls setup for the KSP on each level
336 */
337 #undef __FUNCT__
338 #define __FUNCT__ "PCSetUp_MG"
339 static PetscErrorCode PCSetUp_MG(PC pc)
340 {
341   MG             *mg = (MG*)pc->data;
342   PetscErrorCode ierr;
343   PetscInt       i,n = mg[0]->levels;
344   PC             cpc;
345   PetscTruth     preonly,lu,redundant,monitor = PETSC_FALSE,dump;
346   PetscViewer    ascii;
347   MPI_Comm       comm;
348 
349   PetscFunctionBegin;
350   if (!pc->setupcalled) {
351     ierr = PetscOptionsHasName(0,"-pc_mg_monitor",&monitor);CHKERRQ(ierr);
352 
353     for (i=1; i<n; i++) {
354       if (mg[i]->smoothd) {
355         if (monitor) {
356           ierr = PetscObjectGetComm((PetscObject)mg[i]->smoothd,&comm);CHKERRQ(ierr);
357           ierr = PetscViewerASCIIOpen(comm,"stdout",&ascii);CHKERRQ(ierr);
358           ierr = PetscViewerASCIISetTab(ascii,n-i);CHKERRQ(ierr);
359           ierr = KSPSetMonitor(mg[i]->smoothd,KSPDefaultMonitor,ascii,(PetscErrorCode(*)(void*))PetscViewerDestroy);CHKERRQ(ierr);
360         }
361         ierr = KSPSetFromOptions(mg[i]->smoothd);CHKERRQ(ierr);
362       }
363     }
364     for (i=0; i<n; i++) {
365       if (mg[i]->smoothu && mg[i]->smoothu != mg[i]->smoothd) {
366         if (monitor) {
367           ierr = PetscObjectGetComm((PetscObject)mg[i]->smoothu,&comm);CHKERRQ(ierr);
368           ierr = PetscViewerASCIIOpen(comm,"stdout",&ascii);CHKERRQ(ierr);
369           ierr = PetscViewerASCIISetTab(ascii,n-i);CHKERRQ(ierr);
370           ierr = KSPSetMonitor(mg[i]->smoothu,KSPDefaultMonitor,ascii,(PetscErrorCode(*)(void*))PetscViewerDestroy);CHKERRQ(ierr);
371         }
372         ierr = KSPSetFromOptions(mg[i]->smoothu);CHKERRQ(ierr);
373       }
374     }
375   }
376 
377   for (i=1; i<n; i++) {
378     if (mg[i]->smoothd) {
379       ierr = KSPSetInitialGuessNonzero(mg[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
380       if (mg[i]->eventsetup) {ierr = PetscLogEventBegin(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);}
381       ierr = KSPSetUp(mg[i]->smoothd);CHKERRQ(ierr);
382       if (mg[i]->eventsetup) {ierr = PetscLogEventEnd(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);}
383     }
384   }
385   for (i=0; i<n; i++) {
386     if (mg[i]->smoothu && mg[i]->smoothu != mg[i]->smoothd) {
387         PC           downpc,uppc;
388         Mat          downmat,downpmat,upmat,uppmat;
389         MatStructure matflag;
390 
391       /* check if operators have been set for up, if not use down operators to set them */
392       ierr = KSPGetPC(mg[i]->smoothu,&uppc);CHKERRQ(ierr);
393       ierr = PCGetOperators(uppc,&upmat,&uppmat,PETSC_NULL);CHKERRQ(ierr);
394       if (!upmat) {
395         ierr = KSPGetPC(mg[i]->smoothd,&downpc);CHKERRQ(ierr);
396         ierr = PCGetOperators(downpc,&downmat,&downpmat,&matflag);CHKERRQ(ierr);
397         ierr = KSPSetOperators(mg[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr);
398       }
399 
400       ierr = KSPSetInitialGuessNonzero(mg[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
401       if (mg[i]->eventsetup) {ierr = PetscLogEventBegin(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);}
402       ierr = KSPSetUp(mg[i]->smoothu);CHKERRQ(ierr);
403       if (mg[i]->eventsetup) {ierr = PetscLogEventEnd(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);}
404     }
405   }
406 
407   /*
408       If coarse solver is not direct method then DO NOT USE preonly
409   */
410   ierr = PetscTypeCompare((PetscObject)mg[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr);
411   if (preonly) {
412     ierr = KSPGetPC(mg[0]->smoothd,&cpc);CHKERRQ(ierr);
413     ierr = PetscTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr);
414     ierr = PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr);
415     if (!lu && !redundant) {
416       ierr = KSPSetType(mg[0]->smoothd,KSPGMRES);CHKERRQ(ierr);
417     }
418   }
419 
420   if (!pc->setupcalled) {
421     if (monitor) {
422       ierr = PetscObjectGetComm((PetscObject)mg[0]->smoothd,&comm);CHKERRQ(ierr);
423       ierr = PetscViewerASCIIOpen(comm,"stdout",&ascii);CHKERRQ(ierr);
424       ierr = PetscViewerASCIISetTab(ascii,n);CHKERRQ(ierr);
425       ierr = KSPSetMonitor(mg[0]->smoothd,KSPDefaultMonitor,ascii,(PetscErrorCode(*)(void*))PetscViewerDestroy);CHKERRQ(ierr);
426     }
427     ierr = KSPSetFromOptions(mg[0]->smoothd);CHKERRQ(ierr);
428   }
429 
430   if (mg[0]->eventsetup) {ierr = PetscLogEventBegin(mg[0]->eventsetup,0,0,0,0);CHKERRQ(ierr);}
431   ierr = KSPSetUp(mg[0]->smoothd);CHKERRQ(ierr);
432   if (mg[0]->eventsetup) {ierr = PetscLogEventEnd(mg[0]->eventsetup,0,0,0,0);CHKERRQ(ierr);}
433 
434   /*
435      Dump the interpolation/restriction matrices to matlab plus the
436    Jacobian/stiffness on each level. This allows Matlab users to
437    easily check if the Galerkin condition A_c = R A_f R^T is satisfied */
438   ierr = PetscOptionsHasName(pc->prefix,"-pc_mg_dump_matlab",&dump);CHKERRQ(ierr);
439   if (dump) {
440     for (i=1; i<n; i++) {
441       ierr = MatView(mg[i]->restrct,PETSC_VIEWER_SOCKET_(pc->comm));CHKERRQ(ierr);
442     }
443     for (i=0; i<n; i++) {
444       ierr = KSPGetPC(mg[i]->smoothd,&pc);CHKERRQ(ierr);
445       ierr = MatView(pc->mat,PETSC_VIEWER_SOCKET_(pc->comm));CHKERRQ(ierr);
446     }
447   }
448 
449   PetscFunctionReturn(0);
450 }
451 
452 /* -------------------------------------------------------------------------------------*/
453 
454 #undef __FUNCT__
455 #define __FUNCT__ "MGSetLevels"
456 /*@C
457    MGSetLevels - Sets the number of levels to use with MG.
458    Must be called before any other MG routine.
459 
460    Collective on PC
461 
462    Input Parameters:
463 +  pc - the preconditioner context
464 .  levels - the number of levels
465 -  comms - optional communicators for each level; this is to allow solving the coarser problems
466            on smaller sets of processors. Use PETSC_NULL_OBJECT for default in Fortran
467 
468    Level: intermediate
469 
470    Notes:
471      If the number of levels is one then the multigrid uses the -mg_levels prefix
472   for setting the level options rather than the -mg_coarse prefix.
473 
474 .keywords: MG, set, levels, multigrid
475 
476 .seealso: MGSetType(), MGGetLevels()
477 @*/
478 PetscErrorCode MGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
479 {
480   PetscErrorCode ierr;
481   MG             *mg;
482 
483   PetscFunctionBegin;
484   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
485 
486   if (pc->data) {
487     SETERRQ(PETSC_ERR_ORDER,"Number levels already set for MG\n\
488     make sure that you call MGSetLevels() before KSPSetFromOptions()");
489   }
490   ierr                     = MGCreate_Private(pc->comm,levels,pc,comms,&mg);CHKERRQ(ierr);
491   mg[0]->am                = MGMULTIPLICATIVE;
492   pc->data                 = (void*)mg;
493   pc->ops->applyrichardson = PCApplyRichardson_MG;
494   PetscFunctionReturn(0);
495 }
496 
497 #undef __FUNCT__
498 #define __FUNCT__ "MGGetLevels"
499 /*@
500    MGGetLevels - Gets the number of levels to use with MG.
501 
502    Not Collective
503 
504    Input Parameter:
505 .  pc - the preconditioner context
506 
507    Output parameter:
508 .  levels - the number of levels
509 
510    Level: advanced
511 
512 .keywords: MG, get, levels, multigrid
513 
514 .seealso: MGSetLevels()
515 @*/
516 PetscErrorCode MGGetLevels(PC pc,PetscInt *levels)
517 {
518   MG  *mg;
519 
520   PetscFunctionBegin;
521   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
522   PetscValidIntPointer(levels,2);
523 
524   mg      = (MG*)pc->data;
525   *levels = mg[0]->levels;
526   PetscFunctionReturn(0);
527 }
528 
529 #undef __FUNCT__
530 #define __FUNCT__ "MGSetType"
531 /*@
532    MGSetType - Determines the form of multigrid to use:
533    multiplicative, additive, full, or the Kaskade algorithm.
534 
535    Collective on PC
536 
537    Input Parameters:
538 +  pc - the preconditioner context
539 -  form - multigrid form, one of MGMULTIPLICATIVE, MGADDITIVE,
540    MGFULL, MGKASKADE
541 
542    Options Database Key:
543 .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
544    additive, full, kaskade
545 
546    Level: advanced
547 
548 .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
549 
550 .seealso: MGSetLevels()
551 @*/
552 PetscErrorCode MGSetType(PC pc,MGType form)
553 {
554   MG *mg;
555 
556   PetscFunctionBegin;
557   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
558   mg = (MG*)pc->data;
559 
560   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
561   mg[0]->am = form;
562   if (form == MGMULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
563   else pc->ops->applyrichardson = 0;
564   PetscFunctionReturn(0);
565 }
566 
567 #undef __FUNCT__
568 #define __FUNCT__ "MGSetCycles"
569 /*@
570    MGSetCycles - Sets the type cycles to use.  Use MGSetCyclesOnLevel() for more
571    complicated cycling.
572 
573    Collective on PC
574 
575    Input Parameters:
576 +  mg - the multigrid context
577 -  n - the number of cycles
578 
579    Options Database Key:
580 $  -pc_mg_cycles n - 1 denotes a V-cycle; 2 denotes a W-cycle.
581 
582    Level: advanced
583 
584 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
585 
586 .seealso: MGSetCyclesOnLevel()
587 @*/
588 PetscErrorCode MGSetCycles(PC pc,PetscInt n)
589 {
590   MG       *mg;
591   PetscInt i,levels;
592 
593   PetscFunctionBegin;
594   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
595   mg     = (MG*)pc->data;
596   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
597   levels = mg[0]->levels;
598 
599   for (i=0; i<levels; i++) {
600     mg[i]->cycles  = n;
601   }
602   PetscFunctionReturn(0);
603 }
604 
605 #undef __FUNCT__
606 #define __FUNCT__ "MGCheck"
607 /*@
608    MGCheck - Checks that all components of the MG structure have
609    been set.
610 
611    Collective on PC
612 
613    Input Parameters:
614 .  mg - the MG structure
615 
616    Level: advanced
617 
618 .keywords: MG, check, set, multigrid
619 @*/
620 PetscErrorCode MGCheck(PC pc)
621 {
622   MG       *mg;
623   PetscInt i,n,count = 0;
624 
625   PetscFunctionBegin;
626   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
627   mg = (MG*)pc->data;
628 
629   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
630 
631   n = mg[0]->levels;
632 
633   for (i=1; i<n; i++) {
634     if (!mg[i]->restrct) {
635       (*PetscErrorPrintf)("No restrict set level %D \n",n-i); count++;
636     }
637     if (!mg[i]->interpolate) {
638       (*PetscErrorPrintf)("No interpolate set level %D \n",n-i); count++;
639     }
640     if (!mg[i]->residual) {
641       (*PetscErrorPrintf)("No residual set level %D \n",n-i); count++;
642     }
643     if (!mg[i]->smoothu) {
644       (*PetscErrorPrintf)("No smoothup set level %D \n",n-i); count++;
645     }
646     if (!mg[i]->smoothd) {
647       (*PetscErrorPrintf)("No smoothdown set level %D \n",n-i); count++;
648     }
649     if (!mg[i]->r) {
650       (*PetscErrorPrintf)("No r set level %D \n",n-i); count++;
651     }
652     if (!mg[i-1]->x) {
653       (*PetscErrorPrintf)("No x set level %D \n",n-i); count++;
654     }
655     if (!mg[i-1]->b) {
656       (*PetscErrorPrintf)("No b set level %D \n",n-i); count++;
657     }
658   }
659   PetscFunctionReturn(count);
660 }
661 
662 
663 #undef __FUNCT__
664 #define __FUNCT__ "MGSetNumberSmoothDown"
665 /*@
666    MGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
667    use on all levels. Use MGGetSmootherDown() to set different
668    pre-smoothing steps on different levels.
669 
670    Collective on PC
671 
672    Input Parameters:
673 +  mg - the multigrid context
674 -  n - the number of smoothing steps
675 
676    Options Database Key:
677 .  -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps
678 
679    Level: advanced
680 
681 .keywords: MG, smooth, down, pre-smoothing, steps, multigrid
682 
683 .seealso: MGSetNumberSmoothUp()
684 @*/
685 PetscErrorCode MGSetNumberSmoothDown(PC pc,PetscInt n)
686 {
687   MG             *mg;
688   PetscErrorCode ierr;
689   PetscInt       i,levels;
690 
691   PetscFunctionBegin;
692   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
693   mg     = (MG*)pc->data;
694   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
695   levels = mg[0]->levels;
696 
697   for (i=0; i<levels; i++) {
698     /* make sure smoother up and down are different */
699     ierr = MGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
700     ierr = KSPSetTolerances(mg[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
701     mg[i]->default_smoothd = n;
702   }
703   PetscFunctionReturn(0);
704 }
705 
706 #undef __FUNCT__
707 #define __FUNCT__ "MGSetNumberSmoothUp"
708 /*@
709    MGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
710    on all levels. Use MGGetSmootherUp() to set different numbers of
711    post-smoothing steps on different levels.
712 
713    Collective on PC
714 
715    Input Parameters:
716 +  mg - the multigrid context
717 -  n - the number of smoothing steps
718 
719    Options Database Key:
720 .  -pc_mg_smoothup <n> - Sets number of post-smoothing steps
721 
722    Level: advanced
723 
724    Note: this does not set a value on the coarsest grid, since we assume that
725     there is no seperate smooth up on the coarsest grid. If you want to have a
726     seperate smooth up on the coarsest grid then call MGGetSmoothUp(pc,0,&ksp);
727 
728 .keywords: MG, smooth, up, post-smoothing, steps, multigrid
729 
730 .seealso: MGSetNumberSmoothDown()
731 @*/
732 PetscErrorCode  MGSetNumberSmoothUp(PC pc,PetscInt n)
733 {
734   MG             *mg;
735   PetscErrorCode ierr;
736   PetscInt       i,levels;
737 
738   PetscFunctionBegin;
739   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
740   mg     = (MG*)pc->data;
741   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
742   levels = mg[0]->levels;
743 
744   for (i=1; i<levels; i++) {
745     /* make sure smoother up and down are different */
746     ierr = MGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
747     ierr = KSPSetTolerances(mg[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
748     mg[i]->default_smoothu = n;
749   }
750   PetscFunctionReturn(0);
751 }
752 
753 /* ----------------------------------------------------------------------------------------*/
754 
755 /*MC
756    PCMG - Use geometric multigrid preconditioning. This preconditioner requires you provide additional
757     information about the coarser grid matrices and restriction/interpolation operators.
758 
759    Options Database Keys:
760 +  -pc_mg_levels <nlevels> - number of levels including finest
761 .  -pc_mg_cycles 1 or 2 - for V or W-cycle
762 .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
763 .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
764 .  -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default
765 .  -pc_mg_log - log information about time spent on each level of the solver
766 .  -pc_mg_monitor - print information on the multigrid convergence
767 -  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
768                         to the Socket viewer for reading from Matlab.
769 
770    Notes:
771 
772    Level: intermediate
773 
774    Concepts: multigrid
775 
776 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
777            MGSetLevels(), MGGetLevels(), MGSetType(), MPSetCycles(), MGSetNumberSmoothDown(),
778            MGSetNumberSmoothUp(), MGGetCoarseSolve(), MGSetResidual(), MGSetInterpolation(),
779            MGSetRestriction(), MGGetSmoother(), MGGetSmootherUp(), MGGetSmootherDown(),
780            MGSetCyclesOnLevel(), MGSetRhs(), MGSetX(), MGSetR()
781 M*/
782 
783 EXTERN_C_BEGIN
784 #undef __FUNCT__
785 #define __FUNCT__ "PCCreate_MG"
786 PetscErrorCode PCCreate_MG(PC pc)
787 {
788   PetscFunctionBegin;
789   pc->ops->apply          = PCApply_MG;
790   pc->ops->setup          = PCSetUp_MG;
791   pc->ops->destroy        = PCDestroy_MG;
792   pc->ops->setfromoptions = PCSetFromOptions_MG;
793   pc->ops->view           = PCView_MG;
794 
795   pc->data                = (void*)0;
796   PetscFunctionReturn(0);
797 }
798 EXTERN_C_END
799