xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision ab83eea49124b1d1715a5b5e2ff780d4dd525b5f)
1dba47a55SKris Buschelman 
24b9ad928SBarry Smith /*
34b9ad928SBarry Smith     Defines the multigrid preconditioner interface.
44b9ad928SBarry Smith */
5af0996ceSBarry Smith #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);
225422a814eSBarry Smith     ierr = KSPSetErrorIfNotConverged(mglevels[i]->smoothd,pc->erroriffailure);CHKERRQ(ierr);
226336babb1SJed Brown     ierr = KSPSetType(mglevels[i]->smoothd,KSPCHEBYSHEV);CHKERRQ(ierr);
2270059c7bdSJed Brown     ierr = KSPSetConvergenceTest(mglevels[i]->smoothd,KSPConvergedSkip,NULL,NULL);CHKERRQ(ierr);
228336babb1SJed Brown     ierr = KSPSetNormType(mglevels[i]->smoothd,KSP_NORM_NONE);CHKERRQ(ierr);
229336babb1SJed Brown     ierr = KSPGetPC(mglevels[i]->smoothd,&ipc);CHKERRQ(ierr);
230336babb1SJed Brown     ierr = PCSetType(ipc,PCSOR);CHKERRQ(ierr);
231f3fbd535SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr);
232336babb1SJed Brown     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, i ? mg->default_smoothd : 1);CHKERRQ(ierr);
233f3fbd535SBarry Smith     ierr = KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);CHKERRQ(ierr);
234*ab83eea4SMatthew G. Knepley     ierr = PetscObjectComposedDataSetInt((PetscObject) mglevels[i]->smoothd, PetscMGLevelId, mglevels[i]->level);CHKERRQ(ierr);
235f3fbd535SBarry Smith 
236f3fbd535SBarry Smith     /* do special stuff for coarse grid */
237f3fbd535SBarry Smith     if (!i && levels > 1) {
238f3fbd535SBarry Smith       ierr = KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");CHKERRQ(ierr);
239f3fbd535SBarry Smith 
2409dbfc187SHong Zhang       /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */
241f3fbd535SBarry Smith       ierr = KSPSetType(mglevels[0]->smoothd,KSPPREONLY);CHKERRQ(ierr);
242f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[0]->smoothd,&ipc);CHKERRQ(ierr);
243f3fbd535SBarry Smith       ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
244f3fbd535SBarry Smith       if (size > 1) {
2459dbfc187SHong Zhang         KSP innerksp;
2469dbfc187SHong Zhang         PC  innerpc;
247f3fbd535SBarry Smith         ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr);
2489dbfc187SHong Zhang         ierr = PCRedundantGetKSP(ipc,&innerksp);CHKERRQ(ierr);
2499dbfc187SHong Zhang         ierr = KSPGetPC(innerksp,&innerpc);CHKERRQ(ierr);
2509dbfc187SHong Zhang         ierr = PCFactorSetShiftType(innerpc,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
251f3fbd535SBarry Smith       } else {
252f3fbd535SBarry Smith         ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr);
2539dbfc187SHong Zhang         ierr = PCFactorSetShiftType(ipc,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
254f3fbd535SBarry Smith       }
255f3fbd535SBarry Smith     } else {
256f3fbd535SBarry Smith       char tprefix[128];
257f3fbd535SBarry Smith       sprintf(tprefix,"mg_levels_%d_",(int)i);
258f3fbd535SBarry Smith       ierr = KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);CHKERRQ(ierr);
259f3fbd535SBarry Smith     }
2603bb1ff40SBarry Smith     ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[i]->smoothd);CHKERRQ(ierr);
2612fa5cd67SKarl Rupp 
262f3fbd535SBarry Smith     mglevels[i]->smoothu = mglevels[i]->smoothd;
26331567311SBarry Smith     mg->rtol             = 0.0;
26431567311SBarry Smith     mg->abstol           = 0.0;
26531567311SBarry Smith     mg->dtol             = 0.0;
26631567311SBarry Smith     mg->ttol             = 0.0;
26731567311SBarry Smith     mg->cyclesperpcapply = 1;
268f3fbd535SBarry Smith   }
26931567311SBarry Smith   mg->am                   = PC_MG_MULTIPLICATIVE;
270f3fbd535SBarry Smith   mg->levels               = mglevels;
2714b9ad928SBarry Smith   pc->ops->applyrichardson = PCApplyRichardson_MG;
2724b9ad928SBarry Smith   PetscFunctionReturn(0);
2734b9ad928SBarry Smith }
2744b9ad928SBarry Smith 
275c07bf074SBarry Smith 
276c07bf074SBarry Smith #undef __FUNCT__
277c07bf074SBarry Smith #define __FUNCT__ "PCDestroy_MG"
278c07bf074SBarry Smith PetscErrorCode PCDestroy_MG(PC pc)
279c07bf074SBarry Smith {
280c07bf074SBarry Smith   PetscErrorCode ierr;
281a06653b4SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
282a06653b4SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
283a06653b4SBarry Smith   PetscInt       i,n;
284c07bf074SBarry Smith 
285c07bf074SBarry Smith   PetscFunctionBegin;
286a06653b4SBarry Smith   ierr = PCReset_MG(pc);CHKERRQ(ierr);
287a06653b4SBarry Smith   if (mglevels) {
288a06653b4SBarry Smith     n = mglevels[0]->levels;
289a06653b4SBarry Smith     for (i=0; i<n; i++) {
290a06653b4SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
2916bf464f9SBarry Smith         ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
292a06653b4SBarry Smith       }
2936bf464f9SBarry Smith       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
294a06653b4SBarry Smith       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
295a06653b4SBarry Smith     }
296a06653b4SBarry Smith     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
297a06653b4SBarry Smith   }
298c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
299f3fbd535SBarry Smith   PetscFunctionReturn(0);
300f3fbd535SBarry Smith }
301f3fbd535SBarry Smith 
302f3fbd535SBarry Smith 
303f3fbd535SBarry Smith 
30409573ac7SBarry Smith extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**);
30509573ac7SBarry Smith extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**);
30609573ac7SBarry Smith extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**);
307f3fbd535SBarry Smith 
308f3fbd535SBarry Smith /*
309f3fbd535SBarry Smith    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
310f3fbd535SBarry Smith              or full cycle of multigrid.
311f3fbd535SBarry Smith 
312f3fbd535SBarry Smith   Note:
313f3fbd535SBarry Smith   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
314f3fbd535SBarry Smith */
315f3fbd535SBarry Smith #undef __FUNCT__
316f3fbd535SBarry Smith #define __FUNCT__ "PCApply_MG"
317f3fbd535SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
318f3fbd535SBarry Smith {
319f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
320f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
321f3fbd535SBarry Smith   PetscErrorCode ierr;
322f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
323f3fbd535SBarry Smith 
324f3fbd535SBarry Smith   PetscFunctionBegin;
325a9db87a2SMatthew G Knepley   if (mg->stageApply) {ierr = PetscLogStagePush(mg->stageApply);CHKERRQ(ierr);}
326e1d8e5deSBarry Smith   /* When the DM is supplying the matrix then it will not exist until here */
327298cc208SBarry Smith   for (i=0; i<levels; i++) {
328e1d8e5deSBarry Smith     if (!mglevels[i]->A) {
32923ee1639SBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr);
330298cc208SBarry Smith       ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr);
331e1d8e5deSBarry Smith     }
332e1d8e5deSBarry Smith   }
333e1d8e5deSBarry Smith 
334f3fbd535SBarry Smith   mglevels[levels-1]->b = b;
335f3fbd535SBarry Smith   mglevels[levels-1]->x = x;
33631567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
337f3fbd535SBarry Smith     ierr = VecSet(x,0.0);CHKERRQ(ierr);
33831567311SBarry Smith     for (i=0; i<mg->cyclesperpcapply; i++) {
3390298fd71SBarry Smith       ierr = PCMGMCycle_Private(pc,mglevels+levels-1,NULL);CHKERRQ(ierr);
340f3fbd535SBarry Smith     }
3412fa5cd67SKarl Rupp   } else if (mg->am == PC_MG_ADDITIVE) {
34231567311SBarry Smith     ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr);
3432fa5cd67SKarl Rupp   } else if (mg->am == PC_MG_KASKADE) {
34431567311SBarry Smith     ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr);
3452fa5cd67SKarl Rupp   } else {
346f3fbd535SBarry Smith     ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr);
347f3fbd535SBarry Smith   }
348a9db87a2SMatthew G Knepley   if (mg->stageApply) {ierr = PetscLogStagePop();CHKERRQ(ierr);}
349f3fbd535SBarry Smith   PetscFunctionReturn(0);
350f3fbd535SBarry Smith }
351f3fbd535SBarry Smith 
352f3fbd535SBarry Smith 
353f3fbd535SBarry Smith #undef __FUNCT__
354f3fbd535SBarry Smith #define __FUNCT__ "PCSetFromOptions_MG"
3558c34d3f5SBarry Smith PetscErrorCode PCSetFromOptions_MG(PetscOptions *PetscOptionsObject,PC pc)
356f3fbd535SBarry Smith {
357f3fbd535SBarry Smith   PetscErrorCode ierr;
358f3fbd535SBarry Smith   PetscInt       m,levels = 1,cycles;
359c91913d3SJed Brown   PetscBool      flg,set;
360f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
361f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
362f3fbd535SBarry Smith   PCMGType       mgtype;
363f3fbd535SBarry Smith   PCMGCycleType  mgctype;
364f3fbd535SBarry Smith 
365f3fbd535SBarry Smith   PetscFunctionBegin;
366e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Multigrid options");CHKERRQ(ierr);
3673aeaf226SBarry Smith   if (!mg->levels) {
368f3fbd535SBarry Smith     ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr);
369eb3f98d2SBarry Smith     if (!flg && pc->dm) {
370eb3f98d2SBarry Smith       ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
371eb3f98d2SBarry Smith       levels++;
3723aeaf226SBarry Smith       mg->usedmfornumberoflevels = PETSC_TRUE;
373eb3f98d2SBarry Smith     }
3740298fd71SBarry Smith     ierr = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr);
375f3fbd535SBarry Smith   }
3763aeaf226SBarry Smith   mglevels = mg->levels;
3773aeaf226SBarry Smith 
378f3fbd535SBarry Smith   mgctype = (PCMGCycleType) mglevels[0]->cycles;
379f3fbd535SBarry Smith   ierr    = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr);
380f3fbd535SBarry Smith   if (flg) {
381f3fbd535SBarry Smith     ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr);
3822fa5cd67SKarl Rupp   }
383f3fbd535SBarry Smith   flg  = PETSC_FALSE;
384c91913d3SJed Brown   ierr = PetscOptionsBool("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",flg,&flg,&set);CHKERRQ(ierr);
385c91913d3SJed Brown   if (set) {
386c91913d3SJed Brown     ierr = PCMGSetGalerkin(pc,flg);CHKERRQ(ierr);
387f3fbd535SBarry Smith   }
38894ae4db5SBarry Smith   ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",mg->default_smoothu,&m,&flg);CHKERRQ(ierr);
389f3fbd535SBarry Smith   if (flg) {
390f3fbd535SBarry Smith     ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr);
391f3fbd535SBarry Smith   }
39294ae4db5SBarry Smith   ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",mg->default_smoothd,&m,&flg);CHKERRQ(ierr);
393f3fbd535SBarry Smith   if (flg) {
394f3fbd535SBarry Smith     ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr);
395f3fbd535SBarry Smith   }
39631567311SBarry Smith   mgtype = mg->am;
397f3fbd535SBarry Smith   ierr   = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
398f3fbd535SBarry Smith   if (flg) {
399f3fbd535SBarry Smith     ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
400f3fbd535SBarry Smith   }
40131567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
40231567311SBarry Smith     ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGSetLevels",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr);
403f3fbd535SBarry Smith     if (flg) {
404f3fbd535SBarry Smith       ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr);
405f3fbd535SBarry Smith     }
406f3fbd535SBarry Smith   }
407f3fbd535SBarry Smith   flg  = PETSC_FALSE;
4080298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,NULL);CHKERRQ(ierr);
409f3fbd535SBarry Smith   if (flg) {
410f3fbd535SBarry Smith     PetscInt i;
411f3fbd535SBarry Smith     char     eventname[128];
412ce94432eSBarry Smith     if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
413f3fbd535SBarry Smith     levels = mglevels[0]->levels;
414f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
415f3fbd535SBarry Smith       sprintf(eventname,"MGSetup Level %d",(int)i);
41663e6d426SJed Brown       ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr);
417f3fbd535SBarry Smith       sprintf(eventname,"MGSmooth Level %d",(int)i);
41863e6d426SJed Brown       ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr);
419f3fbd535SBarry Smith       if (i) {
420f3fbd535SBarry Smith         sprintf(eventname,"MGResid Level %d",(int)i);
42163e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr);
422f3fbd535SBarry Smith         sprintf(eventname,"MGInterp Level %d",(int)i);
42363e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr);
424f3fbd535SBarry Smith       }
425f3fbd535SBarry Smith     }
426eec35531SMatthew G Knepley 
427ec60ca73SSatish Balay #if defined(PETSC_USE_LOG)
428eec35531SMatthew G Knepley     {
429eec35531SMatthew G Knepley       const char    *sname = "MG Apply";
430eec35531SMatthew G Knepley       PetscStageLog stageLog;
431eec35531SMatthew G Knepley       PetscInt      st;
432eec35531SMatthew G Knepley 
433eec35531SMatthew G Knepley       PetscFunctionBegin;
434eec35531SMatthew G Knepley       ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
435eec35531SMatthew G Knepley       for (st = 0; st < stageLog->numStages; ++st) {
436eec35531SMatthew G Knepley         PetscBool same;
437eec35531SMatthew G Knepley 
438eec35531SMatthew G Knepley         ierr = PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);CHKERRQ(ierr);
4392fa5cd67SKarl Rupp         if (same) mg->stageApply = st;
440eec35531SMatthew G Knepley       }
441eec35531SMatthew G Knepley       if (!mg->stageApply) {
442eec35531SMatthew G Knepley         ierr = PetscLogStageRegister(sname, &mg->stageApply);CHKERRQ(ierr);
443eec35531SMatthew G Knepley       }
444eec35531SMatthew G Knepley     }
445ec60ca73SSatish Balay #endif
446f3fbd535SBarry Smith   }
447f3fbd535SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
448f3fbd535SBarry Smith   PetscFunctionReturn(0);
449f3fbd535SBarry Smith }
450f3fbd535SBarry Smith 
4516a6fc655SJed Brown const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0};
4526a6fc655SJed Brown const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0};
453f3fbd535SBarry Smith 
4549804daf3SBarry Smith #include <petscdraw.h>
455f3fbd535SBarry Smith #undef __FUNCT__
456f3fbd535SBarry Smith #define __FUNCT__ "PCView_MG"
457f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
458f3fbd535SBarry Smith {
459f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
460f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
461f3fbd535SBarry Smith   PetscErrorCode ierr;
462e3deeeafSJed Brown   PetscInt       levels = mglevels ? mglevels[0]->levels : 0,i;
463219b1045SBarry Smith   PetscBool      iascii,isbinary,isdraw;
464f3fbd535SBarry Smith 
465f3fbd535SBarry Smith   PetscFunctionBegin;
466251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
4675b0b0462SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
468219b1045SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
469f3fbd535SBarry Smith   if (iascii) {
470e3deeeafSJed Brown     const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown";
471e3deeeafSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  MG: type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,cyclename);CHKERRQ(ierr);
47231567311SBarry Smith     if (mg->am == PC_MG_MULTIPLICATIVE) {
47331567311SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr);
474f3fbd535SBarry Smith     }
475218a07d4SBarry Smith     if (mg->galerkin) {
476f3fbd535SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
4774f66f45eSBarry Smith     } else {
4784f66f45eSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
479f3fbd535SBarry Smith     }
4805adeb434SBarry Smith     if (mg->view){
4815adeb434SBarry Smith       ierr = (*mg->view)(pc,viewer);CHKERRQ(ierr);
4825adeb434SBarry Smith     }
483f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
484f3fbd535SBarry Smith       if (!i) {
48589cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);CHKERRQ(ierr);
486f3fbd535SBarry Smith       } else {
48789cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
488f3fbd535SBarry Smith       }
489f3fbd535SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
490f3fbd535SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
491f3fbd535SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
492f3fbd535SBarry Smith       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
493f3fbd535SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
494f3fbd535SBarry Smith       } else if (i) {
49589cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
496f3fbd535SBarry Smith         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
497f3fbd535SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
498f3fbd535SBarry Smith         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
499f3fbd535SBarry Smith       }
500f3fbd535SBarry Smith     }
5015b0b0462SBarry Smith   } else if (isbinary) {
5025b0b0462SBarry Smith     for (i=levels-1; i>=0; i--) {
5035b0b0462SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
5045b0b0462SBarry Smith       if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) {
5055b0b0462SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
5065b0b0462SBarry Smith       }
5075b0b0462SBarry Smith     }
508219b1045SBarry Smith   } else if (isdraw) {
509219b1045SBarry Smith     PetscDraw draw;
510b4375e8dSPeter Brune     PetscReal x,w,y,bottom,th;
511219b1045SBarry Smith     ierr   = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
512219b1045SBarry Smith     ierr   = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
5130298fd71SBarry Smith     ierr   = PetscDrawStringGetSize(draw,NULL,&th);CHKERRQ(ierr);
514219b1045SBarry Smith     bottom = y - th;
515219b1045SBarry Smith     for (i=levels-1; i>=0; i--) {
516b4375e8dSPeter Brune       if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) {
517219b1045SBarry Smith         ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr);
518219b1045SBarry Smith         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
519219b1045SBarry Smith         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
520b4375e8dSPeter Brune       } else {
521b4375e8dSPeter Brune         w    = 0.5*PetscMin(1.0-x,x);
522b4375e8dSPeter Brune         ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr);
523b4375e8dSPeter Brune         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
524b4375e8dSPeter Brune         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
525b4375e8dSPeter Brune         ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr);
526b4375e8dSPeter Brune         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
527b4375e8dSPeter Brune         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
528b4375e8dSPeter Brune       }
5290298fd71SBarry Smith       ierr    = PetscDrawGetBoundingBox(draw,NULL,&bottom,NULL,NULL);CHKERRQ(ierr);
5301cd381d6SBarry Smith       bottom -= th;
531219b1045SBarry Smith     }
5325b0b0462SBarry Smith   }
533f3fbd535SBarry Smith   PetscFunctionReturn(0);
534f3fbd535SBarry Smith }
535f3fbd535SBarry Smith 
536af0996ceSBarry Smith #include <petsc/private/dmimpl.h>
537af0996ceSBarry Smith #include <petsc/private/kspimpl.h>
538cab2e9ccSBarry Smith 
539f3fbd535SBarry Smith /*
540f3fbd535SBarry Smith     Calls setup for the KSP on each level
541f3fbd535SBarry Smith */
542f3fbd535SBarry Smith #undef __FUNCT__
543f3fbd535SBarry Smith #define __FUNCT__ "PCSetUp_MG"
544f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc)
545f3fbd535SBarry Smith {
546f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
547f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
548f3fbd535SBarry Smith   PetscErrorCode ierr;
549f3fbd535SBarry Smith   PetscInt       i,n = mglevels[0]->levels;
55098e04984SBarry Smith   PC             cpc;
5513db492dfSBarry Smith   PetscBool      dump = PETSC_FALSE,opsset,use_amat,missinginterpolate = PETSC_FALSE;
552f3fbd535SBarry Smith   Mat            dA,dB;
553f3fbd535SBarry Smith   Vec            tvec;
554218a07d4SBarry Smith   DM             *dms;
555649052a6SBarry Smith   PetscViewer    viewer = 0;
556f3fbd535SBarry Smith 
557f3fbd535SBarry Smith   PetscFunctionBegin;
55801bc414fSDmitry Karpeev   /* FIX: Move this to PCSetFromOptions_MG? */
5593aeaf226SBarry Smith   if (mg->usedmfornumberoflevels) {
5603aeaf226SBarry Smith     PetscInt levels;
5613aeaf226SBarry Smith     ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
5623aeaf226SBarry Smith     levels++;
5633aeaf226SBarry Smith     if (levels > n) { /* the problem is now being solved on a finer grid */
5640298fd71SBarry Smith       ierr     = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr);
5653aeaf226SBarry Smith       n        = levels;
5663aeaf226SBarry Smith       ierr     = PCSetFromOptions(pc);CHKERRQ(ierr); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
5673aeaf226SBarry Smith       mglevels =  mg->levels;
5683aeaf226SBarry Smith     }
5693aeaf226SBarry Smith   }
57098e04984SBarry Smith   ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr);
5713aeaf226SBarry Smith 
572f3fbd535SBarry Smith 
573f3fbd535SBarry Smith   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
574f3fbd535SBarry Smith   /* so use those from global PC */
575f3fbd535SBarry Smith   /* Is this what we always want? What if user wants to keep old one? */
5760298fd71SBarry Smith   ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset);CHKERRQ(ierr);
57798e04984SBarry Smith   if (opsset) {
57898e04984SBarry Smith     Mat mmat;
57923ee1639SBarry Smith     ierr = KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat);CHKERRQ(ierr);
58098e04984SBarry Smith     if (mmat == pc->pmat) opsset = PETSC_FALSE;
58198e04984SBarry Smith   }
5824a5f07a5SMark F. Adams 
5834a5f07a5SMark F. Adams   if (!opsset) {
58471b23a65SMark F. Adams     ierr = PCGetUseAmat(pc,&use_amat);CHKERRQ(ierr);
5854a5f07a5SMark F. Adams     if(use_amat){
586f3fbd535SBarry 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);
58723ee1639SBarry Smith       ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat);CHKERRQ(ierr);
588f3fbd535SBarry Smith     }
5894a5f07a5SMark F. Adams     else {
5904a5f07a5SMark 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);
59123ee1639SBarry Smith       ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat);CHKERRQ(ierr);
5924a5f07a5SMark F. Adams     }
5934a5f07a5SMark F. Adams   }
594f3fbd535SBarry Smith 
5959c7ffb39SBarry Smith   for (i=n-1; i>0; i--) {
5969c7ffb39SBarry Smith     if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) {
5979c7ffb39SBarry Smith       missinginterpolate = PETSC_TRUE;
5989c7ffb39SBarry Smith       continue;
5999c7ffb39SBarry Smith     }
6009c7ffb39SBarry Smith   }
6019c7ffb39SBarry Smith   /*
6029c7ffb39SBarry Smith    Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS)
6039c7ffb39SBarry Smith    Skipping for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs?
6049c7ffb39SBarry Smith   */
6059c7ffb39SBarry Smith   if (missinginterpolate && pc->dm && mg->galerkin != 2 && !pc->setupcalled) {
6062d2b81a6SBarry Smith     /* construct the interpolation from the DMs */
6072e499ae9SBarry Smith     Mat p;
60873250ac0SBarry Smith     Vec rscale;
609785e854fSJed Brown     ierr     = PetscMalloc1(n,&dms);CHKERRQ(ierr);
610218a07d4SBarry Smith     dms[n-1] = pc->dm;
611ef1267afSMatthew G. Knepley     /* Separately create them so we do not get DMKSP interference between levels */
612ef1267afSMatthew G. Knepley     for (i=n-2; i>-1; i--) {ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr);}
613218a07d4SBarry Smith     for (i=n-2; i>-1; i--) {
614942e3340SBarry Smith       DMKSP kdm;
6153c0c59f3SBarry Smith       ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr);
6163c0c59f3SBarry Smith       if (mg->galerkin) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);}
617942e3340SBarry Smith       ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr);
618d1a3e677SJed Brown       /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take
619d1a3e677SJed Brown        * a bitwise OR of computing the matrix, RHS, and initial iterate. */
6200298fd71SBarry Smith       kdm->ops->computerhs = NULL;
6210298fd71SBarry Smith       kdm->rhsctx          = NULL;
62224c3aa18SBarry Smith       if (!mglevels[i+1]->interpolate) {
623e727c939SJed Brown         ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr);
6242d2b81a6SBarry Smith         ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr);
62500ac8be1SBarry Smith         if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);}
62673250ac0SBarry Smith         ierr = VecDestroy(&rscale);CHKERRQ(ierr);
6276bf464f9SBarry Smith         ierr = MatDestroy(&p);CHKERRQ(ierr);
628218a07d4SBarry Smith       }
62924c3aa18SBarry Smith     }
6302d2b81a6SBarry Smith 
631ef1267afSMatthew G. Knepley     for (i=n-2; i>-1; i--) {ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);}
6322d2b81a6SBarry Smith     ierr = PetscFree(dms);CHKERRQ(ierr);
633ce4cda84SJed Brown   }
634cab2e9ccSBarry Smith 
635ce4cda84SJed Brown   if (pc->dm && !pc->setupcalled) {
636ce4cda84SJed Brown     /* finest smoother also gets DM but it is not active, independent of whether galerkin==2 */
637cab2e9ccSBarry Smith     ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr);
638cab2e9ccSBarry Smith     ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr);
639218a07d4SBarry Smith   }
640218a07d4SBarry Smith 
641c91913d3SJed Brown   if (mg->galerkin == 1) {
642f3fbd535SBarry Smith     Mat B;
643f3fbd535SBarry Smith     /* currently only handle case where mat and pmat are the same on coarser levels */
64423ee1639SBarry Smith     ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB);CHKERRQ(ierr);
645f3fbd535SBarry Smith     if (!pc->setupcalled) {
646f3fbd535SBarry Smith       for (i=n-2; i>-1; i--) {
647b9367914SBarry 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");
648b9367914SBarry Smith         if (!mglevels[i+1]->interpolate) {
649b9367914SBarry Smith           ierr = PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct);CHKERRQ(ierr);
650b9367914SBarry Smith         }
651b9367914SBarry Smith         if (!mglevels[i+1]->restrct) {
652b9367914SBarry Smith           ierr = PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate);CHKERRQ(ierr);
653b9367914SBarry Smith         }
654b9367914SBarry Smith         if (mglevels[i+1]->interpolate == mglevels[i+1]->restrct) {
655f3fbd535SBarry Smith           ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
6567d537d92SJed Brown         } else {
6577d537d92SJed Brown           ierr = MatMatMatMult(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
6587d537d92SJed Brown         }
65923ee1639SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B);CHKERRQ(ierr);
660f3fbd535SBarry Smith         if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
661f3fbd535SBarry Smith         dB = B;
662f3fbd535SBarry Smith       }
663cd9507b2SBarry Smith       if (n > 1) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
664f3fbd535SBarry Smith     } else {
665f3fbd535SBarry Smith       for (i=n-2; i>-1; i--) {
666b9367914SBarry 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");
667b9367914SBarry Smith         if (!mglevels[i+1]->interpolate) {
668b9367914SBarry Smith           ierr = PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct);CHKERRQ(ierr);
669b9367914SBarry Smith         }
670b9367914SBarry Smith         if (!mglevels[i+1]->restrct) {
671b9367914SBarry Smith           ierr = PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate);CHKERRQ(ierr);
672b9367914SBarry Smith         }
67323ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&B);CHKERRQ(ierr);
674b9367914SBarry Smith         if (mglevels[i+1]->interpolate == mglevels[i+1]->restrct) {
675f3fbd535SBarry Smith           ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
6767d537d92SJed Brown         } else {
6777d537d92SJed Brown           ierr = MatMatMatMult(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
6787d537d92SJed Brown         }
67923ee1639SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B);CHKERRQ(ierr);
680f3fbd535SBarry Smith         dB   = B;
681f3fbd535SBarry Smith       }
682f3fbd535SBarry Smith     }
683ce4cda84SJed Brown   } else if (!mg->galerkin && pc->dm && pc->dm->x) {
684cab2e9ccSBarry Smith     /* need to restrict Jacobian location to coarser meshes for evaluation */
685cab2e9ccSBarry Smith     for (i=n-2; i>-1; i--) {
686c88c5224SJed Brown       Mat R;
687c88c5224SJed Brown       Vec rscale;
688cab2e9ccSBarry Smith       if (!mglevels[i]->smoothd->dm->x) {
689cab2e9ccSBarry Smith         Vec *vecs;
6902a7a6963SBarry Smith         ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vecs,0,NULL);CHKERRQ(ierr);
6912fa5cd67SKarl Rupp 
692cab2e9ccSBarry Smith         mglevels[i]->smoothd->dm->x = vecs[0];
6932fa5cd67SKarl Rupp 
694cab2e9ccSBarry Smith         ierr = PetscFree(vecs);CHKERRQ(ierr);
695cab2e9ccSBarry Smith       }
696c88c5224SJed Brown       ierr = PCMGGetRestriction(pc,i+1,&R);CHKERRQ(ierr);
697c88c5224SJed Brown       ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
698c88c5224SJed Brown       ierr = MatRestrict(R,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);CHKERRQ(ierr);
699c88c5224SJed Brown       ierr = VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,rscale);CHKERRQ(ierr);
700cab2e9ccSBarry Smith     }
701f3fbd535SBarry Smith   }
702ccceb50cSJed Brown   if (!mg->galerkin && pc->dm) {
703caa4e7f2SJed Brown     for (i=n-2; i>=0; i--) {
704ccceb50cSJed Brown       DM  dmfine,dmcoarse;
705caa4e7f2SJed Brown       Mat Restrict,Inject;
706caa4e7f2SJed Brown       Vec rscale;
707ccceb50cSJed Brown       ierr   = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr);
708ccceb50cSJed Brown       ierr   = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr);
709caa4e7f2SJed Brown       ierr   = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr);
710caa4e7f2SJed Brown       ierr   = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
7110298fd71SBarry Smith       Inject = NULL;      /* Callback should create it if it needs Injection */
712caa4e7f2SJed Brown       ierr   = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr);
713caa4e7f2SJed Brown     }
714caa4e7f2SJed Brown   }
715f3fbd535SBarry Smith 
716f3fbd535SBarry Smith   if (!pc->setupcalled) {
717f3fbd535SBarry Smith     for (i=0; i<n; i++) {
718f3fbd535SBarry Smith       ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr);
719f3fbd535SBarry Smith     }
720f3fbd535SBarry Smith     for (i=1; i<n; i++) {
721f3fbd535SBarry Smith       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
722f3fbd535SBarry Smith         ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr);
723f3fbd535SBarry Smith       }
724f3fbd535SBarry Smith     }
725f3fbd535SBarry Smith     for (i=1; i<n; i++) {
726c88c5224SJed Brown       ierr = PCMGGetInterpolation(pc,i,&mglevels[i]->interpolate);CHKERRQ(ierr);
727c88c5224SJed Brown       ierr = PCMGGetRestriction(pc,i,&mglevels[i]->restrct);CHKERRQ(ierr);
728f3fbd535SBarry Smith     }
729f3fbd535SBarry Smith     for (i=0; i<n-1; i++) {
730f3fbd535SBarry Smith       if (!mglevels[i]->b) {
731f3fbd535SBarry Smith         Vec *vec;
7322a7a6963SBarry Smith         ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
733f3fbd535SBarry Smith         ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr);
7346bf464f9SBarry Smith         ierr = VecDestroy(vec);CHKERRQ(ierr);
735f3fbd535SBarry Smith         ierr = PetscFree(vec);CHKERRQ(ierr);
736f3fbd535SBarry Smith       }
737f3fbd535SBarry Smith       if (!mglevels[i]->r && i) {
738f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
739f3fbd535SBarry Smith         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
7406bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
741f3fbd535SBarry Smith       }
742f3fbd535SBarry Smith       if (!mglevels[i]->x) {
743f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
744f3fbd535SBarry Smith         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
7456bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
746f3fbd535SBarry Smith       }
747f3fbd535SBarry Smith     }
748f3fbd535SBarry Smith     if (n != 1 && !mglevels[n-1]->r) {
749f3fbd535SBarry Smith       /* PCMGSetR() on the finest level if user did not supply it */
750f3fbd535SBarry Smith       Vec *vec;
7512a7a6963SBarry Smith       ierr = KSPCreateVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
752f3fbd535SBarry Smith       ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr);
7536bf464f9SBarry Smith       ierr = VecDestroy(vec);CHKERRQ(ierr);
754f3fbd535SBarry Smith       ierr = PetscFree(vec);CHKERRQ(ierr);
755f3fbd535SBarry Smith     }
756f3fbd535SBarry Smith   }
757f3fbd535SBarry Smith 
75898e04984SBarry Smith   if (pc->dm) {
75998e04984SBarry Smith     /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */
76098e04984SBarry Smith     for (i=0; i<n-1; i++) {
76198e04984SBarry Smith       if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX;
76298e04984SBarry Smith     }
76398e04984SBarry Smith   }
764f3fbd535SBarry Smith 
765f3fbd535SBarry Smith   for (i=1; i<n; i++) {
7662a21e185SBarry Smith     if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1){
767f3fbd535SBarry Smith       /* if doing only down then initial guess is zero */
768f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
769f3fbd535SBarry Smith     }
77063e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
771f3fbd535SBarry Smith     ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr);
77263e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
773d42688cbSBarry Smith     if (!mglevels[i]->residual) {
774d42688cbSBarry Smith       Mat mat;
77523ee1639SBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&mat);CHKERRQ(ierr);
77654b2cd4bSJed Brown       ierr = PCMGSetResidual(pc,i,PCMGResidualDefault,mat);CHKERRQ(ierr);
777d42688cbSBarry Smith     }
778f3fbd535SBarry Smith   }
779f3fbd535SBarry Smith   for (i=1; i<n; i++) {
780f3fbd535SBarry Smith     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
781f3fbd535SBarry Smith       Mat          downmat,downpmat;
782f3fbd535SBarry Smith 
783f3fbd535SBarry Smith       /* check if operators have been set for up, if not use down operators to set them */
7840298fd71SBarry Smith       ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);CHKERRQ(ierr);
785f3fbd535SBarry Smith       if (!opsset) {
78623ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);CHKERRQ(ierr);
78723ee1639SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat);CHKERRQ(ierr);
788f3fbd535SBarry Smith       }
789f3fbd535SBarry Smith 
790f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
79163e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
792f3fbd535SBarry Smith       ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr);
79363e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
794f3fbd535SBarry Smith     }
795f3fbd535SBarry Smith   }
796f3fbd535SBarry Smith 
79763e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
798f3fbd535SBarry Smith   ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr);
79963e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
800f3fbd535SBarry Smith 
801f3fbd535SBarry Smith   /*
802f3fbd535SBarry Smith      Dump the interpolation/restriction matrices plus the
803e3c5b3baSBarry Smith    Jacobian/stiffness on each level. This allows MATLAB users to
804f3fbd535SBarry Smith    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
805f3fbd535SBarry Smith 
806f3fbd535SBarry Smith    Only support one or the other at the same time.
807f3fbd535SBarry Smith   */
808f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
8090298fd71SBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);CHKERRQ(ierr);
810ce94432eSBarry Smith   if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc));
811f3fbd535SBarry Smith   dump = PETSC_FALSE;
812f3fbd535SBarry Smith #endif
8130298fd71SBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);CHKERRQ(ierr);
814ce94432eSBarry Smith   if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc));
815f3fbd535SBarry Smith 
816f3fbd535SBarry Smith   if (viewer) {
817f3fbd535SBarry Smith     for (i=1; i<n; i++) {
818f3fbd535SBarry Smith       ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr);
819f3fbd535SBarry Smith     }
820f3fbd535SBarry Smith     for (i=0; i<n; i++) {
821f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr);
822f3fbd535SBarry Smith       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
823f3fbd535SBarry Smith     }
824f3fbd535SBarry Smith   }
825f3fbd535SBarry Smith   PetscFunctionReturn(0);
826f3fbd535SBarry Smith }
827f3fbd535SBarry Smith 
828f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/
829f3fbd535SBarry Smith 
830f3fbd535SBarry Smith #undef __FUNCT__
8319dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetLevels"
8324b9ad928SBarry Smith /*@
83397177400SBarry Smith    PCMGGetLevels - Gets the number of levels to use with MG.
8344b9ad928SBarry Smith 
8354b9ad928SBarry Smith    Not Collective
8364b9ad928SBarry Smith 
8374b9ad928SBarry Smith    Input Parameter:
8384b9ad928SBarry Smith .  pc - the preconditioner context
8394b9ad928SBarry Smith 
8404b9ad928SBarry Smith    Output parameter:
8414b9ad928SBarry Smith .  levels - the number of levels
8424b9ad928SBarry Smith 
8434b9ad928SBarry Smith    Level: advanced
8444b9ad928SBarry Smith 
8454b9ad928SBarry Smith .keywords: MG, get, levels, multigrid
8464b9ad928SBarry Smith 
84797177400SBarry Smith .seealso: PCMGSetLevels()
8484b9ad928SBarry Smith @*/
8497087cfbeSBarry Smith PetscErrorCode  PCMGGetLevels(PC pc,PetscInt *levels)
8504b9ad928SBarry Smith {
851f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
8524b9ad928SBarry Smith 
8534b9ad928SBarry Smith   PetscFunctionBegin;
8540700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
8554482741eSBarry Smith   PetscValidIntPointer(levels,2);
856f3fbd535SBarry Smith   *levels = mg->nlevels;
8574b9ad928SBarry Smith   PetscFunctionReturn(0);
8584b9ad928SBarry Smith }
8594b9ad928SBarry Smith 
8604b9ad928SBarry Smith #undef __FUNCT__
8619dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetType"
8624b9ad928SBarry Smith /*@
86397177400SBarry Smith    PCMGSetType - Determines the form of multigrid to use:
8644b9ad928SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
8654b9ad928SBarry Smith 
866ad4df100SBarry Smith    Logically Collective on PC
8674b9ad928SBarry Smith 
8684b9ad928SBarry Smith    Input Parameters:
8694b9ad928SBarry Smith +  pc - the preconditioner context
8709dcbbd2bSBarry Smith -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
8719dcbbd2bSBarry Smith    PC_MG_FULL, PC_MG_KASKADE
8724b9ad928SBarry Smith 
8734b9ad928SBarry Smith    Options Database Key:
8744b9ad928SBarry Smith .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
8754b9ad928SBarry Smith    additive, full, kaskade
8764b9ad928SBarry Smith 
8774b9ad928SBarry Smith    Level: advanced
8784b9ad928SBarry Smith 
8794b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
8804b9ad928SBarry Smith 
88197177400SBarry Smith .seealso: PCMGSetLevels()
8824b9ad928SBarry Smith @*/
8837087cfbeSBarry Smith PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
8844b9ad928SBarry Smith {
885f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
8864b9ad928SBarry Smith 
8874b9ad928SBarry Smith   PetscFunctionBegin;
8880700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
889c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,form,2);
89031567311SBarry Smith   mg->am = form;
8919dcbbd2bSBarry Smith   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
892c60c7ad4SBarry Smith   else pc->ops->applyrichardson = NULL;
893c60c7ad4SBarry Smith   PetscFunctionReturn(0);
894c60c7ad4SBarry Smith }
895c60c7ad4SBarry Smith 
896c60c7ad4SBarry Smith /*@
897c60c7ad4SBarry Smith    PCMGGetType - Determines the form of multigrid to use:
898c60c7ad4SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
899c60c7ad4SBarry Smith 
900c60c7ad4SBarry Smith    Logically Collective on PC
901c60c7ad4SBarry Smith 
902c60c7ad4SBarry Smith    Input Parameter:
903c60c7ad4SBarry Smith .  pc - the preconditioner context
904c60c7ad4SBarry Smith 
905c60c7ad4SBarry Smith    Output Parameter:
906c60c7ad4SBarry Smith .  type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE
907c60c7ad4SBarry Smith 
908c60c7ad4SBarry Smith 
909c60c7ad4SBarry Smith    Level: advanced
910c60c7ad4SBarry Smith 
911c60c7ad4SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
912c60c7ad4SBarry Smith 
913c60c7ad4SBarry Smith .seealso: PCMGSetLevels()
914c60c7ad4SBarry Smith @*/
915c60c7ad4SBarry Smith PetscErrorCode  PCMGGetType(PC pc,PCMGType *type)
916c60c7ad4SBarry Smith {
917c60c7ad4SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
918c60c7ad4SBarry Smith 
919c60c7ad4SBarry Smith   PetscFunctionBegin;
920c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
921c60c7ad4SBarry Smith   *type = mg->am;
9224b9ad928SBarry Smith   PetscFunctionReturn(0);
9234b9ad928SBarry Smith }
9244b9ad928SBarry Smith 
9254b9ad928SBarry Smith #undef __FUNCT__
9260d353602SBarry Smith #define __FUNCT__ "PCMGSetCycleType"
9274b9ad928SBarry Smith /*@
9280d353602SBarry Smith    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
9294b9ad928SBarry Smith    complicated cycling.
9304b9ad928SBarry Smith 
931ad4df100SBarry Smith    Logically Collective on PC
9324b9ad928SBarry Smith 
9334b9ad928SBarry Smith    Input Parameters:
934c2be2410SBarry Smith +  pc - the multigrid context
9350d353602SBarry Smith -  PC_MG_CYCLE_V or PC_MG_CYCLE_W
9364b9ad928SBarry Smith 
9374b9ad928SBarry Smith    Options Database Key:
9384446f3b4SBarry Smith .  -pc_mg_cycle_type <v,w>
9394b9ad928SBarry Smith 
9404b9ad928SBarry Smith    Level: advanced
9414b9ad928SBarry Smith 
9424b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
9434b9ad928SBarry Smith 
9440d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel()
9454b9ad928SBarry Smith @*/
9467087cfbeSBarry Smith PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
9474b9ad928SBarry Smith {
948f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
949f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
95079416396SBarry Smith   PetscInt     i,levels;
9514b9ad928SBarry Smith 
9524b9ad928SBarry Smith   PetscFunctionBegin;
9530700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
954ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
95569fbec6eSBarry Smith   PetscValidLogicalCollectiveEnum(pc,n,2);
956f3fbd535SBarry Smith   levels = mglevels[0]->levels;
9574b9ad928SBarry Smith 
9582fa5cd67SKarl Rupp   for (i=0; i<levels; i++) mglevels[i]->cycles = n;
9594b9ad928SBarry Smith   PetscFunctionReturn(0);
9604b9ad928SBarry Smith }
9614b9ad928SBarry Smith 
9624b9ad928SBarry Smith #undef __FUNCT__
9638cc2d5dfSBarry Smith #define __FUNCT__ "PCMGMultiplicativeSetCycles"
9648cc2d5dfSBarry Smith /*@
9658cc2d5dfSBarry Smith    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
9668cc2d5dfSBarry Smith          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
9678cc2d5dfSBarry Smith 
968ad4df100SBarry Smith    Logically Collective on PC
9698cc2d5dfSBarry Smith 
9708cc2d5dfSBarry Smith    Input Parameters:
9718cc2d5dfSBarry Smith +  pc - the multigrid context
9728cc2d5dfSBarry Smith -  n - number of cycles (default is 1)
9738cc2d5dfSBarry Smith 
9748cc2d5dfSBarry Smith    Options Database Key:
975e1bc860dSBarry Smith .  -pc_mg_multiplicative_cycles n
9768cc2d5dfSBarry Smith 
9778cc2d5dfSBarry Smith    Level: advanced
9788cc2d5dfSBarry Smith 
9798cc2d5dfSBarry Smith    Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
9808cc2d5dfSBarry Smith 
9818cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
9828cc2d5dfSBarry Smith 
9838cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
9848cc2d5dfSBarry Smith @*/
9857087cfbeSBarry Smith PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
9868cc2d5dfSBarry Smith {
987f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
9888cc2d5dfSBarry Smith 
9898cc2d5dfSBarry Smith   PetscFunctionBegin;
9900700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
991c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
9922a21e185SBarry Smith   mg->cyclesperpcapply = n;
9938cc2d5dfSBarry Smith   PetscFunctionReturn(0);
9948cc2d5dfSBarry Smith }
9958cc2d5dfSBarry Smith 
9968cc2d5dfSBarry Smith #undef __FUNCT__
997fb15c04eSBarry Smith #define __FUNCT__ "PCMGSetGalerkin_MG"
998fb15c04eSBarry Smith PetscErrorCode PCMGSetGalerkin_MG(PC pc,PetscBool use)
999fb15c04eSBarry Smith {
1000fb15c04eSBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
1001fb15c04eSBarry Smith 
1002fb15c04eSBarry Smith   PetscFunctionBegin;
1003fb15c04eSBarry Smith   mg->galerkin = use ? 1 : 0;
1004fb15c04eSBarry Smith   PetscFunctionReturn(0);
1005fb15c04eSBarry Smith }
1006fb15c04eSBarry Smith 
1007fb15c04eSBarry Smith #undef __FUNCT__
10089dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetGalerkin"
1009c2be2410SBarry Smith /*@
101097177400SBarry Smith    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
101182b87d87SMatthew G. Knepley       finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i
1012c2be2410SBarry Smith 
1013ad4df100SBarry Smith    Logically Collective on PC
1014c2be2410SBarry Smith 
1015c2be2410SBarry Smith    Input Parameters:
1016c91913d3SJed Brown +  pc - the multigrid context
1017c91913d3SJed Brown -  use - PETSC_TRUE to use the Galerkin process to compute coarse-level operators
1018c2be2410SBarry Smith 
1019c2be2410SBarry Smith    Options Database Key:
10201c1aac46SBarry Smith .  -pc_mg_galerkin <true,false>
1021c2be2410SBarry Smith 
1022c2be2410SBarry Smith    Level: intermediate
1023c2be2410SBarry Smith 
10241c1aac46SBarry Smith    Notes: Some codes that use PCMG such as PCGAMG use Galerkin internally while constructing the hierarchy and thus do not
10251c1aac46SBarry Smith      use the PCMG construction of the coarser grids.
10261c1aac46SBarry Smith 
1027c2be2410SBarry Smith .keywords: MG, set, Galerkin
1028c2be2410SBarry Smith 
10293fc8bf9cSBarry Smith .seealso: PCMGGetGalerkin()
10303fc8bf9cSBarry Smith 
1031c2be2410SBarry Smith @*/
1032c91913d3SJed Brown PetscErrorCode PCMGSetGalerkin(PC pc,PetscBool use)
1033c2be2410SBarry Smith {
1034fb15c04eSBarry Smith   PetscErrorCode ierr;
1035c2be2410SBarry Smith 
1036c2be2410SBarry Smith   PetscFunctionBegin;
10370700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10387a4c0024SBarry Smith   ierr = PetscTryMethod(pc,"PCMGSetGalerkin_C",(PC,PetscBool),(pc,use));CHKERRQ(ierr);
1039c2be2410SBarry Smith   PetscFunctionReturn(0);
1040c2be2410SBarry Smith }
1041c2be2410SBarry Smith 
1042c2be2410SBarry Smith #undef __FUNCT__
10433fc8bf9cSBarry Smith #define __FUNCT__ "PCMGGetGalerkin"
10443fc8bf9cSBarry Smith /*@
10453fc8bf9cSBarry Smith    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
104682b87d87SMatthew G. Knepley       A_i-1 = r_i * A_i * p_i
10473fc8bf9cSBarry Smith 
10483fc8bf9cSBarry Smith    Not Collective
10493fc8bf9cSBarry Smith 
10503fc8bf9cSBarry Smith    Input Parameter:
10513fc8bf9cSBarry Smith .  pc - the multigrid context
10523fc8bf9cSBarry Smith 
10533fc8bf9cSBarry Smith    Output Parameter:
1054e1bc860dSBarry Smith .  galerkin - PETSC_TRUE or PETSC_FALSE
10553fc8bf9cSBarry Smith 
10563fc8bf9cSBarry Smith    Options Database Key:
1057e1bc860dSBarry Smith .  -pc_mg_galerkin
10583fc8bf9cSBarry Smith 
10593fc8bf9cSBarry Smith    Level: intermediate
10603fc8bf9cSBarry Smith 
10613fc8bf9cSBarry Smith .keywords: MG, set, Galerkin
10623fc8bf9cSBarry Smith 
10633fc8bf9cSBarry Smith .seealso: PCMGSetGalerkin()
10643fc8bf9cSBarry Smith 
10653fc8bf9cSBarry Smith @*/
10667087cfbeSBarry Smith PetscErrorCode  PCMGGetGalerkin(PC pc,PetscBool  *galerkin)
10673fc8bf9cSBarry Smith {
1068f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
10693fc8bf9cSBarry Smith 
10703fc8bf9cSBarry Smith   PetscFunctionBegin;
10710700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
107213fdb3acSJose Roman   *galerkin = (PetscBool)mg->galerkin;
10733fc8bf9cSBarry Smith   PetscFunctionReturn(0);
10743fc8bf9cSBarry Smith }
10753fc8bf9cSBarry Smith 
10763fc8bf9cSBarry Smith #undef __FUNCT__
10779dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothDown"
10784b9ad928SBarry Smith /*@
107997177400SBarry Smith    PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
108097177400SBarry Smith    use on all levels. Use PCMGGetSmootherDown() to set different
10814b9ad928SBarry Smith    pre-smoothing steps on different levels.
10824b9ad928SBarry Smith 
1083ad4df100SBarry Smith    Logically Collective on PC
10844b9ad928SBarry Smith 
10854b9ad928SBarry Smith    Input Parameters:
10864b9ad928SBarry Smith +  mg - the multigrid context
10874b9ad928SBarry Smith -  n - the number of smoothing steps
10884b9ad928SBarry Smith 
10894b9ad928SBarry Smith    Options Database Key:
10904b9ad928SBarry Smith .  -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps
10914b9ad928SBarry Smith 
10924b9ad928SBarry Smith    Level: advanced
10934b9ad928SBarry Smith 
10944b9ad928SBarry Smith .keywords: MG, smooth, down, pre-smoothing, steps, multigrid
10954b9ad928SBarry Smith 
109697177400SBarry Smith .seealso: PCMGSetNumberSmoothUp()
10974b9ad928SBarry Smith @*/
10987087cfbeSBarry Smith PetscErrorCode  PCMGSetNumberSmoothDown(PC pc,PetscInt n)
10994b9ad928SBarry Smith {
1100f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
1101f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
11026849ba73SBarry Smith   PetscErrorCode ierr;
110379416396SBarry Smith   PetscInt       i,levels;
11044b9ad928SBarry Smith 
11054b9ad928SBarry Smith   PetscFunctionBegin;
11060700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1107ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1108c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
1109f3fbd535SBarry Smith   levels = mglevels[0]->levels;
11104b9ad928SBarry Smith 
1111b05257ddSBarry Smith   for (i=1; i<levels; i++) {
11124b9ad928SBarry Smith     /* make sure smoother up and down are different */
11130298fd71SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,NULL);CHKERRQ(ierr);
1114f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
11152fa5cd67SKarl Rupp 
111631567311SBarry Smith     mg->default_smoothd = n;
11174b9ad928SBarry Smith   }
11184b9ad928SBarry Smith   PetscFunctionReturn(0);
11194b9ad928SBarry Smith }
11204b9ad928SBarry Smith 
11214b9ad928SBarry Smith #undef __FUNCT__
11229dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothUp"
11234b9ad928SBarry Smith /*@
112497177400SBarry Smith    PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
112597177400SBarry Smith    on all levels. Use PCMGGetSmootherUp() to set different numbers of
11264b9ad928SBarry Smith    post-smoothing steps on different levels.
11274b9ad928SBarry Smith 
1128ad4df100SBarry Smith    Logically Collective on PC
11294b9ad928SBarry Smith 
11304b9ad928SBarry Smith    Input Parameters:
11314b9ad928SBarry Smith +  mg - the multigrid context
11324b9ad928SBarry Smith -  n - the number of smoothing steps
11334b9ad928SBarry Smith 
11344b9ad928SBarry Smith    Options Database Key:
11354b9ad928SBarry Smith .  -pc_mg_smoothup <n> - Sets number of post-smoothing steps
11364b9ad928SBarry Smith 
11374b9ad928SBarry Smith    Level: advanced
11384b9ad928SBarry Smith 
11394b9ad928SBarry Smith    Note: this does not set a value on the coarsest grid, since we assume that
1140a8c7a070SBarry Smith     there is no separate smooth up on the coarsest grid.
11414b9ad928SBarry Smith 
11424b9ad928SBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid
11434b9ad928SBarry Smith 
114497177400SBarry Smith .seealso: PCMGSetNumberSmoothDown()
11454b9ad928SBarry Smith @*/
11467087cfbeSBarry Smith PetscErrorCode  PCMGSetNumberSmoothUp(PC pc,PetscInt n)
11474b9ad928SBarry Smith {
1148f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
1149f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
11506849ba73SBarry Smith   PetscErrorCode ierr;
115179416396SBarry Smith   PetscInt       i,levels;
11524b9ad928SBarry Smith 
11534b9ad928SBarry Smith   PetscFunctionBegin;
11540700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1155ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1156c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
1157f3fbd535SBarry Smith   levels = mglevels[0]->levels;
11584b9ad928SBarry Smith 
11594b9ad928SBarry Smith   for (i=1; i<levels; i++) {
11604b9ad928SBarry Smith     /* make sure smoother up and down are different */
11610298fd71SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,NULL);CHKERRQ(ierr);
1162f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
11632fa5cd67SKarl Rupp 
116431567311SBarry Smith     mg->default_smoothu = n;
11654b9ad928SBarry Smith   }
11664b9ad928SBarry Smith   PetscFunctionReturn(0);
11674b9ad928SBarry Smith }
11684b9ad928SBarry Smith 
11694b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/
11704b9ad928SBarry Smith 
11713b09bd56SBarry Smith /*MC
1172ccb205f8SBarry Smith    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
11733b09bd56SBarry Smith     information about the coarser grid matrices and restriction/interpolation operators.
11743b09bd56SBarry Smith 
11753b09bd56SBarry Smith    Options Database Keys:
11763b09bd56SBarry Smith +  -pc_mg_levels <nlevels> - number of levels including finest
1177b4519db0SPatrick Sanan .  -pc_mg_cycles <v,w> -
117879416396SBarry Smith .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
11793b09bd56SBarry Smith .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
11808c1c2452SJed Brown .  -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default
11813b09bd56SBarry Smith .  -pc_mg_log - log information about time spent on each level of the solver
11828cf6037eSBarry Smith .  -pc_mg_galerkin - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R'
11838cf6037eSBarry Smith .  -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1)
11848cf6037eSBarry Smith .  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1185e3c5b3baSBarry Smith                         to the Socket viewer for reading from MATLAB.
11868cf6037eSBarry Smith -  -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices
11878cf6037eSBarry Smith                         to the binary output file called binaryoutput
11883b09bd56SBarry Smith 
118924c3aa18SBarry 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
11903b09bd56SBarry Smith 
11918cf6037eSBarry Smith        When run with a single level the smoother options are used on that level NOT the coarse grid solver options
11928cf6037eSBarry Smith 
11933b09bd56SBarry Smith    Level: intermediate
11943b09bd56SBarry Smith 
11958f87f92bSBarry Smith    Concepts: multigrid/multilevel
11963b09bd56SBarry Smith 
11978cf6037eSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE
11980d353602SBarry Smith            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
119997177400SBarry Smith            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
120097177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
12010d353602SBarry Smith            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
12023b09bd56SBarry Smith M*/
12033b09bd56SBarry Smith 
12044b9ad928SBarry Smith #undef __FUNCT__
12054b9ad928SBarry Smith #define __FUNCT__ "PCCreate_MG"
12068cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc)
12074b9ad928SBarry Smith {
1208f3fbd535SBarry Smith   PC_MG          *mg;
1209f3fbd535SBarry Smith   PetscErrorCode ierr;
1210f3fbd535SBarry Smith 
12114b9ad928SBarry Smith   PetscFunctionBegin;
1212b00a9115SJed Brown   ierr        = PetscNewLog(pc,&mg);CHKERRQ(ierr);
1213f3fbd535SBarry Smith   pc->data    = (void*)mg;
1214f3fbd535SBarry Smith   mg->nlevels = -1;
1215f3fbd535SBarry Smith 
121637a44384SMark Adams   pc->useAmat = PETSC_TRUE;
121737a44384SMark Adams 
12184b9ad928SBarry Smith   pc->ops->apply          = PCApply_MG;
12194b9ad928SBarry Smith   pc->ops->setup          = PCSetUp_MG;
1220a06653b4SBarry Smith   pc->ops->reset          = PCReset_MG;
12214b9ad928SBarry Smith   pc->ops->destroy        = PCDestroy_MG;
12224b9ad928SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_MG;
12234b9ad928SBarry Smith   pc->ops->view           = PCView_MG;
1224fb15c04eSBarry Smith 
1225fb15c04eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",PCMGSetGalerkin_MG);CHKERRQ(ierr);
12264b9ad928SBarry Smith   PetscFunctionReturn(0);
12274b9ad928SBarry Smith }
1228