xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision e941fc338b6dd16e4abb5ab29938f3746ca3290e)
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;
16a0808db4SHong Zhang   PC             subpc;
17a0808db4SHong Zhang   PCFailedReason pcreason;
184b9ad928SBarry Smith 
194b9ad928SBarry Smith   PetscFunctionBegin;
2063e6d426SJed Brown   if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
2131567311SBarry Smith   ierr = KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x);CHKERRQ(ierr);  /* pre-smooth */
22a0808db4SHong Zhang   ierr = KSPGetPC(mglevels->smoothd,&subpc);CHKERRQ(ierr);
23a0808db4SHong Zhang   ierr = PCGetSetUpFailedReason(subpc,&pcreason);CHKERRQ(ierr);
24a0808db4SHong Zhang   if (pcreason) {
25a0808db4SHong Zhang     pc->failedreason = PC_SUBPC_ERROR;
26a0808db4SHong Zhang   }
2763e6d426SJed Brown   if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
2831567311SBarry Smith   if (mglevels->level) {  /* not the coarsest grid */
2963e6d426SJed Brown     if (mglevels->eventresidual) {ierr = PetscLogEventBegin(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);}
3031567311SBarry Smith     ierr = (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);CHKERRQ(ierr);
3163e6d426SJed Brown     if (mglevels->eventresidual) {ierr = PetscLogEventEnd(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);}
324b9ad928SBarry Smith 
334b9ad928SBarry Smith     /* if on finest level and have convergence criteria set */
3431567311SBarry Smith     if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) {
354b9ad928SBarry Smith       PetscReal rnorm;
3631567311SBarry Smith       ierr = VecNorm(mglevels->r,NORM_2,&rnorm);CHKERRQ(ierr);
374b9ad928SBarry Smith       if (rnorm <= mg->ttol) {
3870441072SBarry Smith         if (rnorm < mg->abstol) {
394d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_ATOL;
4057622a8eSBarry 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);
414b9ad928SBarry Smith         } else {
424d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_RTOL;
4357622a8eSBarry 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);
444b9ad928SBarry Smith         }
454b9ad928SBarry Smith         PetscFunctionReturn(0);
464b9ad928SBarry Smith       }
474b9ad928SBarry Smith     }
484b9ad928SBarry Smith 
4931567311SBarry Smith     mgc = *(mglevelsin - 1);
5063e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
5131567311SBarry Smith     ierr = MatRestrict(mglevels->restrct,mglevels->r,mgc->b);CHKERRQ(ierr);
5263e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
53efb30889SBarry Smith     ierr = VecSet(mgc->x,0.0);CHKERRQ(ierr);
544b9ad928SBarry Smith     while (cycles--) {
5531567311SBarry Smith       ierr = PCMGMCycle_Private(pc,mglevelsin-1,reason);CHKERRQ(ierr);
564b9ad928SBarry Smith     }
5763e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
5831567311SBarry Smith     ierr = MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x);CHKERRQ(ierr);
5963e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
6063e6d426SJed Brown     if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
6131567311SBarry Smith     ierr = KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x);CHKERRQ(ierr);    /* post smooth */
6263e6d426SJed Brown     if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
634b9ad928SBarry Smith   }
644b9ad928SBarry Smith   PetscFunctionReturn(0);
654b9ad928SBarry Smith }
664b9ad928SBarry Smith 
674b9ad928SBarry Smith #undef __FUNCT__
684b9ad928SBarry Smith #define __FUNCT__ "PCApplyRichardson_MG"
69ace3abfcSBarry 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)
704b9ad928SBarry Smith {
71f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
72f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
73dfbe8321SBarry Smith   PetscErrorCode ierr;
74f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
754b9ad928SBarry Smith 
764b9ad928SBarry Smith   PetscFunctionBegin;
77a762d673SBarry Smith   /* When the DM is supplying the matrix then it will not exist until here */
78a762d673SBarry Smith   for (i=0; i<levels; i++) {
79a762d673SBarry Smith     if (!mglevels[i]->A) {
80a762d673SBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr);
81a762d673SBarry Smith       ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr);
82a762d673SBarry Smith     }
83a762d673SBarry Smith   }
84f3fbd535SBarry Smith   mglevels[levels-1]->b = b;
85f3fbd535SBarry Smith   mglevels[levels-1]->x = x;
864b9ad928SBarry Smith 
8731567311SBarry Smith   mg->rtol   = rtol;
8831567311SBarry Smith   mg->abstol = abstol;
8931567311SBarry Smith   mg->dtol   = dtol;
904b9ad928SBarry Smith   if (rtol) {
914b9ad928SBarry Smith     /* compute initial residual norm for relative convergence test */
924b9ad928SBarry Smith     PetscReal rnorm;
937319c654SBarry Smith     if (zeroguess) {
947319c654SBarry Smith       ierr = VecNorm(b,NORM_2,&rnorm);CHKERRQ(ierr);
957319c654SBarry Smith     } else {
96f3fbd535SBarry Smith       ierr = (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);CHKERRQ(ierr);
974b9ad928SBarry Smith       ierr = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr);
987319c654SBarry Smith     }
9931567311SBarry Smith     mg->ttol = PetscMax(rtol*rnorm,abstol);
1002fa5cd67SKarl Rupp   } else if (abstol) mg->ttol = abstol;
1012fa5cd67SKarl Rupp   else mg->ttol = 0.0;
1024b9ad928SBarry Smith 
1034d0a8057SBarry Smith   /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't
104336babb1SJed Brown      stop prematurely due to small residual */
1054d0a8057SBarry Smith   for (i=1; i<levels; i++) {
106f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
107f3fbd535SBarry Smith     if (mglevels[i]->smoothu != mglevels[i]->smoothd) {
108f3fbd535SBarry Smith       ierr = KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
1094b9ad928SBarry Smith     }
1104d0a8057SBarry Smith   }
1114d0a8057SBarry Smith 
1124d0a8057SBarry Smith   *reason = (PCRichardsonConvergedReason)0;
1134d0a8057SBarry Smith   for (i=0; i<its; i++) {
114f3fbd535SBarry Smith     ierr = PCMGMCycle_Private(pc,mglevels+levels-1,reason);CHKERRQ(ierr);
1154d0a8057SBarry Smith     if (*reason) break;
1164d0a8057SBarry Smith   }
1174d0a8057SBarry Smith   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
1184d0a8057SBarry Smith   *outits = i;
1194b9ad928SBarry Smith   PetscFunctionReturn(0);
1204b9ad928SBarry Smith }
1214b9ad928SBarry Smith 
1224b9ad928SBarry Smith #undef __FUNCT__
1233aeaf226SBarry Smith #define __FUNCT__ "PCReset_MG"
1243aeaf226SBarry Smith PetscErrorCode PCReset_MG(PC pc)
1253aeaf226SBarry Smith {
1263aeaf226SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
1273aeaf226SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
1283aeaf226SBarry Smith   PetscErrorCode ierr;
1293aeaf226SBarry Smith   PetscInt       i,n;
1303aeaf226SBarry Smith 
1313aeaf226SBarry Smith   PetscFunctionBegin;
1323aeaf226SBarry Smith   if (mglevels) {
1333aeaf226SBarry Smith     n = mglevels[0]->levels;
1343aeaf226SBarry Smith     for (i=0; i<n-1; i++) {
1353aeaf226SBarry Smith       ierr = VecDestroy(&mglevels[i+1]->r);CHKERRQ(ierr);
1363aeaf226SBarry Smith       ierr = VecDestroy(&mglevels[i]->b);CHKERRQ(ierr);
1373aeaf226SBarry Smith       ierr = VecDestroy(&mglevels[i]->x);CHKERRQ(ierr);
1383aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i+1]->restrct);CHKERRQ(ierr);
1393aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i+1]->interpolate);CHKERRQ(ierr);
14073250ac0SBarry Smith       ierr = VecDestroy(&mglevels[i+1]->rscale);CHKERRQ(ierr);
1413aeaf226SBarry Smith     }
1423aeaf226SBarry Smith 
1433aeaf226SBarry Smith     for (i=0; i<n; i++) {
1443aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i]->A);CHKERRQ(ierr);
1453aeaf226SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
1463aeaf226SBarry Smith         ierr = KSPReset(mglevels[i]->smoothd);CHKERRQ(ierr);
1473aeaf226SBarry Smith       }
1483aeaf226SBarry Smith       ierr = KSPReset(mglevels[i]->smoothu);CHKERRQ(ierr);
1493aeaf226SBarry Smith     }
1503aeaf226SBarry Smith   }
1513aeaf226SBarry Smith   PetscFunctionReturn(0);
1523aeaf226SBarry Smith }
1533aeaf226SBarry Smith 
1543aeaf226SBarry Smith #undef __FUNCT__
1559dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetLevels"
1564b9ad928SBarry Smith /*@C
15797177400SBarry Smith    PCMGSetLevels - Sets the number of levels to use with MG.
1584b9ad928SBarry Smith    Must be called before any other MG routine.
1594b9ad928SBarry Smith 
160ad4df100SBarry Smith    Logically Collective on PC
1614b9ad928SBarry Smith 
1624b9ad928SBarry Smith    Input Parameters:
1634b9ad928SBarry Smith +  pc - the preconditioner context
1644b9ad928SBarry Smith .  levels - the number of levels
1654b9ad928SBarry Smith -  comms - optional communicators for each level; this is to allow solving the coarser problems
1660298fd71SBarry Smith            on smaller sets of processors. Use NULL_OBJECT for default in Fortran
1674b9ad928SBarry Smith 
1684b9ad928SBarry Smith    Level: intermediate
1694b9ad928SBarry Smith 
1704b9ad928SBarry Smith    Notes:
1714b9ad928SBarry Smith      If the number of levels is one then the multigrid uses the -mg_levels prefix
1724b9ad928SBarry Smith   for setting the level options rather than the -mg_coarse prefix.
1734b9ad928SBarry Smith 
1744b9ad928SBarry Smith .keywords: MG, set, levels, multigrid
1754b9ad928SBarry Smith 
17697177400SBarry Smith .seealso: PCMGSetType(), PCMGGetLevels()
1774b9ad928SBarry Smith @*/
1787087cfbeSBarry Smith PetscErrorCode  PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
1794b9ad928SBarry Smith {
180dfbe8321SBarry Smith   PetscErrorCode ierr;
181f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
182ce94432eSBarry Smith   MPI_Comm       comm;
1833aeaf226SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
18410eca3edSLisandro Dalcin   PCMGType       mgtype     = mg->am;
18510167fecSBarry Smith   PetscInt       mgctype    = (PetscInt) PC_MG_CYCLE_V;
186f3fbd535SBarry Smith   PetscInt       i;
187f3fbd535SBarry Smith   PetscMPIInt    size;
188f3fbd535SBarry Smith   const char     *prefix;
189f3fbd535SBarry Smith   PC             ipc;
1903aeaf226SBarry Smith   PetscInt       n;
1914b9ad928SBarry Smith 
1924b9ad928SBarry Smith   PetscFunctionBegin;
1930700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
194548767e0SJed Brown   PetscValidLogicalCollectiveInt(pc,levels,2);
195ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
196548767e0SJed Brown   if (mg->nlevels == levels) PetscFunctionReturn(0);
1973aeaf226SBarry Smith   if (mglevels) {
19810eca3edSLisandro Dalcin     mgctype = mglevels[0]->cycles;
1993aeaf226SBarry Smith     /* changing the number of levels so free up the previous stuff */
2003aeaf226SBarry Smith     ierr = PCReset_MG(pc);CHKERRQ(ierr);
2013aeaf226SBarry Smith     n    = mglevels[0]->levels;
2023aeaf226SBarry Smith     for (i=0; i<n; i++) {
2033aeaf226SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
2043aeaf226SBarry Smith         ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
2053aeaf226SBarry Smith       }
2063aeaf226SBarry Smith       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
2073aeaf226SBarry Smith       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
2083aeaf226SBarry Smith     }
2093aeaf226SBarry Smith     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
2103aeaf226SBarry Smith   }
211f3fbd535SBarry Smith 
212f3fbd535SBarry Smith   mg->nlevels = levels;
213f3fbd535SBarry Smith 
214785e854fSJed Brown   ierr = PetscMalloc1(levels,&mglevels);CHKERRQ(ierr);
2153bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr);
216f3fbd535SBarry Smith 
217f3fbd535SBarry Smith   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
218f3fbd535SBarry Smith 
219a9db87a2SMatthew G Knepley   mg->stageApply = 0;
220f3fbd535SBarry Smith   for (i=0; i<levels; i++) {
221b00a9115SJed Brown     ierr = PetscNewLog(pc,&mglevels[i]);CHKERRQ(ierr);
2222fa5cd67SKarl Rupp 
223f3fbd535SBarry Smith     mglevels[i]->level               = i;
224f3fbd535SBarry Smith     mglevels[i]->levels              = levels;
22510eca3edSLisandro Dalcin     mglevels[i]->cycles              = mgctype;
226336babb1SJed Brown     mg->default_smoothu              = 2;
227336babb1SJed Brown     mg->default_smoothd              = 2;
22863e6d426SJed Brown     mglevels[i]->eventsmoothsetup    = 0;
22963e6d426SJed Brown     mglevels[i]->eventsmoothsolve    = 0;
23063e6d426SJed Brown     mglevels[i]->eventresidual       = 0;
23163e6d426SJed Brown     mglevels[i]->eventinterprestrict = 0;
232f3fbd535SBarry Smith 
233f3fbd535SBarry Smith     if (comms) comm = comms[i];
234f3fbd535SBarry Smith     ierr = KSPCreate(comm,&mglevels[i]->smoothd);CHKERRQ(ierr);
235422a814eSBarry Smith     ierr = KSPSetErrorIfNotConverged(mglevels[i]->smoothd,pc->erroriffailure);CHKERRQ(ierr);
2360ee9a668SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr);
2370ee9a668SBarry Smith     ierr = KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);CHKERRQ(ierr);
2380ee9a668SBarry Smith     ierr = PetscObjectComposedDataSetInt((PetscObject) mglevels[i]->smoothd, PetscMGLevelId, mglevels[i]->level);CHKERRQ(ierr);
2390ee9a668SBarry Smith     if (i || levels == 1) {
2400ee9a668SBarry Smith       char tprefix[128];
2410ee9a668SBarry Smith 
242336babb1SJed Brown       ierr = KSPSetType(mglevels[i]->smoothd,KSPCHEBYSHEV);CHKERRQ(ierr);
2430059c7bdSJed Brown       ierr = KSPSetConvergenceTest(mglevels[i]->smoothd,KSPConvergedSkip,NULL,NULL);CHKERRQ(ierr);
244336babb1SJed Brown       ierr = KSPSetNormType(mglevels[i]->smoothd,KSP_NORM_NONE);CHKERRQ(ierr);
245336babb1SJed Brown       ierr = KSPGetPC(mglevels[i]->smoothd,&ipc);CHKERRQ(ierr);
246336babb1SJed Brown       ierr = PCSetType(ipc,PCSOR);CHKERRQ(ierr);
2470ee9a668SBarry Smith       ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg->default_smoothd);CHKERRQ(ierr);
248f3fbd535SBarry Smith 
2490ee9a668SBarry Smith       sprintf(tprefix,"mg_levels_%d_",(int)i);
2500ee9a668SBarry Smith       ierr = KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);CHKERRQ(ierr);
2510ee9a668SBarry Smith     } else {
252f3fbd535SBarry Smith       ierr = KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");CHKERRQ(ierr);
253f3fbd535SBarry Smith 
2549dbfc187SHong Zhang       /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */
255f3fbd535SBarry Smith       ierr = KSPSetType(mglevels[0]->smoothd,KSPPREONLY);CHKERRQ(ierr);
256f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[0]->smoothd,&ipc);CHKERRQ(ierr);
257f3fbd535SBarry Smith       ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
258f3fbd535SBarry Smith       if (size > 1) {
259f3fbd535SBarry Smith         ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr);
260f3fbd535SBarry Smith       } else {
261f3fbd535SBarry Smith         ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr);
262f3fbd535SBarry Smith       }
263753b7fb9SBarry Smith       ierr = PCFactorSetShiftType(ipc,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
264f3fbd535SBarry Smith     }
2653bb1ff40SBarry Smith     ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[i]->smoothd);CHKERRQ(ierr);
2662fa5cd67SKarl Rupp 
267f3fbd535SBarry Smith     mglevels[i]->smoothu = mglevels[i]->smoothd;
26831567311SBarry Smith     mg->rtol             = 0.0;
26931567311SBarry Smith     mg->abstol           = 0.0;
27031567311SBarry Smith     mg->dtol             = 0.0;
27131567311SBarry Smith     mg->ttol             = 0.0;
27231567311SBarry Smith     mg->cyclesperpcapply = 1;
273f3fbd535SBarry Smith   }
274f3fbd535SBarry Smith   mg->levels = mglevels;
27510eca3edSLisandro Dalcin   ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
2764b9ad928SBarry Smith   PetscFunctionReturn(0);
2774b9ad928SBarry Smith }
2784b9ad928SBarry Smith 
279c07bf074SBarry Smith 
280c07bf074SBarry Smith #undef __FUNCT__
281c07bf074SBarry Smith #define __FUNCT__ "PCDestroy_MG"
282c07bf074SBarry Smith PetscErrorCode PCDestroy_MG(PC pc)
283c07bf074SBarry Smith {
284c07bf074SBarry Smith   PetscErrorCode ierr;
285a06653b4SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
286a06653b4SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
287a06653b4SBarry Smith   PetscInt       i,n;
288c07bf074SBarry Smith 
289c07bf074SBarry Smith   PetscFunctionBegin;
290a06653b4SBarry Smith   ierr = PCReset_MG(pc);CHKERRQ(ierr);
291a06653b4SBarry Smith   if (mglevels) {
292a06653b4SBarry Smith     n = mglevels[0]->levels;
293a06653b4SBarry Smith     for (i=0; i<n; i++) {
294a06653b4SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
2956bf464f9SBarry Smith         ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
296a06653b4SBarry Smith       }
2976bf464f9SBarry Smith       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
298a06653b4SBarry Smith       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
299a06653b4SBarry Smith     }
300a06653b4SBarry Smith     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
301a06653b4SBarry Smith   }
302c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
303f3fbd535SBarry Smith   PetscFunctionReturn(0);
304f3fbd535SBarry Smith }
305f3fbd535SBarry Smith 
306f3fbd535SBarry Smith 
307f3fbd535SBarry Smith 
30809573ac7SBarry Smith extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**);
30909573ac7SBarry Smith extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**);
31009573ac7SBarry Smith extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**);
311f3fbd535SBarry Smith 
312f3fbd535SBarry Smith /*
313f3fbd535SBarry Smith    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
314f3fbd535SBarry Smith              or full cycle of multigrid.
315f3fbd535SBarry Smith 
316f3fbd535SBarry Smith   Note:
317f3fbd535SBarry Smith   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
318f3fbd535SBarry Smith */
319f3fbd535SBarry Smith #undef __FUNCT__
320f3fbd535SBarry Smith #define __FUNCT__ "PCApply_MG"
321f3fbd535SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
322f3fbd535SBarry Smith {
323f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
324f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
325f3fbd535SBarry Smith   PetscErrorCode ierr;
326f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
327f3fbd535SBarry Smith 
328f3fbd535SBarry Smith   PetscFunctionBegin;
329a9db87a2SMatthew G Knepley   if (mg->stageApply) {ierr = PetscLogStagePush(mg->stageApply);CHKERRQ(ierr);}
330e1d8e5deSBarry Smith   /* When the DM is supplying the matrix then it will not exist until here */
331298cc208SBarry Smith   for (i=0; i<levels; i++) {
332e1d8e5deSBarry Smith     if (!mglevels[i]->A) {
33323ee1639SBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr);
334298cc208SBarry Smith       ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr);
335e1d8e5deSBarry Smith     }
336e1d8e5deSBarry Smith   }
337e1d8e5deSBarry Smith 
338f3fbd535SBarry Smith   mglevels[levels-1]->b = b;
339f3fbd535SBarry Smith   mglevels[levels-1]->x = x;
34031567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
341f3fbd535SBarry Smith     ierr = VecSet(x,0.0);CHKERRQ(ierr);
34231567311SBarry Smith     for (i=0; i<mg->cyclesperpcapply; i++) {
3430298fd71SBarry Smith       ierr = PCMGMCycle_Private(pc,mglevels+levels-1,NULL);CHKERRQ(ierr);
344f3fbd535SBarry Smith     }
3452fa5cd67SKarl Rupp   } else if (mg->am == PC_MG_ADDITIVE) {
34631567311SBarry Smith     ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr);
3472fa5cd67SKarl Rupp   } else if (mg->am == PC_MG_KASKADE) {
34831567311SBarry Smith     ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr);
3492fa5cd67SKarl Rupp   } else {
350f3fbd535SBarry Smith     ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr);
351f3fbd535SBarry Smith   }
352a9db87a2SMatthew G Knepley   if (mg->stageApply) {ierr = PetscLogStagePop();CHKERRQ(ierr);}
353f3fbd535SBarry Smith   PetscFunctionReturn(0);
354f3fbd535SBarry Smith }
355f3fbd535SBarry Smith 
356f3fbd535SBarry Smith 
357f3fbd535SBarry Smith #undef __FUNCT__
358f3fbd535SBarry Smith #define __FUNCT__ "PCSetFromOptions_MG"
3594416b707SBarry Smith PetscErrorCode PCSetFromOptions_MG(PetscOptionItems *PetscOptionsObject,PC pc)
360f3fbd535SBarry Smith {
361f3fbd535SBarry Smith   PetscErrorCode ierr;
362f3fbd535SBarry Smith   PetscInt       m,levels = 1,cycles;
363c91913d3SJed Brown   PetscBool      flg,set;
364f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
3653d3eaba7SBarry Smith   PC_MG_Levels   **mglevels;
366f3fbd535SBarry Smith   PCMGType       mgtype;
367f3fbd535SBarry Smith   PCMGCycleType  mgctype;
368f3fbd535SBarry Smith 
369f3fbd535SBarry Smith   PetscFunctionBegin;
370e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Multigrid options");CHKERRQ(ierr);
3713aeaf226SBarry Smith   if (!mg->levels) {
372f3fbd535SBarry Smith     ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr);
373eb3f98d2SBarry Smith     if (!flg && pc->dm) {
374eb3f98d2SBarry Smith       ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
375eb3f98d2SBarry Smith       levels++;
3763aeaf226SBarry Smith       mg->usedmfornumberoflevels = PETSC_TRUE;
377eb3f98d2SBarry Smith     }
3780298fd71SBarry Smith     ierr = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr);
379f3fbd535SBarry Smith   }
3803aeaf226SBarry Smith   mglevels = mg->levels;
3813aeaf226SBarry Smith 
382f3fbd535SBarry Smith   mgctype = (PCMGCycleType) mglevels[0]->cycles;
383f3fbd535SBarry Smith   ierr    = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr);
384f3fbd535SBarry Smith   if (flg) {
385f3fbd535SBarry Smith     ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr);
3862fa5cd67SKarl Rupp   }
387f3fbd535SBarry Smith   flg  = PETSC_FALSE;
388c91913d3SJed Brown   ierr = PetscOptionsBool("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",flg,&flg,&set);CHKERRQ(ierr);
389c91913d3SJed Brown   if (set) {
390c91913d3SJed Brown     ierr = PCMGSetGalerkin(pc,flg);CHKERRQ(ierr);
391f3fbd535SBarry Smith   }
39294ae4db5SBarry Smith   ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",mg->default_smoothu,&m,&flg);CHKERRQ(ierr);
393f3fbd535SBarry Smith   if (flg) {
394f3fbd535SBarry Smith     ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr);
395f3fbd535SBarry Smith   }
39694ae4db5SBarry Smith   ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",mg->default_smoothd,&m,&flg);CHKERRQ(ierr);
397f3fbd535SBarry Smith   if (flg) {
398f3fbd535SBarry Smith     ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr);
399f3fbd535SBarry Smith   }
40031567311SBarry Smith   mgtype = mg->am;
401f3fbd535SBarry Smith   ierr   = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
402f3fbd535SBarry Smith   if (flg) {
403f3fbd535SBarry Smith     ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
404f3fbd535SBarry Smith   }
40531567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
4065363c1e0SLisandro Dalcin     ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGMultiplicativeSetCycles",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr);
407f3fbd535SBarry Smith     if (flg) {
408f3fbd535SBarry Smith       ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr);
409f3fbd535SBarry Smith     }
410f3fbd535SBarry Smith   }
411f3fbd535SBarry Smith   flg  = PETSC_FALSE;
4120298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,NULL);CHKERRQ(ierr);
413f3fbd535SBarry Smith   if (flg) {
414f3fbd535SBarry Smith     PetscInt i;
415f3fbd535SBarry Smith     char     eventname[128];
416ce94432eSBarry Smith     if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
417f3fbd535SBarry Smith     levels = mglevels[0]->levels;
418f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
419f3fbd535SBarry Smith       sprintf(eventname,"MGSetup Level %d",(int)i);
42063e6d426SJed Brown       ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr);
421f3fbd535SBarry Smith       sprintf(eventname,"MGSmooth Level %d",(int)i);
42263e6d426SJed Brown       ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr);
423f3fbd535SBarry Smith       if (i) {
424f3fbd535SBarry Smith         sprintf(eventname,"MGResid Level %d",(int)i);
42563e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr);
426f3fbd535SBarry Smith         sprintf(eventname,"MGInterp Level %d",(int)i);
42763e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr);
428f3fbd535SBarry Smith       }
429f3fbd535SBarry Smith     }
430eec35531SMatthew G Knepley 
431ec60ca73SSatish Balay #if defined(PETSC_USE_LOG)
432eec35531SMatthew G Knepley     {
433eec35531SMatthew G Knepley       const char    *sname = "MG Apply";
434eec35531SMatthew G Knepley       PetscStageLog stageLog;
435eec35531SMatthew G Knepley       PetscInt      st;
436eec35531SMatthew G Knepley 
437eec35531SMatthew G Knepley       PetscFunctionBegin;
438eec35531SMatthew G Knepley       ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
439eec35531SMatthew G Knepley       for (st = 0; st < stageLog->numStages; ++st) {
440eec35531SMatthew G Knepley         PetscBool same;
441eec35531SMatthew G Knepley 
442eec35531SMatthew G Knepley         ierr = PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);CHKERRQ(ierr);
4432fa5cd67SKarl Rupp         if (same) mg->stageApply = st;
444eec35531SMatthew G Knepley       }
445eec35531SMatthew G Knepley       if (!mg->stageApply) {
446eec35531SMatthew G Knepley         ierr = PetscLogStageRegister(sname, &mg->stageApply);CHKERRQ(ierr);
447eec35531SMatthew G Knepley       }
448eec35531SMatthew G Knepley     }
449ec60ca73SSatish Balay #endif
450f3fbd535SBarry Smith   }
451f3fbd535SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
452f3fbd535SBarry Smith   PetscFunctionReturn(0);
453f3fbd535SBarry Smith }
454f3fbd535SBarry Smith 
4556a6fc655SJed Brown const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0};
4566a6fc655SJed Brown const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0};
457f3fbd535SBarry Smith 
4589804daf3SBarry Smith #include <petscdraw.h>
459f3fbd535SBarry Smith #undef __FUNCT__
460f3fbd535SBarry Smith #define __FUNCT__ "PCView_MG"
461f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
462f3fbd535SBarry Smith {
463f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
464f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
465f3fbd535SBarry Smith   PetscErrorCode ierr;
466e3deeeafSJed Brown   PetscInt       levels = mglevels ? mglevels[0]->levels : 0,i;
467219b1045SBarry Smith   PetscBool      iascii,isbinary,isdraw;
468f3fbd535SBarry Smith 
469f3fbd535SBarry Smith   PetscFunctionBegin;
470251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
4715b0b0462SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
472219b1045SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
473f3fbd535SBarry Smith   if (iascii) {
474e3deeeafSJed Brown     const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown";
475e3deeeafSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  MG: type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,cyclename);CHKERRQ(ierr);
47631567311SBarry Smith     if (mg->am == PC_MG_MULTIPLICATIVE) {
47731567311SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr);
478f3fbd535SBarry Smith     }
479218a07d4SBarry Smith     if (mg->galerkin) {
480f3fbd535SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
4814f66f45eSBarry Smith     } else {
4824f66f45eSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
483f3fbd535SBarry Smith     }
4845adeb434SBarry Smith     if (mg->view){
4855adeb434SBarry Smith       ierr = (*mg->view)(pc,viewer);CHKERRQ(ierr);
4865adeb434SBarry Smith     }
487f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
488f3fbd535SBarry Smith       if (!i) {
48989cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);CHKERRQ(ierr);
490f3fbd535SBarry Smith       } else {
49189cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
492f3fbd535SBarry Smith       }
493f3fbd535SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
494f3fbd535SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
495f3fbd535SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
496f3fbd535SBarry Smith       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
497f3fbd535SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
498f3fbd535SBarry Smith       } else if (i) {
49989cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
500f3fbd535SBarry Smith         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
501f3fbd535SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
502f3fbd535SBarry Smith         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
503f3fbd535SBarry Smith       }
504f3fbd535SBarry Smith     }
5055b0b0462SBarry Smith   } else if (isbinary) {
5065b0b0462SBarry Smith     for (i=levels-1; i>=0; i--) {
5075b0b0462SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
5085b0b0462SBarry Smith       if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) {
5095b0b0462SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
5105b0b0462SBarry Smith       }
5115b0b0462SBarry Smith     }
512219b1045SBarry Smith   } else if (isdraw) {
513219b1045SBarry Smith     PetscDraw draw;
514b4375e8dSPeter Brune     PetscReal x,w,y,bottom,th;
515219b1045SBarry Smith     ierr   = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
516219b1045SBarry Smith     ierr   = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
5170298fd71SBarry Smith     ierr   = PetscDrawStringGetSize(draw,NULL,&th);CHKERRQ(ierr);
518219b1045SBarry Smith     bottom = y - th;
519219b1045SBarry Smith     for (i=levels-1; i>=0; i--) {
520b4375e8dSPeter Brune       if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) {
521219b1045SBarry Smith         ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr);
522219b1045SBarry Smith         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
523219b1045SBarry Smith         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
524b4375e8dSPeter Brune       } else {
525b4375e8dSPeter Brune         w    = 0.5*PetscMin(1.0-x,x);
526b4375e8dSPeter Brune         ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr);
527b4375e8dSPeter Brune         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
528b4375e8dSPeter Brune         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
529b4375e8dSPeter Brune         ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr);
530b4375e8dSPeter Brune         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
531b4375e8dSPeter Brune         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
532b4375e8dSPeter Brune       }
5330298fd71SBarry Smith       ierr    = PetscDrawGetBoundingBox(draw,NULL,&bottom,NULL,NULL);CHKERRQ(ierr);
5341cd381d6SBarry Smith       bottom -= th;
535219b1045SBarry Smith     }
5365b0b0462SBarry Smith   }
537f3fbd535SBarry Smith   PetscFunctionReturn(0);
538f3fbd535SBarry Smith }
539f3fbd535SBarry Smith 
540af0996ceSBarry Smith #include <petsc/private/dmimpl.h>
541af0996ceSBarry Smith #include <petsc/private/kspimpl.h>
542cab2e9ccSBarry Smith 
543f3fbd535SBarry Smith /*
544f3fbd535SBarry Smith     Calls setup for the KSP on each level
545f3fbd535SBarry Smith */
546f3fbd535SBarry Smith #undef __FUNCT__
547f3fbd535SBarry Smith #define __FUNCT__ "PCSetUp_MG"
548f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc)
549f3fbd535SBarry Smith {
550f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
551f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
552f3fbd535SBarry Smith   PetscErrorCode ierr;
553f3fbd535SBarry Smith   PetscInt       i,n = mglevels[0]->levels;
55498e04984SBarry Smith   PC             cpc;
5553db492dfSBarry Smith   PetscBool      dump = PETSC_FALSE,opsset,use_amat,missinginterpolate = PETSC_FALSE;
556f3fbd535SBarry Smith   Mat            dA,dB;
557f3fbd535SBarry Smith   Vec            tvec;
558218a07d4SBarry Smith   DM             *dms;
559649052a6SBarry Smith   PetscViewer    viewer = 0;
560f3fbd535SBarry Smith 
561f3fbd535SBarry Smith   PetscFunctionBegin;
56201bc414fSDmitry Karpeev   /* FIX: Move this to PCSetFromOptions_MG? */
5633aeaf226SBarry Smith   if (mg->usedmfornumberoflevels) {
5643aeaf226SBarry Smith     PetscInt levels;
5653aeaf226SBarry Smith     ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
5663aeaf226SBarry Smith     levels++;
5673aeaf226SBarry Smith     if (levels > n) { /* the problem is now being solved on a finer grid */
5680298fd71SBarry Smith       ierr     = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr);
5693aeaf226SBarry Smith       n        = levels;
5703aeaf226SBarry Smith       ierr     = PCSetFromOptions(pc);CHKERRQ(ierr); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
5713aeaf226SBarry Smith       mglevels =  mg->levels;
5723aeaf226SBarry Smith     }
5733aeaf226SBarry Smith   }
57498e04984SBarry Smith   ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr);
5753aeaf226SBarry Smith 
576f3fbd535SBarry Smith 
577f3fbd535SBarry Smith   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
578f3fbd535SBarry Smith   /* so use those from global PC */
579f3fbd535SBarry Smith   /* Is this what we always want? What if user wants to keep old one? */
5800298fd71SBarry Smith   ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset);CHKERRQ(ierr);
58198e04984SBarry Smith   if (opsset) {
58298e04984SBarry Smith     Mat mmat;
58323ee1639SBarry Smith     ierr = KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat);CHKERRQ(ierr);
58498e04984SBarry Smith     if (mmat == pc->pmat) opsset = PETSC_FALSE;
58598e04984SBarry Smith   }
5864a5f07a5SMark F. Adams 
5874a5f07a5SMark F. Adams   if (!opsset) {
58871b23a65SMark F. Adams     ierr = PCGetUseAmat(pc,&use_amat);CHKERRQ(ierr);
5894a5f07a5SMark F. Adams     if(use_amat){
590f3fbd535SBarry 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);
59123ee1639SBarry Smith       ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat);CHKERRQ(ierr);
592f3fbd535SBarry Smith     }
5934a5f07a5SMark F. Adams     else {
5944a5f07a5SMark 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);
59523ee1639SBarry Smith       ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat);CHKERRQ(ierr);
5964a5f07a5SMark F. Adams     }
5974a5f07a5SMark F. Adams   }
598f3fbd535SBarry Smith 
5999c7ffb39SBarry Smith   for (i=n-1; i>0; i--) {
6009c7ffb39SBarry Smith     if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) {
6019c7ffb39SBarry Smith       missinginterpolate = PETSC_TRUE;
6029c7ffb39SBarry Smith       continue;
6039c7ffb39SBarry Smith     }
6049c7ffb39SBarry Smith   }
6059c7ffb39SBarry Smith   /*
6069c7ffb39SBarry Smith    Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS)
6079c7ffb39SBarry Smith    Skipping for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs?
6089c7ffb39SBarry Smith   */
6099c7ffb39SBarry Smith   if (missinginterpolate && pc->dm && mg->galerkin != 2 && !pc->setupcalled) {
6102d2b81a6SBarry Smith     /* construct the interpolation from the DMs */
6112e499ae9SBarry Smith     Mat p;
61273250ac0SBarry Smith     Vec rscale;
613785e854fSJed Brown     ierr     = PetscMalloc1(n,&dms);CHKERRQ(ierr);
614218a07d4SBarry Smith     dms[n-1] = pc->dm;
615ef1267afSMatthew G. Knepley     /* Separately create them so we do not get DMKSP interference between levels */
616ef1267afSMatthew G. Knepley     for (i=n-2; i>-1; i--) {ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr);}
617218a07d4SBarry Smith     for (i=n-2; i>-1; i--) {
618942e3340SBarry Smith       DMKSP     kdm;
6193ad4599aSBarry Smith       PetscBool dmhasrestrict;
6203c0c59f3SBarry Smith       ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr);
6213c0c59f3SBarry Smith       if (mg->galerkin) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);}
622942e3340SBarry Smith       ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr);
623d1a3e677SJed Brown       /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take
624d1a3e677SJed Brown        * a bitwise OR of computing the matrix, RHS, and initial iterate. */
6250298fd71SBarry Smith       kdm->ops->computerhs = NULL;
6260298fd71SBarry Smith       kdm->rhsctx          = NULL;
62724c3aa18SBarry Smith       if (!mglevels[i+1]->interpolate) {
628e727c939SJed Brown         ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr);
6292d2b81a6SBarry Smith         ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr);
63000ac8be1SBarry Smith         if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);}
63173250ac0SBarry Smith         ierr = VecDestroy(&rscale);CHKERRQ(ierr);
6326bf464f9SBarry Smith         ierr = MatDestroy(&p);CHKERRQ(ierr);
633218a07d4SBarry Smith       }
6343ad4599aSBarry Smith       ierr = DMHasCreateRestriction(dms[i],&dmhasrestrict);CHKERRQ(ierr);
6353ad4599aSBarry Smith       if (dmhasrestrict && !mglevels[i+1]->restrct){
6363ad4599aSBarry Smith         ierr = DMCreateRestriction(dms[i],dms[i+1],&p);CHKERRQ(ierr);
6373ad4599aSBarry Smith         ierr = PCMGSetRestriction(pc,i+1,p);CHKERRQ(ierr);
6383ad4599aSBarry Smith         ierr = MatDestroy(&p);CHKERRQ(ierr);
6393ad4599aSBarry Smith       }
64024c3aa18SBarry Smith     }
6412d2b81a6SBarry Smith 
642ef1267afSMatthew G. Knepley     for (i=n-2; i>-1; i--) {ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);}
6432d2b81a6SBarry Smith     ierr = PetscFree(dms);CHKERRQ(ierr);
644ce4cda84SJed Brown   }
645cab2e9ccSBarry Smith 
646ce4cda84SJed Brown   if (pc->dm && !pc->setupcalled) {
647ce4cda84SJed Brown     /* finest smoother also gets DM but it is not active, independent of whether galerkin==2 */
648cab2e9ccSBarry Smith     ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr);
649cab2e9ccSBarry Smith     ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr);
650218a07d4SBarry Smith   }
651218a07d4SBarry Smith 
652c91913d3SJed Brown   if (mg->galerkin == 1) {
653f3fbd535SBarry Smith     Mat B;
654f3fbd535SBarry Smith     /* currently only handle case where mat and pmat are the same on coarser levels */
65523ee1639SBarry Smith     ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB);CHKERRQ(ierr);
656f3fbd535SBarry Smith     if (!pc->setupcalled) {
657f3fbd535SBarry Smith       for (i=n-2; i>-1; i--) {
658b9367914SBarry 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");
659b9367914SBarry Smith         if (!mglevels[i+1]->interpolate) {
660b9367914SBarry Smith           ierr = PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct);CHKERRQ(ierr);
661b9367914SBarry Smith         }
662b9367914SBarry Smith         if (!mglevels[i+1]->restrct) {
663b9367914SBarry Smith           ierr = PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate);CHKERRQ(ierr);
664b9367914SBarry Smith         }
665b9367914SBarry Smith         if (mglevels[i+1]->interpolate == mglevels[i+1]->restrct) {
666f3fbd535SBarry Smith           ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
6677d537d92SJed Brown         } else {
6687d537d92SJed Brown           ierr = MatMatMatMult(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
6697d537d92SJed Brown         }
67023ee1639SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B);CHKERRQ(ierr);
671f3fbd535SBarry Smith         if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
672f3fbd535SBarry Smith         dB = B;
673f3fbd535SBarry Smith       }
674cd9507b2SBarry Smith       if (n > 1) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
675f3fbd535SBarry Smith     } else {
676f3fbd535SBarry Smith       for (i=n-2; i>-1; i--) {
677b9367914SBarry 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");
678b9367914SBarry Smith         if (!mglevels[i+1]->interpolate) {
679b9367914SBarry Smith           ierr = PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct);CHKERRQ(ierr);
680b9367914SBarry Smith         }
681b9367914SBarry Smith         if (!mglevels[i+1]->restrct) {
682b9367914SBarry Smith           ierr = PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate);CHKERRQ(ierr);
683b9367914SBarry Smith         }
68423ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&B);CHKERRQ(ierr);
685b9367914SBarry Smith         if (mglevels[i+1]->interpolate == mglevels[i+1]->restrct) {
686f3fbd535SBarry Smith           ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
6877d537d92SJed Brown         } else {
6887d537d92SJed Brown           ierr = MatMatMatMult(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
6897d537d92SJed Brown         }
69023ee1639SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B);CHKERRQ(ierr);
691f3fbd535SBarry Smith         dB   = B;
692f3fbd535SBarry Smith       }
693f3fbd535SBarry Smith     }
694ce4cda84SJed Brown   } else if (!mg->galerkin && pc->dm && pc->dm->x) {
695cab2e9ccSBarry Smith     /* need to restrict Jacobian location to coarser meshes for evaluation */
696cab2e9ccSBarry Smith     for (i=n-2; i>-1; i--) {
697c88c5224SJed Brown       Mat R;
698c88c5224SJed Brown       Vec rscale;
699cab2e9ccSBarry Smith       if (!mglevels[i]->smoothd->dm->x) {
700cab2e9ccSBarry Smith         Vec *vecs;
7012a7a6963SBarry Smith         ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vecs,0,NULL);CHKERRQ(ierr);
7022fa5cd67SKarl Rupp 
703cab2e9ccSBarry Smith         mglevels[i]->smoothd->dm->x = vecs[0];
7042fa5cd67SKarl Rupp 
705cab2e9ccSBarry Smith         ierr = PetscFree(vecs);CHKERRQ(ierr);
706cab2e9ccSBarry Smith       }
707c88c5224SJed Brown       ierr = PCMGGetRestriction(pc,i+1,&R);CHKERRQ(ierr);
708c88c5224SJed Brown       ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
709c88c5224SJed Brown       ierr = MatRestrict(R,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);CHKERRQ(ierr);
710c88c5224SJed Brown       ierr = VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,rscale);CHKERRQ(ierr);
711cab2e9ccSBarry Smith     }
712f3fbd535SBarry Smith   }
713ccceb50cSJed Brown   if (!mg->galerkin && pc->dm) {
714caa4e7f2SJed Brown     for (i=n-2; i>=0; i--) {
715ccceb50cSJed Brown       DM  dmfine,dmcoarse;
716caa4e7f2SJed Brown       Mat Restrict,Inject;
717caa4e7f2SJed Brown       Vec rscale;
718ccceb50cSJed Brown       ierr   = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr);
719ccceb50cSJed Brown       ierr   = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr);
720caa4e7f2SJed Brown       ierr   = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr);
721caa4e7f2SJed Brown       ierr   = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
7220298fd71SBarry Smith       Inject = NULL;      /* Callback should create it if it needs Injection */
723caa4e7f2SJed Brown       ierr   = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr);
724caa4e7f2SJed Brown     }
725caa4e7f2SJed Brown   }
726f3fbd535SBarry Smith 
727f3fbd535SBarry Smith   if (!pc->setupcalled) {
728f3fbd535SBarry Smith     for (i=0; i<n; i++) {
729f3fbd535SBarry Smith       ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr);
730f3fbd535SBarry Smith     }
731f3fbd535SBarry Smith     for (i=1; i<n; i++) {
732f3fbd535SBarry Smith       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
733f3fbd535SBarry Smith         ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr);
734f3fbd535SBarry Smith       }
735f3fbd535SBarry Smith     }
7363ad4599aSBarry Smith     /* insure that if either interpolation or restriction is set the other other one is set */
737f3fbd535SBarry Smith     for (i=1; i<n; i++) {
7383ad4599aSBarry Smith       ierr = PCMGGetInterpolation(pc,i,NULL);CHKERRQ(ierr);
7393ad4599aSBarry Smith       ierr = PCMGGetRestriction(pc,i,NULL);CHKERRQ(ierr);
740f3fbd535SBarry Smith     }
741f3fbd535SBarry Smith     for (i=0; i<n-1; i++) {
742f3fbd535SBarry Smith       if (!mglevels[i]->b) {
743f3fbd535SBarry Smith         Vec *vec;
7442a7a6963SBarry Smith         ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
745f3fbd535SBarry Smith         ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr);
7466bf464f9SBarry Smith         ierr = VecDestroy(vec);CHKERRQ(ierr);
747f3fbd535SBarry Smith         ierr = PetscFree(vec);CHKERRQ(ierr);
748f3fbd535SBarry Smith       }
749f3fbd535SBarry Smith       if (!mglevels[i]->r && i) {
750f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
751f3fbd535SBarry Smith         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
7526bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
753f3fbd535SBarry Smith       }
754f3fbd535SBarry Smith       if (!mglevels[i]->x) {
755f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
756f3fbd535SBarry Smith         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
7576bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
758f3fbd535SBarry Smith       }
759f3fbd535SBarry Smith     }
760f3fbd535SBarry Smith     if (n != 1 && !mglevels[n-1]->r) {
761f3fbd535SBarry Smith       /* PCMGSetR() on the finest level if user did not supply it */
762f3fbd535SBarry Smith       Vec *vec;
7632a7a6963SBarry Smith       ierr = KSPCreateVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
764f3fbd535SBarry Smith       ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr);
7656bf464f9SBarry Smith       ierr = VecDestroy(vec);CHKERRQ(ierr);
766f3fbd535SBarry Smith       ierr = PetscFree(vec);CHKERRQ(ierr);
767f3fbd535SBarry Smith     }
768f3fbd535SBarry Smith   }
769f3fbd535SBarry Smith 
77098e04984SBarry Smith   if (pc->dm) {
77198e04984SBarry Smith     /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */
77298e04984SBarry Smith     for (i=0; i<n-1; i++) {
77398e04984SBarry Smith       if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX;
77498e04984SBarry Smith     }
77598e04984SBarry Smith   }
776f3fbd535SBarry Smith 
777f3fbd535SBarry Smith   for (i=1; i<n; i++) {
7782a21e185SBarry Smith     if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1){
779f3fbd535SBarry Smith       /* if doing only down then initial guess is zero */
780f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
781f3fbd535SBarry Smith     }
78263e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
783f3fbd535SBarry Smith     ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr);
784899639b0SHong Zhang     if (mglevels[i]->smoothd->reason == KSP_DIVERGED_PCSETUP_FAILED) {
785899639b0SHong Zhang       pc->failedreason = PC_SUBPC_ERROR;
786899639b0SHong Zhang     }
78763e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
788d42688cbSBarry Smith     if (!mglevels[i]->residual) {
789d42688cbSBarry Smith       Mat mat;
79023ee1639SBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&mat);CHKERRQ(ierr);
79154b2cd4bSJed Brown       ierr = PCMGSetResidual(pc,i,PCMGResidualDefault,mat);CHKERRQ(ierr);
792d42688cbSBarry Smith     }
793f3fbd535SBarry Smith   }
794f3fbd535SBarry Smith   for (i=1; i<n; i++) {
795f3fbd535SBarry Smith     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
796f3fbd535SBarry Smith       Mat          downmat,downpmat;
797f3fbd535SBarry Smith 
798f3fbd535SBarry Smith       /* check if operators have been set for up, if not use down operators to set them */
7990298fd71SBarry Smith       ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);CHKERRQ(ierr);
800f3fbd535SBarry Smith       if (!opsset) {
80123ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);CHKERRQ(ierr);
80223ee1639SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat);CHKERRQ(ierr);
803f3fbd535SBarry Smith       }
804f3fbd535SBarry Smith 
805f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
80663e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
807f3fbd535SBarry Smith       ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr);
808899639b0SHong Zhang       if (mglevels[i]->smoothu->reason == KSP_DIVERGED_PCSETUP_FAILED) {
809899639b0SHong Zhang         pc->failedreason = PC_SUBPC_ERROR;
810899639b0SHong Zhang       }
81163e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
812f3fbd535SBarry Smith     }
813f3fbd535SBarry Smith   }
814f3fbd535SBarry Smith 
81563e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
816f3fbd535SBarry Smith   ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr);
817899639b0SHong Zhang   if (mglevels[0]->smoothd->reason == KSP_DIVERGED_PCSETUP_FAILED) {
818899639b0SHong Zhang     pc->failedreason = PC_SUBPC_ERROR;
819899639b0SHong Zhang   }
82063e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
821f3fbd535SBarry Smith 
822f3fbd535SBarry Smith   /*
823f3fbd535SBarry Smith      Dump the interpolation/restriction matrices plus the
824e3c5b3baSBarry Smith    Jacobian/stiffness on each level. This allows MATLAB users to
825f3fbd535SBarry Smith    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
826f3fbd535SBarry Smith 
827f3fbd535SBarry Smith    Only support one or the other at the same time.
828f3fbd535SBarry Smith   */
829f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
830c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);CHKERRQ(ierr);
831ce94432eSBarry Smith   if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc));
832f3fbd535SBarry Smith   dump = PETSC_FALSE;
833f3fbd535SBarry Smith #endif
834c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);CHKERRQ(ierr);
835ce94432eSBarry Smith   if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc));
836f3fbd535SBarry Smith 
837f3fbd535SBarry Smith   if (viewer) {
838f3fbd535SBarry Smith     for (i=1; i<n; i++) {
839f3fbd535SBarry Smith       ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr);
840f3fbd535SBarry Smith     }
841f3fbd535SBarry Smith     for (i=0; i<n; i++) {
842f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr);
843f3fbd535SBarry Smith       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
844f3fbd535SBarry Smith     }
845f3fbd535SBarry Smith   }
846f3fbd535SBarry Smith   PetscFunctionReturn(0);
847f3fbd535SBarry Smith }
848f3fbd535SBarry Smith 
849f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/
850f3fbd535SBarry Smith 
851f3fbd535SBarry Smith #undef __FUNCT__
8529dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetLevels"
8534b9ad928SBarry Smith /*@
85497177400SBarry Smith    PCMGGetLevels - Gets the number of levels to use with MG.
8554b9ad928SBarry Smith 
8564b9ad928SBarry Smith    Not Collective
8574b9ad928SBarry Smith 
8584b9ad928SBarry Smith    Input Parameter:
8594b9ad928SBarry Smith .  pc - the preconditioner context
8604b9ad928SBarry Smith 
8614b9ad928SBarry Smith    Output parameter:
8624b9ad928SBarry Smith .  levels - the number of levels
8634b9ad928SBarry Smith 
8644b9ad928SBarry Smith    Level: advanced
8654b9ad928SBarry Smith 
8664b9ad928SBarry Smith .keywords: MG, get, levels, multigrid
8674b9ad928SBarry Smith 
86897177400SBarry Smith .seealso: PCMGSetLevels()
8694b9ad928SBarry Smith @*/
8707087cfbeSBarry Smith PetscErrorCode  PCMGGetLevels(PC pc,PetscInt *levels)
8714b9ad928SBarry Smith {
872f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
8734b9ad928SBarry Smith 
8744b9ad928SBarry Smith   PetscFunctionBegin;
8750700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
8764482741eSBarry Smith   PetscValidIntPointer(levels,2);
877f3fbd535SBarry Smith   *levels = mg->nlevels;
8784b9ad928SBarry Smith   PetscFunctionReturn(0);
8794b9ad928SBarry Smith }
8804b9ad928SBarry Smith 
8814b9ad928SBarry Smith #undef __FUNCT__
8829dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetType"
8834b9ad928SBarry Smith /*@
88497177400SBarry Smith    PCMGSetType - Determines the form of multigrid to use:
8854b9ad928SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
8864b9ad928SBarry Smith 
887ad4df100SBarry Smith    Logically Collective on PC
8884b9ad928SBarry Smith 
8894b9ad928SBarry Smith    Input Parameters:
8904b9ad928SBarry Smith +  pc - the preconditioner context
8919dcbbd2bSBarry Smith -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
8929dcbbd2bSBarry Smith    PC_MG_FULL, PC_MG_KASKADE
8934b9ad928SBarry Smith 
8944b9ad928SBarry Smith    Options Database Key:
8954b9ad928SBarry Smith .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
8964b9ad928SBarry Smith    additive, full, kaskade
8974b9ad928SBarry Smith 
8984b9ad928SBarry Smith    Level: advanced
8994b9ad928SBarry Smith 
9004b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
9014b9ad928SBarry Smith 
90297177400SBarry Smith .seealso: PCMGSetLevels()
9034b9ad928SBarry Smith @*/
9047087cfbeSBarry Smith PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
9054b9ad928SBarry Smith {
906f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
9074b9ad928SBarry Smith 
9084b9ad928SBarry Smith   PetscFunctionBegin;
9090700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
910c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,form,2);
91131567311SBarry Smith   mg->am = form;
9129dcbbd2bSBarry Smith   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
913c60c7ad4SBarry Smith   else pc->ops->applyrichardson = NULL;
914c60c7ad4SBarry Smith   PetscFunctionReturn(0);
915c60c7ad4SBarry Smith }
916c60c7ad4SBarry Smith 
917f0af28e8SLisandro Dalcin #undef __FUNCT__
918f0af28e8SLisandro Dalcin #define __FUNCT__ "PCMGGetType"
919c60c7ad4SBarry Smith /*@
920c60c7ad4SBarry Smith    PCMGGetType - Determines the form of multigrid to use:
921c60c7ad4SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
922c60c7ad4SBarry Smith 
923c60c7ad4SBarry Smith    Logically Collective on PC
924c60c7ad4SBarry Smith 
925c60c7ad4SBarry Smith    Input Parameter:
926c60c7ad4SBarry Smith .  pc - the preconditioner context
927c60c7ad4SBarry Smith 
928c60c7ad4SBarry Smith    Output Parameter:
929c60c7ad4SBarry Smith .  type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE
930c60c7ad4SBarry Smith 
931c60c7ad4SBarry Smith 
932c60c7ad4SBarry Smith    Level: advanced
933c60c7ad4SBarry Smith 
934c60c7ad4SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
935c60c7ad4SBarry Smith 
936c60c7ad4SBarry Smith .seealso: PCMGSetLevels()
937c60c7ad4SBarry Smith @*/
938c60c7ad4SBarry Smith PetscErrorCode  PCMGGetType(PC pc,PCMGType *type)
939c60c7ad4SBarry Smith {
940c60c7ad4SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
941c60c7ad4SBarry Smith 
942c60c7ad4SBarry Smith   PetscFunctionBegin;
943c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
944c60c7ad4SBarry Smith   *type = mg->am;
9454b9ad928SBarry Smith   PetscFunctionReturn(0);
9464b9ad928SBarry Smith }
9474b9ad928SBarry Smith 
9484b9ad928SBarry Smith #undef __FUNCT__
9490d353602SBarry Smith #define __FUNCT__ "PCMGSetCycleType"
9504b9ad928SBarry Smith /*@
9510d353602SBarry Smith    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
9524b9ad928SBarry Smith    complicated cycling.
9534b9ad928SBarry Smith 
954ad4df100SBarry Smith    Logically Collective on PC
9554b9ad928SBarry Smith 
9564b9ad928SBarry Smith    Input Parameters:
957c2be2410SBarry Smith +  pc - the multigrid context
9580d353602SBarry Smith -  PC_MG_CYCLE_V or PC_MG_CYCLE_W
9594b9ad928SBarry Smith 
9604b9ad928SBarry Smith    Options Database Key:
9614446f3b4SBarry Smith .  -pc_mg_cycle_type <v,w>
9624b9ad928SBarry Smith 
9634b9ad928SBarry Smith    Level: advanced
9644b9ad928SBarry Smith 
9654b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
9664b9ad928SBarry Smith 
9670d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel()
9684b9ad928SBarry Smith @*/
9697087cfbeSBarry Smith PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
9704b9ad928SBarry Smith {
971f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
972f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
97379416396SBarry Smith   PetscInt     i,levels;
9744b9ad928SBarry Smith 
9754b9ad928SBarry Smith   PetscFunctionBegin;
9760700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
977ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
97869fbec6eSBarry Smith   PetscValidLogicalCollectiveEnum(pc,n,2);
979f3fbd535SBarry Smith   levels = mglevels[0]->levels;
9804b9ad928SBarry Smith 
9812fa5cd67SKarl Rupp   for (i=0; i<levels; i++) mglevels[i]->cycles = n;
9824b9ad928SBarry Smith   PetscFunctionReturn(0);
9834b9ad928SBarry Smith }
9844b9ad928SBarry Smith 
9854b9ad928SBarry Smith #undef __FUNCT__
9868cc2d5dfSBarry Smith #define __FUNCT__ "PCMGMultiplicativeSetCycles"
9878cc2d5dfSBarry Smith /*@
9888cc2d5dfSBarry Smith    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
9898cc2d5dfSBarry Smith          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
9908cc2d5dfSBarry Smith 
991ad4df100SBarry Smith    Logically Collective on PC
9928cc2d5dfSBarry Smith 
9938cc2d5dfSBarry Smith    Input Parameters:
9948cc2d5dfSBarry Smith +  pc - the multigrid context
9958cc2d5dfSBarry Smith -  n - number of cycles (default is 1)
9968cc2d5dfSBarry Smith 
9978cc2d5dfSBarry Smith    Options Database Key:
998e1bc860dSBarry Smith .  -pc_mg_multiplicative_cycles n
9998cc2d5dfSBarry Smith 
10008cc2d5dfSBarry Smith    Level: advanced
10018cc2d5dfSBarry Smith 
10028cc2d5dfSBarry Smith    Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
10038cc2d5dfSBarry Smith 
10048cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
10058cc2d5dfSBarry Smith 
10068cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
10078cc2d5dfSBarry Smith @*/
10087087cfbeSBarry Smith PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
10098cc2d5dfSBarry Smith {
1010f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
10118cc2d5dfSBarry Smith 
10128cc2d5dfSBarry Smith   PetscFunctionBegin;
10130700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1014c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
10152a21e185SBarry Smith   mg->cyclesperpcapply = n;
10168cc2d5dfSBarry Smith   PetscFunctionReturn(0);
10178cc2d5dfSBarry Smith }
10188cc2d5dfSBarry Smith 
10198cc2d5dfSBarry Smith #undef __FUNCT__
1020fb15c04eSBarry Smith #define __FUNCT__ "PCMGSetGalerkin_MG"
1021fb15c04eSBarry Smith PetscErrorCode PCMGSetGalerkin_MG(PC pc,PetscBool use)
1022fb15c04eSBarry Smith {
1023fb15c04eSBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
1024fb15c04eSBarry Smith 
1025fb15c04eSBarry Smith   PetscFunctionBegin;
1026fb15c04eSBarry Smith   mg->galerkin = use ? 1 : 0;
1027fb15c04eSBarry Smith   PetscFunctionReturn(0);
1028fb15c04eSBarry Smith }
1029fb15c04eSBarry Smith 
1030fb15c04eSBarry Smith #undef __FUNCT__
10319dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetGalerkin"
1032c2be2410SBarry Smith /*@
103397177400SBarry Smith    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
103482b87d87SMatthew G. Knepley       finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i
1035c2be2410SBarry Smith 
1036ad4df100SBarry Smith    Logically Collective on PC
1037c2be2410SBarry Smith 
1038c2be2410SBarry Smith    Input Parameters:
1039c91913d3SJed Brown +  pc - the multigrid context
1040c91913d3SJed Brown -  use - PETSC_TRUE to use the Galerkin process to compute coarse-level operators
1041c2be2410SBarry Smith 
1042c2be2410SBarry Smith    Options Database Key:
10431c1aac46SBarry Smith .  -pc_mg_galerkin <true,false>
1044c2be2410SBarry Smith 
1045c2be2410SBarry Smith    Level: intermediate
1046c2be2410SBarry Smith 
10471c1aac46SBarry Smith    Notes: Some codes that use PCMG such as PCGAMG use Galerkin internally while constructing the hierarchy and thus do not
10481c1aac46SBarry Smith      use the PCMG construction of the coarser grids.
10491c1aac46SBarry Smith 
1050c2be2410SBarry Smith .keywords: MG, set, Galerkin
1051c2be2410SBarry Smith 
10523fc8bf9cSBarry Smith .seealso: PCMGGetGalerkin()
10533fc8bf9cSBarry Smith 
1054c2be2410SBarry Smith @*/
1055c91913d3SJed Brown PetscErrorCode PCMGSetGalerkin(PC pc,PetscBool use)
1056c2be2410SBarry Smith {
1057fb15c04eSBarry Smith   PetscErrorCode ierr;
1058c2be2410SBarry Smith 
1059c2be2410SBarry Smith   PetscFunctionBegin;
10600700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10617a4c0024SBarry Smith   ierr = PetscTryMethod(pc,"PCMGSetGalerkin_C",(PC,PetscBool),(pc,use));CHKERRQ(ierr);
1062c2be2410SBarry Smith   PetscFunctionReturn(0);
1063c2be2410SBarry Smith }
1064c2be2410SBarry Smith 
1065c2be2410SBarry Smith #undef __FUNCT__
10663fc8bf9cSBarry Smith #define __FUNCT__ "PCMGGetGalerkin"
10673fc8bf9cSBarry Smith /*@
10683fc8bf9cSBarry Smith    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
106982b87d87SMatthew G. Knepley       A_i-1 = r_i * A_i * p_i
10703fc8bf9cSBarry Smith 
10713fc8bf9cSBarry Smith    Not Collective
10723fc8bf9cSBarry Smith 
10733fc8bf9cSBarry Smith    Input Parameter:
10743fc8bf9cSBarry Smith .  pc - the multigrid context
10753fc8bf9cSBarry Smith 
10763fc8bf9cSBarry Smith    Output Parameter:
1077e1bc860dSBarry Smith .  galerkin - PETSC_TRUE or PETSC_FALSE
10783fc8bf9cSBarry Smith 
10793fc8bf9cSBarry Smith    Options Database Key:
1080e1bc860dSBarry Smith .  -pc_mg_galerkin
10813fc8bf9cSBarry Smith 
10823fc8bf9cSBarry Smith    Level: intermediate
10833fc8bf9cSBarry Smith 
10843fc8bf9cSBarry Smith .keywords: MG, set, Galerkin
10853fc8bf9cSBarry Smith 
10863fc8bf9cSBarry Smith .seealso: PCMGSetGalerkin()
10873fc8bf9cSBarry Smith 
10883fc8bf9cSBarry Smith @*/
10897087cfbeSBarry Smith PetscErrorCode  PCMGGetGalerkin(PC pc,PetscBool  *galerkin)
10903fc8bf9cSBarry Smith {
1091f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
10923fc8bf9cSBarry Smith 
10933fc8bf9cSBarry Smith   PetscFunctionBegin;
10940700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
109513fdb3acSJose Roman   *galerkin = (PetscBool)mg->galerkin;
10963fc8bf9cSBarry Smith   PetscFunctionReturn(0);
10973fc8bf9cSBarry Smith }
10983fc8bf9cSBarry Smith 
10993fc8bf9cSBarry Smith #undef __FUNCT__
11009dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothDown"
11014b9ad928SBarry Smith /*@
110297177400SBarry Smith    PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
110397177400SBarry Smith    use on all levels. Use PCMGGetSmootherDown() to set different
11044b9ad928SBarry Smith    pre-smoothing steps on different levels.
11054b9ad928SBarry Smith 
1106ad4df100SBarry Smith    Logically Collective on PC
11074b9ad928SBarry Smith 
11084b9ad928SBarry Smith    Input Parameters:
11094b9ad928SBarry Smith +  mg - the multigrid context
11104b9ad928SBarry Smith -  n - the number of smoothing steps
11114b9ad928SBarry Smith 
11124b9ad928SBarry Smith    Options Database Key:
11134b9ad928SBarry Smith .  -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps
11144b9ad928SBarry Smith 
11154b9ad928SBarry Smith    Level: advanced
11164b9ad928SBarry Smith 
11174b9ad928SBarry Smith .keywords: MG, smooth, down, pre-smoothing, steps, multigrid
11184b9ad928SBarry Smith 
111997177400SBarry Smith .seealso: PCMGSetNumberSmoothUp()
11204b9ad928SBarry Smith @*/
11217087cfbeSBarry Smith PetscErrorCode  PCMGSetNumberSmoothDown(PC pc,PetscInt n)
11224b9ad928SBarry Smith {
1123f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
1124f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
11256849ba73SBarry Smith   PetscErrorCode ierr;
112679416396SBarry Smith   PetscInt       i,levels;
11274b9ad928SBarry Smith 
11284b9ad928SBarry Smith   PetscFunctionBegin;
11290700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1130ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1131c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
1132f3fbd535SBarry Smith   levels = mglevels[0]->levels;
11334b9ad928SBarry Smith 
1134b05257ddSBarry Smith   for (i=1; i<levels; i++) {
11354b9ad928SBarry Smith     /* make sure smoother up and down are different */
11360298fd71SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,NULL);CHKERRQ(ierr);
1137f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
11382fa5cd67SKarl Rupp 
113931567311SBarry Smith     mg->default_smoothd = n;
11404b9ad928SBarry Smith   }
11414b9ad928SBarry Smith   PetscFunctionReturn(0);
11424b9ad928SBarry Smith }
11434b9ad928SBarry Smith 
11444b9ad928SBarry Smith #undef __FUNCT__
11459dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothUp"
11464b9ad928SBarry Smith /*@
114797177400SBarry Smith    PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
114897177400SBarry Smith    on all levels. Use PCMGGetSmootherUp() to set different numbers of
11494b9ad928SBarry Smith    post-smoothing steps on different levels.
11504b9ad928SBarry Smith 
1151ad4df100SBarry Smith    Logically Collective on PC
11524b9ad928SBarry Smith 
11534b9ad928SBarry Smith    Input Parameters:
11544b9ad928SBarry Smith +  mg - the multigrid context
11554b9ad928SBarry Smith -  n - the number of smoothing steps
11564b9ad928SBarry Smith 
11574b9ad928SBarry Smith    Options Database Key:
11584b9ad928SBarry Smith .  -pc_mg_smoothup <n> - Sets number of post-smoothing steps
11594b9ad928SBarry Smith 
11604b9ad928SBarry Smith    Level: advanced
11614b9ad928SBarry Smith 
11624b9ad928SBarry Smith    Note: this does not set a value on the coarsest grid, since we assume that
1163a8c7a070SBarry Smith     there is no separate smooth up on the coarsest grid.
11644b9ad928SBarry Smith 
11654b9ad928SBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid
11664b9ad928SBarry Smith 
116797177400SBarry Smith .seealso: PCMGSetNumberSmoothDown()
11684b9ad928SBarry Smith @*/
11697087cfbeSBarry Smith PetscErrorCode  PCMGSetNumberSmoothUp(PC pc,PetscInt n)
11704b9ad928SBarry Smith {
1171f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
1172f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
11736849ba73SBarry Smith   PetscErrorCode ierr;
117479416396SBarry Smith   PetscInt       i,levels;
11754b9ad928SBarry Smith 
11764b9ad928SBarry Smith   PetscFunctionBegin;
11770700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1178ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1179c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
1180f3fbd535SBarry Smith   levels = mglevels[0]->levels;
11814b9ad928SBarry Smith 
11824b9ad928SBarry Smith   for (i=1; i<levels; i++) {
11834b9ad928SBarry Smith     /* make sure smoother up and down are different */
11840298fd71SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,NULL);CHKERRQ(ierr);
1185f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
11862fa5cd67SKarl Rupp 
118731567311SBarry Smith     mg->default_smoothu = n;
11884b9ad928SBarry Smith   }
11894b9ad928SBarry Smith   PetscFunctionReturn(0);
11904b9ad928SBarry Smith }
11914b9ad928SBarry Smith 
11924b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/
11934b9ad928SBarry Smith 
11943b09bd56SBarry Smith /*MC
1195ccb205f8SBarry Smith    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
11963b09bd56SBarry Smith     information about the coarser grid matrices and restriction/interpolation operators.
11973b09bd56SBarry Smith 
11983b09bd56SBarry Smith    Options Database Keys:
11993b09bd56SBarry Smith +  -pc_mg_levels <nlevels> - number of levels including finest
12004afba7b0SPatrick Sanan .  -pc_mg_cycle_type <v,w> -
120179416396SBarry Smith .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
12023b09bd56SBarry Smith .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
12038c1c2452SJed Brown .  -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default
12043b09bd56SBarry Smith .  -pc_mg_log - log information about time spent on each level of the solver
12058cf6037eSBarry Smith .  -pc_mg_galerkin - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R'
12068cf6037eSBarry Smith .  -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1)
12078cf6037eSBarry Smith .  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1208e3c5b3baSBarry Smith                         to the Socket viewer for reading from MATLAB.
12098cf6037eSBarry Smith -  -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices
12108cf6037eSBarry Smith                         to the binary output file called binaryoutput
12113b09bd56SBarry Smith 
1212*e941fc33SBarry Smith    Notes: If one uses a Krylov method such GMRES or CG as the smoother than one must use KSPFGMRES, KSPGCG, or KSPRICHARDSON as the outer Krylov method
12133b09bd56SBarry Smith 
12148cf6037eSBarry Smith        When run with a single level the smoother options are used on that level NOT the coarse grid solver options
12158cf6037eSBarry Smith 
12163b09bd56SBarry Smith    Level: intermediate
12173b09bd56SBarry Smith 
12188f87f92bSBarry Smith    Concepts: multigrid/multilevel
12193b09bd56SBarry Smith 
12208cf6037eSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE
12210d353602SBarry Smith            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
122297177400SBarry Smith            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
122397177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
12240d353602SBarry Smith            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
12253b09bd56SBarry Smith M*/
12263b09bd56SBarry Smith 
12274b9ad928SBarry Smith #undef __FUNCT__
12284b9ad928SBarry Smith #define __FUNCT__ "PCCreate_MG"
12298cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc)
12304b9ad928SBarry Smith {
1231f3fbd535SBarry Smith   PC_MG          *mg;
1232f3fbd535SBarry Smith   PetscErrorCode ierr;
1233f3fbd535SBarry Smith 
12344b9ad928SBarry Smith   PetscFunctionBegin;
1235b00a9115SJed Brown   ierr        = PetscNewLog(pc,&mg);CHKERRQ(ierr);
1236f3fbd535SBarry Smith   pc->data    = (void*)mg;
1237f3fbd535SBarry Smith   mg->nlevels = -1;
123810eca3edSLisandro Dalcin   mg->am      = PC_MG_MULTIPLICATIVE;
1239f3fbd535SBarry Smith 
124037a44384SMark Adams   pc->useAmat = PETSC_TRUE;
124137a44384SMark Adams 
12424b9ad928SBarry Smith   pc->ops->apply          = PCApply_MG;
12434b9ad928SBarry Smith   pc->ops->setup          = PCSetUp_MG;
1244a06653b4SBarry Smith   pc->ops->reset          = PCReset_MG;
12454b9ad928SBarry Smith   pc->ops->destroy        = PCDestroy_MG;
12464b9ad928SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_MG;
12474b9ad928SBarry Smith   pc->ops->view           = PCView_MG;
1248fb15c04eSBarry Smith 
1249fb15c04eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",PCMGSetGalerkin_MG);CHKERRQ(ierr);
12504b9ad928SBarry Smith   PetscFunctionReturn(0);
12514b9ad928SBarry Smith }
1252