xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision a58c3bc3391eee32bc3fd94ac7edeea38fe57aae)
1 #define PETSCKSP_DLL
2 
3 /*
4     Defines the multigrid preconditioner interface.
5 */
6 #include "src/ksp/pc/impls/mg/mgimpl.h"                    /*I "petscmg.h" I*/
7 
8 
9 #undef __FUNCT__
10 #define __FUNCT__ "PCMGMCycle_Private"
11 PetscErrorCode PCMGMCycle_Private(PC_MG **mglevels,PetscTruth *converged)
12 {
13   PC_MG          *mg = *mglevels,*mgc;
14   PetscErrorCode ierr;
15   PetscInt       cycles = (PetscInt) mg->cycles;
16 
17   PetscFunctionBegin;
18   if (converged) *converged = PETSC_FALSE;
19 
20   if (mg->eventsmoothsolve) {ierr = PetscLogEventBegin(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
21   ierr = KSPSolve(mg->smoothd,mg->b,mg->x);CHKERRQ(ierr);  /* pre-smooth */
22   if (mg->eventsmoothsolve) {ierr = PetscLogEventEnd(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
23   if (mg->level) {  /* not the coarsest grid */
24     if (mg->eventresidual) {ierr = PetscLogEventBegin(mg->eventresidual,0,0,0,0);CHKERRQ(ierr);}
25     ierr = (*mg->residual)(mg->A,mg->b,mg->x,mg->r);CHKERRQ(ierr);
26     if (mg->eventresidual) {ierr = PetscLogEventEnd(mg->eventresidual,0,0,0,0);CHKERRQ(ierr);}
27 
28     /* if on finest level and have convergence criteria set */
29     if (mg->level == mg->levels-1 && mg->ttol) {
30       PetscReal rnorm;
31       ierr = VecNorm(mg->r,NORM_2,&rnorm);CHKERRQ(ierr);
32       if (rnorm <= mg->ttol) {
33         *converged = PETSC_TRUE;
34         if (rnorm < mg->abstol) {
35           ierr = PetscInfo2(0,"Linear solver has converged. Residual norm %G is less than absolute tolerance %G\n",rnorm,mg->abstol);CHKERRQ(ierr);
36         } else {
37           ierr = PetscInfo2(0,"Linear solver has converged. Residual norm %G is less than relative tolerance times initial residual norm %G\n",rnorm,mg->ttol);CHKERRQ(ierr);
38         }
39         PetscFunctionReturn(0);
40       }
41     }
42 
43     mgc = *(mglevels - 1);
44     if (mg->eventinterprestrict) {ierr = PetscLogEventBegin(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
45     ierr = MatRestrict(mg->restrct,mg->r,mgc->b);CHKERRQ(ierr);
46     if (mg->eventinterprestrict) {ierr = PetscLogEventEnd(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
47     ierr = VecSet(mgc->x,0.0);CHKERRQ(ierr);
48     while (cycles--) {
49       ierr = PCMGMCycle_Private(mglevels-1,converged);CHKERRQ(ierr);
50     }
51     if (mg->eventinterprestrict) {ierr = PetscLogEventBegin(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
52     ierr = MatInterpolateAdd(mg->interpolate,mgc->x,mg->x,mg->x);CHKERRQ(ierr);
53     if (mg->eventinterprestrict) {ierr = PetscLogEventEnd(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
54     if (mg->eventsmoothsolve) {ierr = PetscLogEventBegin(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
55     ierr = KSPSolve(mg->smoothu,mg->b,mg->x);CHKERRQ(ierr);    /* post smooth */
56     if (mg->eventsmoothsolve) {ierr = PetscLogEventEnd(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
57   }
58   PetscFunctionReturn(0);
59 }
60 
61 /*
62        PCMGCreate_Private - Creates a PC_MG structure for use with the
63                multigrid code. Level 0 is the coarsest. (But the
64                finest level is stored first in the array).
65 
66 */
67 #undef __FUNCT__
68 #define __FUNCT__ "PCMGCreate_Private"
69 static PetscErrorCode PCMGCreate_Private(MPI_Comm comm,PetscInt levels,PC pc,MPI_Comm *comms,PC_MG ***result)
70 {
71   PC_MG          **mg;
72   PetscErrorCode ierr;
73   PetscInt       i;
74   PetscMPIInt    size;
75   const char     *prefix;
76   PC             ipc;
77 
78   PetscFunctionBegin;
79   ierr = PetscMalloc(levels*sizeof(PC_MG*),&mg);CHKERRQ(ierr);
80   ierr = PetscLogObjectMemory(pc,levels*(sizeof(PC_MG*)+sizeof(PC_MG)));CHKERRQ(ierr);
81 
82   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
83 
84   for (i=0; i<levels; i++) {
85     ierr = PetscNew(PC_MG,&mg[i]);CHKERRQ(ierr);
86     mg[i]->level           = i;
87     mg[i]->levels          = levels;
88     mg[i]->cycles          = PC_MG_CYCLE_V;
89     mg[i]->galerkin        = PETSC_FALSE;
90     mg[i]->galerkinused    = PETSC_FALSE;
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       } else {
110         ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr);
111       }
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     ierr = PetscLogObjectParent(pc,mg[i]->smoothd);CHKERRQ(ierr);
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]->eventsmoothsetup    = 0;
125     mg[i]->eventsmoothsolve    = 0;
126     mg[i]->eventresidual       = 0;
127     mg[i]->eventinterprestrict = 0;
128     mg[i]->cyclesperpcapply    = 1;
129   }
130   *result = mg;
131   PetscFunctionReturn(0);
132 }
133 
134 #undef __FUNCT__
135 #define __FUNCT__ "PCDestroy_MG"
136 static PetscErrorCode PCDestroy_MG(PC pc)
137 {
138   PC_MG          **mg = (PC_MG**)pc->data;
139   PetscErrorCode ierr;
140   PetscInt       i,n = mg[0]->levels;
141 
142   PetscFunctionBegin;
143   for (i=0; i<n-1; i++) {
144     if (mg[i+1]->r) {ierr = VecDestroy(mg[i+1]->r);CHKERRQ(ierr);}
145     if (mg[i]->b) {ierr = VecDestroy(mg[i]->b);CHKERRQ(ierr);}
146     if (mg[i]->x) {ierr = VecDestroy(mg[i]->x);CHKERRQ(ierr);}
147     if (mg[i+1]->restrct) {ierr = MatDestroy(mg[i+1]->restrct);CHKERRQ(ierr);}
148     if (mg[i+1]->interpolate) {ierr = MatDestroy(mg[i+1]->interpolate);CHKERRQ(ierr);}
149   }
150 
151   for (i=0; i<n; i++) {
152     if (mg[i]->smoothd != mg[i]->smoothu) {
153       ierr = KSPDestroy(mg[i]->smoothd);CHKERRQ(ierr);
154     }
155     ierr = KSPDestroy(mg[i]->smoothu);CHKERRQ(ierr);
156     ierr = PetscFree(mg[i]);CHKERRQ(ierr);
157   }
158   ierr = PetscFree(mg);CHKERRQ(ierr);
159   PetscFunctionReturn(0);
160 }
161 
162 
163 
164 EXTERN PetscErrorCode PCMGACycle_Private(PC_MG**);
165 EXTERN PetscErrorCode PCMGFCycle_Private(PC_MG**);
166 EXTERN PetscErrorCode PCMGKCycle_Private(PC_MG**);
167 
168 /*
169    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
170              or full cycle of multigrid.
171 
172   Note:
173   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
174 */
175 #undef __FUNCT__
176 #define __FUNCT__ "PCApply_MG"
177 static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
178 {
179   PC_MG          **mg = (PC_MG**)pc->data;
180   PetscErrorCode ierr;
181   PetscInt       levels = mg[0]->levels,i;
182 
183   PetscFunctionBegin;
184   mg[levels-1]->b = b;
185   mg[levels-1]->x = x;
186   if (!mg[levels-1]->r && mg[0]->am != PC_MG_ADDITIVE && levels > 1) {
187     Vec tvec;
188     ierr = VecDuplicate(mg[levels-1]->b,&tvec);CHKERRQ(ierr);
189     ierr = PCMGSetR(pc,levels-1,tvec);CHKERRQ(ierr);
190     ierr = VecDestroy(tvec);CHKERRQ(ierr);
191   }
192   if (mg[0]->am == PC_MG_MULTIPLICATIVE) {
193     ierr = VecSet(x,0.0);CHKERRQ(ierr);
194     for (i=0; i<mg[0]->cyclesperpcapply; i++) {
195       ierr = PCMGMCycle_Private(mg+levels-1,PETSC_NULL);CHKERRQ(ierr);
196     }
197   }
198   else if (mg[0]->am == PC_MG_ADDITIVE) {
199     ierr = PCMGACycle_Private(mg);CHKERRQ(ierr);
200   }
201   else if (mg[0]->am == PC_MG_KASKADE) {
202     ierr = PCMGKCycle_Private(mg);CHKERRQ(ierr);
203   }
204   else {
205     ierr = PCMGFCycle_Private(mg);CHKERRQ(ierr);
206   }
207   PetscFunctionReturn(0);
208 }
209 
210 #undef __FUNCT__
211 #define __FUNCT__ "PCApplyRichardson_MG"
212 static PetscErrorCode PCApplyRichardson_MG(PC pc,Vec b,Vec x,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its)
213 {
214   PC_MG          **mg = (PC_MG**)pc->data;
215   PetscErrorCode ierr;
216   PetscInt       levels = mg[0]->levels;
217   PetscTruth     converged = PETSC_FALSE;
218 
219   PetscFunctionBegin;
220   mg[levels-1]->b    = b;
221   mg[levels-1]->x    = x;
222 
223   mg[levels-1]->rtol = rtol;
224   mg[levels-1]->abstol = abstol;
225   mg[levels-1]->dtol = dtol;
226   if (rtol) {
227     /* compute initial residual norm for relative convergence test */
228     PetscReal rnorm;
229     ierr               = (*mg[levels-1]->residual)(mg[levels-1]->A,b,x,w);CHKERRQ(ierr);
230     ierr               = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr);
231     mg[levels-1]->ttol = PetscMax(rtol*rnorm,abstol);
232   } else if (abstol) {
233     mg[levels-1]->ttol = abstol;
234   } else {
235     mg[levels-1]->ttol = 0.0;
236   }
237 
238   while (its-- && !converged) {
239     ierr = PCMGMCycle_Private(mg+levels-1,&converged);CHKERRQ(ierr);
240   }
241   PetscFunctionReturn(0);
242 }
243 
244 #undef __FUNCT__
245 #define __FUNCT__ "PCSetFromOptions_MG"
246 PetscErrorCode PCSetFromOptions_MG(PC pc)
247 {
248   PetscErrorCode ierr;
249   PetscInt       m,levels = 1,cycles;
250   PetscTruth     flg;
251   PC_MG          **mg = (PC_MG**)pc->data;
252   PCMGType       mgtype = PC_MG_ADDITIVE;
253   PCMGCycleType  mgctype;
254 
255   PetscFunctionBegin;
256   ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr);
257     if (!pc->data) {
258       ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr);
259       ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr);
260       mg = (PC_MG**)pc->data;
261     }
262     mgctype = (PCMGCycleType) mg[0]->cycles;
263     ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr);
264     if (flg) {
265       ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr);
266     };
267     ierr = PetscOptionsName("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",&flg);CHKERRQ(ierr);
268     if (flg) {
269       ierr = PCMGSetGalerkin(pc);CHKERRQ(ierr);
270     }
271     ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr);
272     if (flg) {
273       ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr);
274     }
275     ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr);
276     if (flg) {
277       ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr);
278     }
279     ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
280     if (flg) {
281       ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
282     }
283     if (mg[0]->am == PC_MG_MULTIPLICATIVE) {
284       ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGSetLevels",mg[0]->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr);
285       if (flg) {
286 	ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr);
287       }
288     }
289     ierr = PetscOptionsName("-pc_mg_log","Log times for each multigrid level","None",&flg);CHKERRQ(ierr);
290     if (flg) {
291       PetscInt i;
292       char     eventname[128];
293       if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
294       levels = mg[0]->levels;
295       for (i=0; i<levels; i++) {
296         sprintf(eventname,"MGSetup Level %d",(int)i);
297         ierr = PetscLogEventRegister(&mg[i]->eventsmoothsetup,eventname,pc->cookie);CHKERRQ(ierr);
298         sprintf(eventname,"MGSmooth Level %d",(int)i);
299         ierr = PetscLogEventRegister(&mg[i]->eventsmoothsolve,eventname,pc->cookie);CHKERRQ(ierr);
300         if (i) {
301           sprintf(eventname,"MGResid Level %d",(int)i);
302           ierr = PetscLogEventRegister(&mg[i]->eventresidual,eventname,pc->cookie);CHKERRQ(ierr);
303           sprintf(eventname,"MGInterp Level %d",(int)i);
304           ierr = PetscLogEventRegister(&mg[i]->eventinterprestrict,eventname,pc->cookie);CHKERRQ(ierr);
305         }
306       }
307     }
308   ierr = PetscOptionsTail();CHKERRQ(ierr);
309   PetscFunctionReturn(0);
310 }
311 
312 const char *PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0};
313 const char *PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0};
314 
315 #undef __FUNCT__
316 #define __FUNCT__ "PCView_MG"
317 static PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
318 {
319   PC_MG          **mg = (PC_MG**)pc->data;
320   PetscErrorCode ierr;
321   PetscInt       levels = mg[0]->levels,i;
322   PetscTruth     iascii;
323 
324   PetscFunctionBegin;
325   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
326   if (iascii) {
327     ierr = PetscViewerASCIIPrintf(viewer,"  MG: type is %s, levels=%D cycles=%s, pre-smooths=%D, post-smooths=%D\n",
328 				  PCMGTypes[mg[0]->am],levels,(mg[0]->cycles == PC_MG_CYCLE_V) ? "v" : "w",
329                                   mg[0]->default_smoothd,mg[0]->default_smoothu);CHKERRQ(ierr);
330     if (mg[0]->galerkin) {
331       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
332     }
333     for (i=0; i<levels; i++) {
334       if (!i) {
335         ierr = PetscViewerASCIIPrintf(viewer,"Coarse gride solver -- level %D -------------------------------\n",i);CHKERRQ(ierr);
336       } else {
337         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
338       }
339       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
340       ierr = KSPView(mg[i]->smoothd,viewer);CHKERRQ(ierr);
341       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
342       if (i && mg[i]->smoothd == mg[i]->smoothu) {
343         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
344       } else if (i){
345         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
346         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
347         ierr = KSPView(mg[i]->smoothu,viewer);CHKERRQ(ierr);
348         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
349       }
350     }
351   } else {
352     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCMG",((PetscObject)viewer)->type_name);
353   }
354   PetscFunctionReturn(0);
355 }
356 
357 /*
358     Calls setup for the KSP on each level
359 */
360 #undef __FUNCT__
361 #define __FUNCT__ "PCSetUp_MG"
362 static PetscErrorCode PCSetUp_MG(PC pc)
363 {
364   PC_MG                   **mg = (PC_MG**)pc->data;
365   PetscErrorCode          ierr;
366   PetscInt                i,n = mg[0]->levels;
367   PC                      cpc,mpc;
368   PetscTruth              preonly,lu,redundant,cholesky,monitor = PETSC_FALSE,dump,opsset;
369   PetscViewerASCIIMonitor ascii;
370   PetscViewer             viewer = PETSC_NULL;
371   MPI_Comm                comm;
372   Mat                     dA,dB;
373   MatStructure            uflag;
374   Vec                     tvec;
375 
376   PetscFunctionBegin;
377 
378   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
379   /* so use those from global PC */
380   /* Is this what we always want? What if user wants to keep old one? */
381   ierr = KSPGetOperatorsSet(mg[n-1]->smoothd,PETSC_NULL,&opsset);CHKERRQ(ierr);
382   ierr = KSPGetPC(mg[0]->smoothd,&cpc);CHKERRQ(ierr);
383   ierr = KSPGetPC(mg[n-1]->smoothd,&mpc);CHKERRQ(ierr);
384   if (!opsset || ((cpc->setupcalled == 1) && (mpc->setupcalled == 2))) {
385     ierr = PetscInfo(pc,"Using outer operators to define finest grid operator \n  because PCMGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n");CHKERRQ(ierr);
386     ierr = KSPSetOperators(mg[n-1]->smoothd,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr);
387   }
388 
389   if (mg[0]->galerkin) {
390     Mat B;
391     mg[0]->galerkinused = PETSC_TRUE;
392     /* currently only handle case where mat and pmat are the same on coarser levels */
393     ierr = KSPGetOperators(mg[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr);
394     if (!pc->setupcalled) {
395       for (i=n-2; i>-1; i--) {
396         ierr = MatPtAP(dB,mg[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
397         ierr = KSPSetOperators(mg[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
398 	if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
399         dB   = B;
400       }
401       ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);
402     } else {
403       for (i=n-2; i>-1; i--) {
404         ierr = KSPGetOperators(mg[i]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
405         ierr = MatPtAP(dB,mg[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
406         ierr = KSPSetOperators(mg[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
407         dB   = B;
408       }
409     }
410   }
411 
412   if (!pc->setupcalled) {
413     ierr = PetscOptionsHasName(0,"-pc_mg_monitor",&monitor);CHKERRQ(ierr);
414 
415     for (i=0; i<n; i++) {
416       if (monitor) {
417         ierr = PetscObjectGetComm((PetscObject)mg[i]->smoothd,&comm);CHKERRQ(ierr);
418         ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr);
419         ierr = KSPMonitorSet(mg[i]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
420       }
421       ierr = KSPSetFromOptions(mg[i]->smoothd);CHKERRQ(ierr);
422     }
423     for (i=1; i<n; i++) {
424       if (mg[i]->smoothu && (mg[i]->smoothu != mg[i]->smoothd)) {
425         if (monitor) {
426           ierr = PetscObjectGetComm((PetscObject)mg[i]->smoothu,&comm);CHKERRQ(ierr);
427           ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr);
428           ierr = KSPMonitorSet(mg[i]->smoothu,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
429         }
430         ierr = KSPSetFromOptions(mg[i]->smoothu);CHKERRQ(ierr);
431       }
432     }
433     for (i=1; i<n; i++) {
434       if (!mg[i]->residual) {
435         Mat mat;
436         ierr = KSPGetOperators(mg[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr);
437         ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr);
438       }
439       if (mg[i]->restrct && !mg[i]->interpolate) {
440         ierr = PCMGSetInterpolation(pc,i,mg[i]->restrct);CHKERRQ(ierr);
441       }
442       if (!mg[i]->restrct && mg[i]->interpolate) {
443         ierr = PCMGSetRestriction(pc,i,mg[i]->interpolate);CHKERRQ(ierr);
444       }
445 #if defined(PETSC_USE_DEBUG)
446       if (!mg[i]->restrct || !mg[i]->interpolate) {
447         SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Need to set restriction or interpolation on level %d",(int)i);
448       }
449 #endif
450     }
451     for (i=0; i<n-1; i++) {
452       if (!mg[i]->b) {
453         Vec *vec;
454         ierr = KSPGetVecs(mg[i]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
455         ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr);
456         ierr = PetscFree(vec);CHKERRQ(ierr);
457       }
458       if (!mg[i]->r && i) {
459         ierr = VecDuplicate(mg[i]->b,&tvec);CHKERRQ(ierr);
460         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
461         ierr = VecDestroy(tvec);CHKERRQ(ierr);
462       }
463       if (!mg[i]->x) {
464         ierr = VecDuplicate(mg[i]->b,&tvec);CHKERRQ(ierr);
465         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
466         ierr = VecDestroy(tvec);CHKERRQ(ierr);
467       }
468     }
469   }
470 
471 
472   for (i=1; i<n; i++) {
473     if (mg[i]->smoothu == mg[i]->smoothd) {
474       /* if doing only down then initial guess is zero */
475       ierr = KSPSetInitialGuessNonzero(mg[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
476     }
477     if (mg[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mg[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
478     ierr = KSPSetUp(mg[i]->smoothd);CHKERRQ(ierr);
479     if (mg[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mg[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
480   }
481   for (i=1; i<n; i++) {
482     if (mg[i]->smoothu && mg[i]->smoothu != mg[i]->smoothd) {
483       Mat          downmat,downpmat;
484       MatStructure matflag;
485       PetscTruth   opsset;
486 
487       /* check if operators have been set for up, if not use down operators to set them */
488       ierr = KSPGetOperatorsSet(mg[i]->smoothu,&opsset,PETSC_NULL);CHKERRQ(ierr);
489       if (!opsset) {
490         ierr = KSPGetOperators(mg[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr);
491         ierr = KSPSetOperators(mg[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr);
492       }
493 
494       ierr = KSPSetInitialGuessNonzero(mg[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
495       if (mg[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mg[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
496       ierr = KSPSetUp(mg[i]->smoothu);CHKERRQ(ierr);
497       if (mg[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mg[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
498     }
499   }
500 
501   /*
502       If coarse solver is not direct method then DO NOT USE preonly
503   */
504   ierr = PetscTypeCompare((PetscObject)mg[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr);
505   if (preonly) {
506     ierr = PetscTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr);
507     ierr = PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr);
508     ierr = PetscTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr);
509     if (!lu && !redundant && !cholesky) {
510       ierr = KSPSetType(mg[0]->smoothd,KSPGMRES);CHKERRQ(ierr);
511     }
512   }
513 
514   if (!pc->setupcalled) {
515     if (monitor) {
516       ierr = PetscObjectGetComm((PetscObject)mg[0]->smoothd,&comm);CHKERRQ(ierr);
517       ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n,&ascii);CHKERRQ(ierr);
518       ierr = KSPMonitorSet(mg[0]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
519     }
520     ierr = KSPSetFromOptions(mg[0]->smoothd);CHKERRQ(ierr);
521   }
522 
523   if (mg[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mg[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
524   ierr = KSPSetUp(mg[0]->smoothd);CHKERRQ(ierr);
525   if (mg[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mg[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
526 
527   /*
528      Dump the interpolation/restriction matrices plus the
529    Jacobian/stiffness on each level. This allows Matlab users to
530    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
531 
532    Only support one or the other at the same time.
533   */
534 #if defined(PETSC_USE_SOCKET_VIEWER)
535   ierr = PetscOptionsHasName(pc->prefix,"-pc_mg_dump_matlab",&dump);CHKERRQ(ierr);
536   if (dump) {
537     viewer = PETSC_VIEWER_SOCKET_(pc->comm);
538   }
539 #endif
540   ierr = PetscOptionsHasName(pc->prefix,"-pc_mg_dump_binary",&dump);CHKERRQ(ierr);
541   if (dump) {
542     viewer = PETSC_VIEWER_BINARY_(pc->comm);
543   }
544 
545   if (viewer) {
546     for (i=1; i<n; i++) {
547       ierr = MatView(mg[i]->restrct,viewer);CHKERRQ(ierr);
548     }
549     for (i=0; i<n; i++) {
550       ierr = KSPGetPC(mg[i]->smoothd,&pc);CHKERRQ(ierr);
551       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
552     }
553   }
554   PetscFunctionReturn(0);
555 }
556 
557 /* -------------------------------------------------------------------------------------*/
558 
559 #undef __FUNCT__
560 #define __FUNCT__ "PCMGSetLevels"
561 /*@C
562    PCMGSetLevels - Sets the number of levels to use with MG.
563    Must be called before any other MG routine.
564 
565    Collective on PC
566 
567    Input Parameters:
568 +  pc - the preconditioner context
569 .  levels - the number of levels
570 -  comms - optional communicators for each level; this is to allow solving the coarser problems
571            on smaller sets of processors. Use PETSC_NULL_OBJECT for default in Fortran
572 
573    Level: intermediate
574 
575    Notes:
576      If the number of levels is one then the multigrid uses the -mg_levels prefix
577   for setting the level options rather than the -mg_coarse prefix.
578 
579 .keywords: MG, set, levels, multigrid
580 
581 .seealso: PCMGSetType(), PCMGGetLevels()
582 @*/
583 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
584 {
585   PetscErrorCode ierr;
586   PC_MG          **mg=0;
587 
588   PetscFunctionBegin;
589   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
590 
591   if (pc->data) {
592     SETERRQ(PETSC_ERR_ORDER,"Number levels already set for MG\n\
593     make sure that you call PCMGSetLevels() before KSPSetFromOptions()");
594   }
595   ierr                     = PCMGCreate_Private(pc->comm,levels,pc,comms,&mg);CHKERRQ(ierr);
596   mg[0]->am                = PC_MG_MULTIPLICATIVE;
597   pc->data                 = (void*)mg;
598   pc->ops->applyrichardson = PCApplyRichardson_MG;
599   PetscFunctionReturn(0);
600 }
601 
602 #undef __FUNCT__
603 #define __FUNCT__ "PCMGGetLevels"
604 /*@
605    PCMGGetLevels - Gets the number of levels to use with MG.
606 
607    Not Collective
608 
609    Input Parameter:
610 .  pc - the preconditioner context
611 
612    Output parameter:
613 .  levels - the number of levels
614 
615    Level: advanced
616 
617 .keywords: MG, get, levels, multigrid
618 
619 .seealso: PCMGSetLevels()
620 @*/
621 PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetLevels(PC pc,PetscInt *levels)
622 {
623   PC_MG  **mg;
624 
625   PetscFunctionBegin;
626   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
627   PetscValidIntPointer(levels,2);
628 
629   mg      = (PC_MG**)pc->data;
630   *levels = mg[0]->levels;
631   PetscFunctionReturn(0);
632 }
633 
634 #undef __FUNCT__
635 #define __FUNCT__ "PCMGSetType"
636 /*@
637    PCMGSetType - Determines the form of multigrid to use:
638    multiplicative, additive, full, or the Kaskade algorithm.
639 
640    Collective on PC
641 
642    Input Parameters:
643 +  pc - the preconditioner context
644 -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
645    PC_MG_FULL, PC_MG_KASKADE
646 
647    Options Database Key:
648 .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
649    additive, full, kaskade
650 
651    Level: advanced
652 
653 .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
654 
655 .seealso: PCMGSetLevels()
656 @*/
657 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetType(PC pc,PCMGType form)
658 {
659   PC_MG **mg;
660 
661   PetscFunctionBegin;
662   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
663   mg = (PC_MG**)pc->data;
664 
665   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
666   mg[0]->am = form;
667   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
668   else pc->ops->applyrichardson = 0;
669   PetscFunctionReturn(0);
670 }
671 
672 #undef __FUNCT__
673 #define __FUNCT__ "PCMGSetCycleType"
674 /*@
675    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
676    complicated cycling.
677 
678    Collective on PC
679 
680    Input Parameters:
681 +  pc - the multigrid context
682 -  PC_MG_CYCLE_V or PC_MG_CYCLE_W
683 
684    Options Database Key:
685 $  -pc_mg_cycle_type v or w
686 
687    Level: advanced
688 
689 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
690 
691 .seealso: PCMGSetCycleTypeOnLevel()
692 @*/
693 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetCycleType(PC pc,PCMGCycleType n)
694 {
695   PC_MG    **mg;
696   PetscInt i,levels;
697 
698   PetscFunctionBegin;
699   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
700   mg     = (PC_MG**)pc->data;
701   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
702   levels = mg[0]->levels;
703 
704   for (i=0; i<levels; i++) {
705     mg[i]->cycles  = n;
706   }
707   PetscFunctionReturn(0);
708 }
709 
710 #undef __FUNCT__
711 #define __FUNCT__ "PCMGMultiplicativeSetCycles"
712 /*@
713    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
714          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
715 
716    Collective on PC
717 
718    Input Parameters:
719 +  pc - the multigrid context
720 -  n - number of cycles (default is 1)
721 
722    Options Database Key:
723 $  -pc_mg_multiplicative_cycles n
724 
725    Level: advanced
726 
727    Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
728 
729 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
730 
731 .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
732 @*/
733 PetscErrorCode PETSCKSP_DLLEXPORT PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
734 {
735   PC_MG    **mg;
736   PetscInt i,levels;
737 
738   PetscFunctionBegin;
739   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
740   mg     = (PC_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=0; i<levels; i++) {
745     mg[i]->cyclesperpcapply  = n;
746   }
747   PetscFunctionReturn(0);
748 }
749 
750 #undef __FUNCT__
751 #define __FUNCT__ "PCMGSetGalerkin"
752 /*@
753    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
754       finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t
755 
756    Collective on PC
757 
758    Input Parameters:
759 .  pc - the multigrid context
760 
761    Options Database Key:
762 $  -pc_mg_galerkin
763 
764    Level: intermediate
765 
766 .keywords: MG, set, Galerkin
767 
768 .seealso: PCMGGetGalerkin()
769 
770 @*/
771 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetGalerkin(PC pc)
772 {
773   PC_MG    **mg;
774   PetscInt i,levels;
775 
776   PetscFunctionBegin;
777   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
778   mg     = (PC_MG**)pc->data;
779   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
780   levels = mg[0]->levels;
781 
782   for (i=0; i<levels; i++) {
783     mg[i]->galerkin = PETSC_TRUE;
784   }
785   PetscFunctionReturn(0);
786 }
787 
788 #undef __FUNCT__
789 #define __FUNCT__ "PCMGGetGalerkin"
790 /*@
791    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
792       A_i-1 = r_i * A_i * r_i^t
793 
794    Not Collective
795 
796    Input Parameter:
797 .  pc - the multigrid context
798 
799    Output Parameter:
800 .  gelerkin - PETSC_TRUE or PETSC_FALSE
801 
802    Options Database Key:
803 $  -pc_mg_galerkin
804 
805    Level: intermediate
806 
807 .keywords: MG, set, Galerkin
808 
809 .seealso: PCMGSetGalerkin()
810 
811 @*/
812 PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetGalerkin(PC pc,PetscTruth *galerkin)
813 {
814   PC_MG    **mg;
815 
816   PetscFunctionBegin;
817   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
818   mg     = (PC_MG**)pc->data;
819   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
820   *galerkin = mg[0]->galerkin;
821   PetscFunctionReturn(0);
822 }
823 
824 #undef __FUNCT__
825 #define __FUNCT__ "PCMGSetNumberSmoothDown"
826 /*@
827    PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
828    use on all levels. Use PCMGGetSmootherDown() to set different
829    pre-smoothing steps on different levels.
830 
831    Collective on PC
832 
833    Input Parameters:
834 +  mg - the multigrid context
835 -  n - the number of smoothing steps
836 
837    Options Database Key:
838 .  -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps
839 
840    Level: advanced
841 
842 .keywords: MG, smooth, down, pre-smoothing, steps, multigrid
843 
844 .seealso: PCMGSetNumberSmoothUp()
845 @*/
846 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothDown(PC pc,PetscInt n)
847 {
848   PC_MG          **mg;
849   PetscErrorCode ierr;
850   PetscInt       i,levels;
851 
852   PetscFunctionBegin;
853   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
854   mg     = (PC_MG**)pc->data;
855   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
856   levels = mg[0]->levels;
857 
858   for (i=1; i<levels; i++) {
859     /* make sure smoother up and down are different */
860     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
861     ierr = KSPSetTolerances(mg[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
862     mg[i]->default_smoothd = n;
863   }
864   PetscFunctionReturn(0);
865 }
866 
867 #undef __FUNCT__
868 #define __FUNCT__ "PCMGSetNumberSmoothUp"
869 /*@
870    PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
871    on all levels. Use PCMGGetSmootherUp() to set different numbers of
872    post-smoothing steps on different levels.
873 
874    Collective on PC
875 
876    Input Parameters:
877 +  mg - the multigrid context
878 -  n - the number of smoothing steps
879 
880    Options Database Key:
881 .  -pc_mg_smoothup <n> - Sets number of post-smoothing steps
882 
883    Level: advanced
884 
885    Note: this does not set a value on the coarsest grid, since we assume that
886     there is no separate smooth up on the coarsest grid.
887 
888 .keywords: MG, smooth, up, post-smoothing, steps, multigrid
889 
890 .seealso: PCMGSetNumberSmoothDown()
891 @*/
892 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothUp(PC pc,PetscInt n)
893 {
894   PC_MG          **mg;
895   PetscErrorCode ierr;
896   PetscInt       i,levels;
897 
898   PetscFunctionBegin;
899   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
900   mg     = (PC_MG**)pc->data;
901   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
902   levels = mg[0]->levels;
903 
904   for (i=1; i<levels; i++) {
905     /* make sure smoother up and down are different */
906     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
907     ierr = KSPSetTolerances(mg[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
908     mg[i]->default_smoothu = n;
909   }
910   PetscFunctionReturn(0);
911 }
912 
913 /* ----------------------------------------------------------------------------------------*/
914 
915 /*MC
916    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
917     information about the coarser grid matrices and restriction/interpolation operators.
918 
919    Options Database Keys:
920 +  -pc_mg_levels <nlevels> - number of levels including finest
921 .  -pc_mg_cycles v or w
922 .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
923 .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
924 .  -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default
925 .  -pc_mg_log - log information about time spent on each level of the solver
926 .  -pc_mg_monitor - print information on the multigrid convergence
927 .  -pc_mg_galerkin - use Galerkin process to compute coarser operators
928 -  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
929                         to the Socket viewer for reading from Matlab.
930 
931    Notes:
932 
933    Level: intermediate
934 
935    Concepts: multigrid
936 
937 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
938            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
939            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
940            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
941            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
942 M*/
943 
944 EXTERN_C_BEGIN
945 #undef __FUNCT__
946 #define __FUNCT__ "PCCreate_MG"
947 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_MG(PC pc)
948 {
949   PetscFunctionBegin;
950   pc->ops->apply          = PCApply_MG;
951   pc->ops->setup          = PCSetUp_MG;
952   pc->ops->destroy        = PCDestroy_MG;
953   pc->ops->setfromoptions = PCSetFromOptions_MG;
954   pc->ops->view           = PCView_MG;
955 
956   pc->data                = (void*)0;
957   PetscFunctionReturn(0);
958 }
959 EXTERN_C_END
960