xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision c60c7ad4eb3baa38c2a8deee7f23b741b9653612)
1dba47a55SKris Buschelman 
24b9ad928SBarry Smith /*
34b9ad928SBarry Smith     Defines the multigrid preconditioner interface.
44b9ad928SBarry Smith */
5b4876bcbSBarry Smith #include <../src/ksp/pc/impls/mg/mgimpl.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"
353f3fbd535SBarry Smith PetscErrorCode PCSetFromOptions_MG(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;
364f3fbd535SBarry Smith   ierr = PetscOptionsHead("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     }
478f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
479f3fbd535SBarry Smith       if (!i) {
48089cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);CHKERRQ(ierr);
481f3fbd535SBarry Smith       } else {
48289cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
483f3fbd535SBarry Smith       }
484f3fbd535SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
485f3fbd535SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
486f3fbd535SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
487f3fbd535SBarry Smith       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
488f3fbd535SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
489f3fbd535SBarry Smith       } else if (i) {
49089cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
491f3fbd535SBarry Smith         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
492f3fbd535SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
493f3fbd535SBarry Smith         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
494f3fbd535SBarry Smith       }
495f3fbd535SBarry Smith     }
4965b0b0462SBarry Smith   } else if (isbinary) {
4975b0b0462SBarry Smith     for (i=levels-1; i>=0; i--) {
4985b0b0462SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
4995b0b0462SBarry Smith       if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) {
5005b0b0462SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
5015b0b0462SBarry Smith       }
5025b0b0462SBarry Smith     }
503219b1045SBarry Smith   } else if (isdraw) {
504219b1045SBarry Smith     PetscDraw draw;
505b4375e8dSPeter Brune     PetscReal x,w,y,bottom,th;
506219b1045SBarry Smith     ierr   = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
507219b1045SBarry Smith     ierr   = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
5080298fd71SBarry Smith     ierr   = PetscDrawStringGetSize(draw,NULL,&th);CHKERRQ(ierr);
509219b1045SBarry Smith     bottom = y - th;
510219b1045SBarry Smith     for (i=levels-1; i>=0; i--) {
511b4375e8dSPeter Brune       if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) {
512219b1045SBarry Smith         ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr);
513219b1045SBarry Smith         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
514219b1045SBarry Smith         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
515b4375e8dSPeter Brune       } else {
516b4375e8dSPeter Brune         w    = 0.5*PetscMin(1.0-x,x);
517b4375e8dSPeter Brune         ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr);
518b4375e8dSPeter Brune         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
519b4375e8dSPeter Brune         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
520b4375e8dSPeter Brune         ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr);
521b4375e8dSPeter Brune         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
522b4375e8dSPeter Brune         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
523b4375e8dSPeter Brune       }
5240298fd71SBarry Smith       ierr    = PetscDrawGetBoundingBox(draw,NULL,&bottom,NULL,NULL);CHKERRQ(ierr);
5251cd381d6SBarry Smith       bottom -= th;
526219b1045SBarry Smith     }
5275b0b0462SBarry Smith   }
528f3fbd535SBarry Smith   PetscFunctionReturn(0);
529f3fbd535SBarry Smith }
530f3fbd535SBarry Smith 
531b45d2f2cSJed Brown #include <petsc-private/dmimpl.h>
532b45d2f2cSJed Brown #include <petsc-private/kspimpl.h>
533cab2e9ccSBarry Smith 
534f3fbd535SBarry Smith /*
535f3fbd535SBarry Smith     Calls setup for the KSP on each level
536f3fbd535SBarry Smith */
537f3fbd535SBarry Smith #undef __FUNCT__
538f3fbd535SBarry Smith #define __FUNCT__ "PCSetUp_MG"
539f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc)
540f3fbd535SBarry Smith {
541f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
542f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
543f3fbd535SBarry Smith   PetscErrorCode ierr;
544f3fbd535SBarry Smith   PetscInt       i,n = mglevels[0]->levels;
54598e04984SBarry Smith   PC             cpc;
54671b23a65SMark F. Adams   PetscBool      preonly,lu,redundant,cholesky,svd,dump = PETSC_FALSE,opsset,use_amat;
547f3fbd535SBarry Smith   Mat            dA,dB;
548f3fbd535SBarry Smith   Vec            tvec;
549218a07d4SBarry Smith   DM             *dms;
550649052a6SBarry Smith   PetscViewer    viewer = 0;
551f3fbd535SBarry Smith 
552f3fbd535SBarry Smith   PetscFunctionBegin;
55301bc414fSDmitry Karpeev   /* FIX: Move this to PCSetFromOptions_MG? */
5543aeaf226SBarry Smith   if (mg->usedmfornumberoflevels) {
5553aeaf226SBarry Smith     PetscInt levels;
5563aeaf226SBarry Smith     ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
5573aeaf226SBarry Smith     levels++;
5583aeaf226SBarry Smith     if (levels > n) { /* the problem is now being solved on a finer grid */
5590298fd71SBarry Smith       ierr     = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr);
5603aeaf226SBarry Smith       n        = levels;
5613aeaf226SBarry Smith       ierr     = PCSetFromOptions(pc);CHKERRQ(ierr); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
5623aeaf226SBarry Smith       mglevels =  mg->levels;
5633aeaf226SBarry Smith     }
5643aeaf226SBarry Smith   }
56598e04984SBarry Smith   ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr);
5663aeaf226SBarry Smith 
567f3fbd535SBarry Smith 
568f3fbd535SBarry Smith   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
569f3fbd535SBarry Smith   /* so use those from global PC */
570f3fbd535SBarry Smith   /* Is this what we always want? What if user wants to keep old one? */
5710298fd71SBarry Smith   ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset);CHKERRQ(ierr);
57298e04984SBarry Smith   if (opsset) {
57398e04984SBarry Smith     Mat mmat;
57423ee1639SBarry Smith     ierr = KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat);CHKERRQ(ierr);
57598e04984SBarry Smith     if (mmat == pc->pmat) opsset = PETSC_FALSE;
57698e04984SBarry Smith   }
5774a5f07a5SMark F. Adams 
5784a5f07a5SMark F. Adams   if (!opsset) {
57971b23a65SMark F. Adams     ierr = PCGetUseAmat(pc,&use_amat);CHKERRQ(ierr);
5804a5f07a5SMark F. Adams     if(use_amat){
581f3fbd535SBarry 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);
58223ee1639SBarry Smith       ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat);CHKERRQ(ierr);
583f3fbd535SBarry Smith     }
5844a5f07a5SMark F. Adams     else {
5854a5f07a5SMark 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);
58623ee1639SBarry Smith       ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat);CHKERRQ(ierr);
5874a5f07a5SMark F. Adams     }
5884a5f07a5SMark F. Adams   }
589f3fbd535SBarry Smith 
59001bc414fSDmitry 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? */
591ce4cda84SJed Brown   if (pc->dm && mg->galerkin != 2 && !pc->setupcalled) {
5922d2b81a6SBarry Smith     /* construct the interpolation from the DMs */
5932e499ae9SBarry Smith     Mat p;
59473250ac0SBarry Smith     Vec rscale;
595785e854fSJed Brown     ierr     = PetscMalloc1(n,&dms);CHKERRQ(ierr);
596218a07d4SBarry Smith     dms[n-1] = pc->dm;
597218a07d4SBarry Smith     for (i=n-2; i>-1; i--) {
598942e3340SBarry Smith       DMKSP kdm;
5992ee06e3bSJed Brown       ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr);
6003c0c59f3SBarry Smith       ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr);
6013c0c59f3SBarry Smith       if (mg->galerkin) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);}
602942e3340SBarry Smith       ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr);
603d1a3e677SJed Brown       /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take
604d1a3e677SJed Brown        * a bitwise OR of computing the matrix, RHS, and initial iterate. */
6050298fd71SBarry Smith       kdm->ops->computerhs = NULL;
6060298fd71SBarry Smith       kdm->rhsctx          = NULL;
60724c3aa18SBarry Smith       if (!mglevels[i+1]->interpolate) {
608e727c939SJed Brown         ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr);
6092d2b81a6SBarry Smith         ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr);
61000ac8be1SBarry Smith         if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);}
61173250ac0SBarry Smith         ierr = VecDestroy(&rscale);CHKERRQ(ierr);
6126bf464f9SBarry Smith         ierr = MatDestroy(&p);CHKERRQ(ierr);
613218a07d4SBarry Smith       }
61424c3aa18SBarry Smith     }
6152d2b81a6SBarry Smith 
616218a07d4SBarry Smith     for (i=n-2; i>-1; i--) {
6176bf464f9SBarry Smith       ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);
618218a07d4SBarry Smith     }
6192d2b81a6SBarry Smith     ierr = PetscFree(dms);CHKERRQ(ierr);
620ce4cda84SJed Brown   }
621cab2e9ccSBarry Smith 
622ce4cda84SJed Brown   if (pc->dm && !pc->setupcalled) {
623ce4cda84SJed Brown     /* finest smoother also gets DM but it is not active, independent of whether galerkin==2 */
624cab2e9ccSBarry Smith     ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr);
625cab2e9ccSBarry Smith     ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr);
626218a07d4SBarry Smith   }
627218a07d4SBarry Smith 
628c91913d3SJed Brown   if (mg->galerkin == 1) {
629f3fbd535SBarry Smith     Mat B;
630f3fbd535SBarry Smith     /* currently only handle case where mat and pmat are the same on coarser levels */
63123ee1639SBarry Smith     ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB);CHKERRQ(ierr);
632f3fbd535SBarry Smith     if (!pc->setupcalled) {
633f3fbd535SBarry Smith       for (i=n-2; i>-1; i--) {
6349f3d9e6bSJed Brown         if (mglevels[i+1]->interpolate == mglevels[i+1]->restrct || !mglevels[i+1]->restrct) {
635f3fbd535SBarry Smith           ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
6367d537d92SJed Brown         } else {
6377d537d92SJed Brown           ierr = MatMatMatMult(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
6387d537d92SJed Brown         }
63923ee1639SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B);CHKERRQ(ierr);
640f3fbd535SBarry Smith         if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
641f3fbd535SBarry Smith         dB = B;
642f3fbd535SBarry Smith       }
643cd9507b2SBarry Smith       if (n > 1) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
644f3fbd535SBarry Smith     } else {
645f3fbd535SBarry Smith       for (i=n-2; i>-1; i--) {
64623ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&B);CHKERRQ(ierr);
6479f3d9e6bSJed Brown         if (mglevels[i+1]->interpolate == mglevels[i+1]->restrct || !mglevels[i+1]->restrct) {
648f3fbd535SBarry Smith           ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
6497d537d92SJed Brown         } else {
6507d537d92SJed Brown           ierr = MatMatMatMult(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
6517d537d92SJed Brown         }
65223ee1639SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B);CHKERRQ(ierr);
653f3fbd535SBarry Smith         dB   = B;
654f3fbd535SBarry Smith       }
655f3fbd535SBarry Smith     }
656ce4cda84SJed Brown   } else if (!mg->galerkin && pc->dm && pc->dm->x) {
657cab2e9ccSBarry Smith     /* need to restrict Jacobian location to coarser meshes for evaluation */
658cab2e9ccSBarry Smith     for (i=n-2; i>-1; i--) {
659c88c5224SJed Brown       Mat R;
660c88c5224SJed Brown       Vec rscale;
661cab2e9ccSBarry Smith       if (!mglevels[i]->smoothd->dm->x) {
662cab2e9ccSBarry Smith         Vec *vecs;
6632a7a6963SBarry Smith         ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vecs,0,NULL);CHKERRQ(ierr);
6642fa5cd67SKarl Rupp 
665cab2e9ccSBarry Smith         mglevels[i]->smoothd->dm->x = vecs[0];
6662fa5cd67SKarl Rupp 
667cab2e9ccSBarry Smith         ierr = PetscFree(vecs);CHKERRQ(ierr);
668cab2e9ccSBarry Smith       }
669c88c5224SJed Brown       ierr = PCMGGetRestriction(pc,i+1,&R);CHKERRQ(ierr);
670c88c5224SJed Brown       ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
671c88c5224SJed Brown       ierr = MatRestrict(R,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);CHKERRQ(ierr);
672c88c5224SJed Brown       ierr = VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,rscale);CHKERRQ(ierr);
673cab2e9ccSBarry Smith     }
674f3fbd535SBarry Smith   }
675ccceb50cSJed Brown   if (!mg->galerkin && pc->dm) {
676caa4e7f2SJed Brown     for (i=n-2; i>=0; i--) {
677ccceb50cSJed Brown       DM  dmfine,dmcoarse;
678caa4e7f2SJed Brown       Mat Restrict,Inject;
679caa4e7f2SJed Brown       Vec rscale;
680ccceb50cSJed Brown       ierr   = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr);
681ccceb50cSJed Brown       ierr   = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr);
682caa4e7f2SJed Brown       ierr   = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr);
683caa4e7f2SJed Brown       ierr   = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
6840298fd71SBarry Smith       Inject = NULL;      /* Callback should create it if it needs Injection */
685caa4e7f2SJed Brown       ierr   = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr);
686caa4e7f2SJed Brown     }
687caa4e7f2SJed Brown   }
688f3fbd535SBarry Smith 
689f3fbd535SBarry Smith   if (!pc->setupcalled) {
690f3fbd535SBarry Smith     for (i=0; i<n; i++) {
691f3fbd535SBarry Smith       ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr);
692f3fbd535SBarry Smith     }
693f3fbd535SBarry Smith     for (i=1; i<n; i++) {
694f3fbd535SBarry Smith       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
695f3fbd535SBarry Smith         ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr);
696f3fbd535SBarry Smith       }
697f3fbd535SBarry Smith     }
698f3fbd535SBarry Smith     for (i=1; i<n; i++) {
699c88c5224SJed Brown       ierr = PCMGGetInterpolation(pc,i,&mglevels[i]->interpolate);CHKERRQ(ierr);
700c88c5224SJed Brown       ierr = PCMGGetRestriction(pc,i,&mglevels[i]->restrct);CHKERRQ(ierr);
701f3fbd535SBarry Smith     }
702f3fbd535SBarry Smith     for (i=0; i<n-1; i++) {
703f3fbd535SBarry Smith       if (!mglevels[i]->b) {
704f3fbd535SBarry Smith         Vec *vec;
7052a7a6963SBarry Smith         ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
706f3fbd535SBarry Smith         ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr);
7076bf464f9SBarry Smith         ierr = VecDestroy(vec);CHKERRQ(ierr);
708f3fbd535SBarry Smith         ierr = PetscFree(vec);CHKERRQ(ierr);
709f3fbd535SBarry Smith       }
710f3fbd535SBarry Smith       if (!mglevels[i]->r && i) {
711f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
712f3fbd535SBarry Smith         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
7136bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
714f3fbd535SBarry Smith       }
715f3fbd535SBarry Smith       if (!mglevels[i]->x) {
716f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
717f3fbd535SBarry Smith         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
7186bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
719f3fbd535SBarry Smith       }
720f3fbd535SBarry Smith     }
721f3fbd535SBarry Smith     if (n != 1 && !mglevels[n-1]->r) {
722f3fbd535SBarry Smith       /* PCMGSetR() on the finest level if user did not supply it */
723f3fbd535SBarry Smith       Vec *vec;
7242a7a6963SBarry Smith       ierr = KSPCreateVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
725f3fbd535SBarry Smith       ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr);
7266bf464f9SBarry Smith       ierr = VecDestroy(vec);CHKERRQ(ierr);
727f3fbd535SBarry Smith       ierr = PetscFree(vec);CHKERRQ(ierr);
728f3fbd535SBarry Smith     }
729f3fbd535SBarry Smith   }
730f3fbd535SBarry Smith 
73198e04984SBarry Smith   if (pc->dm) {
73298e04984SBarry Smith     /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */
73398e04984SBarry Smith     for (i=0; i<n-1; i++) {
73498e04984SBarry Smith       if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX;
73598e04984SBarry Smith     }
73698e04984SBarry Smith   }
737f3fbd535SBarry Smith 
738f3fbd535SBarry Smith   for (i=1; i<n; i++) {
739f3fbd535SBarry Smith     if (mglevels[i]->smoothu == mglevels[i]->smoothd) {
740f3fbd535SBarry Smith       /* if doing only down then initial guess is zero */
741f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
742f3fbd535SBarry Smith     }
74363e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
744f3fbd535SBarry Smith     ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr);
74563e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
746d42688cbSBarry Smith     if (!mglevels[i]->residual) {
747d42688cbSBarry Smith       Mat mat;
74823ee1639SBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&mat);CHKERRQ(ierr);
74954b2cd4bSJed Brown       ierr = PCMGSetResidual(pc,i,PCMGResidualDefault,mat);CHKERRQ(ierr);
750d42688cbSBarry Smith     }
751f3fbd535SBarry Smith   }
752f3fbd535SBarry Smith   for (i=1; i<n; i++) {
753f3fbd535SBarry Smith     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
754f3fbd535SBarry Smith       Mat          downmat,downpmat;
755f3fbd535SBarry Smith 
756f3fbd535SBarry Smith       /* check if operators have been set for up, if not use down operators to set them */
7570298fd71SBarry Smith       ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);CHKERRQ(ierr);
758f3fbd535SBarry Smith       if (!opsset) {
75923ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);CHKERRQ(ierr);
76023ee1639SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat);CHKERRQ(ierr);
761f3fbd535SBarry Smith       }
762f3fbd535SBarry Smith 
763f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
76463e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
765f3fbd535SBarry Smith       ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr);
76663e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
767f3fbd535SBarry Smith     }
768f3fbd535SBarry Smith   }
769f3fbd535SBarry Smith 
770f3fbd535SBarry Smith   /*
771f3fbd535SBarry Smith       If coarse solver is not direct method then DO NOT USE preonly
772f3fbd535SBarry Smith   */
773251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mglevels[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr);
774f3fbd535SBarry Smith   if (preonly) {
775251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr);
776251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr);
777251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr);
778251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)cpc,PCSVD,&svd);CHKERRQ(ierr);
779fd1303e9SJungho Lee     if (!lu && !redundant && !cholesky && !svd) {
780f3fbd535SBarry Smith       ierr = KSPSetType(mglevels[0]->smoothd,KSPGMRES);CHKERRQ(ierr);
781f3fbd535SBarry Smith     }
782f3fbd535SBarry Smith   }
783f3fbd535SBarry Smith 
784f3fbd535SBarry Smith   if (!pc->setupcalled) {
785f3fbd535SBarry Smith     ierr = KSPSetFromOptions(mglevels[0]->smoothd);CHKERRQ(ierr);
786f3fbd535SBarry Smith   }
787f3fbd535SBarry Smith 
78863e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
789f3fbd535SBarry Smith   ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr);
79063e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
791f3fbd535SBarry Smith 
792f3fbd535SBarry Smith   /*
793f3fbd535SBarry Smith      Dump the interpolation/restriction matrices plus the
794e3c5b3baSBarry Smith    Jacobian/stiffness on each level. This allows MATLAB users to
795f3fbd535SBarry Smith    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
796f3fbd535SBarry Smith 
797f3fbd535SBarry Smith    Only support one or the other at the same time.
798f3fbd535SBarry Smith   */
799f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
8000298fd71SBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);CHKERRQ(ierr);
801ce94432eSBarry Smith   if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc));
802f3fbd535SBarry Smith   dump = PETSC_FALSE;
803f3fbd535SBarry Smith #endif
8040298fd71SBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);CHKERRQ(ierr);
805ce94432eSBarry Smith   if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc));
806f3fbd535SBarry Smith 
807f3fbd535SBarry Smith   if (viewer) {
808f3fbd535SBarry Smith     for (i=1; i<n; i++) {
809f3fbd535SBarry Smith       ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr);
810f3fbd535SBarry Smith     }
811f3fbd535SBarry Smith     for (i=0; i<n; i++) {
812f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr);
813f3fbd535SBarry Smith       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
814f3fbd535SBarry Smith     }
815f3fbd535SBarry Smith   }
816f3fbd535SBarry Smith   PetscFunctionReturn(0);
817f3fbd535SBarry Smith }
818f3fbd535SBarry Smith 
819f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/
820f3fbd535SBarry Smith 
821f3fbd535SBarry Smith #undef __FUNCT__
8229dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetLevels"
8234b9ad928SBarry Smith /*@
82497177400SBarry Smith    PCMGGetLevels - Gets the number of levels to use with MG.
8254b9ad928SBarry Smith 
8264b9ad928SBarry Smith    Not Collective
8274b9ad928SBarry Smith 
8284b9ad928SBarry Smith    Input Parameter:
8294b9ad928SBarry Smith .  pc - the preconditioner context
8304b9ad928SBarry Smith 
8314b9ad928SBarry Smith    Output parameter:
8324b9ad928SBarry Smith .  levels - the number of levels
8334b9ad928SBarry Smith 
8344b9ad928SBarry Smith    Level: advanced
8354b9ad928SBarry Smith 
8364b9ad928SBarry Smith .keywords: MG, get, levels, multigrid
8374b9ad928SBarry Smith 
83897177400SBarry Smith .seealso: PCMGSetLevels()
8394b9ad928SBarry Smith @*/
8407087cfbeSBarry Smith PetscErrorCode  PCMGGetLevels(PC pc,PetscInt *levels)
8414b9ad928SBarry Smith {
842f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
8434b9ad928SBarry Smith 
8444b9ad928SBarry Smith   PetscFunctionBegin;
8450700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
8464482741eSBarry Smith   PetscValidIntPointer(levels,2);
847f3fbd535SBarry Smith   *levels = mg->nlevels;
8484b9ad928SBarry Smith   PetscFunctionReturn(0);
8494b9ad928SBarry Smith }
8504b9ad928SBarry Smith 
8514b9ad928SBarry Smith #undef __FUNCT__
8529dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetType"
8534b9ad928SBarry Smith /*@
85497177400SBarry Smith    PCMGSetType - Determines the form of multigrid to use:
8554b9ad928SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
8564b9ad928SBarry Smith 
857ad4df100SBarry Smith    Logically Collective on PC
8584b9ad928SBarry Smith 
8594b9ad928SBarry Smith    Input Parameters:
8604b9ad928SBarry Smith +  pc - the preconditioner context
8619dcbbd2bSBarry Smith -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
8629dcbbd2bSBarry Smith    PC_MG_FULL, PC_MG_KASKADE
8634b9ad928SBarry Smith 
8644b9ad928SBarry Smith    Options Database Key:
8654b9ad928SBarry Smith .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
8664b9ad928SBarry Smith    additive, full, kaskade
8674b9ad928SBarry Smith 
8684b9ad928SBarry Smith    Level: advanced
8694b9ad928SBarry Smith 
8704b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
8714b9ad928SBarry Smith 
87297177400SBarry Smith .seealso: PCMGSetLevels()
8734b9ad928SBarry Smith @*/
8747087cfbeSBarry Smith PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
8754b9ad928SBarry Smith {
876f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
8774b9ad928SBarry Smith 
8784b9ad928SBarry Smith   PetscFunctionBegin;
8790700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
880c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,form,2);
88131567311SBarry Smith   mg->am = form;
8829dcbbd2bSBarry Smith   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
883*c60c7ad4SBarry Smith   else pc->ops->applyrichardson = NULL;
884*c60c7ad4SBarry Smith   PetscFunctionReturn(0);
885*c60c7ad4SBarry Smith }
886*c60c7ad4SBarry Smith 
887*c60c7ad4SBarry Smith /*@
888*c60c7ad4SBarry Smith    PCMGGetType - Determines the form of multigrid to use:
889*c60c7ad4SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
890*c60c7ad4SBarry Smith 
891*c60c7ad4SBarry Smith    Logically Collective on PC
892*c60c7ad4SBarry Smith 
893*c60c7ad4SBarry Smith    Input Parameter:
894*c60c7ad4SBarry Smith .  pc - the preconditioner context
895*c60c7ad4SBarry Smith 
896*c60c7ad4SBarry Smith    Output Parameter:
897*c60c7ad4SBarry Smith .  type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE
898*c60c7ad4SBarry Smith 
899*c60c7ad4SBarry Smith 
900*c60c7ad4SBarry Smith    Level: advanced
901*c60c7ad4SBarry Smith 
902*c60c7ad4SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
903*c60c7ad4SBarry Smith 
904*c60c7ad4SBarry Smith .seealso: PCMGSetLevels()
905*c60c7ad4SBarry Smith @*/
906*c60c7ad4SBarry Smith PetscErrorCode  PCMGGetType(PC pc,PCMGType *type)
907*c60c7ad4SBarry Smith {
908*c60c7ad4SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
909*c60c7ad4SBarry Smith 
910*c60c7ad4SBarry Smith   PetscFunctionBegin;
911*c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
912*c60c7ad4SBarry Smith   *type = mg->am;
9134b9ad928SBarry Smith   PetscFunctionReturn(0);
9144b9ad928SBarry Smith }
9154b9ad928SBarry Smith 
9164b9ad928SBarry Smith #undef __FUNCT__
9170d353602SBarry Smith #define __FUNCT__ "PCMGSetCycleType"
9184b9ad928SBarry Smith /*@
9190d353602SBarry Smith    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
9204b9ad928SBarry Smith    complicated cycling.
9214b9ad928SBarry Smith 
922ad4df100SBarry Smith    Logically Collective on PC
9234b9ad928SBarry Smith 
9244b9ad928SBarry Smith    Input Parameters:
925c2be2410SBarry Smith +  pc - the multigrid context
9260d353602SBarry Smith -  PC_MG_CYCLE_V or PC_MG_CYCLE_W
9274b9ad928SBarry Smith 
9284b9ad928SBarry Smith    Options Database Key:
929b4519db0SPatrick Sanan $  -pc_mg_cycle_type <v,w>
9304b9ad928SBarry Smith 
9314b9ad928SBarry Smith    Level: advanced
9324b9ad928SBarry Smith 
9334b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
9344b9ad928SBarry Smith 
9350d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel()
9364b9ad928SBarry Smith @*/
9377087cfbeSBarry Smith PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
9384b9ad928SBarry Smith {
939f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
940f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
94179416396SBarry Smith   PetscInt     i,levels;
9424b9ad928SBarry Smith 
9434b9ad928SBarry Smith   PetscFunctionBegin;
9440700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
945ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
946c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
947f3fbd535SBarry Smith   levels = mglevels[0]->levels;
9484b9ad928SBarry Smith 
9492fa5cd67SKarl Rupp   for (i=0; i<levels; i++) mglevels[i]->cycles = n;
9504b9ad928SBarry Smith   PetscFunctionReturn(0);
9514b9ad928SBarry Smith }
9524b9ad928SBarry Smith 
9534b9ad928SBarry Smith #undef __FUNCT__
9548cc2d5dfSBarry Smith #define __FUNCT__ "PCMGMultiplicativeSetCycles"
9558cc2d5dfSBarry Smith /*@
9568cc2d5dfSBarry Smith    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
9578cc2d5dfSBarry Smith          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
9588cc2d5dfSBarry Smith 
959ad4df100SBarry Smith    Logically Collective on PC
9608cc2d5dfSBarry Smith 
9618cc2d5dfSBarry Smith    Input Parameters:
9628cc2d5dfSBarry Smith +  pc - the multigrid context
9638cc2d5dfSBarry Smith -  n - number of cycles (default is 1)
9648cc2d5dfSBarry Smith 
9658cc2d5dfSBarry Smith    Options Database Key:
9668cc2d5dfSBarry Smith $  -pc_mg_multiplicative_cycles n
9678cc2d5dfSBarry Smith 
9688cc2d5dfSBarry Smith    Level: advanced
9698cc2d5dfSBarry Smith 
9708cc2d5dfSBarry Smith    Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
9718cc2d5dfSBarry Smith 
9728cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
9738cc2d5dfSBarry Smith 
9748cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
9758cc2d5dfSBarry Smith @*/
9767087cfbeSBarry Smith PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
9778cc2d5dfSBarry Smith {
978f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
979f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
9808cc2d5dfSBarry Smith   PetscInt     i,levels;
9818cc2d5dfSBarry Smith 
9828cc2d5dfSBarry Smith   PetscFunctionBegin;
9830700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
984ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
985c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
986f3fbd535SBarry Smith   levels = mglevels[0]->levels;
9878cc2d5dfSBarry Smith 
9882fa5cd67SKarl Rupp   for (i=0; i<levels; i++) mg->cyclesperpcapply = n;
9898cc2d5dfSBarry Smith   PetscFunctionReturn(0);
9908cc2d5dfSBarry Smith }
9918cc2d5dfSBarry Smith 
9928cc2d5dfSBarry Smith #undef __FUNCT__
9939dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetGalerkin"
994c2be2410SBarry Smith /*@
99597177400SBarry Smith    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
99682b87d87SMatthew G. Knepley       finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i
997c2be2410SBarry Smith 
998ad4df100SBarry Smith    Logically Collective on PC
999c2be2410SBarry Smith 
1000c2be2410SBarry Smith    Input Parameters:
1001c91913d3SJed Brown +  pc - the multigrid context
1002c91913d3SJed Brown -  use - PETSC_TRUE to use the Galerkin process to compute coarse-level operators
1003c2be2410SBarry Smith 
1004c2be2410SBarry Smith    Options Database Key:
1005c2be2410SBarry Smith $  -pc_mg_galerkin
1006c2be2410SBarry Smith 
1007c2be2410SBarry Smith    Level: intermediate
1008c2be2410SBarry Smith 
1009c2be2410SBarry Smith .keywords: MG, set, Galerkin
1010c2be2410SBarry Smith 
10113fc8bf9cSBarry Smith .seealso: PCMGGetGalerkin()
10123fc8bf9cSBarry Smith 
1013c2be2410SBarry Smith @*/
1014c91913d3SJed Brown PetscErrorCode PCMGSetGalerkin(PC pc,PetscBool use)
1015c2be2410SBarry Smith {
1016f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
1017c2be2410SBarry Smith 
1018c2be2410SBarry Smith   PetscFunctionBegin;
10190700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1020789726d7SBarry Smith   mg->galerkin = use ? 1 : 0;
1021c2be2410SBarry Smith   PetscFunctionReturn(0);
1022c2be2410SBarry Smith }
1023c2be2410SBarry Smith 
1024c2be2410SBarry Smith #undef __FUNCT__
10253fc8bf9cSBarry Smith #define __FUNCT__ "PCMGGetGalerkin"
10263fc8bf9cSBarry Smith /*@
10273fc8bf9cSBarry Smith    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
102882b87d87SMatthew G. Knepley       A_i-1 = r_i * A_i * p_i
10293fc8bf9cSBarry Smith 
10303fc8bf9cSBarry Smith    Not Collective
10313fc8bf9cSBarry Smith 
10323fc8bf9cSBarry Smith    Input Parameter:
10333fc8bf9cSBarry Smith .  pc - the multigrid context
10343fc8bf9cSBarry Smith 
10353fc8bf9cSBarry Smith    Output Parameter:
10363fc8bf9cSBarry Smith .  gelerkin - PETSC_TRUE or PETSC_FALSE
10373fc8bf9cSBarry Smith 
10383fc8bf9cSBarry Smith    Options Database Key:
10393fc8bf9cSBarry Smith $  -pc_mg_galerkin
10403fc8bf9cSBarry Smith 
10413fc8bf9cSBarry Smith    Level: intermediate
10423fc8bf9cSBarry Smith 
10433fc8bf9cSBarry Smith .keywords: MG, set, Galerkin
10443fc8bf9cSBarry Smith 
10453fc8bf9cSBarry Smith .seealso: PCMGSetGalerkin()
10463fc8bf9cSBarry Smith 
10473fc8bf9cSBarry Smith @*/
10487087cfbeSBarry Smith PetscErrorCode  PCMGGetGalerkin(PC pc,PetscBool  *galerkin)
10493fc8bf9cSBarry Smith {
1050f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
10513fc8bf9cSBarry Smith 
10523fc8bf9cSBarry Smith   PetscFunctionBegin;
10530700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
105413fdb3acSJose Roman   *galerkin = (PetscBool)mg->galerkin;
10553fc8bf9cSBarry Smith   PetscFunctionReturn(0);
10563fc8bf9cSBarry Smith }
10573fc8bf9cSBarry Smith 
10583fc8bf9cSBarry Smith #undef __FUNCT__
10599dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothDown"
10604b9ad928SBarry Smith /*@
106197177400SBarry Smith    PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
106297177400SBarry Smith    use on all levels. Use PCMGGetSmootherDown() to set different
10634b9ad928SBarry Smith    pre-smoothing steps on different levels.
10644b9ad928SBarry Smith 
1065ad4df100SBarry Smith    Logically Collective on PC
10664b9ad928SBarry Smith 
10674b9ad928SBarry Smith    Input Parameters:
10684b9ad928SBarry Smith +  mg - the multigrid context
10694b9ad928SBarry Smith -  n - the number of smoothing steps
10704b9ad928SBarry Smith 
10714b9ad928SBarry Smith    Options Database Key:
10724b9ad928SBarry Smith .  -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps
10734b9ad928SBarry Smith 
10744b9ad928SBarry Smith    Level: advanced
10754b9ad928SBarry Smith 
10764b9ad928SBarry Smith .keywords: MG, smooth, down, pre-smoothing, steps, multigrid
10774b9ad928SBarry Smith 
107897177400SBarry Smith .seealso: PCMGSetNumberSmoothUp()
10794b9ad928SBarry Smith @*/
10807087cfbeSBarry Smith PetscErrorCode  PCMGSetNumberSmoothDown(PC pc,PetscInt n)
10814b9ad928SBarry Smith {
1082f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
1083f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
10846849ba73SBarry Smith   PetscErrorCode ierr;
108579416396SBarry Smith   PetscInt       i,levels;
10864b9ad928SBarry Smith 
10874b9ad928SBarry Smith   PetscFunctionBegin;
10880700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1089ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1090c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
1091f3fbd535SBarry Smith   levels = mglevels[0]->levels;
10924b9ad928SBarry Smith 
1093b05257ddSBarry Smith   for (i=1; i<levels; i++) {
10944b9ad928SBarry Smith     /* make sure smoother up and down are different */
10950298fd71SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,NULL);CHKERRQ(ierr);
1096f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
10972fa5cd67SKarl Rupp 
109831567311SBarry Smith     mg->default_smoothd = n;
10994b9ad928SBarry Smith   }
11004b9ad928SBarry Smith   PetscFunctionReturn(0);
11014b9ad928SBarry Smith }
11024b9ad928SBarry Smith 
11034b9ad928SBarry Smith #undef __FUNCT__
11049dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothUp"
11054b9ad928SBarry Smith /*@
110697177400SBarry Smith    PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
110797177400SBarry Smith    on all levels. Use PCMGGetSmootherUp() to set different numbers of
11084b9ad928SBarry Smith    post-smoothing steps on different levels.
11094b9ad928SBarry Smith 
1110ad4df100SBarry Smith    Logically Collective on PC
11114b9ad928SBarry Smith 
11124b9ad928SBarry Smith    Input Parameters:
11134b9ad928SBarry Smith +  mg - the multigrid context
11144b9ad928SBarry Smith -  n - the number of smoothing steps
11154b9ad928SBarry Smith 
11164b9ad928SBarry Smith    Options Database Key:
11174b9ad928SBarry Smith .  -pc_mg_smoothup <n> - Sets number of post-smoothing steps
11184b9ad928SBarry Smith 
11194b9ad928SBarry Smith    Level: advanced
11204b9ad928SBarry Smith 
11214b9ad928SBarry Smith    Note: this does not set a value on the coarsest grid, since we assume that
1122a8c7a070SBarry Smith     there is no separate smooth up on the coarsest grid.
11234b9ad928SBarry Smith 
11244b9ad928SBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid
11254b9ad928SBarry Smith 
112697177400SBarry Smith .seealso: PCMGSetNumberSmoothDown()
11274b9ad928SBarry Smith @*/
11287087cfbeSBarry Smith PetscErrorCode  PCMGSetNumberSmoothUp(PC pc,PetscInt n)
11294b9ad928SBarry Smith {
1130f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
1131f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
11326849ba73SBarry Smith   PetscErrorCode ierr;
113379416396SBarry Smith   PetscInt       i,levels;
11344b9ad928SBarry Smith 
11354b9ad928SBarry Smith   PetscFunctionBegin;
11360700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1137ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1138c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
1139f3fbd535SBarry Smith   levels = mglevels[0]->levels;
11404b9ad928SBarry Smith 
11414b9ad928SBarry Smith   for (i=1; i<levels; i++) {
11424b9ad928SBarry Smith     /* make sure smoother up and down are different */
11430298fd71SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,NULL);CHKERRQ(ierr);
1144f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
11452fa5cd67SKarl Rupp 
114631567311SBarry Smith     mg->default_smoothu = n;
11474b9ad928SBarry Smith   }
11484b9ad928SBarry Smith   PetscFunctionReturn(0);
11494b9ad928SBarry Smith }
11504b9ad928SBarry Smith 
11514b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/
11524b9ad928SBarry Smith 
11533b09bd56SBarry Smith /*MC
1154ccb205f8SBarry Smith    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
11553b09bd56SBarry Smith     information about the coarser grid matrices and restriction/interpolation operators.
11563b09bd56SBarry Smith 
11573b09bd56SBarry Smith    Options Database Keys:
11583b09bd56SBarry Smith +  -pc_mg_levels <nlevels> - number of levels including finest
1159b4519db0SPatrick Sanan .  -pc_mg_cycles <v,w> -
116079416396SBarry Smith .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
11613b09bd56SBarry Smith .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
11628c1c2452SJed Brown .  -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default
11633b09bd56SBarry Smith .  -pc_mg_log - log information about time spent on each level of the solver
11643b09bd56SBarry Smith .  -pc_mg_monitor - print information on the multigrid convergence
11658cf6037eSBarry Smith .  -pc_mg_galerkin - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R'
11668cf6037eSBarry Smith .  -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1)
11678cf6037eSBarry Smith .  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1168e3c5b3baSBarry Smith                         to the Socket viewer for reading from MATLAB.
11698cf6037eSBarry Smith -  -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices
11708cf6037eSBarry Smith                         to the binary output file called binaryoutput
11713b09bd56SBarry Smith 
117224c3aa18SBarry 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
11733b09bd56SBarry Smith 
11748cf6037eSBarry Smith        When run with a single level the smoother options are used on that level NOT the coarse grid solver options
11758cf6037eSBarry Smith 
11763b09bd56SBarry Smith    Level: intermediate
11773b09bd56SBarry Smith 
11788f87f92bSBarry Smith    Concepts: multigrid/multilevel
11793b09bd56SBarry Smith 
11808cf6037eSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE
11810d353602SBarry Smith            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
118297177400SBarry Smith            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
118397177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
11840d353602SBarry Smith            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
11853b09bd56SBarry Smith M*/
11863b09bd56SBarry Smith 
11874b9ad928SBarry Smith #undef __FUNCT__
11884b9ad928SBarry Smith #define __FUNCT__ "PCCreate_MG"
11898cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc)
11904b9ad928SBarry Smith {
1191f3fbd535SBarry Smith   PC_MG          *mg;
1192f3fbd535SBarry Smith   PetscErrorCode ierr;
1193f3fbd535SBarry Smith 
11944b9ad928SBarry Smith   PetscFunctionBegin;
1195b00a9115SJed Brown   ierr        = PetscNewLog(pc,&mg);CHKERRQ(ierr);
1196f3fbd535SBarry Smith   pc->data    = (void*)mg;
1197f3fbd535SBarry Smith   mg->nlevels = -1;
1198f3fbd535SBarry Smith 
119937a44384SMark Adams   pc->useAmat = PETSC_TRUE;
120037a44384SMark Adams 
12014b9ad928SBarry Smith   pc->ops->apply          = PCApply_MG;
12024b9ad928SBarry Smith   pc->ops->setup          = PCSetUp_MG;
1203a06653b4SBarry Smith   pc->ops->reset          = PCReset_MG;
12044b9ad928SBarry Smith   pc->ops->destroy        = PCDestroy_MG;
12054b9ad928SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_MG;
12064b9ad928SBarry Smith   pc->ops->view           = PCView_MG;
12074b9ad928SBarry Smith   PetscFunctionReturn(0);
12084b9ad928SBarry Smith }
1209