xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision cb497662fae1ebe98e095979598fa3ebb2f149fa)
1 
2 /*
3     Defines the multigrid preconditioner interface.
4 */
5 #include <petsc/private/pcmgimpl.h>                    /*I "petscksp.h" I*/
6 #include <petscdm.h>
7 PETSC_INTERN PetscErrorCode PCPreSolveChangeRHS(PC,PetscBool*);
8 
9 /*
10    Contains the list of registered coarse space construction routines
11 */
12 PetscFunctionList PCMGCoarseList = NULL;
13 
14 PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG_Levels **mglevelsin,PCRichardsonConvergedReason *reason)
15 {
16   PC_MG          *mg = (PC_MG*)pc->data;
17   PC_MG_Levels   *mgc,*mglevels = *mglevelsin;
18   PetscErrorCode ierr;
19   PetscInt       cycles = (mglevels->level == 1) ? 1 : (PetscInt) mglevels->cycles;
20 
21   PetscFunctionBegin;
22   if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
23   ierr = KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x);CHKERRQ(ierr);  /* pre-smooth */
24   ierr = KSPCheckSolve(mglevels->smoothd,pc,mglevels->x);CHKERRQ(ierr);
25   if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
26   if (mglevels->level) {  /* not the coarsest grid */
27     if (mglevels->eventresidual) {ierr = PetscLogEventBegin(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);}
28     ierr = (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);CHKERRQ(ierr);
29     if (mglevels->eventresidual) {ierr = PetscLogEventEnd(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);}
30 
31     /* if on finest level and have convergence criteria set */
32     if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) {
33       PetscReal rnorm;
34       ierr = VecNorm(mglevels->r,NORM_2,&rnorm);CHKERRQ(ierr);
35       if (rnorm <= mg->ttol) {
36         if (rnorm < mg->abstol) {
37           *reason = PCRICHARDSON_CONVERGED_ATOL;
38           ierr    = PetscInfo2(pc,"Linear solver has converged. Residual norm %g is less than absolute tolerance %g\n",(double)rnorm,(double)mg->abstol);CHKERRQ(ierr);
39         } else {
40           *reason = PCRICHARDSON_CONVERGED_RTOL;
41           ierr    = PetscInfo2(pc,"Linear solver has converged. Residual norm %g is less than relative tolerance times initial residual norm %g\n",(double)rnorm,(double)mg->ttol);CHKERRQ(ierr);
42         }
43         PetscFunctionReturn(0);
44       }
45     }
46 
47     mgc = *(mglevelsin - 1);
48     if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
49     ierr = MatRestrict(mglevels->restrct,mglevels->r,mgc->b);CHKERRQ(ierr);
50     if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
51     ierr = VecSet(mgc->x,0.0);CHKERRQ(ierr);
52     while (cycles--) {
53       ierr = PCMGMCycle_Private(pc,mglevelsin-1,reason);CHKERRQ(ierr);
54     }
55     if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
56     ierr = MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x);CHKERRQ(ierr);
57     if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
58     if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
59     ierr = KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x);CHKERRQ(ierr);    /* post smooth */
60     ierr = KSPCheckSolve(mglevels->smoothu,pc,mglevels->x);CHKERRQ(ierr);
61     if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
62   }
63   PetscFunctionReturn(0);
64 }
65 
66 static PetscErrorCode PCApplyRichardson_MG(PC pc,Vec b,Vec x,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscBool zeroguess,PetscInt *outits,PCRichardsonConvergedReason *reason)
67 {
68   PC_MG          *mg        = (PC_MG*)pc->data;
69   PC_MG_Levels   **mglevels = mg->levels;
70   PetscErrorCode ierr;
71   PC             tpc;
72   PetscBool      changeu,changed;
73   PetscInt       levels = mglevels[0]->levels,i;
74 
75   PetscFunctionBegin;
76   /* When the DM is supplying the matrix then it will not exist until here */
77   for (i=0; i<levels; i++) {
78     if (!mglevels[i]->A) {
79       ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr);
80       ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr);
81     }
82   }
83 
84   ierr = KSPGetPC(mglevels[levels-1]->smoothd,&tpc);CHKERRQ(ierr);
85   ierr = PCPreSolveChangeRHS(tpc,&changed);CHKERRQ(ierr);
86   ierr = KSPGetPC(mglevels[levels-1]->smoothu,&tpc);CHKERRQ(ierr);
87   ierr = PCPreSolveChangeRHS(tpc,&changeu);CHKERRQ(ierr);
88   if (!changed && !changeu) {
89     ierr = VecDestroy(&mglevels[levels-1]->b);CHKERRQ(ierr);
90     mglevels[levels-1]->b = b;
91   } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */
92     if (!mglevels[levels-1]->b) {
93       Vec *vec;
94 
95       ierr = KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
96       mglevels[levels-1]->b = *vec;
97       ierr = PetscFree(vec);CHKERRQ(ierr);
98     }
99     ierr = VecCopy(b,mglevels[levels-1]->b);CHKERRQ(ierr);
100   }
101   mglevels[levels-1]->x = x;
102 
103   mg->rtol   = rtol;
104   mg->abstol = abstol;
105   mg->dtol   = dtol;
106   if (rtol) {
107     /* compute initial residual norm for relative convergence test */
108     PetscReal rnorm;
109     if (zeroguess) {
110       ierr = VecNorm(b,NORM_2,&rnorm);CHKERRQ(ierr);
111     } else {
112       ierr = (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);CHKERRQ(ierr);
113       ierr = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr);
114     }
115     mg->ttol = PetscMax(rtol*rnorm,abstol);
116   } else if (abstol) mg->ttol = abstol;
117   else mg->ttol = 0.0;
118 
119   /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't
120      stop prematurely due to small residual */
121   for (i=1; i<levels; i++) {
122     ierr = KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
123     if (mglevels[i]->smoothu != mglevels[i]->smoothd) {
124       /* For Richardson the initial guess is nonzero since it is solving in each cycle the original system not just applying as a preconditioner */
125       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
126       ierr = KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
127     }
128   }
129 
130   *reason = (PCRichardsonConvergedReason)0;
131   for (i=0; i<its; i++) {
132     ierr = PCMGMCycle_Private(pc,mglevels+levels-1,reason);CHKERRQ(ierr);
133     if (*reason) break;
134   }
135   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
136   *outits = i;
137   if (!changed && !changeu) mglevels[levels-1]->b = NULL;
138   PetscFunctionReturn(0);
139 }
140 
141 PetscErrorCode PCReset_MG(PC pc)
142 {
143   PC_MG          *mg        = (PC_MG*)pc->data;
144   PC_MG_Levels   **mglevels = mg->levels;
145   PetscErrorCode ierr;
146   PetscInt       i,c,n;
147 
148   PetscFunctionBegin;
149   if (mglevels) {
150     n = mglevels[0]->levels;
151     for (i=0; i<n-1; i++) {
152       ierr = VecDestroy(&mglevels[i+1]->r);CHKERRQ(ierr);
153       ierr = VecDestroy(&mglevels[i]->b);CHKERRQ(ierr);
154       ierr = VecDestroy(&mglevels[i]->x);CHKERRQ(ierr);
155       ierr = MatDestroy(&mglevels[i+1]->restrct);CHKERRQ(ierr);
156       ierr = MatDestroy(&mglevels[i+1]->interpolate);CHKERRQ(ierr);
157       ierr = MatDestroy(&mglevels[i+1]->inject);CHKERRQ(ierr);
158       ierr = VecDestroy(&mglevels[i+1]->rscale);CHKERRQ(ierr);
159     }
160     /* this is not null only if the smoother on the finest level
161        changes the rhs during PreSolve */
162     ierr = VecDestroy(&mglevels[n-1]->b);CHKERRQ(ierr);
163 
164     for (i=0; i<n; i++) {
165       if (mglevels[i]->coarseSpace) for (c = 0; c < mg->Nc; ++c) {ierr = VecDestroy(&mglevels[i]->coarseSpace[c]);CHKERRQ(ierr);}
166       ierr = PetscFree(mglevels[i]->coarseSpace);CHKERRQ(ierr);
167       mglevels[i]->coarseSpace = NULL;
168       ierr = MatDestroy(&mglevels[i]->A);CHKERRQ(ierr);
169       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
170         ierr = KSPReset(mglevels[i]->smoothd);CHKERRQ(ierr);
171       }
172       ierr = KSPReset(mglevels[i]->smoothu);CHKERRQ(ierr);
173     }
174     mg->Nc = 0;
175   }
176   PetscFunctionReturn(0);
177 }
178 
179 PetscErrorCode PCMGSetLevels_MG(PC pc,PetscInt levels,MPI_Comm *comms)
180 {
181   PetscErrorCode ierr;
182   PC_MG          *mg        = (PC_MG*)pc->data;
183   MPI_Comm       comm;
184   PC_MG_Levels   **mglevels = mg->levels;
185   PCMGType       mgtype     = mg->am;
186   PetscInt       mgctype    = (PetscInt) PC_MG_CYCLE_V;
187   PetscInt       i;
188   PetscMPIInt    size;
189   const char     *prefix;
190   PC             ipc;
191   PetscInt       n;
192 
193   PetscFunctionBegin;
194   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
195   PetscValidLogicalCollectiveInt(pc,levels,2);
196   if (mg->nlevels == levels) PetscFunctionReturn(0);
197   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
198   if (mglevels) {
199     mgctype = mglevels[0]->cycles;
200     /* changing the number of levels so free up the previous stuff */
201     ierr = PCReset_MG(pc);CHKERRQ(ierr);
202     n    = mglevels[0]->levels;
203     for (i=0; i<n; i++) {
204       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
205         ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
206       }
207       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
208       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
209     }
210     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
211   }
212 
213   mg->nlevels = levels;
214 
215   ierr = PetscMalloc1(levels,&mglevels);CHKERRQ(ierr);
216   ierr = PetscLogObjectMemory((PetscObject)pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr);
217 
218   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
219 
220   mg->stageApply = 0;
221   for (i=0; i<levels; i++) {
222     ierr = PetscNewLog(pc,&mglevels[i]);CHKERRQ(ierr);
223 
224     mglevels[i]->level               = i;
225     mglevels[i]->levels              = levels;
226     mglevels[i]->cycles              = mgctype;
227     mg->default_smoothu              = 2;
228     mg->default_smoothd              = 2;
229     mglevels[i]->eventsmoothsetup    = 0;
230     mglevels[i]->eventsmoothsolve    = 0;
231     mglevels[i]->eventresidual       = 0;
232     mglevels[i]->eventinterprestrict = 0;
233 
234     if (comms) comm = comms[i];
235     if (comm != MPI_COMM_NULL) {
236       ierr = KSPCreate(comm,&mglevels[i]->smoothd);CHKERRQ(ierr);
237       ierr = KSPSetErrorIfNotConverged(mglevels[i]->smoothd,pc->erroriffailure);CHKERRQ(ierr);
238       ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr);
239       ierr = KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);CHKERRQ(ierr);
240       ierr = PetscObjectComposedDataSetInt((PetscObject) mglevels[i]->smoothd, PetscMGLevelId, mglevels[i]->level);CHKERRQ(ierr);
241       if (i || levels == 1) {
242         char tprefix[128];
243 
244         ierr = KSPSetType(mglevels[i]->smoothd,KSPCHEBYSHEV);CHKERRQ(ierr);
245         ierr = KSPSetConvergenceTest(mglevels[i]->smoothd,KSPConvergedSkip,NULL,NULL);CHKERRQ(ierr);
246         ierr = KSPSetNormType(mglevels[i]->smoothd,KSP_NORM_NONE);CHKERRQ(ierr);
247         ierr = KSPGetPC(mglevels[i]->smoothd,&ipc);CHKERRQ(ierr);
248         ierr = PCSetType(ipc,PCSOR);CHKERRQ(ierr);
249         ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg->default_smoothd);CHKERRQ(ierr);
250 
251         ierr = PetscSNPrintf(tprefix,128,"mg_levels_%d_",(int)i);CHKERRQ(ierr);
252         ierr = KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);CHKERRQ(ierr);
253       } else {
254         ierr = KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");CHKERRQ(ierr);
255 
256         /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */
257         ierr = KSPSetType(mglevels[0]->smoothd,KSPPREONLY);CHKERRQ(ierr);
258         ierr = KSPGetPC(mglevels[0]->smoothd,&ipc);CHKERRQ(ierr);
259         ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
260         if (size > 1) {
261           ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr);
262         } else {
263           ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr);
264         }
265         ierr = PCFactorSetShiftType(ipc,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
266       }
267       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[i]->smoothd);CHKERRQ(ierr);
268     }
269     mglevels[i]->smoothu = mglevels[i]->smoothd;
270     mg->rtol             = 0.0;
271     mg->abstol           = 0.0;
272     mg->dtol             = 0.0;
273     mg->ttol             = 0.0;
274     mg->cyclesperpcapply = 1;
275   }
276   mg->levels = mglevels;
277   ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
278   PetscFunctionReturn(0);
279 }
280 
281 /*@C
282    PCMGSetLevels - Sets the number of levels to use with MG.
283    Must be called before any other MG routine.
284 
285    Logically Collective on PC
286 
287    Input Parameters:
288 +  pc - the preconditioner context
289 .  levels - the number of levels
290 -  comms - optional communicators for each level; this is to allow solving the coarser problems
291            on smaller sets of processes. For processes that are not included in the computation
292            you must pass MPI_COMM_NULL.
293 
294    Level: intermediate
295 
296    Notes:
297      If the number of levels is one then the multigrid uses the -mg_levels prefix
298      for setting the level options rather than the -mg_coarse prefix.
299 
300      You can free the information in comms after this routine is called.
301 
302      The array of MPI communicators must contain MPI_COMM_NULL for those ranks that at each level
303      are not participating in the coarser solve. For example, with 2 levels and 1 and 2 ranks on
304      the two levels, rank 0 in the original communicator will pass in an array of 2 communicators
305      of size 2 and 1, while rank 1 in the original communicator will pass in array of 2 communicators
306      the first of size 2 and the second of value MPI_COMM_NULL since the rank 1 does not participate
307      in the coarse grid solve.
308 
309      Since each coarser level may have a new MPI_Comm with fewer ranks than the previous, one
310      must take special care in providing the restriction and interpolation operation. We recommend
311      providing these as two step operations; first perform a standard restriction or interpolation on
312      the full number of ranks for that level and then use an MPI call to copy the resulting vector
313      array entries (after calls to VecGetArray()) to the smaller or larger number of ranks, not in both
314      cases the MPI calls must be made on the larger of the two communicators. Traditional MPI send and
315      recieves or MPI_AlltoAllv() could be used to do the reshuffling of the vector entries.
316 
317 
318 .seealso: PCMGSetType(), PCMGGetLevels()
319 @*/
320 PetscErrorCode PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
321 {
322   PetscErrorCode ierr;
323 
324   PetscFunctionBegin;
325   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
326   if (comms) PetscValidPointer(comms,3);
327   ierr = PetscTryMethod(pc,"PCMGSetLevels_C",(PC,PetscInt,MPI_Comm*),(pc,levels,comms));CHKERRQ(ierr);
328   PetscFunctionReturn(0);
329 }
330 
331 
332 PetscErrorCode PCDestroy_MG(PC pc)
333 {
334   PetscErrorCode ierr;
335   PC_MG          *mg        = (PC_MG*)pc->data;
336   PC_MG_Levels   **mglevels = mg->levels;
337   PetscInt       i,n;
338 
339   PetscFunctionBegin;
340   ierr = PCReset_MG(pc);CHKERRQ(ierr);
341   if (mglevels) {
342     n = mglevels[0]->levels;
343     for (i=0; i<n; i++) {
344       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
345         ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
346       }
347       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
348       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
349     }
350     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
351   }
352   ierr = PetscFree(pc->data);CHKERRQ(ierr);
353   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGetInterpolations_C",NULL);CHKERRQ(ierr);
354   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGetCoarseOperators_C",NULL);CHKERRQ(ierr);
355   PetscFunctionReturn(0);
356 }
357 
358 
359 
360 extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**);
361 extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**);
362 extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**);
363 
364 /*
365    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
366              or full cycle of multigrid.
367 
368   Note:
369   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
370 */
371 static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
372 {
373   PC_MG          *mg        = (PC_MG*)pc->data;
374   PC_MG_Levels   **mglevels = mg->levels;
375   PetscErrorCode ierr;
376   PC             tpc;
377   PetscInt       levels = mglevels[0]->levels,i;
378   PetscBool      changeu,changed;
379 
380   PetscFunctionBegin;
381   if (mg->stageApply) {ierr = PetscLogStagePush(mg->stageApply);CHKERRQ(ierr);}
382   /* When the DM is supplying the matrix then it will not exist until here */
383   for (i=0; i<levels; i++) {
384     if (!mglevels[i]->A) {
385       ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr);
386       ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr);
387     }
388   }
389 
390   ierr = KSPGetPC(mglevels[levels-1]->smoothd,&tpc);CHKERRQ(ierr);
391   ierr = PCPreSolveChangeRHS(tpc,&changed);CHKERRQ(ierr);
392   ierr = KSPGetPC(mglevels[levels-1]->smoothu,&tpc);CHKERRQ(ierr);
393   ierr = PCPreSolveChangeRHS(tpc,&changeu);CHKERRQ(ierr);
394   if (!changeu && !changed) {
395     ierr = VecDestroy(&mglevels[levels-1]->b);CHKERRQ(ierr);
396     mglevels[levels-1]->b = b;
397   } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */
398     if (!mglevels[levels-1]->b) {
399       Vec *vec;
400 
401       ierr = KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
402       mglevels[levels-1]->b = *vec;
403       ierr = PetscFree(vec);CHKERRQ(ierr);
404     }
405     ierr = VecCopy(b,mglevels[levels-1]->b);CHKERRQ(ierr);
406   }
407   mglevels[levels-1]->x = x;
408 
409   if (mg->am == PC_MG_MULTIPLICATIVE) {
410     ierr = VecSet(x,0.0);CHKERRQ(ierr);
411     for (i=0; i<mg->cyclesperpcapply; i++) {
412       ierr = PCMGMCycle_Private(pc,mglevels+levels-1,NULL);CHKERRQ(ierr);
413     }
414   } else if (mg->am == PC_MG_ADDITIVE) {
415     ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr);
416   } else if (mg->am == PC_MG_KASKADE) {
417     ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr);
418   } else {
419     ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr);
420   }
421   if (mg->stageApply) {ierr = PetscLogStagePop();CHKERRQ(ierr);}
422   if (!changeu && !changed) mglevels[levels-1]->b = NULL;
423   PetscFunctionReturn(0);
424 }
425 
426 
427 PetscErrorCode PCSetFromOptions_MG(PetscOptionItems *PetscOptionsObject,PC pc)
428 {
429   PetscErrorCode   ierr;
430   PetscInt         levels,cycles;
431   PetscBool        flg, flg2;
432   PC_MG            *mg = (PC_MG*)pc->data;
433   PC_MG_Levels     **mglevels;
434   PCMGType         mgtype;
435   PCMGCycleType    mgctype;
436   PCMGGalerkinType gtype;
437 
438   PetscFunctionBegin;
439   levels = PetscMax(mg->nlevels,1);
440   ierr = PetscOptionsHead(PetscOptionsObject,"Multigrid options");CHKERRQ(ierr);
441   ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr);
442   if (!flg && !mg->levels && pc->dm) {
443     ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
444     levels++;
445     mg->usedmfornumberoflevels = PETSC_TRUE;
446   }
447   ierr = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr);
448   mglevels = mg->levels;
449 
450   mgctype = (PCMGCycleType) mglevels[0]->cycles;
451   ierr    = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr);
452   if (flg) {
453     ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr);
454   }
455   gtype = mg->galerkin;
456   ierr = PetscOptionsEnum("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",PCMGGalerkinTypes,(PetscEnum)gtype,(PetscEnum*)&gtype,&flg);CHKERRQ(ierr);
457   if (flg) {
458     ierr = PCMGSetGalerkin(pc,gtype);CHKERRQ(ierr);
459   }
460   flg2 = PETSC_FALSE;
461   ierr = PetscOptionsBool("-pc_mg_adapt_interp","Adapt interpolation using some coarse space","PCMGSetAdaptInterpolation",PETSC_FALSE,&flg2,&flg);CHKERRQ(ierr);
462   if (flg) {ierr = PCMGSetAdaptInterpolation(pc, flg2);CHKERRQ(ierr);}
463   ierr = PetscOptionsInt("-pc_mg_adapt_interp_n","Size of the coarse space for adaptive interpolation","PCMGSetCoarseSpace",mg->Nc,&mg->Nc,&flg);CHKERRQ(ierr);
464   ierr = PetscOptionsEnum("-pc_mg_adapt_interp_coarse_space","Type of coarse space: polynomial, harmonic, eigenvector, generalized_eigenvector","PCMGSetAdaptCoarseSpaceType",PCMGCoarseSpaceTypes,(PetscEnum)mg->coarseSpaceType,(PetscEnum*)&mg->coarseSpaceType,&flg);CHKERRQ(ierr);
465   ierr = PetscOptionsBool("-pc_mg_mesp_monitor","Monitor the multilevel eigensolver","PCMGSetAdaptInterpolation",PETSC_FALSE,&mg->mespMonitor,&flg);CHKERRQ(ierr);
466   flg = PETSC_FALSE;
467   ierr = PetscOptionsBool("-pc_mg_distinct_smoothup","Create separate smoothup KSP and append the prefix _up","PCMGSetDistinctSmoothUp",PETSC_FALSE,&flg,NULL);CHKERRQ(ierr);
468   if (flg) {
469     ierr = PCMGSetDistinctSmoothUp(pc);CHKERRQ(ierr);
470   }
471   mgtype = mg->am;
472   ierr   = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
473   if (flg) {
474     ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
475   }
476   if (mg->am == PC_MG_MULTIPLICATIVE) {
477     ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGMultiplicativeSetCycles",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr);
478     if (flg) {
479       ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr);
480     }
481   }
482   flg  = PETSC_FALSE;
483   ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,NULL);CHKERRQ(ierr);
484   if (flg) {
485     PetscInt i;
486     char     eventname[128];
487 
488     levels = mglevels[0]->levels;
489     for (i=0; i<levels; i++) {
490       sprintf(eventname,"MGSetup Level %d",(int)i);
491       ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr);
492       sprintf(eventname,"MGSmooth Level %d",(int)i);
493       ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr);
494       if (i) {
495         sprintf(eventname,"MGResid Level %d",(int)i);
496         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr);
497         sprintf(eventname,"MGInterp Level %d",(int)i);
498         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr);
499       }
500     }
501 
502 #if defined(PETSC_USE_LOG)
503     {
504       const char    *sname = "MG Apply";
505       PetscStageLog stageLog;
506       PetscInt      st;
507 
508       ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
509       for (st = 0; st < stageLog->numStages; ++st) {
510         PetscBool same;
511 
512         ierr = PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);CHKERRQ(ierr);
513         if (same) mg->stageApply = st;
514       }
515       if (!mg->stageApply) {
516         ierr = PetscLogStageRegister(sname, &mg->stageApply);CHKERRQ(ierr);
517       }
518     }
519 #endif
520   }
521   ierr = PetscOptionsTail();CHKERRQ(ierr);
522   /* Check option consistency */
523   ierr = PCMGGetGalerkin(pc, &gtype);CHKERRQ(ierr);
524   ierr = PCMGGetAdaptInterpolation(pc, &flg);CHKERRQ(ierr);
525   if (flg && (gtype >= PC_MG_GALERKIN_NONE)) SETERRQ(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_INCOMP, "Must use Galerkin coarse operators when adapting the interpolator");
526   PetscFunctionReturn(0);
527 }
528 
529 const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",NULL};
530 const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",NULL};
531 const char *const PCMGGalerkinTypes[] = {"both","pmat","mat","none","external","PCMGGalerkinType","PC_MG_GALERKIN",NULL};
532 const char *const PCMGCoarseSpaceTypes[] = {"polynomial","harmonic","eigenvector","generalized_eigenvector","PCMGCoarseSpaceType","PCMG_POLYNOMIAL",NULL};
533 
534 #include <petscdraw.h>
535 PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
536 {
537   PC_MG          *mg        = (PC_MG*)pc->data;
538   PC_MG_Levels   **mglevels = mg->levels;
539   PetscErrorCode ierr;
540   PetscInt       levels = mglevels ? mglevels[0]->levels : 0,i;
541   PetscBool      iascii,isbinary,isdraw;
542 
543   PetscFunctionBegin;
544   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
545   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
546   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
547   if (iascii) {
548     const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown";
549     ierr = PetscViewerASCIIPrintf(viewer,"  type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,cyclename);CHKERRQ(ierr);
550     if (mg->am == PC_MG_MULTIPLICATIVE) {
551       ierr = PetscViewerASCIIPrintf(viewer,"    Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr);
552     }
553     if (mg->galerkin == PC_MG_GALERKIN_BOTH) {
554       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
555     } else if (mg->galerkin == PC_MG_GALERKIN_PMAT) {
556       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices for pmat\n");CHKERRQ(ierr);
557     } else if (mg->galerkin == PC_MG_GALERKIN_MAT) {
558       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices for mat\n");CHKERRQ(ierr);
559     } else if (mg->galerkin == PC_MG_GALERKIN_EXTERNAL) {
560       ierr = PetscViewerASCIIPrintf(viewer,"    Using externally compute Galerkin coarse grid matrices\n");CHKERRQ(ierr);
561     } else {
562       ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
563     }
564     if (mg->view){
565       ierr = (*mg->view)(pc,viewer);CHKERRQ(ierr);
566     }
567     for (i=0; i<levels; i++) {
568       if (!i) {
569         ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);CHKERRQ(ierr);
570       } else {
571         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
572       }
573       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
574       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
575       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
576       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
577         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
578       } else if (i) {
579         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
580         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
581         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
582         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
583       }
584     }
585   } else if (isbinary) {
586     for (i=levels-1; i>=0; i--) {
587       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
588       if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) {
589         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
590       }
591     }
592   } else if (isdraw) {
593     PetscDraw draw;
594     PetscReal x,w,y,bottom,th;
595     ierr   = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
596     ierr   = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
597     ierr   = PetscDrawStringGetSize(draw,NULL,&th);CHKERRQ(ierr);
598     bottom = y - th;
599     for (i=levels-1; i>=0; i--) {
600       if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) {
601         ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr);
602         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
603         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
604       } else {
605         w    = 0.5*PetscMin(1.0-x,x);
606         ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr);
607         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
608         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
609         ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr);
610         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
611         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
612       }
613       ierr    = PetscDrawGetBoundingBox(draw,NULL,&bottom,NULL,NULL);CHKERRQ(ierr);
614       bottom -= th;
615     }
616   }
617   PetscFunctionReturn(0);
618 }
619 
620 #include <petsc/private/kspimpl.h>
621 
622 /*
623     Calls setup for the KSP on each level
624 */
625 PetscErrorCode PCSetUp_MG(PC pc)
626 {
627   PC_MG          *mg        = (PC_MG*)pc->data;
628   PC_MG_Levels   **mglevels = mg->levels;
629   PetscErrorCode ierr;
630   PetscInt       i,n;
631   PC             cpc;
632   PetscBool      dump = PETSC_FALSE,opsset,use_amat,missinginterpolate = PETSC_FALSE;
633   Mat            dA,dB;
634   Vec            tvec;
635   DM             *dms;
636   PetscViewer    viewer = NULL;
637   PetscBool      dAeqdB = PETSC_FALSE, needRestricts = PETSC_FALSE;
638 
639   PetscFunctionBegin;
640   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels with PCMGSetLevels() before setting up");
641   n = mglevels[0]->levels;
642   /* FIX: Move this to PCSetFromOptions_MG? */
643   if (mg->usedmfornumberoflevels) {
644     PetscInt levels;
645     ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
646     levels++;
647     if (levels > n) { /* the problem is now being solved on a finer grid */
648       ierr     = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr);
649       n        = levels;
650       ierr     = PCSetFromOptions(pc);CHKERRQ(ierr); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
651       mglevels = mg->levels;
652     }
653   }
654   ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr);
655 
656 
657   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
658   /* so use those from global PC */
659   /* Is this what we always want? What if user wants to keep old one? */
660   ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset);CHKERRQ(ierr);
661   if (opsset) {
662     Mat mmat;
663     ierr = KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat);CHKERRQ(ierr);
664     if (mmat == pc->pmat) opsset = PETSC_FALSE;
665   }
666 
667   if (!opsset) {
668     ierr = PCGetUseAmat(pc,&use_amat);CHKERRQ(ierr);
669     if (use_amat) {
670       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);
671       ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat);CHKERRQ(ierr);
672     } else {
673       ierr = PetscInfo(pc,"Using matrix (pmat) operators to define finest grid operator \n  because PCMGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n");CHKERRQ(ierr);
674       ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat);CHKERRQ(ierr);
675     }
676   }
677 
678   for (i=n-1; i>0; i--) {
679     if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) {
680       missinginterpolate = PETSC_TRUE;
681       continue;
682     }
683   }
684 
685   ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB);CHKERRQ(ierr);
686   if (dA == dB) dAeqdB = PETSC_TRUE;
687   if ((mg->galerkin == PC_MG_GALERKIN_NONE) || (((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_MAT)) && !dAeqdB)) {
688     needRestricts = PETSC_TRUE;  /* user must compute either mat, pmat, or both so must restrict x to coarser levels */
689   }
690 
691 
692   /*
693    Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS)
694    Skipping for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs?
695   */
696   if (missinginterpolate && pc->dm && mg->galerkin != PC_MG_GALERKIN_EXTERNAL && !pc->setupcalled) {
697         /* construct the interpolation from the DMs */
698     Mat p;
699     Vec rscale;
700     ierr     = PetscMalloc1(n,&dms);CHKERRQ(ierr);
701     dms[n-1] = pc->dm;
702     /* Separately create them so we do not get DMKSP interference between levels */
703     for (i=n-2; i>-1; i--) {ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr);}
704         /*
705            Force the mat type of coarse level operator to be AIJ because usually we want to use LU for coarse level.
706            Notice that it can be overwritten by -mat_type because KSPSetUp() reads command line options.
707            But it is safe to use -dm_mat_type.
708 
709            The mat type should not be hardcoded like this, we need to find a better way.
710     ierr = DMSetMatType(dms[0],MATAIJ);CHKERRQ(ierr);
711     */
712     for (i=n-2; i>-1; i--) {
713       DMKSP     kdm;
714       PetscBool dmhasrestrict, dmhasinject;
715       ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr);
716       if (!needRestricts) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);}
717       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
718         ierr = KSPSetDM(mglevels[i]->smoothu,dms[i]);CHKERRQ(ierr);
719         if (!needRestricts) {ierr = KSPSetDMActive(mglevels[i]->smoothu,PETSC_FALSE);CHKERRQ(ierr);}
720       }
721       ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr);
722       /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take
723        * a bitwise OR of computing the matrix, RHS, and initial iterate. */
724       kdm->ops->computerhs = NULL;
725       kdm->rhsctx          = NULL;
726       if (!mglevels[i+1]->interpolate) {
727         ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr);
728         ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr);
729         if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);}
730         ierr = VecDestroy(&rscale);CHKERRQ(ierr);
731         ierr = MatDestroy(&p);CHKERRQ(ierr);
732       }
733       ierr = DMHasCreateRestriction(dms[i],&dmhasrestrict);CHKERRQ(ierr);
734       if (dmhasrestrict && !mglevels[i+1]->restrct){
735         ierr = DMCreateRestriction(dms[i],dms[i+1],&p);CHKERRQ(ierr);
736         ierr = PCMGSetRestriction(pc,i+1,p);CHKERRQ(ierr);
737         ierr = MatDestroy(&p);CHKERRQ(ierr);
738       }
739       ierr = DMHasCreateInjection(dms[i],&dmhasinject);CHKERRQ(ierr);
740       if (dmhasinject && !mglevels[i+1]->inject){
741         ierr = DMCreateInjection(dms[i],dms[i+1],&p);CHKERRQ(ierr);
742         ierr = PCMGSetInjection(pc,i+1,p);CHKERRQ(ierr);
743         ierr = MatDestroy(&p);CHKERRQ(ierr);
744       }
745     }
746 
747     for (i=n-2; i>-1; i--) {ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);}
748     ierr = PetscFree(dms);CHKERRQ(ierr);
749   }
750 
751   if (pc->dm && !pc->setupcalled) {
752     /* finest smoother also gets DM but it is not active, independent of whether galerkin==PC_MG_GALERKIN_EXTERNAL */
753     ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr);
754     ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr);
755     if (mglevels[n-1]->smoothd != mglevels[n-1]->smoothu) {
756       ierr = KSPSetDM(mglevels[n-1]->smoothu,pc->dm);CHKERRQ(ierr);
757       ierr = KSPSetDMActive(mglevels[n-1]->smoothu,PETSC_FALSE);CHKERRQ(ierr);
758     }
759   }
760 
761   if (mg->galerkin < PC_MG_GALERKIN_NONE) {
762     Mat       A,B;
763     PetscBool doA = PETSC_FALSE,doB = PETSC_FALSE;
764     MatReuse  reuse = MAT_INITIAL_MATRIX;
765 
766     if ((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_BOTH)) doB = PETSC_TRUE;
767     if ((mg->galerkin == PC_MG_GALERKIN_MAT) || ((mg->galerkin == PC_MG_GALERKIN_BOTH) && (dA != dB))) doA = PETSC_TRUE;
768     if (pc->setupcalled) reuse = MAT_REUSE_MATRIX;
769     for (i=n-2; i>-1; i--) {
770       if (!mglevels[i+1]->restrct && !mglevels[i+1]->interpolate) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must provide interpolation or restriction for each MG level except level 0");
771       if (!mglevels[i+1]->interpolate) {
772         ierr = PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct);CHKERRQ(ierr);
773       }
774       if (!mglevels[i+1]->restrct) {
775         ierr = PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate);CHKERRQ(ierr);
776       }
777       if (reuse == MAT_REUSE_MATRIX) {
778         ierr = KSPGetOperators(mglevels[i]->smoothd,&A,&B);CHKERRQ(ierr);
779       }
780       if (doA) {
781         ierr = MatGalerkin(mglevels[i+1]->restrct,dA,mglevels[i+1]->interpolate,reuse,1.0,&A);CHKERRQ(ierr);
782       }
783       if (doB) {
784         ierr = MatGalerkin(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,reuse,1.0,&B);CHKERRQ(ierr);
785       }
786       /* the management of the PetscObjectReference() and PetscObjecDereference() below is rather delicate */
787       if (!doA && dAeqdB) {
788         if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);}
789         A = B;
790       } else if (!doA && reuse == MAT_INITIAL_MATRIX) {
791         ierr = KSPGetOperators(mglevels[i]->smoothd,&A,NULL);CHKERRQ(ierr);
792         ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
793       }
794       if (!doB && dAeqdB) {
795         if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);}
796         B = A;
797       } else if (!doB && reuse == MAT_INITIAL_MATRIX) {
798         ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&B);CHKERRQ(ierr);
799         ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
800       }
801       if (reuse == MAT_INITIAL_MATRIX) {
802         ierr = KSPSetOperators(mglevels[i]->smoothd,A,B);CHKERRQ(ierr);
803         ierr = PetscObjectDereference((PetscObject)A);CHKERRQ(ierr);
804         ierr = PetscObjectDereference((PetscObject)B);CHKERRQ(ierr);
805       }
806       dA = A;
807       dB = B;
808     }
809   }
810 
811 
812   /* Adapt interpolation matrices */
813   if (mg->adaptInterpolation) {
814     mg->Nc = mg->Nc < 0 ? 6 : mg->Nc; /* Default to 6 modes */
815     for (i = 0; i < n; ++i) {
816       ierr = PCMGComputeCoarseSpace_Internal(pc, i, mg->coarseSpaceType, mg->Nc, !i ? NULL : mglevels[i-1]->coarseSpace, &mglevels[i]->coarseSpace);CHKERRQ(ierr);
817       if (i) {ierr = PCMGAdaptInterpolator_Internal(pc, i, mglevels[i-1]->smoothu, mglevels[i]->smoothu, mg->Nc, mglevels[i-1]->coarseSpace, mglevels[i]->coarseSpace);CHKERRQ(ierr);}
818     }
819     for (i = n-2; i > -1; --i) {
820       ierr = PCMGRecomputeLevelOperators_Internal(pc, i);CHKERRQ(ierr);
821     }
822   }
823 
824   if (needRestricts && pc->dm) {
825     for (i=n-2; i>=0; i--) {
826       DM  dmfine,dmcoarse;
827       Mat Restrict,Inject;
828       Vec rscale;
829       ierr = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr);
830       ierr = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr);
831       ierr = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr);
832       ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
833       ierr = PCMGGetInjection(pc,i+1,&Inject);CHKERRQ(ierr);
834       ierr = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr);
835     }
836   }
837 
838   if (!pc->setupcalled) {
839     for (i=0; i<n; i++) {
840       ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr);
841     }
842     for (i=1; i<n; i++) {
843       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
844         ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr);
845       }
846     }
847     /* insure that if either interpolation or restriction is set the other other one is set */
848     for (i=1; i<n; i++) {
849       ierr = PCMGGetInterpolation(pc,i,NULL);CHKERRQ(ierr);
850       ierr = PCMGGetRestriction(pc,i,NULL);CHKERRQ(ierr);
851     }
852     for (i=0; i<n-1; i++) {
853       if (!mglevels[i]->b) {
854         Vec *vec;
855         ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
856         ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr);
857         ierr = VecDestroy(vec);CHKERRQ(ierr);
858         ierr = PetscFree(vec);CHKERRQ(ierr);
859       }
860       if (!mglevels[i]->r && i) {
861         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
862         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
863         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
864       }
865       if (!mglevels[i]->x) {
866         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
867         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
868         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
869       }
870     }
871     if (n != 1 && !mglevels[n-1]->r) {
872       /* PCMGSetR() on the finest level if user did not supply it */
873       Vec *vec;
874       ierr = KSPCreateVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
875       ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr);
876       ierr = VecDestroy(vec);CHKERRQ(ierr);
877       ierr = PetscFree(vec);CHKERRQ(ierr);
878     }
879   }
880 
881   if (pc->dm) {
882     /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */
883     for (i=0; i<n-1; i++) {
884       if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX;
885     }
886   }
887 
888   for (i=1; i<n; i++) {
889     if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1){
890       /* if doing only down then initial guess is zero */
891       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
892     }
893     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
894     ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr);
895     if (mglevels[i]->smoothd->reason == KSP_DIVERGED_PC_FAILED) {
896       pc->failedreason = PC_SUBPC_ERROR;
897     }
898     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
899     if (!mglevels[i]->residual) {
900       Mat mat;
901       ierr = KSPGetOperators(mglevels[i]->smoothd,&mat,NULL);CHKERRQ(ierr);
902       ierr = PCMGSetResidual(pc,i,PCMGResidualDefault,mat);CHKERRQ(ierr);
903     }
904   }
905   for (i=1; i<n; i++) {
906     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
907       Mat downmat,downpmat;
908 
909       /* check if operators have been set for up, if not use down operators to set them */
910       ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);CHKERRQ(ierr);
911       if (!opsset) {
912         ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);CHKERRQ(ierr);
913         ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat);CHKERRQ(ierr);
914       }
915 
916       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
917       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
918       ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr);
919       if (mglevels[i]->smoothu->reason == KSP_DIVERGED_PC_FAILED) {
920         pc->failedreason = PC_SUBPC_ERROR;
921       }
922       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
923     }
924   }
925 
926   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
927   ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr);
928   if (mglevels[0]->smoothd->reason == KSP_DIVERGED_PC_FAILED) {
929     pc->failedreason = PC_SUBPC_ERROR;
930   }
931   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
932 
933   /*
934      Dump the interpolation/restriction matrices plus the
935    Jacobian/stiffness on each level. This allows MATLAB users to
936    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
937 
938    Only support one or the other at the same time.
939   */
940 #if defined(PETSC_USE_SOCKET_VIEWER)
941   ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);CHKERRQ(ierr);
942   if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc));
943   dump = PETSC_FALSE;
944 #endif
945   ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);CHKERRQ(ierr);
946   if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc));
947 
948   if (viewer) {
949     for (i=1; i<n; i++) {
950       ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr);
951     }
952     for (i=0; i<n; i++) {
953       ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr);
954       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
955     }
956   }
957   PetscFunctionReturn(0);
958 }
959 
960 /* -------------------------------------------------------------------------------------*/
961 
962 PetscErrorCode PCMGGetLevels_MG(PC pc, PetscInt *levels)
963 {
964   PC_MG *mg = (PC_MG *) pc->data;
965 
966   PetscFunctionBegin;
967   *levels = mg->nlevels;
968   PetscFunctionReturn(0);
969 }
970 
971 /*@
972    PCMGGetLevels - Gets the number of levels to use with MG.
973 
974    Not Collective
975 
976    Input Parameter:
977 .  pc - the preconditioner context
978 
979    Output parameter:
980 .  levels - the number of levels
981 
982    Level: advanced
983 
984 .seealso: PCMGSetLevels()
985 @*/
986 PetscErrorCode PCMGGetLevels(PC pc,PetscInt *levels)
987 {
988   PetscErrorCode ierr;
989 
990   PetscFunctionBegin;
991   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
992   PetscValidIntPointer(levels,2);
993   *levels = 0;
994   ierr = PetscTryMethod(pc,"PCMGGetLevels_C",(PC,PetscInt*),(pc,levels));CHKERRQ(ierr);
995   PetscFunctionReturn(0);
996 }
997 
998 /*@
999    PCMGSetType - Determines the form of multigrid to use:
1000    multiplicative, additive, full, or the Kaskade algorithm.
1001 
1002    Logically Collective on PC
1003 
1004    Input Parameters:
1005 +  pc - the preconditioner context
1006 -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
1007    PC_MG_FULL, PC_MG_KASKADE
1008 
1009    Options Database Key:
1010 .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
1011    additive, full, kaskade
1012 
1013    Level: advanced
1014 
1015 .seealso: PCMGSetLevels()
1016 @*/
1017 PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
1018 {
1019   PC_MG *mg = (PC_MG*)pc->data;
1020 
1021   PetscFunctionBegin;
1022   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1023   PetscValidLogicalCollectiveEnum(pc,form,2);
1024   mg->am = form;
1025   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
1026   else pc->ops->applyrichardson = NULL;
1027   PetscFunctionReturn(0);
1028 }
1029 
1030 /*@
1031    PCMGGetType - Determines the form of multigrid to use:
1032    multiplicative, additive, full, or the Kaskade algorithm.
1033 
1034    Logically Collective on PC
1035 
1036    Input Parameter:
1037 .  pc - the preconditioner context
1038 
1039    Output Parameter:
1040 .  type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE
1041 
1042 
1043    Level: advanced
1044 
1045 .seealso: PCMGSetLevels()
1046 @*/
1047 PetscErrorCode  PCMGGetType(PC pc,PCMGType *type)
1048 {
1049   PC_MG *mg = (PC_MG*)pc->data;
1050 
1051   PetscFunctionBegin;
1052   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1053   *type = mg->am;
1054   PetscFunctionReturn(0);
1055 }
1056 
1057 /*@
1058    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
1059    complicated cycling.
1060 
1061    Logically Collective on PC
1062 
1063    Input Parameters:
1064 +  pc - the multigrid context
1065 -  n - either PC_MG_CYCLE_V or PC_MG_CYCLE_W
1066 
1067    Options Database Key:
1068 .  -pc_mg_cycle_type <v,w> - provide the cycle desired
1069 
1070    Level: advanced
1071 
1072 .seealso: PCMGSetCycleTypeOnLevel()
1073 @*/
1074 PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
1075 {
1076   PC_MG        *mg        = (PC_MG*)pc->data;
1077   PC_MG_Levels **mglevels = mg->levels;
1078   PetscInt     i,levels;
1079 
1080   PetscFunctionBegin;
1081   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1082   PetscValidLogicalCollectiveEnum(pc,n,2);
1083   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
1084   levels = mglevels[0]->levels;
1085   for (i=0; i<levels; i++) mglevels[i]->cycles = n;
1086   PetscFunctionReturn(0);
1087 }
1088 
1089 /*@
1090    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
1091          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
1092 
1093    Logically Collective on PC
1094 
1095    Input Parameters:
1096 +  pc - the multigrid context
1097 -  n - number of cycles (default is 1)
1098 
1099    Options Database Key:
1100 .  -pc_mg_multiplicative_cycles n
1101 
1102    Level: advanced
1103 
1104    Notes:
1105     This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
1106 
1107 .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
1108 @*/
1109 PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
1110 {
1111   PC_MG        *mg        = (PC_MG*)pc->data;
1112 
1113   PetscFunctionBegin;
1114   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1115   PetscValidLogicalCollectiveInt(pc,n,2);
1116   mg->cyclesperpcapply = n;
1117   PetscFunctionReturn(0);
1118 }
1119 
1120 PetscErrorCode PCMGSetGalerkin_MG(PC pc,PCMGGalerkinType use)
1121 {
1122   PC_MG *mg = (PC_MG*)pc->data;
1123 
1124   PetscFunctionBegin;
1125   mg->galerkin = use;
1126   PetscFunctionReturn(0);
1127 }
1128 
1129 /*@
1130    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
1131       finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i
1132 
1133    Logically Collective on PC
1134 
1135    Input Parameters:
1136 +  pc - the multigrid context
1137 -  use - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, or PC_MG_GALERKIN_NONE
1138 
1139    Options Database Key:
1140 .  -pc_mg_galerkin <both,pmat,mat,none>
1141 
1142    Level: intermediate
1143 
1144    Notes:
1145     Some codes that use PCMG such as PCGAMG use Galerkin internally while constructing the hierarchy and thus do not
1146      use the PCMG construction of the coarser grids.
1147 
1148 .seealso: PCMGGetGalerkin(), PCMGGalerkinType
1149 
1150 @*/
1151 PetscErrorCode PCMGSetGalerkin(PC pc,PCMGGalerkinType use)
1152 {
1153   PetscErrorCode ierr;
1154 
1155   PetscFunctionBegin;
1156   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1157   ierr = PetscTryMethod(pc,"PCMGSetGalerkin_C",(PC,PCMGGalerkinType),(pc,use));CHKERRQ(ierr);
1158   PetscFunctionReturn(0);
1159 }
1160 
1161 /*@
1162    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
1163       A_i-1 = r_i * A_i * p_i
1164 
1165    Not Collective
1166 
1167    Input Parameter:
1168 .  pc - the multigrid context
1169 
1170    Output Parameter:
1171 .  galerkin - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, PC_MG_GALERKIN_NONE, or PC_MG_GALERKIN_EXTERNAL
1172 
1173    Level: intermediate
1174 
1175 .seealso: PCMGSetGalerkin(), PCMGGalerkinType
1176 
1177 @*/
1178 PetscErrorCode  PCMGGetGalerkin(PC pc,PCMGGalerkinType  *galerkin)
1179 {
1180   PC_MG *mg = (PC_MG*)pc->data;
1181 
1182   PetscFunctionBegin;
1183   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1184   *galerkin = mg->galerkin;
1185   PetscFunctionReturn(0);
1186 }
1187 
1188 PetscErrorCode PCMGSetAdaptInterpolation_MG(PC pc, PetscBool adapt)
1189 {
1190   PC_MG *mg = (PC_MG *) pc->data;
1191 
1192   PetscFunctionBegin;
1193   mg->adaptInterpolation = adapt;
1194   PetscFunctionReturn(0);
1195 }
1196 
1197 PetscErrorCode PCMGGetAdaptInterpolation_MG(PC pc, PetscBool *adapt)
1198 {
1199   PC_MG *mg = (PC_MG *) pc->data;
1200 
1201   PetscFunctionBegin;
1202   *adapt = mg->adaptInterpolation;
1203   PetscFunctionReturn(0);
1204 }
1205 
1206 /*@
1207   PCMGSetAdaptInterpolation - Adapt the interpolator based upon a vector space which should be accurately captured by the next coarser mesh, and thus accurately interpolated.
1208 
1209   Logically Collective on PC
1210 
1211   Input Parameters:
1212 + pc    - the multigrid context
1213 - adapt - flag for adaptation of the interpolator
1214 
1215   Options Database Keys:
1216 + -pc_mg_adapt_interp                     - Turn on adaptation
1217 . -pc_mg_adapt_interp_n <int>             - The number of modes to use, should be divisible by dimension
1218 - -pc_mg_adapt_interp_coarse_space <type> - The type of coarse space: polynomial, harmonic, eigenvector, generalized_eigenvector
1219 
1220   Level: intermediate
1221 
1222 .keywords: MG, set, Galerkin
1223 .seealso: PCMGGetAdaptInterpolation(), PCMGSetGalerkin()
1224 @*/
1225 PetscErrorCode PCMGSetAdaptInterpolation(PC pc, PetscBool adapt)
1226 {
1227   PetscErrorCode ierr;
1228 
1229   PetscFunctionBegin;
1230   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1231   ierr = PetscTryMethod(pc,"PCMGSetAdaptInterpolation_C",(PC,PetscBool),(pc,adapt));CHKERRQ(ierr);
1232   PetscFunctionReturn(0);
1233 }
1234 
1235 /*@
1236   PCMGGetAdaptInterpolation - Get the flag to adapt the interpolator based upon a vector space which should be accurately captured by the next coarser mesh, and thus accurately interpolated.
1237 
1238   Not collective
1239 
1240   Input Parameter:
1241 . pc    - the multigrid context
1242 
1243   Output Parameter:
1244 . adapt - flag for adaptation of the interpolator
1245 
1246   Level: intermediate
1247 
1248 .keywords: MG, set, Galerkin
1249 .seealso: PCMGSetAdaptInterpolation(), PCMGSetGalerkin()
1250 @*/
1251 PetscErrorCode PCMGGetAdaptInterpolation(PC pc, PetscBool *adapt)
1252 {
1253   PetscErrorCode ierr;
1254 
1255   PetscFunctionBegin;
1256   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1257   PetscValidPointer(adapt, 2);
1258   ierr = PetscUseMethod(pc,"PCMGGetAdaptInterpolation_C",(PC,PetscBool*),(pc,adapt));CHKERRQ(ierr);
1259   PetscFunctionReturn(0);
1260 }
1261 
1262 /*@
1263    PCMGSetNumberSmooth - Sets the number of pre and post-smoothing steps to use
1264    on all levels.  Use PCMGDistinctSmoothUp() to create separate up and down smoothers if you want different numbers of
1265    pre- and post-smoothing steps.
1266 
1267    Logically Collective on PC
1268 
1269    Input Parameters:
1270 +  mg - the multigrid context
1271 -  n - the number of smoothing steps
1272 
1273    Options Database Key:
1274 .  -mg_levels_ksp_max_it <n> - Sets number of pre and post-smoothing steps
1275 
1276    Level: advanced
1277 
1278    Notes:
1279     this does not set a value on the coarsest grid, since we assume that
1280     there is no separate smooth up on the coarsest grid.
1281 
1282 .seealso: PCMGSetDistinctSmoothUp()
1283 @*/
1284 PetscErrorCode  PCMGSetNumberSmooth(PC pc,PetscInt n)
1285 {
1286   PC_MG          *mg        = (PC_MG*)pc->data;
1287   PC_MG_Levels   **mglevels = mg->levels;
1288   PetscErrorCode ierr;
1289   PetscInt       i,levels;
1290 
1291   PetscFunctionBegin;
1292   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1293   PetscValidLogicalCollectiveInt(pc,n,2);
1294   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
1295   levels = mglevels[0]->levels;
1296 
1297   for (i=1; i<levels; i++) {
1298     ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
1299     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
1300     mg->default_smoothu = n;
1301     mg->default_smoothd = n;
1302   }
1303   PetscFunctionReturn(0);
1304 }
1305 
1306 /*@
1307    PCMGSetDistinctSmoothUp - sets the up (post) smoother to be a separate KSP from the down (pre) smoother on all levels
1308        and adds the suffix _up to the options name
1309 
1310    Logically Collective on PC
1311 
1312    Input Parameters:
1313 .  pc - the preconditioner context
1314 
1315    Options Database Key:
1316 .  -pc_mg_distinct_smoothup
1317 
1318    Level: advanced
1319 
1320    Notes:
1321     this does not set a value on the coarsest grid, since we assume that
1322     there is no separate smooth up on the coarsest grid.
1323 
1324 .seealso: PCMGSetNumberSmooth()
1325 @*/
1326 PetscErrorCode  PCMGSetDistinctSmoothUp(PC pc)
1327 {
1328   PC_MG          *mg        = (PC_MG*)pc->data;
1329   PC_MG_Levels   **mglevels = mg->levels;
1330   PetscErrorCode ierr;
1331   PetscInt       i,levels;
1332   KSP            subksp;
1333 
1334   PetscFunctionBegin;
1335   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1336   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
1337   levels = mglevels[0]->levels;
1338 
1339   for (i=1; i<levels; i++) {
1340     const char *prefix = NULL;
1341     /* make sure smoother up and down are different */
1342     ierr = PCMGGetSmootherUp(pc,i,&subksp);CHKERRQ(ierr);
1343     ierr = KSPGetOptionsPrefix(mglevels[i]->smoothd,&prefix);CHKERRQ(ierr);
1344     ierr = KSPSetOptionsPrefix(subksp,prefix);CHKERRQ(ierr);
1345     ierr = KSPAppendOptionsPrefix(subksp,"up_");CHKERRQ(ierr);
1346   }
1347   PetscFunctionReturn(0);
1348 }
1349 
1350 /* No new matrices are created, and the coarse operator matrices are the references to the original ones */
1351 PetscErrorCode  PCGetInterpolations_MG(PC pc,PetscInt *num_levels,Mat *interpolations[])
1352 {
1353   PC_MG          *mg        = (PC_MG*)pc->data;
1354   PC_MG_Levels   **mglevels = mg->levels;
1355   Mat            *mat;
1356   PetscInt       l;
1357   PetscErrorCode ierr;
1358 
1359   PetscFunctionBegin;
1360   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1361   ierr = PetscMalloc1(mg->nlevels,&mat);CHKERRQ(ierr);
1362   for (l=1; l< mg->nlevels; l++) {
1363     mat[l-1] = mglevels[l]->interpolate;
1364     ierr = PetscObjectReference((PetscObject)mat[l-1]);CHKERRQ(ierr);
1365   }
1366   *num_levels = mg->nlevels;
1367   *interpolations = mat;
1368   PetscFunctionReturn(0);
1369 }
1370 
1371 /* No new matrices are created, and the coarse operator matrices are the references to the original ones */
1372 PetscErrorCode  PCGetCoarseOperators_MG(PC pc,PetscInt *num_levels,Mat *coarseOperators[])
1373 {
1374   PC_MG          *mg        = (PC_MG*)pc->data;
1375   PC_MG_Levels   **mglevels = mg->levels;
1376   PetscInt       l;
1377   Mat            *mat;
1378   PetscErrorCode ierr;
1379 
1380   PetscFunctionBegin;
1381   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1382   ierr = PetscMalloc1(mg->nlevels,&mat);CHKERRQ(ierr);
1383   for (l=0; l<mg->nlevels-1; l++) {
1384     ierr = KSPGetOperators(mglevels[l]->smoothd,NULL,&(mat[l]));CHKERRQ(ierr);
1385     ierr = PetscObjectReference((PetscObject)mat[l]);CHKERRQ(ierr);
1386   }
1387   *num_levels = mg->nlevels;
1388   *coarseOperators = mat;
1389   PetscFunctionReturn(0);
1390 }
1391 
1392 /*@C
1393   PCMGRegisterCoarseSpaceConstructor -  Adds a method to the PCMG package for coarse space construction.
1394 
1395   Not collective
1396 
1397   Input Parameters:
1398 + name     - name of the constructor
1399 - function - constructor routine
1400 
1401   Notes:
1402   Calling sequence for the routine:
1403 $ my_csp(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, const Vec initGuess[], Vec **coarseSp)
1404 $   pc        - The PC object
1405 $   l         - The multigrid level, 0 is the coarse level
1406 $   dm        - The DM for this level
1407 $   smooth    - The level smoother
1408 $   Nc        - The size of the coarse space
1409 $   initGuess - Basis for an initial guess for the space
1410 $   coarseSp  - A basis for the computed coarse space
1411 
1412   Level: advanced
1413 
1414 .seealso: PCMGGetCoarseSpaceConstructor(), PCRegister()
1415 @*/
1416 PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char name[], PetscErrorCode (*function)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **))
1417 {
1418   PetscErrorCode ierr;
1419 
1420   PetscFunctionBegin;
1421   ierr = PCInitializePackage();CHKERRQ(ierr);
1422   ierr = PetscFunctionListAdd(&PCMGCoarseList,name,function);CHKERRQ(ierr);
1423   PetscFunctionReturn(0);
1424 }
1425 
1426 /*@C
1427   PCMGGetCoarseSpaceConstructor -  Returns the given coarse space construction method.
1428 
1429   Not collective
1430 
1431   Input Parameter:
1432 . name     - name of the constructor
1433 
1434   Output Parameter:
1435 . function - constructor routine
1436 
1437   Notes:
1438   Calling sequence for the routine:
1439 $ my_csp(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, const Vec initGuess[], Vec **coarseSp)
1440 $   pc        - The PC object
1441 $   l         - The multigrid level, 0 is the coarse level
1442 $   dm        - The DM for this level
1443 $   smooth    - The level smoother
1444 $   Nc        - The size of the coarse space
1445 $   initGuess - Basis for an initial guess for the space
1446 $   coarseSp  - A basis for the computed coarse space
1447 
1448   Level: advanced
1449 
1450 .seealso: PCMGRegisterCoarseSpaceConstructor(), PCRegister()
1451 @*/
1452 PetscErrorCode PCMGGetCoarseSpaceConstructor(const char name[], PetscErrorCode (**function)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **))
1453 {
1454   PetscErrorCode ierr;
1455 
1456   PetscFunctionBegin;
1457   ierr = PetscFunctionListFind(PCMGCoarseList,name,function);CHKERRQ(ierr);
1458   PetscFunctionReturn(0);
1459 }
1460 
1461 /* ----------------------------------------------------------------------------------------*/
1462 
1463 /*MC
1464    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
1465     information about the coarser grid matrices and restriction/interpolation operators.
1466 
1467    Options Database Keys:
1468 +  -pc_mg_levels <nlevels> - number of levels including finest
1469 .  -pc_mg_cycle_type <v,w> - provide the cycle desired
1470 .  -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default
1471 .  -pc_mg_log - log information about time spent on each level of the solver
1472 .  -pc_mg_distinct_smoothup - configure up (after interpolation) and down (before restriction) smoothers separately (with different options prefixes)
1473 .  -pc_mg_galerkin <both,pmat,mat,none> - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R'
1474 .  -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1)
1475 .  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1476                         to the Socket viewer for reading from MATLAB.
1477 -  -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices
1478                         to the binary output file called binaryoutput
1479 
1480    Notes:
1481     If one uses a Krylov method such GMRES or CG as the smoother then one must use KSPFGMRES, KSPGCR, or KSPRICHARDSON as the outer Krylov method
1482 
1483        When run with a single level the smoother options are used on that level NOT the coarse grid solver options
1484 
1485        When run with KSPRICHARDSON the convergence test changes slightly if monitor is turned on. The iteration count may change slightly. This
1486        is because without monitoring the residual norm is computed WITHIN each multigrid cycle on the finest level after the pre-smoothing
1487        (because the residual has just been computed for the multigrid algorithm and is hence available for free) while with monitoring the
1488        residual is computed at the end of each cycle.
1489 
1490    Level: intermediate
1491 
1492 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE
1493            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(),
1494            PCMGSetDistinctSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
1495            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
1496            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
1497 M*/
1498 
1499 PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc)
1500 {
1501   PC_MG          *mg;
1502   PetscErrorCode ierr;
1503 
1504   PetscFunctionBegin;
1505   ierr         = PetscNewLog(pc,&mg);CHKERRQ(ierr);
1506   pc->data     = (void*)mg;
1507   mg->nlevels  = -1;
1508   mg->am       = PC_MG_MULTIPLICATIVE;
1509   mg->galerkin = PC_MG_GALERKIN_NONE;
1510   mg->adaptInterpolation = PETSC_FALSE;
1511   mg->Nc                 = -1;
1512   mg->eigenvalue         = -1;
1513 
1514   pc->useAmat = PETSC_TRUE;
1515 
1516   pc->ops->apply          = PCApply_MG;
1517   pc->ops->setup          = PCSetUp_MG;
1518   pc->ops->reset          = PCReset_MG;
1519   pc->ops->destroy        = PCDestroy_MG;
1520   pc->ops->setfromoptions = PCSetFromOptions_MG;
1521   pc->ops->view           = PCView_MG;
1522 
1523   ierr = PetscObjectComposedDataRegister(&mg->eigenvalue);CHKERRQ(ierr);
1524   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",PCMGSetGalerkin_MG);CHKERRQ(ierr);
1525   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGGetLevels_C",PCMGGetLevels_MG);CHKERRQ(ierr);
1526   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetLevels_C",PCMGSetLevels_MG);CHKERRQ(ierr);
1527   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGetInterpolations_C",PCGetInterpolations_MG);CHKERRQ(ierr);
1528   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGetCoarseOperators_C",PCGetCoarseOperators_MG);CHKERRQ(ierr);
1529   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetAdaptInterpolation_C",PCMGSetAdaptInterpolation_MG);CHKERRQ(ierr);
1530   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGGetAdaptInterpolation_C",PCMGGetAdaptInterpolation_MG);CHKERRQ(ierr);
1531   PetscFunctionReturn(0);
1532 }
1533