xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision 5adeb434db1cbecfaacaa61eb22541e5543bf6b2)
1dba47a55SKris Buschelman 
24b9ad928SBarry Smith /*
34b9ad928SBarry Smith     Defines the multigrid preconditioner interface.
44b9ad928SBarry Smith */
51996caf8SToby Isaac #include <petsc-private/pcmgimpl.h>                    /*I "petscksp.h" I*/
61e25c274SJed Brown #include <petscdm.h>
74b9ad928SBarry Smith 
84b9ad928SBarry Smith #undef __FUNCT__
99dcbbd2bSBarry Smith #define __FUNCT__ "PCMGMCycle_Private"
1031567311SBarry Smith PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG_Levels **mglevelsin,PCRichardsonConvergedReason *reason)
114b9ad928SBarry Smith {
1231567311SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
1331567311SBarry Smith   PC_MG_Levels   *mgc,*mglevels = *mglevelsin;
146849ba73SBarry Smith   PetscErrorCode ierr;
1531567311SBarry Smith   PetscInt       cycles = (mglevels->level == 1) ? 1 : (PetscInt) mglevels->cycles;
164b9ad928SBarry Smith 
174b9ad928SBarry Smith   PetscFunctionBegin;
1863e6d426SJed Brown   if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
1931567311SBarry Smith   ierr = KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x);CHKERRQ(ierr);  /* pre-smooth */
2063e6d426SJed Brown   if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
2131567311SBarry Smith   if (mglevels->level) {  /* not the coarsest grid */
2263e6d426SJed Brown     if (mglevels->eventresidual) {ierr = PetscLogEventBegin(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);}
2331567311SBarry Smith     ierr = (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);CHKERRQ(ierr);
2463e6d426SJed Brown     if (mglevels->eventresidual) {ierr = PetscLogEventEnd(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);}
254b9ad928SBarry Smith 
264b9ad928SBarry Smith     /* if on finest level and have convergence criteria set */
2731567311SBarry Smith     if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) {
284b9ad928SBarry Smith       PetscReal rnorm;
2931567311SBarry Smith       ierr = VecNorm(mglevels->r,NORM_2,&rnorm);CHKERRQ(ierr);
304b9ad928SBarry Smith       if (rnorm <= mg->ttol) {
3170441072SBarry Smith         if (rnorm < mg->abstol) {
324d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_ATOL;
3357622a8eSBarry Smith           ierr    = PetscInfo2(pc,"Linear solver has converged. Residual norm %g is less than absolute tolerance %g\n",(double)rnorm,(double)mg->abstol);CHKERRQ(ierr);
344b9ad928SBarry Smith         } else {
354d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_RTOL;
3657622a8eSBarry Smith           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);
374b9ad928SBarry Smith         }
384b9ad928SBarry Smith         PetscFunctionReturn(0);
394b9ad928SBarry Smith       }
404b9ad928SBarry Smith     }
414b9ad928SBarry Smith 
4231567311SBarry Smith     mgc = *(mglevelsin - 1);
4363e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
4431567311SBarry Smith     ierr = MatRestrict(mglevels->restrct,mglevels->r,mgc->b);CHKERRQ(ierr);
4563e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
46efb30889SBarry Smith     ierr = VecSet(mgc->x,0.0);CHKERRQ(ierr);
474b9ad928SBarry Smith     while (cycles--) {
4831567311SBarry Smith       ierr = PCMGMCycle_Private(pc,mglevelsin-1,reason);CHKERRQ(ierr);
494b9ad928SBarry Smith     }
5063e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
5131567311SBarry Smith     ierr = MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x);CHKERRQ(ierr);
5263e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
5363e6d426SJed Brown     if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
5431567311SBarry Smith     ierr = KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x);CHKERRQ(ierr);    /* post smooth */
5563e6d426SJed Brown     if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
564b9ad928SBarry Smith   }
574b9ad928SBarry Smith   PetscFunctionReturn(0);
584b9ad928SBarry Smith }
594b9ad928SBarry Smith 
604b9ad928SBarry Smith #undef __FUNCT__
614b9ad928SBarry Smith #define __FUNCT__ "PCApplyRichardson_MG"
62ace3abfcSBarry Smith 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)
634b9ad928SBarry Smith {
64f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
65f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
66dfbe8321SBarry Smith   PetscErrorCode ierr;
67f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
684b9ad928SBarry Smith 
694b9ad928SBarry Smith   PetscFunctionBegin;
70a762d673SBarry Smith   /* When the DM is supplying the matrix then it will not exist until here */
71a762d673SBarry Smith   for (i=0; i<levels; i++) {
72a762d673SBarry Smith     if (!mglevels[i]->A) {
73a762d673SBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr);
74a762d673SBarry Smith       ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr);
75a762d673SBarry Smith     }
76a762d673SBarry Smith   }
77f3fbd535SBarry Smith   mglevels[levels-1]->b = b;
78f3fbd535SBarry Smith   mglevels[levels-1]->x = x;
794b9ad928SBarry Smith 
8031567311SBarry Smith   mg->rtol   = rtol;
8131567311SBarry Smith   mg->abstol = abstol;
8231567311SBarry Smith   mg->dtol   = dtol;
834b9ad928SBarry Smith   if (rtol) {
844b9ad928SBarry Smith     /* compute initial residual norm for relative convergence test */
854b9ad928SBarry Smith     PetscReal rnorm;
867319c654SBarry Smith     if (zeroguess) {
877319c654SBarry Smith       ierr = VecNorm(b,NORM_2,&rnorm);CHKERRQ(ierr);
887319c654SBarry Smith     } else {
89f3fbd535SBarry Smith       ierr = (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);CHKERRQ(ierr);
904b9ad928SBarry Smith       ierr = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr);
917319c654SBarry Smith     }
9231567311SBarry Smith     mg->ttol = PetscMax(rtol*rnorm,abstol);
932fa5cd67SKarl Rupp   } else if (abstol) mg->ttol = abstol;
942fa5cd67SKarl Rupp   else mg->ttol = 0.0;
954b9ad928SBarry Smith 
964d0a8057SBarry Smith   /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't
97336babb1SJed Brown      stop prematurely due to small residual */
984d0a8057SBarry Smith   for (i=1; i<levels; i++) {
99f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
100f3fbd535SBarry Smith     if (mglevels[i]->smoothu != mglevels[i]->smoothd) {
101f3fbd535SBarry Smith       ierr = KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
1024b9ad928SBarry Smith     }
1034d0a8057SBarry Smith   }
1044d0a8057SBarry Smith 
1054d0a8057SBarry Smith   *reason = (PCRichardsonConvergedReason)0;
1064d0a8057SBarry Smith   for (i=0; i<its; i++) {
107f3fbd535SBarry Smith     ierr = PCMGMCycle_Private(pc,mglevels+levels-1,reason);CHKERRQ(ierr);
1084d0a8057SBarry Smith     if (*reason) break;
1094d0a8057SBarry Smith   }
1104d0a8057SBarry Smith   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
1114d0a8057SBarry Smith   *outits = i;
1124b9ad928SBarry Smith   PetscFunctionReturn(0);
1134b9ad928SBarry Smith }
1144b9ad928SBarry Smith 
1154b9ad928SBarry Smith #undef __FUNCT__
1163aeaf226SBarry Smith #define __FUNCT__ "PCReset_MG"
1173aeaf226SBarry Smith PetscErrorCode PCReset_MG(PC pc)
1183aeaf226SBarry Smith {
1193aeaf226SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
1203aeaf226SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
1213aeaf226SBarry Smith   PetscErrorCode ierr;
1223aeaf226SBarry Smith   PetscInt       i,n;
1233aeaf226SBarry Smith 
1243aeaf226SBarry Smith   PetscFunctionBegin;
1253aeaf226SBarry Smith   if (mglevels) {
1263aeaf226SBarry Smith     n = mglevels[0]->levels;
1273aeaf226SBarry Smith     for (i=0; i<n-1; i++) {
1283aeaf226SBarry Smith       ierr = VecDestroy(&mglevels[i+1]->r);CHKERRQ(ierr);
1293aeaf226SBarry Smith       ierr = VecDestroy(&mglevels[i]->b);CHKERRQ(ierr);
1303aeaf226SBarry Smith       ierr = VecDestroy(&mglevels[i]->x);CHKERRQ(ierr);
1313aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i+1]->restrct);CHKERRQ(ierr);
1323aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i+1]->interpolate);CHKERRQ(ierr);
13373250ac0SBarry Smith       ierr = VecDestroy(&mglevels[i+1]->rscale);CHKERRQ(ierr);
1343aeaf226SBarry Smith     }
1353aeaf226SBarry Smith 
1363aeaf226SBarry Smith     for (i=0; i<n; i++) {
1373aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i]->A);CHKERRQ(ierr);
1383aeaf226SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
1393aeaf226SBarry Smith         ierr = KSPReset(mglevels[i]->smoothd);CHKERRQ(ierr);
1403aeaf226SBarry Smith       }
1413aeaf226SBarry Smith       ierr = KSPReset(mglevels[i]->smoothu);CHKERRQ(ierr);
1423aeaf226SBarry Smith     }
1433aeaf226SBarry Smith   }
1443aeaf226SBarry Smith   PetscFunctionReturn(0);
1453aeaf226SBarry Smith }
1463aeaf226SBarry Smith 
1473aeaf226SBarry Smith #undef __FUNCT__
1489dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetLevels"
1494b9ad928SBarry Smith /*@C
15097177400SBarry Smith    PCMGSetLevels - Sets the number of levels to use with MG.
1514b9ad928SBarry Smith    Must be called before any other MG routine.
1524b9ad928SBarry Smith 
153ad4df100SBarry Smith    Logically Collective on PC
1544b9ad928SBarry Smith 
1554b9ad928SBarry Smith    Input Parameters:
1564b9ad928SBarry Smith +  pc - the preconditioner context
1574b9ad928SBarry Smith .  levels - the number of levels
1584b9ad928SBarry Smith -  comms - optional communicators for each level; this is to allow solving the coarser problems
1590298fd71SBarry Smith            on smaller sets of processors. Use NULL_OBJECT for default in Fortran
1604b9ad928SBarry Smith 
1614b9ad928SBarry Smith    Level: intermediate
1624b9ad928SBarry Smith 
1634b9ad928SBarry Smith    Notes:
1644b9ad928SBarry Smith      If the number of levels is one then the multigrid uses the -mg_levels prefix
1654b9ad928SBarry Smith   for setting the level options rather than the -mg_coarse prefix.
1664b9ad928SBarry Smith 
1674b9ad928SBarry Smith .keywords: MG, set, levels, multigrid
1684b9ad928SBarry Smith 
16997177400SBarry Smith .seealso: PCMGSetType(), PCMGGetLevels()
1704b9ad928SBarry Smith @*/
1717087cfbeSBarry Smith PetscErrorCode  PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
1724b9ad928SBarry Smith {
173dfbe8321SBarry Smith   PetscErrorCode ierr;
174f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
175ce94432eSBarry Smith   MPI_Comm       comm;
1763aeaf226SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
177f3fbd535SBarry Smith   PetscInt       i;
178f3fbd535SBarry Smith   PetscMPIInt    size;
179f3fbd535SBarry Smith   const char     *prefix;
180f3fbd535SBarry Smith   PC             ipc;
1813aeaf226SBarry Smith   PetscInt       n;
1824b9ad928SBarry Smith 
1834b9ad928SBarry Smith   PetscFunctionBegin;
1840700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
185548767e0SJed Brown   PetscValidLogicalCollectiveInt(pc,levels,2);
186ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
187548767e0SJed Brown   if (mg->nlevels == levels) PetscFunctionReturn(0);
1883aeaf226SBarry Smith   if (mglevels) {
1893aeaf226SBarry Smith     /* changing the number of levels so free up the previous stuff */
1903aeaf226SBarry Smith     ierr = PCReset_MG(pc);CHKERRQ(ierr);
1913aeaf226SBarry Smith     n    = mglevels[0]->levels;
1923aeaf226SBarry Smith     for (i=0; i<n; i++) {
1933aeaf226SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
1943aeaf226SBarry Smith         ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
1953aeaf226SBarry Smith       }
1963aeaf226SBarry Smith       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
1973aeaf226SBarry Smith       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
1983aeaf226SBarry Smith     }
1993aeaf226SBarry Smith     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
2003aeaf226SBarry Smith   }
201f3fbd535SBarry Smith 
202f3fbd535SBarry Smith   mg->nlevels = levels;
203f3fbd535SBarry Smith 
204785e854fSJed Brown   ierr = PetscMalloc1(levels,&mglevels);CHKERRQ(ierr);
2053bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr);
206f3fbd535SBarry Smith 
207f3fbd535SBarry Smith   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
208f3fbd535SBarry Smith 
209a9db87a2SMatthew G Knepley   mg->stageApply = 0;
210f3fbd535SBarry Smith   for (i=0; i<levels; i++) {
211b00a9115SJed Brown     ierr = PetscNewLog(pc,&mglevels[i]);CHKERRQ(ierr);
2122fa5cd67SKarl Rupp 
213f3fbd535SBarry Smith     mglevels[i]->level               = i;
214f3fbd535SBarry Smith     mglevels[i]->levels              = levels;
215f3fbd535SBarry Smith     mglevels[i]->cycles              = PC_MG_CYCLE_V;
216336babb1SJed Brown     mg->default_smoothu              = 2;
217336babb1SJed Brown     mg->default_smoothd              = 2;
21863e6d426SJed Brown     mglevels[i]->eventsmoothsetup    = 0;
21963e6d426SJed Brown     mglevels[i]->eventsmoothsolve    = 0;
22063e6d426SJed Brown     mglevels[i]->eventresidual       = 0;
22163e6d426SJed Brown     mglevels[i]->eventinterprestrict = 0;
222f3fbd535SBarry Smith 
223f3fbd535SBarry Smith     if (comms) comm = comms[i];
224f3fbd535SBarry Smith     ierr = KSPCreate(comm,&mglevels[i]->smoothd);CHKERRQ(ierr);
225336babb1SJed Brown     ierr = KSPSetType(mglevels[i]->smoothd,KSPCHEBYSHEV);CHKERRQ(ierr);
2260059c7bdSJed Brown     ierr = KSPSetConvergenceTest(mglevels[i]->smoothd,KSPConvergedSkip,NULL,NULL);CHKERRQ(ierr);
227336babb1SJed Brown     ierr = KSPSetNormType(mglevels[i]->smoothd,KSP_NORM_NONE);CHKERRQ(ierr);
228336babb1SJed Brown     ierr = KSPGetPC(mglevels[i]->smoothd,&ipc);CHKERRQ(ierr);
229336babb1SJed Brown     ierr = PCSetType(ipc,PCSOR);CHKERRQ(ierr);
230f3fbd535SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr);
231336babb1SJed Brown     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, i ? mg->default_smoothd : 1);CHKERRQ(ierr);
232f3fbd535SBarry Smith     ierr = KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);CHKERRQ(ierr);
233f3fbd535SBarry Smith 
234f3fbd535SBarry Smith     /* do special stuff for coarse grid */
235f3fbd535SBarry Smith     if (!i && levels > 1) {
236f3fbd535SBarry Smith       ierr = KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");CHKERRQ(ierr);
237f3fbd535SBarry Smith 
2389dbfc187SHong Zhang       /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */
239f3fbd535SBarry Smith       ierr = KSPSetType(mglevels[0]->smoothd,KSPPREONLY);CHKERRQ(ierr);
240f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[0]->smoothd,&ipc);CHKERRQ(ierr);
241f3fbd535SBarry Smith       ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
242f3fbd535SBarry Smith       if (size > 1) {
2439dbfc187SHong Zhang         KSP innerksp;
2449dbfc187SHong Zhang         PC  innerpc;
245f3fbd535SBarry Smith         ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr);
2469dbfc187SHong Zhang         ierr = PCRedundantGetKSP(ipc,&innerksp);CHKERRQ(ierr);
2479dbfc187SHong Zhang         ierr = KSPGetPC(innerksp,&innerpc);CHKERRQ(ierr);
2489dbfc187SHong Zhang         ierr = PCFactorSetShiftType(innerpc,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
249f3fbd535SBarry Smith       } else {
250f3fbd535SBarry Smith         ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr);
2519dbfc187SHong Zhang         ierr = PCFactorSetShiftType(ipc,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
252f3fbd535SBarry Smith       }
253f3fbd535SBarry Smith     } else {
254f3fbd535SBarry Smith       char tprefix[128];
255f3fbd535SBarry Smith       sprintf(tprefix,"mg_levels_%d_",(int)i);
256f3fbd535SBarry Smith       ierr = KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);CHKERRQ(ierr);
257f3fbd535SBarry Smith     }
2583bb1ff40SBarry Smith     ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[i]->smoothd);CHKERRQ(ierr);
2592fa5cd67SKarl Rupp 
260f3fbd535SBarry Smith     mglevels[i]->smoothu = mglevels[i]->smoothd;
26131567311SBarry Smith     mg->rtol             = 0.0;
26231567311SBarry Smith     mg->abstol           = 0.0;
26331567311SBarry Smith     mg->dtol             = 0.0;
26431567311SBarry Smith     mg->ttol             = 0.0;
26531567311SBarry Smith     mg->cyclesperpcapply = 1;
266f3fbd535SBarry Smith   }
26731567311SBarry Smith   mg->am                   = PC_MG_MULTIPLICATIVE;
268f3fbd535SBarry Smith   mg->levels               = mglevels;
2694b9ad928SBarry Smith   pc->ops->applyrichardson = PCApplyRichardson_MG;
2704b9ad928SBarry Smith   PetscFunctionReturn(0);
2714b9ad928SBarry Smith }
2724b9ad928SBarry Smith 
273c07bf074SBarry Smith 
274c07bf074SBarry Smith #undef __FUNCT__
275c07bf074SBarry Smith #define __FUNCT__ "PCDestroy_MG"
276c07bf074SBarry Smith PetscErrorCode PCDestroy_MG(PC pc)
277c07bf074SBarry Smith {
278c07bf074SBarry Smith   PetscErrorCode ierr;
279a06653b4SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
280a06653b4SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
281a06653b4SBarry Smith   PetscInt       i,n;
282c07bf074SBarry Smith 
283c07bf074SBarry Smith   PetscFunctionBegin;
284a06653b4SBarry Smith   ierr = PCReset_MG(pc);CHKERRQ(ierr);
285a06653b4SBarry Smith   if (mglevels) {
286a06653b4SBarry Smith     n = mglevels[0]->levels;
287a06653b4SBarry Smith     for (i=0; i<n; i++) {
288a06653b4SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
2896bf464f9SBarry Smith         ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
290a06653b4SBarry Smith       }
2916bf464f9SBarry Smith       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
292a06653b4SBarry Smith       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
293a06653b4SBarry Smith     }
294a06653b4SBarry Smith     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
295a06653b4SBarry Smith   }
296c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
297f3fbd535SBarry Smith   PetscFunctionReturn(0);
298f3fbd535SBarry Smith }
299f3fbd535SBarry Smith 
300f3fbd535SBarry Smith 
301f3fbd535SBarry Smith 
30209573ac7SBarry Smith extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**);
30309573ac7SBarry Smith extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**);
30409573ac7SBarry Smith extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**);
305f3fbd535SBarry Smith 
306f3fbd535SBarry Smith /*
307f3fbd535SBarry Smith    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
308f3fbd535SBarry Smith              or full cycle of multigrid.
309f3fbd535SBarry Smith 
310f3fbd535SBarry Smith   Note:
311f3fbd535SBarry Smith   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
312f3fbd535SBarry Smith */
313f3fbd535SBarry Smith #undef __FUNCT__
314f3fbd535SBarry Smith #define __FUNCT__ "PCApply_MG"
315f3fbd535SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
316f3fbd535SBarry Smith {
317f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
318f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
319f3fbd535SBarry Smith   PetscErrorCode ierr;
320f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
321f3fbd535SBarry Smith 
322f3fbd535SBarry Smith   PetscFunctionBegin;
323a9db87a2SMatthew G Knepley   if (mg->stageApply) {ierr = PetscLogStagePush(mg->stageApply);CHKERRQ(ierr);}
324e1d8e5deSBarry Smith   /* When the DM is supplying the matrix then it will not exist until here */
325298cc208SBarry Smith   for (i=0; i<levels; i++) {
326e1d8e5deSBarry Smith     if (!mglevels[i]->A) {
32723ee1639SBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr);
328298cc208SBarry Smith       ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr);
329e1d8e5deSBarry Smith     }
330e1d8e5deSBarry Smith   }
331e1d8e5deSBarry Smith 
332f3fbd535SBarry Smith   mglevels[levels-1]->b = b;
333f3fbd535SBarry Smith   mglevels[levels-1]->x = x;
33431567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
335f3fbd535SBarry Smith     ierr = VecSet(x,0.0);CHKERRQ(ierr);
33631567311SBarry Smith     for (i=0; i<mg->cyclesperpcapply; i++) {
3370298fd71SBarry Smith       ierr = PCMGMCycle_Private(pc,mglevels+levels-1,NULL);CHKERRQ(ierr);
338f3fbd535SBarry Smith     }
3392fa5cd67SKarl Rupp   } else if (mg->am == PC_MG_ADDITIVE) {
34031567311SBarry Smith     ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr);
3412fa5cd67SKarl Rupp   } else if (mg->am == PC_MG_KASKADE) {
34231567311SBarry Smith     ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr);
3432fa5cd67SKarl Rupp   } else {
344f3fbd535SBarry Smith     ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr);
345f3fbd535SBarry Smith   }
346a9db87a2SMatthew G Knepley   if (mg->stageApply) {ierr = PetscLogStagePop();CHKERRQ(ierr);}
347f3fbd535SBarry Smith   PetscFunctionReturn(0);
348f3fbd535SBarry Smith }
349f3fbd535SBarry Smith 
350f3fbd535SBarry Smith 
351f3fbd535SBarry Smith #undef __FUNCT__
352f3fbd535SBarry Smith #define __FUNCT__ "PCSetFromOptions_MG"
3538c34d3f5SBarry Smith PetscErrorCode PCSetFromOptions_MG(PetscOptions *PetscOptionsObject,PC pc)
354f3fbd535SBarry Smith {
355f3fbd535SBarry Smith   PetscErrorCode ierr;
356f3fbd535SBarry Smith   PetscInt       m,levels = 1,cycles;
357c91913d3SJed Brown   PetscBool      flg,set;
358f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
359f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
360f3fbd535SBarry Smith   PCMGType       mgtype;
361f3fbd535SBarry Smith   PCMGCycleType  mgctype;
362f3fbd535SBarry Smith 
363f3fbd535SBarry Smith   PetscFunctionBegin;
364e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Multigrid options");CHKERRQ(ierr);
3653aeaf226SBarry Smith   if (!mg->levels) {
366f3fbd535SBarry Smith     ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr);
367eb3f98d2SBarry Smith     if (!flg && pc->dm) {
368eb3f98d2SBarry Smith       ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
369eb3f98d2SBarry Smith       levels++;
3703aeaf226SBarry Smith       mg->usedmfornumberoflevels = PETSC_TRUE;
371eb3f98d2SBarry Smith     }
3720298fd71SBarry Smith     ierr = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr);
373f3fbd535SBarry Smith   }
3743aeaf226SBarry Smith   mglevels = mg->levels;
3753aeaf226SBarry Smith 
376f3fbd535SBarry Smith   mgctype = (PCMGCycleType) mglevels[0]->cycles;
377f3fbd535SBarry Smith   ierr    = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr);
378f3fbd535SBarry Smith   if (flg) {
379f3fbd535SBarry Smith     ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr);
3802fa5cd67SKarl Rupp   }
381f3fbd535SBarry Smith   flg  = PETSC_FALSE;
382c91913d3SJed Brown   ierr = PetscOptionsBool("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",flg,&flg,&set);CHKERRQ(ierr);
383c91913d3SJed Brown   if (set) {
384c91913d3SJed Brown     ierr = PCMGSetGalerkin(pc,flg);CHKERRQ(ierr);
385f3fbd535SBarry Smith   }
38694ae4db5SBarry Smith   ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",mg->default_smoothu,&m,&flg);CHKERRQ(ierr);
387f3fbd535SBarry Smith   if (flg) {
388f3fbd535SBarry Smith     ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr);
389f3fbd535SBarry Smith   }
39094ae4db5SBarry Smith   ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",mg->default_smoothd,&m,&flg);CHKERRQ(ierr);
391f3fbd535SBarry Smith   if (flg) {
392f3fbd535SBarry Smith     ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr);
393f3fbd535SBarry Smith   }
39431567311SBarry Smith   mgtype = mg->am;
395f3fbd535SBarry Smith   ierr   = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
396f3fbd535SBarry Smith   if (flg) {
397f3fbd535SBarry Smith     ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
398f3fbd535SBarry Smith   }
39931567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
40031567311SBarry Smith     ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGSetLevels",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr);
401f3fbd535SBarry Smith     if (flg) {
402f3fbd535SBarry Smith       ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr);
403f3fbd535SBarry Smith     }
404f3fbd535SBarry Smith   }
405f3fbd535SBarry Smith   flg  = PETSC_FALSE;
4060298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,NULL);CHKERRQ(ierr);
407f3fbd535SBarry Smith   if (flg) {
408f3fbd535SBarry Smith     PetscInt i;
409f3fbd535SBarry Smith     char     eventname[128];
410ce94432eSBarry Smith     if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
411f3fbd535SBarry Smith     levels = mglevels[0]->levels;
412f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
413f3fbd535SBarry Smith       sprintf(eventname,"MGSetup Level %d",(int)i);
41463e6d426SJed Brown       ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr);
415f3fbd535SBarry Smith       sprintf(eventname,"MGSmooth Level %d",(int)i);
41663e6d426SJed Brown       ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr);
417f3fbd535SBarry Smith       if (i) {
418f3fbd535SBarry Smith         sprintf(eventname,"MGResid Level %d",(int)i);
41963e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr);
420f3fbd535SBarry Smith         sprintf(eventname,"MGInterp Level %d",(int)i);
42163e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr);
422f3fbd535SBarry Smith       }
423f3fbd535SBarry Smith     }
424eec35531SMatthew G Knepley 
425ec60ca73SSatish Balay #if defined(PETSC_USE_LOG)
426eec35531SMatthew G Knepley     {
427eec35531SMatthew G Knepley       const char    *sname = "MG Apply";
428eec35531SMatthew G Knepley       PetscStageLog stageLog;
429eec35531SMatthew G Knepley       PetscInt      st;
430eec35531SMatthew G Knepley 
431eec35531SMatthew G Knepley       PetscFunctionBegin;
432eec35531SMatthew G Knepley       ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
433eec35531SMatthew G Knepley       for (st = 0; st < stageLog->numStages; ++st) {
434eec35531SMatthew G Knepley         PetscBool same;
435eec35531SMatthew G Knepley 
436eec35531SMatthew G Knepley         ierr = PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);CHKERRQ(ierr);
4372fa5cd67SKarl Rupp         if (same) mg->stageApply = st;
438eec35531SMatthew G Knepley       }
439eec35531SMatthew G Knepley       if (!mg->stageApply) {
440eec35531SMatthew G Knepley         ierr = PetscLogStageRegister(sname, &mg->stageApply);CHKERRQ(ierr);
441eec35531SMatthew G Knepley       }
442eec35531SMatthew G Knepley     }
443ec60ca73SSatish Balay #endif
444f3fbd535SBarry Smith   }
445f3fbd535SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
446f3fbd535SBarry Smith   PetscFunctionReturn(0);
447f3fbd535SBarry Smith }
448f3fbd535SBarry Smith 
4496a6fc655SJed Brown const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0};
4506a6fc655SJed Brown const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0};
451f3fbd535SBarry Smith 
4529804daf3SBarry Smith #include <petscdraw.h>
453f3fbd535SBarry Smith #undef __FUNCT__
454f3fbd535SBarry Smith #define __FUNCT__ "PCView_MG"
455f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
456f3fbd535SBarry Smith {
457f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
458f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
459f3fbd535SBarry Smith   PetscErrorCode ierr;
460e3deeeafSJed Brown   PetscInt       levels = mglevels ? mglevels[0]->levels : 0,i;
461219b1045SBarry Smith   PetscBool      iascii,isbinary,isdraw;
462f3fbd535SBarry Smith 
463f3fbd535SBarry Smith   PetscFunctionBegin;
464251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
4655b0b0462SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
466219b1045SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
467f3fbd535SBarry Smith   if (iascii) {
468e3deeeafSJed Brown     const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown";
469e3deeeafSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  MG: type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,cyclename);CHKERRQ(ierr);
47031567311SBarry Smith     if (mg->am == PC_MG_MULTIPLICATIVE) {
47131567311SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr);
472f3fbd535SBarry Smith     }
473218a07d4SBarry Smith     if (mg->galerkin) {
474f3fbd535SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
4754f66f45eSBarry Smith     } else {
4764f66f45eSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
477f3fbd535SBarry Smith     }
478*5adeb434SBarry Smith     if (mg->view){
479*5adeb434SBarry Smith       ierr = (*mg->view)(pc,viewer);CHKERRQ(ierr);
480*5adeb434SBarry Smith     }
481f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
482f3fbd535SBarry Smith       if (!i) {
48389cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);CHKERRQ(ierr);
484f3fbd535SBarry Smith       } else {
48589cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
486f3fbd535SBarry Smith       }
487f3fbd535SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
488f3fbd535SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
489f3fbd535SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
490f3fbd535SBarry Smith       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
491f3fbd535SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
492f3fbd535SBarry Smith       } else if (i) {
49389cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
494f3fbd535SBarry Smith         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
495f3fbd535SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
496f3fbd535SBarry Smith         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
497f3fbd535SBarry Smith       }
498f3fbd535SBarry Smith     }
4995b0b0462SBarry Smith   } else if (isbinary) {
5005b0b0462SBarry Smith     for (i=levels-1; i>=0; i--) {
5015b0b0462SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
5025b0b0462SBarry Smith       if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) {
5035b0b0462SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
5045b0b0462SBarry Smith       }
5055b0b0462SBarry Smith     }
506219b1045SBarry Smith   } else if (isdraw) {
507219b1045SBarry Smith     PetscDraw draw;
508b4375e8dSPeter Brune     PetscReal x,w,y,bottom,th;
509219b1045SBarry Smith     ierr   = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
510219b1045SBarry Smith     ierr   = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
5110298fd71SBarry Smith     ierr   = PetscDrawStringGetSize(draw,NULL,&th);CHKERRQ(ierr);
512219b1045SBarry Smith     bottom = y - th;
513219b1045SBarry Smith     for (i=levels-1; i>=0; i--) {
514b4375e8dSPeter Brune       if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) {
515219b1045SBarry Smith         ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr);
516219b1045SBarry Smith         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
517219b1045SBarry Smith         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
518b4375e8dSPeter Brune       } else {
519b4375e8dSPeter Brune         w    = 0.5*PetscMin(1.0-x,x);
520b4375e8dSPeter Brune         ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr);
521b4375e8dSPeter Brune         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
522b4375e8dSPeter Brune         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
523b4375e8dSPeter Brune         ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr);
524b4375e8dSPeter Brune         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
525b4375e8dSPeter Brune         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
526b4375e8dSPeter Brune       }
5270298fd71SBarry Smith       ierr    = PetscDrawGetBoundingBox(draw,NULL,&bottom,NULL,NULL);CHKERRQ(ierr);
5281cd381d6SBarry Smith       bottom -= th;
529219b1045SBarry Smith     }
5305b0b0462SBarry Smith   }
531f3fbd535SBarry Smith   PetscFunctionReturn(0);
532f3fbd535SBarry Smith }
533f3fbd535SBarry Smith 
534b45d2f2cSJed Brown #include <petsc-private/dmimpl.h>
535b45d2f2cSJed Brown #include <petsc-private/kspimpl.h>
536cab2e9ccSBarry Smith 
537f3fbd535SBarry Smith /*
538f3fbd535SBarry Smith     Calls setup for the KSP on each level
539f3fbd535SBarry Smith */
540f3fbd535SBarry Smith #undef __FUNCT__
541f3fbd535SBarry Smith #define __FUNCT__ "PCSetUp_MG"
542f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc)
543f3fbd535SBarry Smith {
544f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
545f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
546f3fbd535SBarry Smith   PetscErrorCode ierr;
547f3fbd535SBarry Smith   PetscInt       i,n = mglevels[0]->levels;
54898e04984SBarry Smith   PC             cpc;
54971b23a65SMark F. Adams   PetscBool      preonly,lu,redundant,cholesky,svd,dump = PETSC_FALSE,opsset,use_amat;
550f3fbd535SBarry Smith   Mat            dA,dB;
551f3fbd535SBarry Smith   Vec            tvec;
552218a07d4SBarry Smith   DM             *dms;
553649052a6SBarry Smith   PetscViewer    viewer = 0;
554f3fbd535SBarry Smith 
555f3fbd535SBarry Smith   PetscFunctionBegin;
55601bc414fSDmitry Karpeev   /* FIX: Move this to PCSetFromOptions_MG? */
5573aeaf226SBarry Smith   if (mg->usedmfornumberoflevels) {
5583aeaf226SBarry Smith     PetscInt levels;
5593aeaf226SBarry Smith     ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
5603aeaf226SBarry Smith     levels++;
5613aeaf226SBarry Smith     if (levels > n) { /* the problem is now being solved on a finer grid */
5620298fd71SBarry Smith       ierr     = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr);
5633aeaf226SBarry Smith       n        = levels;
5643aeaf226SBarry Smith       ierr     = PCSetFromOptions(pc);CHKERRQ(ierr); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
5653aeaf226SBarry Smith       mglevels =  mg->levels;
5663aeaf226SBarry Smith     }
5673aeaf226SBarry Smith   }
56898e04984SBarry Smith   ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr);
5693aeaf226SBarry Smith 
570f3fbd535SBarry Smith 
571f3fbd535SBarry Smith   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
572f3fbd535SBarry Smith   /* so use those from global PC */
573f3fbd535SBarry Smith   /* Is this what we always want? What if user wants to keep old one? */
5740298fd71SBarry Smith   ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset);CHKERRQ(ierr);
57598e04984SBarry Smith   if (opsset) {
57698e04984SBarry Smith     Mat mmat;
57723ee1639SBarry Smith     ierr = KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat);CHKERRQ(ierr);
57898e04984SBarry Smith     if (mmat == pc->pmat) opsset = PETSC_FALSE;
57998e04984SBarry Smith   }
5804a5f07a5SMark F. Adams 
5814a5f07a5SMark F. Adams   if (!opsset) {
58271b23a65SMark F. Adams     ierr = PCGetUseAmat(pc,&use_amat);CHKERRQ(ierr);
5834a5f07a5SMark F. Adams     if(use_amat){
584f3fbd535SBarry Smith       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);
58523ee1639SBarry Smith       ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat);CHKERRQ(ierr);
586f3fbd535SBarry Smith     }
5874a5f07a5SMark F. Adams     else {
5884a5f07a5SMark F. Adams       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);
58923ee1639SBarry Smith       ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat);CHKERRQ(ierr);
5904a5f07a5SMark F. Adams     }
5914a5f07a5SMark F. Adams   }
592f3fbd535SBarry Smith 
59301bc414fSDmitry Karpeev   /* Skipping this for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs? */
594ce4cda84SJed Brown   if (pc->dm && mg->galerkin != 2 && !pc->setupcalled) {
5952d2b81a6SBarry Smith     /* construct the interpolation from the DMs */
5962e499ae9SBarry Smith     Mat p;
59773250ac0SBarry Smith     Vec rscale;
598785e854fSJed Brown     ierr     = PetscMalloc1(n,&dms);CHKERRQ(ierr);
599218a07d4SBarry Smith     dms[n-1] = pc->dm;
600ef1267afSMatthew G. Knepley     /* Separately create them so we do not get DMKSP interference between levels */
601ef1267afSMatthew G. Knepley     for (i=n-2; i>-1; i--) {ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr);}
602218a07d4SBarry Smith     for (i=n-2; i>-1; i--) {
603942e3340SBarry Smith       DMKSP kdm;
6043c0c59f3SBarry Smith       ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr);
6053c0c59f3SBarry Smith       if (mg->galerkin) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);}
606942e3340SBarry Smith       ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr);
607d1a3e677SJed Brown       /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take
608d1a3e677SJed Brown        * a bitwise OR of computing the matrix, RHS, and initial iterate. */
6090298fd71SBarry Smith       kdm->ops->computerhs = NULL;
6100298fd71SBarry Smith       kdm->rhsctx          = NULL;
61124c3aa18SBarry Smith       if (!mglevels[i+1]->interpolate) {
612e727c939SJed Brown         ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr);
6132d2b81a6SBarry Smith         ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr);
61400ac8be1SBarry Smith         if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);}
61573250ac0SBarry Smith         ierr = VecDestroy(&rscale);CHKERRQ(ierr);
6166bf464f9SBarry Smith         ierr = MatDestroy(&p);CHKERRQ(ierr);
617218a07d4SBarry Smith       }
61824c3aa18SBarry Smith     }
6192d2b81a6SBarry Smith 
620ef1267afSMatthew G. Knepley     for (i=n-2; i>-1; i--) {ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);}
6212d2b81a6SBarry Smith     ierr = PetscFree(dms);CHKERRQ(ierr);
622ce4cda84SJed Brown   }
623cab2e9ccSBarry Smith 
624ce4cda84SJed Brown   if (pc->dm && !pc->setupcalled) {
625ce4cda84SJed Brown     /* finest smoother also gets DM but it is not active, independent of whether galerkin==2 */
626cab2e9ccSBarry Smith     ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr);
627cab2e9ccSBarry Smith     ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr);
628218a07d4SBarry Smith   }
629218a07d4SBarry Smith 
630c91913d3SJed Brown   if (mg->galerkin == 1) {
631f3fbd535SBarry Smith     Mat B;
632f3fbd535SBarry Smith     /* currently only handle case where mat and pmat are the same on coarser levels */
63323ee1639SBarry Smith     ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB);CHKERRQ(ierr);
634f3fbd535SBarry Smith     if (!pc->setupcalled) {
635f3fbd535SBarry Smith       for (i=n-2; i>-1; i--) {
636b9367914SBarry Smith         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");
637b9367914SBarry Smith         if (!mglevels[i+1]->interpolate) {
638b9367914SBarry Smith           ierr = PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct);CHKERRQ(ierr);
639b9367914SBarry Smith         }
640b9367914SBarry Smith         if (!mglevels[i+1]->restrct) {
641b9367914SBarry Smith           ierr = PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate);CHKERRQ(ierr);
642b9367914SBarry Smith         }
643b9367914SBarry Smith         if (mglevels[i+1]->interpolate == mglevels[i+1]->restrct) {
644f3fbd535SBarry Smith           ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
6457d537d92SJed Brown         } else {
6467d537d92SJed Brown           ierr = MatMatMatMult(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
6477d537d92SJed Brown         }
64823ee1639SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B);CHKERRQ(ierr);
649f3fbd535SBarry Smith         if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
650f3fbd535SBarry Smith         dB = B;
651f3fbd535SBarry Smith       }
652cd9507b2SBarry Smith       if (n > 1) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
653f3fbd535SBarry Smith     } else {
654f3fbd535SBarry Smith       for (i=n-2; i>-1; i--) {
655b9367914SBarry Smith         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");
656b9367914SBarry Smith         if (!mglevels[i+1]->interpolate) {
657b9367914SBarry Smith           ierr = PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct);CHKERRQ(ierr);
658b9367914SBarry Smith         }
659b9367914SBarry Smith         if (!mglevels[i+1]->restrct) {
660b9367914SBarry Smith           ierr = PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate);CHKERRQ(ierr);
661b9367914SBarry Smith         }
66223ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&B);CHKERRQ(ierr);
663b9367914SBarry Smith         if (mglevels[i+1]->interpolate == mglevels[i+1]->restrct) {
664f3fbd535SBarry Smith           ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
6657d537d92SJed Brown         } else {
6667d537d92SJed Brown           ierr = MatMatMatMult(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
6677d537d92SJed Brown         }
66823ee1639SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B);CHKERRQ(ierr);
669f3fbd535SBarry Smith         dB   = B;
670f3fbd535SBarry Smith       }
671f3fbd535SBarry Smith     }
672ce4cda84SJed Brown   } else if (!mg->galerkin && pc->dm && pc->dm->x) {
673cab2e9ccSBarry Smith     /* need to restrict Jacobian location to coarser meshes for evaluation */
674cab2e9ccSBarry Smith     for (i=n-2; i>-1; i--) {
675c88c5224SJed Brown       Mat R;
676c88c5224SJed Brown       Vec rscale;
677cab2e9ccSBarry Smith       if (!mglevels[i]->smoothd->dm->x) {
678cab2e9ccSBarry Smith         Vec *vecs;
6792a7a6963SBarry Smith         ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vecs,0,NULL);CHKERRQ(ierr);
6802fa5cd67SKarl Rupp 
681cab2e9ccSBarry Smith         mglevels[i]->smoothd->dm->x = vecs[0];
6822fa5cd67SKarl Rupp 
683cab2e9ccSBarry Smith         ierr = PetscFree(vecs);CHKERRQ(ierr);
684cab2e9ccSBarry Smith       }
685c88c5224SJed Brown       ierr = PCMGGetRestriction(pc,i+1,&R);CHKERRQ(ierr);
686c88c5224SJed Brown       ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
687c88c5224SJed Brown       ierr = MatRestrict(R,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);CHKERRQ(ierr);
688c88c5224SJed Brown       ierr = VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,rscale);CHKERRQ(ierr);
689cab2e9ccSBarry Smith     }
690f3fbd535SBarry Smith   }
691ccceb50cSJed Brown   if (!mg->galerkin && pc->dm) {
692caa4e7f2SJed Brown     for (i=n-2; i>=0; i--) {
693ccceb50cSJed Brown       DM  dmfine,dmcoarse;
694caa4e7f2SJed Brown       Mat Restrict,Inject;
695caa4e7f2SJed Brown       Vec rscale;
696ccceb50cSJed Brown       ierr   = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr);
697ccceb50cSJed Brown       ierr   = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr);
698caa4e7f2SJed Brown       ierr   = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr);
699caa4e7f2SJed Brown       ierr   = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
7000298fd71SBarry Smith       Inject = NULL;      /* Callback should create it if it needs Injection */
701caa4e7f2SJed Brown       ierr   = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr);
702caa4e7f2SJed Brown     }
703caa4e7f2SJed Brown   }
704f3fbd535SBarry Smith 
705f3fbd535SBarry Smith   if (!pc->setupcalled) {
706f3fbd535SBarry Smith     for (i=0; i<n; i++) {
707f3fbd535SBarry Smith       ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr);
708f3fbd535SBarry Smith     }
709f3fbd535SBarry Smith     for (i=1; i<n; i++) {
710f3fbd535SBarry Smith       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
711f3fbd535SBarry Smith         ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr);
712f3fbd535SBarry Smith       }
713f3fbd535SBarry Smith     }
714f3fbd535SBarry Smith     for (i=1; i<n; i++) {
715c88c5224SJed Brown       ierr = PCMGGetInterpolation(pc,i,&mglevels[i]->interpolate);CHKERRQ(ierr);
716c88c5224SJed Brown       ierr = PCMGGetRestriction(pc,i,&mglevels[i]->restrct);CHKERRQ(ierr);
717f3fbd535SBarry Smith     }
718f3fbd535SBarry Smith     for (i=0; i<n-1; i++) {
719f3fbd535SBarry Smith       if (!mglevels[i]->b) {
720f3fbd535SBarry Smith         Vec *vec;
7212a7a6963SBarry Smith         ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
722f3fbd535SBarry Smith         ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr);
7236bf464f9SBarry Smith         ierr = VecDestroy(vec);CHKERRQ(ierr);
724f3fbd535SBarry Smith         ierr = PetscFree(vec);CHKERRQ(ierr);
725f3fbd535SBarry Smith       }
726f3fbd535SBarry Smith       if (!mglevels[i]->r && i) {
727f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
728f3fbd535SBarry Smith         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
7296bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
730f3fbd535SBarry Smith       }
731f3fbd535SBarry Smith       if (!mglevels[i]->x) {
732f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
733f3fbd535SBarry Smith         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
7346bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
735f3fbd535SBarry Smith       }
736f3fbd535SBarry Smith     }
737f3fbd535SBarry Smith     if (n != 1 && !mglevels[n-1]->r) {
738f3fbd535SBarry Smith       /* PCMGSetR() on the finest level if user did not supply it */
739f3fbd535SBarry Smith       Vec *vec;
7402a7a6963SBarry Smith       ierr = KSPCreateVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
741f3fbd535SBarry Smith       ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr);
7426bf464f9SBarry Smith       ierr = VecDestroy(vec);CHKERRQ(ierr);
743f3fbd535SBarry Smith       ierr = PetscFree(vec);CHKERRQ(ierr);
744f3fbd535SBarry Smith     }
745f3fbd535SBarry Smith   }
746f3fbd535SBarry Smith 
74798e04984SBarry Smith   if (pc->dm) {
74898e04984SBarry Smith     /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */
74998e04984SBarry Smith     for (i=0; i<n-1; i++) {
75098e04984SBarry Smith       if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX;
75198e04984SBarry Smith     }
75298e04984SBarry Smith   }
753f3fbd535SBarry Smith 
754f3fbd535SBarry Smith   for (i=1; i<n; i++) {
7552a21e185SBarry Smith     if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1){
756f3fbd535SBarry Smith       /* if doing only down then initial guess is zero */
757f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
758f3fbd535SBarry Smith     }
75963e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
760f3fbd535SBarry Smith     ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr);
76163e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
762d42688cbSBarry Smith     if (!mglevels[i]->residual) {
763d42688cbSBarry Smith       Mat mat;
76423ee1639SBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&mat);CHKERRQ(ierr);
76554b2cd4bSJed Brown       ierr = PCMGSetResidual(pc,i,PCMGResidualDefault,mat);CHKERRQ(ierr);
766d42688cbSBarry Smith     }
767f3fbd535SBarry Smith   }
768f3fbd535SBarry Smith   for (i=1; i<n; i++) {
769f3fbd535SBarry Smith     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
770f3fbd535SBarry Smith       Mat          downmat,downpmat;
771f3fbd535SBarry Smith 
772f3fbd535SBarry Smith       /* check if operators have been set for up, if not use down operators to set them */
7730298fd71SBarry Smith       ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);CHKERRQ(ierr);
774f3fbd535SBarry Smith       if (!opsset) {
77523ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);CHKERRQ(ierr);
77623ee1639SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat);CHKERRQ(ierr);
777f3fbd535SBarry Smith       }
778f3fbd535SBarry Smith 
779f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
78063e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
781f3fbd535SBarry Smith       ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr);
78263e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
783f3fbd535SBarry Smith     }
784f3fbd535SBarry Smith   }
785f3fbd535SBarry Smith 
786f3fbd535SBarry Smith   /*
787f3fbd535SBarry Smith       If coarse solver is not direct method then DO NOT USE preonly
788f3fbd535SBarry Smith   */
789251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mglevels[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr);
790f3fbd535SBarry Smith   if (preonly) {
791251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr);
792251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr);
793251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr);
794251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)cpc,PCSVD,&svd);CHKERRQ(ierr);
795fd1303e9SJungho Lee     if (!lu && !redundant && !cholesky && !svd) {
796f3fbd535SBarry Smith       ierr = KSPSetType(mglevels[0]->smoothd,KSPGMRES);CHKERRQ(ierr);
797f3fbd535SBarry Smith     }
798f3fbd535SBarry Smith   }
799f3fbd535SBarry Smith 
80063e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
801f3fbd535SBarry Smith   ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr);
80263e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
803f3fbd535SBarry Smith 
804f3fbd535SBarry Smith   /*
805f3fbd535SBarry Smith      Dump the interpolation/restriction matrices plus the
806e3c5b3baSBarry Smith    Jacobian/stiffness on each level. This allows MATLAB users to
807f3fbd535SBarry Smith    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
808f3fbd535SBarry Smith 
809f3fbd535SBarry Smith    Only support one or the other at the same time.
810f3fbd535SBarry Smith   */
811f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
8120298fd71SBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);CHKERRQ(ierr);
813ce94432eSBarry Smith   if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc));
814f3fbd535SBarry Smith   dump = PETSC_FALSE;
815f3fbd535SBarry Smith #endif
8160298fd71SBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);CHKERRQ(ierr);
817ce94432eSBarry Smith   if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc));
818f3fbd535SBarry Smith 
819f3fbd535SBarry Smith   if (viewer) {
820f3fbd535SBarry Smith     for (i=1; i<n; i++) {
821f3fbd535SBarry Smith       ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr);
822f3fbd535SBarry Smith     }
823f3fbd535SBarry Smith     for (i=0; i<n; i++) {
824f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr);
825f3fbd535SBarry Smith       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
826f3fbd535SBarry Smith     }
827f3fbd535SBarry Smith   }
828f3fbd535SBarry Smith   PetscFunctionReturn(0);
829f3fbd535SBarry Smith }
830f3fbd535SBarry Smith 
831f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/
832f3fbd535SBarry Smith 
833f3fbd535SBarry Smith #undef __FUNCT__
8349dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetLevels"
8354b9ad928SBarry Smith /*@
83697177400SBarry Smith    PCMGGetLevels - Gets the number of levels to use with MG.
8374b9ad928SBarry Smith 
8384b9ad928SBarry Smith    Not Collective
8394b9ad928SBarry Smith 
8404b9ad928SBarry Smith    Input Parameter:
8414b9ad928SBarry Smith .  pc - the preconditioner context
8424b9ad928SBarry Smith 
8434b9ad928SBarry Smith    Output parameter:
8444b9ad928SBarry Smith .  levels - the number of levels
8454b9ad928SBarry Smith 
8464b9ad928SBarry Smith    Level: advanced
8474b9ad928SBarry Smith 
8484b9ad928SBarry Smith .keywords: MG, get, levels, multigrid
8494b9ad928SBarry Smith 
85097177400SBarry Smith .seealso: PCMGSetLevels()
8514b9ad928SBarry Smith @*/
8527087cfbeSBarry Smith PetscErrorCode  PCMGGetLevels(PC pc,PetscInt *levels)
8534b9ad928SBarry Smith {
854f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
8554b9ad928SBarry Smith 
8564b9ad928SBarry Smith   PetscFunctionBegin;
8570700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
8584482741eSBarry Smith   PetscValidIntPointer(levels,2);
859f3fbd535SBarry Smith   *levels = mg->nlevels;
8604b9ad928SBarry Smith   PetscFunctionReturn(0);
8614b9ad928SBarry Smith }
8624b9ad928SBarry Smith 
8634b9ad928SBarry Smith #undef __FUNCT__
8649dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetType"
8654b9ad928SBarry Smith /*@
86697177400SBarry Smith    PCMGSetType - Determines the form of multigrid to use:
8674b9ad928SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
8684b9ad928SBarry Smith 
869ad4df100SBarry Smith    Logically Collective on PC
8704b9ad928SBarry Smith 
8714b9ad928SBarry Smith    Input Parameters:
8724b9ad928SBarry Smith +  pc - the preconditioner context
8739dcbbd2bSBarry Smith -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
8749dcbbd2bSBarry Smith    PC_MG_FULL, PC_MG_KASKADE
8754b9ad928SBarry Smith 
8764b9ad928SBarry Smith    Options Database Key:
8774b9ad928SBarry Smith .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
8784b9ad928SBarry Smith    additive, full, kaskade
8794b9ad928SBarry Smith 
8804b9ad928SBarry Smith    Level: advanced
8814b9ad928SBarry Smith 
8824b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
8834b9ad928SBarry Smith 
88497177400SBarry Smith .seealso: PCMGSetLevels()
8854b9ad928SBarry Smith @*/
8867087cfbeSBarry Smith PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
8874b9ad928SBarry Smith {
888f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
8894b9ad928SBarry Smith 
8904b9ad928SBarry Smith   PetscFunctionBegin;
8910700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
892c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,form,2);
89331567311SBarry Smith   mg->am = form;
8949dcbbd2bSBarry Smith   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
895c60c7ad4SBarry Smith   else pc->ops->applyrichardson = NULL;
896c60c7ad4SBarry Smith   PetscFunctionReturn(0);
897c60c7ad4SBarry Smith }
898c60c7ad4SBarry Smith 
899c60c7ad4SBarry Smith /*@
900c60c7ad4SBarry Smith    PCMGGetType - Determines the form of multigrid to use:
901c60c7ad4SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
902c60c7ad4SBarry Smith 
903c60c7ad4SBarry Smith    Logically Collective on PC
904c60c7ad4SBarry Smith 
905c60c7ad4SBarry Smith    Input Parameter:
906c60c7ad4SBarry Smith .  pc - the preconditioner context
907c60c7ad4SBarry Smith 
908c60c7ad4SBarry Smith    Output Parameter:
909c60c7ad4SBarry Smith .  type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE
910c60c7ad4SBarry Smith 
911c60c7ad4SBarry Smith 
912c60c7ad4SBarry Smith    Level: advanced
913c60c7ad4SBarry Smith 
914c60c7ad4SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
915c60c7ad4SBarry Smith 
916c60c7ad4SBarry Smith .seealso: PCMGSetLevels()
917c60c7ad4SBarry Smith @*/
918c60c7ad4SBarry Smith PetscErrorCode  PCMGGetType(PC pc,PCMGType *type)
919c60c7ad4SBarry Smith {
920c60c7ad4SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
921c60c7ad4SBarry Smith 
922c60c7ad4SBarry Smith   PetscFunctionBegin;
923c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
924c60c7ad4SBarry Smith   *type = mg->am;
9254b9ad928SBarry Smith   PetscFunctionReturn(0);
9264b9ad928SBarry Smith }
9274b9ad928SBarry Smith 
9284b9ad928SBarry Smith #undef __FUNCT__
9290d353602SBarry Smith #define __FUNCT__ "PCMGSetCycleType"
9304b9ad928SBarry Smith /*@
9310d353602SBarry Smith    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
9324b9ad928SBarry Smith    complicated cycling.
9334b9ad928SBarry Smith 
934ad4df100SBarry Smith    Logically Collective on PC
9354b9ad928SBarry Smith 
9364b9ad928SBarry Smith    Input Parameters:
937c2be2410SBarry Smith +  pc - the multigrid context
9380d353602SBarry Smith -  PC_MG_CYCLE_V or PC_MG_CYCLE_W
9394b9ad928SBarry Smith 
9404b9ad928SBarry Smith    Options Database Key:
9414446f3b4SBarry Smith .  -pc_mg_cycle_type <v,w>
9424b9ad928SBarry Smith 
9434b9ad928SBarry Smith    Level: advanced
9444b9ad928SBarry Smith 
9454b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
9464b9ad928SBarry Smith 
9470d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel()
9484b9ad928SBarry Smith @*/
9497087cfbeSBarry Smith PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
9504b9ad928SBarry Smith {
951f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
952f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
95379416396SBarry Smith   PetscInt     i,levels;
9544b9ad928SBarry Smith 
9554b9ad928SBarry Smith   PetscFunctionBegin;
9560700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
957ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
958c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
959f3fbd535SBarry Smith   levels = mglevels[0]->levels;
9604b9ad928SBarry Smith 
9612fa5cd67SKarl Rupp   for (i=0; i<levels; i++) mglevels[i]->cycles = n;
9624b9ad928SBarry Smith   PetscFunctionReturn(0);
9634b9ad928SBarry Smith }
9644b9ad928SBarry Smith 
9654b9ad928SBarry Smith #undef __FUNCT__
9668cc2d5dfSBarry Smith #define __FUNCT__ "PCMGMultiplicativeSetCycles"
9678cc2d5dfSBarry Smith /*@
9688cc2d5dfSBarry Smith    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
9698cc2d5dfSBarry Smith          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
9708cc2d5dfSBarry Smith 
971ad4df100SBarry Smith    Logically Collective on PC
9728cc2d5dfSBarry Smith 
9738cc2d5dfSBarry Smith    Input Parameters:
9748cc2d5dfSBarry Smith +  pc - the multigrid context
9758cc2d5dfSBarry Smith -  n - number of cycles (default is 1)
9768cc2d5dfSBarry Smith 
9778cc2d5dfSBarry Smith    Options Database Key:
978e1bc860dSBarry Smith .  -pc_mg_multiplicative_cycles n
9798cc2d5dfSBarry Smith 
9808cc2d5dfSBarry Smith    Level: advanced
9818cc2d5dfSBarry Smith 
9828cc2d5dfSBarry Smith    Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
9838cc2d5dfSBarry Smith 
9848cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
9858cc2d5dfSBarry Smith 
9868cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
9878cc2d5dfSBarry Smith @*/
9887087cfbeSBarry Smith PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
9898cc2d5dfSBarry Smith {
990f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
9918cc2d5dfSBarry Smith 
9928cc2d5dfSBarry Smith   PetscFunctionBegin;
9930700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
994c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
9952a21e185SBarry Smith   mg->cyclesperpcapply = n;
9968cc2d5dfSBarry Smith   PetscFunctionReturn(0);
9978cc2d5dfSBarry Smith }
9988cc2d5dfSBarry Smith 
9998cc2d5dfSBarry Smith #undef __FUNCT__
10009dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetGalerkin"
1001c2be2410SBarry Smith /*@
100297177400SBarry Smith    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
100382b87d87SMatthew G. Knepley       finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i
1004c2be2410SBarry Smith 
1005ad4df100SBarry Smith    Logically Collective on PC
1006c2be2410SBarry Smith 
1007c2be2410SBarry Smith    Input Parameters:
1008c91913d3SJed Brown +  pc - the multigrid context
1009c91913d3SJed Brown -  use - PETSC_TRUE to use the Galerkin process to compute coarse-level operators
1010c2be2410SBarry Smith 
1011c2be2410SBarry Smith    Options Database Key:
1012e1bc860dSBarry Smith .  -pc_mg_galerkin
1013c2be2410SBarry Smith 
1014c2be2410SBarry Smith    Level: intermediate
1015c2be2410SBarry Smith 
1016c2be2410SBarry Smith .keywords: MG, set, Galerkin
1017c2be2410SBarry Smith 
10183fc8bf9cSBarry Smith .seealso: PCMGGetGalerkin()
10193fc8bf9cSBarry Smith 
1020c2be2410SBarry Smith @*/
1021c91913d3SJed Brown PetscErrorCode PCMGSetGalerkin(PC pc,PetscBool use)
1022c2be2410SBarry Smith {
1023f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
1024c2be2410SBarry Smith 
1025c2be2410SBarry Smith   PetscFunctionBegin;
10260700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1027789726d7SBarry Smith   mg->galerkin = use ? 1 : 0;
1028c2be2410SBarry Smith   PetscFunctionReturn(0);
1029c2be2410SBarry Smith }
1030c2be2410SBarry Smith 
1031c2be2410SBarry Smith #undef __FUNCT__
10323fc8bf9cSBarry Smith #define __FUNCT__ "PCMGGetGalerkin"
10333fc8bf9cSBarry Smith /*@
10343fc8bf9cSBarry Smith    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
103582b87d87SMatthew G. Knepley       A_i-1 = r_i * A_i * p_i
10363fc8bf9cSBarry Smith 
10373fc8bf9cSBarry Smith    Not Collective
10383fc8bf9cSBarry Smith 
10393fc8bf9cSBarry Smith    Input Parameter:
10403fc8bf9cSBarry Smith .  pc - the multigrid context
10413fc8bf9cSBarry Smith 
10423fc8bf9cSBarry Smith    Output Parameter:
1043e1bc860dSBarry Smith .  galerkin - PETSC_TRUE or PETSC_FALSE
10443fc8bf9cSBarry Smith 
10453fc8bf9cSBarry Smith    Options Database Key:
1046e1bc860dSBarry Smith .  -pc_mg_galerkin
10473fc8bf9cSBarry Smith 
10483fc8bf9cSBarry Smith    Level: intermediate
10493fc8bf9cSBarry Smith 
10503fc8bf9cSBarry Smith .keywords: MG, set, Galerkin
10513fc8bf9cSBarry Smith 
10523fc8bf9cSBarry Smith .seealso: PCMGSetGalerkin()
10533fc8bf9cSBarry Smith 
10543fc8bf9cSBarry Smith @*/
10557087cfbeSBarry Smith PetscErrorCode  PCMGGetGalerkin(PC pc,PetscBool  *galerkin)
10563fc8bf9cSBarry Smith {
1057f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
10583fc8bf9cSBarry Smith 
10593fc8bf9cSBarry Smith   PetscFunctionBegin;
10600700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
106113fdb3acSJose Roman   *galerkin = (PetscBool)mg->galerkin;
10623fc8bf9cSBarry Smith   PetscFunctionReturn(0);
10633fc8bf9cSBarry Smith }
10643fc8bf9cSBarry Smith 
10653fc8bf9cSBarry Smith #undef __FUNCT__
10669dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothDown"
10674b9ad928SBarry Smith /*@
106897177400SBarry Smith    PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
106997177400SBarry Smith    use on all levels. Use PCMGGetSmootherDown() to set different
10704b9ad928SBarry Smith    pre-smoothing steps on different levels.
10714b9ad928SBarry Smith 
1072ad4df100SBarry Smith    Logically Collective on PC
10734b9ad928SBarry Smith 
10744b9ad928SBarry Smith    Input Parameters:
10754b9ad928SBarry Smith +  mg - the multigrid context
10764b9ad928SBarry Smith -  n - the number of smoothing steps
10774b9ad928SBarry Smith 
10784b9ad928SBarry Smith    Options Database Key:
10794b9ad928SBarry Smith .  -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps
10804b9ad928SBarry Smith 
10814b9ad928SBarry Smith    Level: advanced
10824b9ad928SBarry Smith 
10834b9ad928SBarry Smith .keywords: MG, smooth, down, pre-smoothing, steps, multigrid
10844b9ad928SBarry Smith 
108597177400SBarry Smith .seealso: PCMGSetNumberSmoothUp()
10864b9ad928SBarry Smith @*/
10877087cfbeSBarry Smith PetscErrorCode  PCMGSetNumberSmoothDown(PC pc,PetscInt n)
10884b9ad928SBarry Smith {
1089f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
1090f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
10916849ba73SBarry Smith   PetscErrorCode ierr;
109279416396SBarry Smith   PetscInt       i,levels;
10934b9ad928SBarry Smith 
10944b9ad928SBarry Smith   PetscFunctionBegin;
10950700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1096ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1097c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
1098f3fbd535SBarry Smith   levels = mglevels[0]->levels;
10994b9ad928SBarry Smith 
1100b05257ddSBarry Smith   for (i=1; i<levels; i++) {
11014b9ad928SBarry Smith     /* make sure smoother up and down are different */
11020298fd71SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,NULL);CHKERRQ(ierr);
1103f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
11042fa5cd67SKarl Rupp 
110531567311SBarry Smith     mg->default_smoothd = n;
11064b9ad928SBarry Smith   }
11074b9ad928SBarry Smith   PetscFunctionReturn(0);
11084b9ad928SBarry Smith }
11094b9ad928SBarry Smith 
11104b9ad928SBarry Smith #undef __FUNCT__
11119dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothUp"
11124b9ad928SBarry Smith /*@
111397177400SBarry Smith    PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
111497177400SBarry Smith    on all levels. Use PCMGGetSmootherUp() to set different numbers of
11154b9ad928SBarry Smith    post-smoothing steps on different levels.
11164b9ad928SBarry Smith 
1117ad4df100SBarry Smith    Logically Collective on PC
11184b9ad928SBarry Smith 
11194b9ad928SBarry Smith    Input Parameters:
11204b9ad928SBarry Smith +  mg - the multigrid context
11214b9ad928SBarry Smith -  n - the number of smoothing steps
11224b9ad928SBarry Smith 
11234b9ad928SBarry Smith    Options Database Key:
11244b9ad928SBarry Smith .  -pc_mg_smoothup <n> - Sets number of post-smoothing steps
11254b9ad928SBarry Smith 
11264b9ad928SBarry Smith    Level: advanced
11274b9ad928SBarry Smith 
11284b9ad928SBarry Smith    Note: this does not set a value on the coarsest grid, since we assume that
1129a8c7a070SBarry Smith     there is no separate smooth up on the coarsest grid.
11304b9ad928SBarry Smith 
11314b9ad928SBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid
11324b9ad928SBarry Smith 
113397177400SBarry Smith .seealso: PCMGSetNumberSmoothDown()
11344b9ad928SBarry Smith @*/
11357087cfbeSBarry Smith PetscErrorCode  PCMGSetNumberSmoothUp(PC pc,PetscInt n)
11364b9ad928SBarry Smith {
1137f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
1138f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
11396849ba73SBarry Smith   PetscErrorCode ierr;
114079416396SBarry Smith   PetscInt       i,levels;
11414b9ad928SBarry Smith 
11424b9ad928SBarry Smith   PetscFunctionBegin;
11430700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1144ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1145c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
1146f3fbd535SBarry Smith   levels = mglevels[0]->levels;
11474b9ad928SBarry Smith 
11484b9ad928SBarry Smith   for (i=1; i<levels; i++) {
11494b9ad928SBarry Smith     /* make sure smoother up and down are different */
11500298fd71SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,NULL);CHKERRQ(ierr);
1151f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
11522fa5cd67SKarl Rupp 
115331567311SBarry Smith     mg->default_smoothu = n;
11544b9ad928SBarry Smith   }
11554b9ad928SBarry Smith   PetscFunctionReturn(0);
11564b9ad928SBarry Smith }
11574b9ad928SBarry Smith 
11584b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/
11594b9ad928SBarry Smith 
11603b09bd56SBarry Smith /*MC
1161ccb205f8SBarry Smith    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
11623b09bd56SBarry Smith     information about the coarser grid matrices and restriction/interpolation operators.
11633b09bd56SBarry Smith 
11643b09bd56SBarry Smith    Options Database Keys:
11653b09bd56SBarry Smith +  -pc_mg_levels <nlevels> - number of levels including finest
1166b4519db0SPatrick Sanan .  -pc_mg_cycles <v,w> -
116779416396SBarry Smith .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
11683b09bd56SBarry Smith .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
11698c1c2452SJed Brown .  -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default
11703b09bd56SBarry Smith .  -pc_mg_log - log information about time spent on each level of the solver
11713b09bd56SBarry Smith .  -pc_mg_monitor - print information on the multigrid convergence
11728cf6037eSBarry Smith .  -pc_mg_galerkin - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R'
11738cf6037eSBarry Smith .  -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1)
11748cf6037eSBarry Smith .  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1175e3c5b3baSBarry Smith                         to the Socket viewer for reading from MATLAB.
11768cf6037eSBarry Smith -  -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices
11778cf6037eSBarry Smith                         to the binary output file called binaryoutput
11783b09bd56SBarry Smith 
117924c3aa18SBarry Smith    Notes: By default this uses GMRES on the fine grid smoother so this should be used with KSPFGMRES or the smoother changed to not use GMRES
11803b09bd56SBarry Smith 
11818cf6037eSBarry Smith        When run with a single level the smoother options are used on that level NOT the coarse grid solver options
11828cf6037eSBarry Smith 
11833b09bd56SBarry Smith    Level: intermediate
11843b09bd56SBarry Smith 
11858f87f92bSBarry Smith    Concepts: multigrid/multilevel
11863b09bd56SBarry Smith 
11878cf6037eSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE
11880d353602SBarry Smith            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
118997177400SBarry Smith            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
119097177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
11910d353602SBarry Smith            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
11923b09bd56SBarry Smith M*/
11933b09bd56SBarry Smith 
11944b9ad928SBarry Smith #undef __FUNCT__
11954b9ad928SBarry Smith #define __FUNCT__ "PCCreate_MG"
11968cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc)
11974b9ad928SBarry Smith {
1198f3fbd535SBarry Smith   PC_MG          *mg;
1199f3fbd535SBarry Smith   PetscErrorCode ierr;
1200f3fbd535SBarry Smith 
12014b9ad928SBarry Smith   PetscFunctionBegin;
1202b00a9115SJed Brown   ierr        = PetscNewLog(pc,&mg);CHKERRQ(ierr);
1203f3fbd535SBarry Smith   pc->data    = (void*)mg;
1204f3fbd535SBarry Smith   mg->nlevels = -1;
1205f3fbd535SBarry Smith 
120637a44384SMark Adams   pc->useAmat = PETSC_TRUE;
120737a44384SMark Adams 
12084b9ad928SBarry Smith   pc->ops->apply          = PCApply_MG;
12094b9ad928SBarry Smith   pc->ops->setup          = PCSetUp_MG;
1210a06653b4SBarry Smith   pc->ops->reset          = PCReset_MG;
12114b9ad928SBarry Smith   pc->ops->destroy        = PCDestroy_MG;
12124b9ad928SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_MG;
12134b9ad928SBarry Smith   pc->ops->view           = PCView_MG;
12144b9ad928SBarry Smith   PetscFunctionReturn(0);
12154b9ad928SBarry Smith }
1216