xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision 10167fec99cacc94e6fe858d292b158e9be1fe37)
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;
17710eca3edSLisandro Dalcin   PCMGType       mgtype     = mg->am;
178*10167fecSBarry Smith   PetscInt       mgctype    = (PetscInt) PC_MG_CYCLE_V;
179f3fbd535SBarry Smith   PetscInt       i;
180f3fbd535SBarry Smith   PetscMPIInt    size;
181f3fbd535SBarry Smith   const char     *prefix;
182f3fbd535SBarry Smith   PC             ipc;
1833aeaf226SBarry Smith   PetscInt       n;
1844b9ad928SBarry Smith 
1854b9ad928SBarry Smith   PetscFunctionBegin;
1860700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
187548767e0SJed Brown   PetscValidLogicalCollectiveInt(pc,levels,2);
188ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
189548767e0SJed Brown   if (mg->nlevels == levels) PetscFunctionReturn(0);
1903aeaf226SBarry Smith   if (mglevels) {
19110eca3edSLisandro Dalcin     mgctype = mglevels[0]->cycles;
1923aeaf226SBarry Smith     /* changing the number of levels so free up the previous stuff */
1933aeaf226SBarry Smith     ierr = PCReset_MG(pc);CHKERRQ(ierr);
1943aeaf226SBarry Smith     n    = mglevels[0]->levels;
1953aeaf226SBarry Smith     for (i=0; i<n; i++) {
1963aeaf226SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
1973aeaf226SBarry Smith         ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
1983aeaf226SBarry Smith       }
1993aeaf226SBarry Smith       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
2003aeaf226SBarry Smith       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
2013aeaf226SBarry Smith     }
2023aeaf226SBarry Smith     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
2033aeaf226SBarry Smith   }
204f3fbd535SBarry Smith 
205f3fbd535SBarry Smith   mg->nlevels = levels;
206f3fbd535SBarry Smith 
207785e854fSJed Brown   ierr = PetscMalloc1(levels,&mglevels);CHKERRQ(ierr);
2083bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr);
209f3fbd535SBarry Smith 
210f3fbd535SBarry Smith   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
211f3fbd535SBarry Smith 
212a9db87a2SMatthew G Knepley   mg->stageApply = 0;
213f3fbd535SBarry Smith   for (i=0; i<levels; i++) {
214b00a9115SJed Brown     ierr = PetscNewLog(pc,&mglevels[i]);CHKERRQ(ierr);
2152fa5cd67SKarl Rupp 
216f3fbd535SBarry Smith     mglevels[i]->level               = i;
217f3fbd535SBarry Smith     mglevels[i]->levels              = levels;
21810eca3edSLisandro Dalcin     mglevels[i]->cycles              = mgctype;
219336babb1SJed Brown     mg->default_smoothu              = 2;
220336babb1SJed Brown     mg->default_smoothd              = 2;
22163e6d426SJed Brown     mglevels[i]->eventsmoothsetup    = 0;
22263e6d426SJed Brown     mglevels[i]->eventsmoothsolve    = 0;
22363e6d426SJed Brown     mglevels[i]->eventresidual       = 0;
22463e6d426SJed Brown     mglevels[i]->eventinterprestrict = 0;
225f3fbd535SBarry Smith 
226f3fbd535SBarry Smith     if (comms) comm = comms[i];
227f3fbd535SBarry Smith     ierr = KSPCreate(comm,&mglevels[i]->smoothd);CHKERRQ(ierr);
228422a814eSBarry Smith     ierr = KSPSetErrorIfNotConverged(mglevels[i]->smoothd,pc->erroriffailure);CHKERRQ(ierr);
229336babb1SJed Brown     ierr = KSPSetType(mglevels[i]->smoothd,KSPCHEBYSHEV);CHKERRQ(ierr);
2300059c7bdSJed Brown     ierr = KSPSetConvergenceTest(mglevels[i]->smoothd,KSPConvergedSkip,NULL,NULL);CHKERRQ(ierr);
231336babb1SJed Brown     ierr = KSPSetNormType(mglevels[i]->smoothd,KSP_NORM_NONE);CHKERRQ(ierr);
232336babb1SJed Brown     ierr = KSPGetPC(mglevels[i]->smoothd,&ipc);CHKERRQ(ierr);
233336babb1SJed Brown     ierr = PCSetType(ipc,PCSOR);CHKERRQ(ierr);
234f3fbd535SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr);
235336babb1SJed Brown     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, i ? mg->default_smoothd : 1);CHKERRQ(ierr);
236f3fbd535SBarry Smith     ierr = KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);CHKERRQ(ierr);
237f3fbd535SBarry Smith 
238f3fbd535SBarry Smith     /* do special stuff for coarse grid */
239f3fbd535SBarry Smith     if (!i && levels > 1) {
240f3fbd535SBarry Smith       ierr = KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");CHKERRQ(ierr);
241f3fbd535SBarry Smith 
2429dbfc187SHong Zhang       /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */
243f3fbd535SBarry Smith       ierr = KSPSetType(mglevels[0]->smoothd,KSPPREONLY);CHKERRQ(ierr);
244f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[0]->smoothd,&ipc);CHKERRQ(ierr);
245f3fbd535SBarry Smith       ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
246f3fbd535SBarry Smith       if (size > 1) {
2479dbfc187SHong Zhang         KSP innerksp;
2489dbfc187SHong Zhang         PC  innerpc;
249f3fbd535SBarry Smith         ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr);
2509dbfc187SHong Zhang         ierr = PCRedundantGetKSP(ipc,&innerksp);CHKERRQ(ierr);
2519dbfc187SHong Zhang         ierr = KSPGetPC(innerksp,&innerpc);CHKERRQ(ierr);
2529dbfc187SHong Zhang         ierr = PCFactorSetShiftType(innerpc,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
253f3fbd535SBarry Smith       } else {
254f3fbd535SBarry Smith         ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr);
2559dbfc187SHong Zhang         ierr = PCFactorSetShiftType(ipc,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
256f3fbd535SBarry Smith       }
257f3fbd535SBarry Smith     } else {
258f3fbd535SBarry Smith       char tprefix[128];
259f3fbd535SBarry Smith       sprintf(tprefix,"mg_levels_%d_",(int)i);
260f3fbd535SBarry Smith       ierr = KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);CHKERRQ(ierr);
261f3fbd535SBarry Smith     }
2623bb1ff40SBarry Smith     ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[i]->smoothd);CHKERRQ(ierr);
2632fa5cd67SKarl Rupp 
264f3fbd535SBarry Smith     mglevels[i]->smoothu = mglevels[i]->smoothd;
26531567311SBarry Smith     mg->rtol             = 0.0;
26631567311SBarry Smith     mg->abstol           = 0.0;
26731567311SBarry Smith     mg->dtol             = 0.0;
26831567311SBarry Smith     mg->ttol             = 0.0;
26931567311SBarry Smith     mg->cyclesperpcapply = 1;
270f3fbd535SBarry Smith   }
271f3fbd535SBarry Smith   mg->levels = mglevels;
27210eca3edSLisandro Dalcin   ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
2734b9ad928SBarry Smith   PetscFunctionReturn(0);
2744b9ad928SBarry Smith }
2754b9ad928SBarry Smith 
276c07bf074SBarry Smith 
277c07bf074SBarry Smith #undef __FUNCT__
278c07bf074SBarry Smith #define __FUNCT__ "PCDestroy_MG"
279c07bf074SBarry Smith PetscErrorCode PCDestroy_MG(PC pc)
280c07bf074SBarry Smith {
281c07bf074SBarry Smith   PetscErrorCode ierr;
282a06653b4SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
283a06653b4SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
284a06653b4SBarry Smith   PetscInt       i,n;
285c07bf074SBarry Smith 
286c07bf074SBarry Smith   PetscFunctionBegin;
287a06653b4SBarry Smith   ierr = PCReset_MG(pc);CHKERRQ(ierr);
288a06653b4SBarry Smith   if (mglevels) {
289a06653b4SBarry Smith     n = mglevels[0]->levels;
290a06653b4SBarry Smith     for (i=0; i<n; i++) {
291a06653b4SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
2926bf464f9SBarry Smith         ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
293a06653b4SBarry Smith       }
2946bf464f9SBarry Smith       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
295a06653b4SBarry Smith       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
296a06653b4SBarry Smith     }
297a06653b4SBarry Smith     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
298a06653b4SBarry Smith   }
299c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
300f3fbd535SBarry Smith   PetscFunctionReturn(0);
301f3fbd535SBarry Smith }
302f3fbd535SBarry Smith 
303f3fbd535SBarry Smith 
304f3fbd535SBarry Smith 
30509573ac7SBarry Smith extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**);
30609573ac7SBarry Smith extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**);
30709573ac7SBarry Smith extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**);
308f3fbd535SBarry Smith 
309f3fbd535SBarry Smith /*
310f3fbd535SBarry Smith    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
311f3fbd535SBarry Smith              or full cycle of multigrid.
312f3fbd535SBarry Smith 
313f3fbd535SBarry Smith   Note:
314f3fbd535SBarry Smith   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
315f3fbd535SBarry Smith */
316f3fbd535SBarry Smith #undef __FUNCT__
317f3fbd535SBarry Smith #define __FUNCT__ "PCApply_MG"
318f3fbd535SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
319f3fbd535SBarry Smith {
320f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
321f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
322f3fbd535SBarry Smith   PetscErrorCode ierr;
323f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
324f3fbd535SBarry Smith 
325f3fbd535SBarry Smith   PetscFunctionBegin;
326a9db87a2SMatthew G Knepley   if (mg->stageApply) {ierr = PetscLogStagePush(mg->stageApply);CHKERRQ(ierr);}
327e1d8e5deSBarry Smith   /* When the DM is supplying the matrix then it will not exist until here */
328298cc208SBarry Smith   for (i=0; i<levels; i++) {
329e1d8e5deSBarry Smith     if (!mglevels[i]->A) {
33023ee1639SBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr);
331298cc208SBarry Smith       ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr);
332e1d8e5deSBarry Smith     }
333e1d8e5deSBarry Smith   }
334e1d8e5deSBarry Smith 
335f3fbd535SBarry Smith   mglevels[levels-1]->b = b;
336f3fbd535SBarry Smith   mglevels[levels-1]->x = x;
33731567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
338f3fbd535SBarry Smith     ierr = VecSet(x,0.0);CHKERRQ(ierr);
33931567311SBarry Smith     for (i=0; i<mg->cyclesperpcapply; i++) {
3400298fd71SBarry Smith       ierr = PCMGMCycle_Private(pc,mglevels+levels-1,NULL);CHKERRQ(ierr);
341f3fbd535SBarry Smith     }
3422fa5cd67SKarl Rupp   } else if (mg->am == PC_MG_ADDITIVE) {
34331567311SBarry Smith     ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr);
3442fa5cd67SKarl Rupp   } else if (mg->am == PC_MG_KASKADE) {
34531567311SBarry Smith     ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr);
3462fa5cd67SKarl Rupp   } else {
347f3fbd535SBarry Smith     ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr);
348f3fbd535SBarry Smith   }
349a9db87a2SMatthew G Knepley   if (mg->stageApply) {ierr = PetscLogStagePop();CHKERRQ(ierr);}
350f3fbd535SBarry Smith   PetscFunctionReturn(0);
351f3fbd535SBarry Smith }
352f3fbd535SBarry Smith 
353f3fbd535SBarry Smith 
354f3fbd535SBarry Smith #undef __FUNCT__
355f3fbd535SBarry Smith #define __FUNCT__ "PCSetFromOptions_MG"
3568c34d3f5SBarry Smith PetscErrorCode PCSetFromOptions_MG(PetscOptions *PetscOptionsObject,PC pc)
357f3fbd535SBarry Smith {
358f3fbd535SBarry Smith   PetscErrorCode ierr;
359f3fbd535SBarry Smith   PetscInt       m,levels = 1,cycles;
360c91913d3SJed Brown   PetscBool      flg,set;
361f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
362f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
363f3fbd535SBarry Smith   PCMGType       mgtype;
364f3fbd535SBarry Smith   PCMGCycleType  mgctype;
365f3fbd535SBarry Smith 
366f3fbd535SBarry Smith   PetscFunctionBegin;
367e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Multigrid options");CHKERRQ(ierr);
3683aeaf226SBarry Smith   if (!mg->levels) {
369f3fbd535SBarry Smith     ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr);
370eb3f98d2SBarry Smith     if (!flg && pc->dm) {
371eb3f98d2SBarry Smith       ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
372eb3f98d2SBarry Smith       levels++;
3733aeaf226SBarry Smith       mg->usedmfornumberoflevels = PETSC_TRUE;
374eb3f98d2SBarry Smith     }
3750298fd71SBarry Smith     ierr = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr);
376f3fbd535SBarry Smith   }
3773aeaf226SBarry Smith   mglevels = mg->levels;
3783aeaf226SBarry Smith 
379f3fbd535SBarry Smith   mgctype = (PCMGCycleType) mglevels[0]->cycles;
380f3fbd535SBarry Smith   ierr    = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr);
381f3fbd535SBarry Smith   if (flg) {
382f3fbd535SBarry Smith     ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr);
3832fa5cd67SKarl Rupp   }
384f3fbd535SBarry Smith   flg  = PETSC_FALSE;
385c91913d3SJed Brown   ierr = PetscOptionsBool("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",flg,&flg,&set);CHKERRQ(ierr);
386c91913d3SJed Brown   if (set) {
387c91913d3SJed Brown     ierr = PCMGSetGalerkin(pc,flg);CHKERRQ(ierr);
388f3fbd535SBarry Smith   }
38994ae4db5SBarry Smith   ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",mg->default_smoothu,&m,&flg);CHKERRQ(ierr);
390f3fbd535SBarry Smith   if (flg) {
391f3fbd535SBarry Smith     ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr);
392f3fbd535SBarry Smith   }
39394ae4db5SBarry Smith   ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",mg->default_smoothd,&m,&flg);CHKERRQ(ierr);
394f3fbd535SBarry Smith   if (flg) {
395f3fbd535SBarry Smith     ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr);
396f3fbd535SBarry Smith   }
39731567311SBarry Smith   mgtype = mg->am;
398f3fbd535SBarry Smith   ierr   = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
399f3fbd535SBarry Smith   if (flg) {
400f3fbd535SBarry Smith     ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
401f3fbd535SBarry Smith   }
40231567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
4035363c1e0SLisandro Dalcin     ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGMultiplicativeSetCycles",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr);
404f3fbd535SBarry Smith     if (flg) {
405f3fbd535SBarry Smith       ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr);
406f3fbd535SBarry Smith     }
407f3fbd535SBarry Smith   }
408f3fbd535SBarry Smith   flg  = PETSC_FALSE;
4090298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,NULL);CHKERRQ(ierr);
410f3fbd535SBarry Smith   if (flg) {
411f3fbd535SBarry Smith     PetscInt i;
412f3fbd535SBarry Smith     char     eventname[128];
413ce94432eSBarry Smith     if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
414f3fbd535SBarry Smith     levels = mglevels[0]->levels;
415f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
416f3fbd535SBarry Smith       sprintf(eventname,"MGSetup Level %d",(int)i);
41763e6d426SJed Brown       ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr);
418f3fbd535SBarry Smith       sprintf(eventname,"MGSmooth Level %d",(int)i);
41963e6d426SJed Brown       ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr);
420f3fbd535SBarry Smith       if (i) {
421f3fbd535SBarry Smith         sprintf(eventname,"MGResid Level %d",(int)i);
42263e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr);
423f3fbd535SBarry Smith         sprintf(eventname,"MGInterp Level %d",(int)i);
42463e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr);
425f3fbd535SBarry Smith       }
426f3fbd535SBarry Smith     }
427eec35531SMatthew G Knepley 
428ec60ca73SSatish Balay #if defined(PETSC_USE_LOG)
429eec35531SMatthew G Knepley     {
430eec35531SMatthew G Knepley       const char    *sname = "MG Apply";
431eec35531SMatthew G Knepley       PetscStageLog stageLog;
432eec35531SMatthew G Knepley       PetscInt      st;
433eec35531SMatthew G Knepley 
434eec35531SMatthew G Knepley       PetscFunctionBegin;
435eec35531SMatthew G Knepley       ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
436eec35531SMatthew G Knepley       for (st = 0; st < stageLog->numStages; ++st) {
437eec35531SMatthew G Knepley         PetscBool same;
438eec35531SMatthew G Knepley 
439eec35531SMatthew G Knepley         ierr = PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);CHKERRQ(ierr);
4402fa5cd67SKarl Rupp         if (same) mg->stageApply = st;
441eec35531SMatthew G Knepley       }
442eec35531SMatthew G Knepley       if (!mg->stageApply) {
443eec35531SMatthew G Knepley         ierr = PetscLogStageRegister(sname, &mg->stageApply);CHKERRQ(ierr);
444eec35531SMatthew G Knepley       }
445eec35531SMatthew G Knepley     }
446ec60ca73SSatish Balay #endif
447f3fbd535SBarry Smith   }
448f3fbd535SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
449f3fbd535SBarry Smith   PetscFunctionReturn(0);
450f3fbd535SBarry Smith }
451f3fbd535SBarry Smith 
4526a6fc655SJed Brown const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0};
4536a6fc655SJed Brown const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0};
454f3fbd535SBarry Smith 
4559804daf3SBarry Smith #include <petscdraw.h>
456f3fbd535SBarry Smith #undef __FUNCT__
457f3fbd535SBarry Smith #define __FUNCT__ "PCView_MG"
458f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
459f3fbd535SBarry Smith {
460f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
461f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
462f3fbd535SBarry Smith   PetscErrorCode ierr;
463e3deeeafSJed Brown   PetscInt       levels = mglevels ? mglevels[0]->levels : 0,i;
464219b1045SBarry Smith   PetscBool      iascii,isbinary,isdraw;
465f3fbd535SBarry Smith 
466f3fbd535SBarry Smith   PetscFunctionBegin;
467251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
4685b0b0462SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
469219b1045SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
470f3fbd535SBarry Smith   if (iascii) {
471e3deeeafSJed Brown     const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown";
472e3deeeafSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  MG: type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,cyclename);CHKERRQ(ierr);
47331567311SBarry Smith     if (mg->am == PC_MG_MULTIPLICATIVE) {
47431567311SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr);
475f3fbd535SBarry Smith     }
476218a07d4SBarry Smith     if (mg->galerkin) {
477f3fbd535SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
4784f66f45eSBarry Smith     } else {
4794f66f45eSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
480f3fbd535SBarry Smith     }
4815adeb434SBarry Smith     if (mg->view){
4825adeb434SBarry Smith       ierr = (*mg->view)(pc,viewer);CHKERRQ(ierr);
4835adeb434SBarry Smith     }
484f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
485f3fbd535SBarry Smith       if (!i) {
48689cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);CHKERRQ(ierr);
487f3fbd535SBarry Smith       } else {
48889cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
489f3fbd535SBarry Smith       }
490f3fbd535SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
491f3fbd535SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
492f3fbd535SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
493f3fbd535SBarry Smith       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
494f3fbd535SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
495f3fbd535SBarry Smith       } else if (i) {
49689cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
497f3fbd535SBarry Smith         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
498f3fbd535SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
499f3fbd535SBarry Smith         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
500f3fbd535SBarry Smith       }
501f3fbd535SBarry Smith     }
5025b0b0462SBarry Smith   } else if (isbinary) {
5035b0b0462SBarry Smith     for (i=levels-1; i>=0; i--) {
5045b0b0462SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
5055b0b0462SBarry Smith       if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) {
5065b0b0462SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
5075b0b0462SBarry Smith       }
5085b0b0462SBarry Smith     }
509219b1045SBarry Smith   } else if (isdraw) {
510219b1045SBarry Smith     PetscDraw draw;
511b4375e8dSPeter Brune     PetscReal x,w,y,bottom,th;
512219b1045SBarry Smith     ierr   = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
513219b1045SBarry Smith     ierr   = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
5140298fd71SBarry Smith     ierr   = PetscDrawStringGetSize(draw,NULL,&th);CHKERRQ(ierr);
515219b1045SBarry Smith     bottom = y - th;
516219b1045SBarry Smith     for (i=levels-1; i>=0; i--) {
517b4375e8dSPeter Brune       if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) {
518219b1045SBarry Smith         ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr);
519219b1045SBarry Smith         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
520219b1045SBarry Smith         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
521b4375e8dSPeter Brune       } else {
522b4375e8dSPeter Brune         w    = 0.5*PetscMin(1.0-x,x);
523b4375e8dSPeter Brune         ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr);
524b4375e8dSPeter Brune         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
525b4375e8dSPeter Brune         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
526b4375e8dSPeter Brune         ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr);
527b4375e8dSPeter Brune         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
528b4375e8dSPeter Brune         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
529b4375e8dSPeter Brune       }
5300298fd71SBarry Smith       ierr    = PetscDrawGetBoundingBox(draw,NULL,&bottom,NULL,NULL);CHKERRQ(ierr);
5311cd381d6SBarry Smith       bottom -= th;
532219b1045SBarry Smith     }
5335b0b0462SBarry Smith   }
534f3fbd535SBarry Smith   PetscFunctionReturn(0);
535f3fbd535SBarry Smith }
536f3fbd535SBarry Smith 
537af0996ceSBarry Smith #include <petsc/private/dmimpl.h>
538af0996ceSBarry Smith #include <petsc/private/kspimpl.h>
539cab2e9ccSBarry Smith 
540f3fbd535SBarry Smith /*
541f3fbd535SBarry Smith     Calls setup for the KSP on each level
542f3fbd535SBarry Smith */
543f3fbd535SBarry Smith #undef __FUNCT__
544f3fbd535SBarry Smith #define __FUNCT__ "PCSetUp_MG"
545f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc)
546f3fbd535SBarry Smith {
547f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
548f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
549f3fbd535SBarry Smith   PetscErrorCode ierr;
550f3fbd535SBarry Smith   PetscInt       i,n = mglevels[0]->levels;
55198e04984SBarry Smith   PC             cpc;
5529c7ffb39SBarry Smith   PetscBool      preonly,lu,redundant,cholesky,svd,dump = PETSC_FALSE,opsset,use_amat,missinginterpolate = PETSC_FALSE;
553f3fbd535SBarry Smith   Mat            dA,dB;
554f3fbd535SBarry Smith   Vec            tvec;
555218a07d4SBarry Smith   DM             *dms;
556649052a6SBarry Smith   PetscViewer    viewer = 0;
557f3fbd535SBarry Smith 
558f3fbd535SBarry Smith   PetscFunctionBegin;
55901bc414fSDmitry Karpeev   /* FIX: Move this to PCSetFromOptions_MG? */
5603aeaf226SBarry Smith   if (mg->usedmfornumberoflevels) {
5613aeaf226SBarry Smith     PetscInt levels;
5623aeaf226SBarry Smith     ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
5633aeaf226SBarry Smith     levels++;
5643aeaf226SBarry Smith     if (levels > n) { /* the problem is now being solved on a finer grid */
5650298fd71SBarry Smith       ierr     = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr);
5663aeaf226SBarry Smith       n        = levels;
5673aeaf226SBarry Smith       ierr     = PCSetFromOptions(pc);CHKERRQ(ierr); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
5683aeaf226SBarry Smith       mglevels =  mg->levels;
5693aeaf226SBarry Smith     }
5703aeaf226SBarry Smith   }
57198e04984SBarry Smith   ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr);
5723aeaf226SBarry Smith 
573f3fbd535SBarry Smith 
574f3fbd535SBarry Smith   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
575f3fbd535SBarry Smith   /* so use those from global PC */
576f3fbd535SBarry Smith   /* Is this what we always want? What if user wants to keep old one? */
5770298fd71SBarry Smith   ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset);CHKERRQ(ierr);
57898e04984SBarry Smith   if (opsset) {
57998e04984SBarry Smith     Mat mmat;
58023ee1639SBarry Smith     ierr = KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat);CHKERRQ(ierr);
58198e04984SBarry Smith     if (mmat == pc->pmat) opsset = PETSC_FALSE;
58298e04984SBarry Smith   }
5834a5f07a5SMark F. Adams 
5844a5f07a5SMark F. Adams   if (!opsset) {
58571b23a65SMark F. Adams     ierr = PCGetUseAmat(pc,&use_amat);CHKERRQ(ierr);
5864a5f07a5SMark F. Adams     if(use_amat){
587f3fbd535SBarry 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);
58823ee1639SBarry Smith       ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat);CHKERRQ(ierr);
589f3fbd535SBarry Smith     }
5904a5f07a5SMark F. Adams     else {
5914a5f07a5SMark 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);
59223ee1639SBarry Smith       ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat);CHKERRQ(ierr);
5934a5f07a5SMark F. Adams     }
5944a5f07a5SMark F. Adams   }
595f3fbd535SBarry Smith 
5969c7ffb39SBarry Smith   for (i=n-1; i>0; i--) {
5979c7ffb39SBarry Smith     if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) {
5989c7ffb39SBarry Smith       missinginterpolate = PETSC_TRUE;
5999c7ffb39SBarry Smith       continue;
6009c7ffb39SBarry Smith     }
6019c7ffb39SBarry Smith   }
6029c7ffb39SBarry Smith   /*
6039c7ffb39SBarry Smith    Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS)
6049c7ffb39SBarry Smith    Skipping for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs?
6059c7ffb39SBarry Smith   */
6069c7ffb39SBarry Smith   if (missinginterpolate && pc->dm && mg->galerkin != 2 && !pc->setupcalled) {
6072d2b81a6SBarry Smith     /* construct the interpolation from the DMs */
6082e499ae9SBarry Smith     Mat p;
60973250ac0SBarry Smith     Vec rscale;
610785e854fSJed Brown     ierr     = PetscMalloc1(n,&dms);CHKERRQ(ierr);
611218a07d4SBarry Smith     dms[n-1] = pc->dm;
612ef1267afSMatthew G. Knepley     /* Separately create them so we do not get DMKSP interference between levels */
613ef1267afSMatthew G. Knepley     for (i=n-2; i>-1; i--) {ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr);}
614218a07d4SBarry Smith     for (i=n-2; i>-1; i--) {
615942e3340SBarry Smith       DMKSP kdm;
6163c0c59f3SBarry Smith       ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr);
6173c0c59f3SBarry Smith       if (mg->galerkin) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);}
618942e3340SBarry Smith       ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr);
619d1a3e677SJed Brown       /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take
620d1a3e677SJed Brown        * a bitwise OR of computing the matrix, RHS, and initial iterate. */
6210298fd71SBarry Smith       kdm->ops->computerhs = NULL;
6220298fd71SBarry Smith       kdm->rhsctx          = NULL;
62324c3aa18SBarry Smith       if (!mglevels[i+1]->interpolate) {
624e727c939SJed Brown         ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr);
6252d2b81a6SBarry Smith         ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr);
62600ac8be1SBarry Smith         if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);}
62773250ac0SBarry Smith         ierr = VecDestroy(&rscale);CHKERRQ(ierr);
6286bf464f9SBarry Smith         ierr = MatDestroy(&p);CHKERRQ(ierr);
629218a07d4SBarry Smith       }
63024c3aa18SBarry Smith     }
6312d2b81a6SBarry Smith 
632ef1267afSMatthew G. Knepley     for (i=n-2; i>-1; i--) {ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);}
6332d2b81a6SBarry Smith     ierr = PetscFree(dms);CHKERRQ(ierr);
634ce4cda84SJed Brown   }
635cab2e9ccSBarry Smith 
636ce4cda84SJed Brown   if (pc->dm && !pc->setupcalled) {
637ce4cda84SJed Brown     /* finest smoother also gets DM but it is not active, independent of whether galerkin==2 */
638cab2e9ccSBarry Smith     ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr);
639cab2e9ccSBarry Smith     ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr);
640218a07d4SBarry Smith   }
641218a07d4SBarry Smith 
642c91913d3SJed Brown   if (mg->galerkin == 1) {
643f3fbd535SBarry Smith     Mat B;
644f3fbd535SBarry Smith     /* currently only handle case where mat and pmat are the same on coarser levels */
64523ee1639SBarry Smith     ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB);CHKERRQ(ierr);
646f3fbd535SBarry Smith     if (!pc->setupcalled) {
647f3fbd535SBarry Smith       for (i=n-2; i>-1; i--) {
648b9367914SBarry 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");
649b9367914SBarry Smith         if (!mglevels[i+1]->interpolate) {
650b9367914SBarry Smith           ierr = PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct);CHKERRQ(ierr);
651b9367914SBarry Smith         }
652b9367914SBarry Smith         if (!mglevels[i+1]->restrct) {
653b9367914SBarry Smith           ierr = PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate);CHKERRQ(ierr);
654b9367914SBarry Smith         }
655b9367914SBarry Smith         if (mglevels[i+1]->interpolate == mglevels[i+1]->restrct) {
656f3fbd535SBarry Smith           ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
6577d537d92SJed Brown         } else {
6587d537d92SJed Brown           ierr = MatMatMatMult(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
6597d537d92SJed Brown         }
66023ee1639SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B);CHKERRQ(ierr);
661f3fbd535SBarry Smith         if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
662f3fbd535SBarry Smith         dB = B;
663f3fbd535SBarry Smith       }
664cd9507b2SBarry Smith       if (n > 1) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
665f3fbd535SBarry Smith     } else {
666f3fbd535SBarry Smith       for (i=n-2; i>-1; i--) {
667b9367914SBarry 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");
668b9367914SBarry Smith         if (!mglevels[i+1]->interpolate) {
669b9367914SBarry Smith           ierr = PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct);CHKERRQ(ierr);
670b9367914SBarry Smith         }
671b9367914SBarry Smith         if (!mglevels[i+1]->restrct) {
672b9367914SBarry Smith           ierr = PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate);CHKERRQ(ierr);
673b9367914SBarry Smith         }
67423ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&B);CHKERRQ(ierr);
675b9367914SBarry Smith         if (mglevels[i+1]->interpolate == mglevels[i+1]->restrct) {
676f3fbd535SBarry Smith           ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
6777d537d92SJed Brown         } else {
6787d537d92SJed Brown           ierr = MatMatMatMult(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
6797d537d92SJed Brown         }
68023ee1639SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B);CHKERRQ(ierr);
681f3fbd535SBarry Smith         dB   = B;
682f3fbd535SBarry Smith       }
683f3fbd535SBarry Smith     }
684ce4cda84SJed Brown   } else if (!mg->galerkin && pc->dm && pc->dm->x) {
685cab2e9ccSBarry Smith     /* need to restrict Jacobian location to coarser meshes for evaluation */
686cab2e9ccSBarry Smith     for (i=n-2; i>-1; i--) {
687c88c5224SJed Brown       Mat R;
688c88c5224SJed Brown       Vec rscale;
689cab2e9ccSBarry Smith       if (!mglevels[i]->smoothd->dm->x) {
690cab2e9ccSBarry Smith         Vec *vecs;
6912a7a6963SBarry Smith         ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vecs,0,NULL);CHKERRQ(ierr);
6922fa5cd67SKarl Rupp 
693cab2e9ccSBarry Smith         mglevels[i]->smoothd->dm->x = vecs[0];
6942fa5cd67SKarl Rupp 
695cab2e9ccSBarry Smith         ierr = PetscFree(vecs);CHKERRQ(ierr);
696cab2e9ccSBarry Smith       }
697c88c5224SJed Brown       ierr = PCMGGetRestriction(pc,i+1,&R);CHKERRQ(ierr);
698c88c5224SJed Brown       ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
699c88c5224SJed Brown       ierr = MatRestrict(R,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);CHKERRQ(ierr);
700c88c5224SJed Brown       ierr = VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,rscale);CHKERRQ(ierr);
701cab2e9ccSBarry Smith     }
702f3fbd535SBarry Smith   }
703ccceb50cSJed Brown   if (!mg->galerkin && pc->dm) {
704caa4e7f2SJed Brown     for (i=n-2; i>=0; i--) {
705ccceb50cSJed Brown       DM  dmfine,dmcoarse;
706caa4e7f2SJed Brown       Mat Restrict,Inject;
707caa4e7f2SJed Brown       Vec rscale;
708ccceb50cSJed Brown       ierr   = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr);
709ccceb50cSJed Brown       ierr   = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr);
710caa4e7f2SJed Brown       ierr   = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr);
711caa4e7f2SJed Brown       ierr   = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
7120298fd71SBarry Smith       Inject = NULL;      /* Callback should create it if it needs Injection */
713caa4e7f2SJed Brown       ierr   = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr);
714caa4e7f2SJed Brown     }
715caa4e7f2SJed Brown   }
716f3fbd535SBarry Smith 
717f3fbd535SBarry Smith   if (!pc->setupcalled) {
718f3fbd535SBarry Smith     for (i=0; i<n; i++) {
719f3fbd535SBarry Smith       ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr);
720f3fbd535SBarry Smith     }
721f3fbd535SBarry Smith     for (i=1; i<n; i++) {
722f3fbd535SBarry Smith       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
723f3fbd535SBarry Smith         ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr);
724f3fbd535SBarry Smith       }
725f3fbd535SBarry Smith     }
726f3fbd535SBarry Smith     for (i=1; i<n; i++) {
727c88c5224SJed Brown       ierr = PCMGGetInterpolation(pc,i,&mglevels[i]->interpolate);CHKERRQ(ierr);
728c88c5224SJed Brown       ierr = PCMGGetRestriction(pc,i,&mglevels[i]->restrct);CHKERRQ(ierr);
729f3fbd535SBarry Smith     }
730f3fbd535SBarry Smith     for (i=0; i<n-1; i++) {
731f3fbd535SBarry Smith       if (!mglevels[i]->b) {
732f3fbd535SBarry Smith         Vec *vec;
7332a7a6963SBarry Smith         ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
734f3fbd535SBarry Smith         ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr);
7356bf464f9SBarry Smith         ierr = VecDestroy(vec);CHKERRQ(ierr);
736f3fbd535SBarry Smith         ierr = PetscFree(vec);CHKERRQ(ierr);
737f3fbd535SBarry Smith       }
738f3fbd535SBarry Smith       if (!mglevels[i]->r && i) {
739f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
740f3fbd535SBarry Smith         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
7416bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
742f3fbd535SBarry Smith       }
743f3fbd535SBarry Smith       if (!mglevels[i]->x) {
744f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
745f3fbd535SBarry Smith         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
7466bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
747f3fbd535SBarry Smith       }
748f3fbd535SBarry Smith     }
749f3fbd535SBarry Smith     if (n != 1 && !mglevels[n-1]->r) {
750f3fbd535SBarry Smith       /* PCMGSetR() on the finest level if user did not supply it */
751f3fbd535SBarry Smith       Vec *vec;
7522a7a6963SBarry Smith       ierr = KSPCreateVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
753f3fbd535SBarry Smith       ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr);
7546bf464f9SBarry Smith       ierr = VecDestroy(vec);CHKERRQ(ierr);
755f3fbd535SBarry Smith       ierr = PetscFree(vec);CHKERRQ(ierr);
756f3fbd535SBarry Smith     }
757f3fbd535SBarry Smith   }
758f3fbd535SBarry Smith 
75998e04984SBarry Smith   if (pc->dm) {
76098e04984SBarry Smith     /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */
76198e04984SBarry Smith     for (i=0; i<n-1; i++) {
76298e04984SBarry Smith       if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX;
76398e04984SBarry Smith     }
76498e04984SBarry Smith   }
765f3fbd535SBarry Smith 
766f3fbd535SBarry Smith   for (i=1; i<n; i++) {
7672a21e185SBarry Smith     if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1){
768f3fbd535SBarry Smith       /* if doing only down then initial guess is zero */
769f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
770f3fbd535SBarry Smith     }
77163e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
772f3fbd535SBarry Smith     ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr);
77363e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
774d42688cbSBarry Smith     if (!mglevels[i]->residual) {
775d42688cbSBarry Smith       Mat mat;
77623ee1639SBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&mat);CHKERRQ(ierr);
77754b2cd4bSJed Brown       ierr = PCMGSetResidual(pc,i,PCMGResidualDefault,mat);CHKERRQ(ierr);
778d42688cbSBarry Smith     }
779f3fbd535SBarry Smith   }
780f3fbd535SBarry Smith   for (i=1; i<n; i++) {
781f3fbd535SBarry Smith     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
782f3fbd535SBarry Smith       Mat          downmat,downpmat;
783f3fbd535SBarry Smith 
784f3fbd535SBarry Smith       /* check if operators have been set for up, if not use down operators to set them */
7850298fd71SBarry Smith       ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);CHKERRQ(ierr);
786f3fbd535SBarry Smith       if (!opsset) {
78723ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);CHKERRQ(ierr);
78823ee1639SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat);CHKERRQ(ierr);
789f3fbd535SBarry Smith       }
790f3fbd535SBarry Smith 
791f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
79263e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
793f3fbd535SBarry Smith       ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr);
79463e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
795f3fbd535SBarry Smith     }
796f3fbd535SBarry Smith   }
797f3fbd535SBarry Smith 
798f3fbd535SBarry Smith   /*
799f3fbd535SBarry Smith       If coarse solver is not direct method then DO NOT USE preonly
800f3fbd535SBarry Smith   */
801251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mglevels[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr);
802f3fbd535SBarry Smith   if (preonly) {
803251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr);
804251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr);
805251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr);
806251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)cpc,PCSVD,&svd);CHKERRQ(ierr);
807fd1303e9SJungho Lee     if (!lu && !redundant && !cholesky && !svd) {
808f3fbd535SBarry Smith       ierr = KSPSetType(mglevels[0]->smoothd,KSPGMRES);CHKERRQ(ierr);
809f3fbd535SBarry Smith     }
810f3fbd535SBarry Smith   }
811f3fbd535SBarry Smith 
81263e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
813f3fbd535SBarry Smith   ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr);
81463e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
815f3fbd535SBarry Smith 
816f3fbd535SBarry Smith   /*
817f3fbd535SBarry Smith      Dump the interpolation/restriction matrices plus the
818e3c5b3baSBarry Smith    Jacobian/stiffness on each level. This allows MATLAB users to
819f3fbd535SBarry Smith    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
820f3fbd535SBarry Smith 
821f3fbd535SBarry Smith    Only support one or the other at the same time.
822f3fbd535SBarry Smith   */
823f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
8240298fd71SBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);CHKERRQ(ierr);
825ce94432eSBarry Smith   if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc));
826f3fbd535SBarry Smith   dump = PETSC_FALSE;
827f3fbd535SBarry Smith #endif
8280298fd71SBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);CHKERRQ(ierr);
829ce94432eSBarry Smith   if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc));
830f3fbd535SBarry Smith 
831f3fbd535SBarry Smith   if (viewer) {
832f3fbd535SBarry Smith     for (i=1; i<n; i++) {
833f3fbd535SBarry Smith       ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr);
834f3fbd535SBarry Smith     }
835f3fbd535SBarry Smith     for (i=0; i<n; i++) {
836f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr);
837f3fbd535SBarry Smith       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
838f3fbd535SBarry Smith     }
839f3fbd535SBarry Smith   }
840f3fbd535SBarry Smith   PetscFunctionReturn(0);
841f3fbd535SBarry Smith }
842f3fbd535SBarry Smith 
843f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/
844f3fbd535SBarry Smith 
845f3fbd535SBarry Smith #undef __FUNCT__
8469dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetLevels"
8474b9ad928SBarry Smith /*@
84897177400SBarry Smith    PCMGGetLevels - Gets the number of levels to use with MG.
8494b9ad928SBarry Smith 
8504b9ad928SBarry Smith    Not Collective
8514b9ad928SBarry Smith 
8524b9ad928SBarry Smith    Input Parameter:
8534b9ad928SBarry Smith .  pc - the preconditioner context
8544b9ad928SBarry Smith 
8554b9ad928SBarry Smith    Output parameter:
8564b9ad928SBarry Smith .  levels - the number of levels
8574b9ad928SBarry Smith 
8584b9ad928SBarry Smith    Level: advanced
8594b9ad928SBarry Smith 
8604b9ad928SBarry Smith .keywords: MG, get, levels, multigrid
8614b9ad928SBarry Smith 
86297177400SBarry Smith .seealso: PCMGSetLevels()
8634b9ad928SBarry Smith @*/
8647087cfbeSBarry Smith PetscErrorCode  PCMGGetLevels(PC pc,PetscInt *levels)
8654b9ad928SBarry Smith {
866f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
8674b9ad928SBarry Smith 
8684b9ad928SBarry Smith   PetscFunctionBegin;
8690700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
8704482741eSBarry Smith   PetscValidIntPointer(levels,2);
871f3fbd535SBarry Smith   *levels = mg->nlevels;
8724b9ad928SBarry Smith   PetscFunctionReturn(0);
8734b9ad928SBarry Smith }
8744b9ad928SBarry Smith 
8754b9ad928SBarry Smith #undef __FUNCT__
8769dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetType"
8774b9ad928SBarry Smith /*@
87897177400SBarry Smith    PCMGSetType - Determines the form of multigrid to use:
8794b9ad928SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
8804b9ad928SBarry Smith 
881ad4df100SBarry Smith    Logically Collective on PC
8824b9ad928SBarry Smith 
8834b9ad928SBarry Smith    Input Parameters:
8844b9ad928SBarry Smith +  pc - the preconditioner context
8859dcbbd2bSBarry Smith -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
8869dcbbd2bSBarry Smith    PC_MG_FULL, PC_MG_KASKADE
8874b9ad928SBarry Smith 
8884b9ad928SBarry Smith    Options Database Key:
8894b9ad928SBarry Smith .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
8904b9ad928SBarry Smith    additive, full, kaskade
8914b9ad928SBarry Smith 
8924b9ad928SBarry Smith    Level: advanced
8934b9ad928SBarry Smith 
8944b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
8954b9ad928SBarry Smith 
89697177400SBarry Smith .seealso: PCMGSetLevels()
8974b9ad928SBarry Smith @*/
8987087cfbeSBarry Smith PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
8994b9ad928SBarry Smith {
900f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
9014b9ad928SBarry Smith 
9024b9ad928SBarry Smith   PetscFunctionBegin;
9030700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
904c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,form,2);
90531567311SBarry Smith   mg->am = form;
9069dcbbd2bSBarry Smith   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
907c60c7ad4SBarry Smith   else pc->ops->applyrichardson = NULL;
908c60c7ad4SBarry Smith   PetscFunctionReturn(0);
909c60c7ad4SBarry Smith }
910c60c7ad4SBarry Smith 
911f0af28e8SLisandro Dalcin #undef __FUNCT__
912f0af28e8SLisandro Dalcin #define __FUNCT__ "PCMGGetType"
913c60c7ad4SBarry Smith /*@
914c60c7ad4SBarry Smith    PCMGGetType - Determines the form of multigrid to use:
915c60c7ad4SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
916c60c7ad4SBarry Smith 
917c60c7ad4SBarry Smith    Logically Collective on PC
918c60c7ad4SBarry Smith 
919c60c7ad4SBarry Smith    Input Parameter:
920c60c7ad4SBarry Smith .  pc - the preconditioner context
921c60c7ad4SBarry Smith 
922c60c7ad4SBarry Smith    Output Parameter:
923c60c7ad4SBarry Smith .  type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE
924c60c7ad4SBarry Smith 
925c60c7ad4SBarry Smith 
926c60c7ad4SBarry Smith    Level: advanced
927c60c7ad4SBarry Smith 
928c60c7ad4SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
929c60c7ad4SBarry Smith 
930c60c7ad4SBarry Smith .seealso: PCMGSetLevels()
931c60c7ad4SBarry Smith @*/
932c60c7ad4SBarry Smith PetscErrorCode  PCMGGetType(PC pc,PCMGType *type)
933c60c7ad4SBarry Smith {
934c60c7ad4SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
935c60c7ad4SBarry Smith 
936c60c7ad4SBarry Smith   PetscFunctionBegin;
937c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
938c60c7ad4SBarry Smith   *type = mg->am;
9394b9ad928SBarry Smith   PetscFunctionReturn(0);
9404b9ad928SBarry Smith }
9414b9ad928SBarry Smith 
9424b9ad928SBarry Smith #undef __FUNCT__
9430d353602SBarry Smith #define __FUNCT__ "PCMGSetCycleType"
9444b9ad928SBarry Smith /*@
9450d353602SBarry Smith    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
9464b9ad928SBarry Smith    complicated cycling.
9474b9ad928SBarry Smith 
948ad4df100SBarry Smith    Logically Collective on PC
9494b9ad928SBarry Smith 
9504b9ad928SBarry Smith    Input Parameters:
951c2be2410SBarry Smith +  pc - the multigrid context
9520d353602SBarry Smith -  PC_MG_CYCLE_V or PC_MG_CYCLE_W
9534b9ad928SBarry Smith 
9544b9ad928SBarry Smith    Options Database Key:
9554446f3b4SBarry Smith .  -pc_mg_cycle_type <v,w>
9564b9ad928SBarry Smith 
9574b9ad928SBarry Smith    Level: advanced
9584b9ad928SBarry Smith 
9594b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
9604b9ad928SBarry Smith 
9610d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel()
9624b9ad928SBarry Smith @*/
9637087cfbeSBarry Smith PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
9644b9ad928SBarry Smith {
965f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
966f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
96779416396SBarry Smith   PetscInt     i,levels;
9684b9ad928SBarry Smith 
9694b9ad928SBarry Smith   PetscFunctionBegin;
9700700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
971ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
97269fbec6eSBarry Smith   PetscValidLogicalCollectiveEnum(pc,n,2);
973f3fbd535SBarry Smith   levels = mglevels[0]->levels;
9744b9ad928SBarry Smith 
9752fa5cd67SKarl Rupp   for (i=0; i<levels; i++) mglevels[i]->cycles = n;
9764b9ad928SBarry Smith   PetscFunctionReturn(0);
9774b9ad928SBarry Smith }
9784b9ad928SBarry Smith 
9794b9ad928SBarry Smith #undef __FUNCT__
9808cc2d5dfSBarry Smith #define __FUNCT__ "PCMGMultiplicativeSetCycles"
9818cc2d5dfSBarry Smith /*@
9828cc2d5dfSBarry Smith    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
9838cc2d5dfSBarry Smith          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
9848cc2d5dfSBarry Smith 
985ad4df100SBarry Smith    Logically Collective on PC
9868cc2d5dfSBarry Smith 
9878cc2d5dfSBarry Smith    Input Parameters:
9888cc2d5dfSBarry Smith +  pc - the multigrid context
9898cc2d5dfSBarry Smith -  n - number of cycles (default is 1)
9908cc2d5dfSBarry Smith 
9918cc2d5dfSBarry Smith    Options Database Key:
992e1bc860dSBarry Smith .  -pc_mg_multiplicative_cycles n
9938cc2d5dfSBarry Smith 
9948cc2d5dfSBarry Smith    Level: advanced
9958cc2d5dfSBarry Smith 
9968cc2d5dfSBarry Smith    Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
9978cc2d5dfSBarry Smith 
9988cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
9998cc2d5dfSBarry Smith 
10008cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
10018cc2d5dfSBarry Smith @*/
10027087cfbeSBarry Smith PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
10038cc2d5dfSBarry Smith {
1004f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
10058cc2d5dfSBarry Smith 
10068cc2d5dfSBarry Smith   PetscFunctionBegin;
10070700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1008c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
10092a21e185SBarry Smith   mg->cyclesperpcapply = n;
10108cc2d5dfSBarry Smith   PetscFunctionReturn(0);
10118cc2d5dfSBarry Smith }
10128cc2d5dfSBarry Smith 
10138cc2d5dfSBarry Smith #undef __FUNCT__
1014fb15c04eSBarry Smith #define __FUNCT__ "PCMGSetGalerkin_MG"
1015fb15c04eSBarry Smith PetscErrorCode PCMGSetGalerkin_MG(PC pc,PetscBool use)
1016fb15c04eSBarry Smith {
1017fb15c04eSBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
1018fb15c04eSBarry Smith 
1019fb15c04eSBarry Smith   PetscFunctionBegin;
1020fb15c04eSBarry Smith   mg->galerkin = use ? 1 : 0;
1021fb15c04eSBarry Smith   PetscFunctionReturn(0);
1022fb15c04eSBarry Smith }
1023fb15c04eSBarry Smith 
1024fb15c04eSBarry Smith #undef __FUNCT__
10259dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetGalerkin"
1026c2be2410SBarry Smith /*@
102797177400SBarry Smith    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
102882b87d87SMatthew G. Knepley       finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i
1029c2be2410SBarry Smith 
1030ad4df100SBarry Smith    Logically Collective on PC
1031c2be2410SBarry Smith 
1032c2be2410SBarry Smith    Input Parameters:
1033c91913d3SJed Brown +  pc - the multigrid context
1034c91913d3SJed Brown -  use - PETSC_TRUE to use the Galerkin process to compute coarse-level operators
1035c2be2410SBarry Smith 
1036c2be2410SBarry Smith    Options Database Key:
10371c1aac46SBarry Smith .  -pc_mg_galerkin <true,false>
1038c2be2410SBarry Smith 
1039c2be2410SBarry Smith    Level: intermediate
1040c2be2410SBarry Smith 
10411c1aac46SBarry Smith    Notes: Some codes that use PCMG such as PCGAMG use Galerkin internally while constructing the hierarchy and thus do not
10421c1aac46SBarry Smith      use the PCMG construction of the coarser grids.
10431c1aac46SBarry Smith 
1044c2be2410SBarry Smith .keywords: MG, set, Galerkin
1045c2be2410SBarry Smith 
10463fc8bf9cSBarry Smith .seealso: PCMGGetGalerkin()
10473fc8bf9cSBarry Smith 
1048c2be2410SBarry Smith @*/
1049c91913d3SJed Brown PetscErrorCode PCMGSetGalerkin(PC pc,PetscBool use)
1050c2be2410SBarry Smith {
1051fb15c04eSBarry Smith   PetscErrorCode ierr;
1052c2be2410SBarry Smith 
1053c2be2410SBarry Smith   PetscFunctionBegin;
10540700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10557a4c0024SBarry Smith   ierr = PetscTryMethod(pc,"PCMGSetGalerkin_C",(PC,PetscBool),(pc,use));CHKERRQ(ierr);
1056c2be2410SBarry Smith   PetscFunctionReturn(0);
1057c2be2410SBarry Smith }
1058c2be2410SBarry Smith 
1059c2be2410SBarry Smith #undef __FUNCT__
10603fc8bf9cSBarry Smith #define __FUNCT__ "PCMGGetGalerkin"
10613fc8bf9cSBarry Smith /*@
10623fc8bf9cSBarry Smith    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
106382b87d87SMatthew G. Knepley       A_i-1 = r_i * A_i * p_i
10643fc8bf9cSBarry Smith 
10653fc8bf9cSBarry Smith    Not Collective
10663fc8bf9cSBarry Smith 
10673fc8bf9cSBarry Smith    Input Parameter:
10683fc8bf9cSBarry Smith .  pc - the multigrid context
10693fc8bf9cSBarry Smith 
10703fc8bf9cSBarry Smith    Output Parameter:
1071e1bc860dSBarry Smith .  galerkin - PETSC_TRUE or PETSC_FALSE
10723fc8bf9cSBarry Smith 
10733fc8bf9cSBarry Smith    Options Database Key:
1074e1bc860dSBarry Smith .  -pc_mg_galerkin
10753fc8bf9cSBarry Smith 
10763fc8bf9cSBarry Smith    Level: intermediate
10773fc8bf9cSBarry Smith 
10783fc8bf9cSBarry Smith .keywords: MG, set, Galerkin
10793fc8bf9cSBarry Smith 
10803fc8bf9cSBarry Smith .seealso: PCMGSetGalerkin()
10813fc8bf9cSBarry Smith 
10823fc8bf9cSBarry Smith @*/
10837087cfbeSBarry Smith PetscErrorCode  PCMGGetGalerkin(PC pc,PetscBool  *galerkin)
10843fc8bf9cSBarry Smith {
1085f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
10863fc8bf9cSBarry Smith 
10873fc8bf9cSBarry Smith   PetscFunctionBegin;
10880700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
108913fdb3acSJose Roman   *galerkin = (PetscBool)mg->galerkin;
10903fc8bf9cSBarry Smith   PetscFunctionReturn(0);
10913fc8bf9cSBarry Smith }
10923fc8bf9cSBarry Smith 
10933fc8bf9cSBarry Smith #undef __FUNCT__
10949dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothDown"
10954b9ad928SBarry Smith /*@
109697177400SBarry Smith    PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
109797177400SBarry Smith    use on all levels. Use PCMGGetSmootherDown() to set different
10984b9ad928SBarry Smith    pre-smoothing steps on different levels.
10994b9ad928SBarry Smith 
1100ad4df100SBarry Smith    Logically Collective on PC
11014b9ad928SBarry Smith 
11024b9ad928SBarry Smith    Input Parameters:
11034b9ad928SBarry Smith +  mg - the multigrid context
11044b9ad928SBarry Smith -  n - the number of smoothing steps
11054b9ad928SBarry Smith 
11064b9ad928SBarry Smith    Options Database Key:
11074b9ad928SBarry Smith .  -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps
11084b9ad928SBarry Smith 
11094b9ad928SBarry Smith    Level: advanced
11104b9ad928SBarry Smith 
11114b9ad928SBarry Smith .keywords: MG, smooth, down, pre-smoothing, steps, multigrid
11124b9ad928SBarry Smith 
111397177400SBarry Smith .seealso: PCMGSetNumberSmoothUp()
11144b9ad928SBarry Smith @*/
11157087cfbeSBarry Smith PetscErrorCode  PCMGSetNumberSmoothDown(PC pc,PetscInt n)
11164b9ad928SBarry Smith {
1117f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
1118f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
11196849ba73SBarry Smith   PetscErrorCode ierr;
112079416396SBarry Smith   PetscInt       i,levels;
11214b9ad928SBarry Smith 
11224b9ad928SBarry Smith   PetscFunctionBegin;
11230700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1124ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1125c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
1126f3fbd535SBarry Smith   levels = mglevels[0]->levels;
11274b9ad928SBarry Smith 
1128b05257ddSBarry Smith   for (i=1; i<levels; i++) {
11294b9ad928SBarry Smith     /* make sure smoother up and down are different */
11300298fd71SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,NULL);CHKERRQ(ierr);
1131f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
11322fa5cd67SKarl Rupp 
113331567311SBarry Smith     mg->default_smoothd = n;
11344b9ad928SBarry Smith   }
11354b9ad928SBarry Smith   PetscFunctionReturn(0);
11364b9ad928SBarry Smith }
11374b9ad928SBarry Smith 
11384b9ad928SBarry Smith #undef __FUNCT__
11399dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothUp"
11404b9ad928SBarry Smith /*@
114197177400SBarry Smith    PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
114297177400SBarry Smith    on all levels. Use PCMGGetSmootherUp() to set different numbers of
11434b9ad928SBarry Smith    post-smoothing steps on different levels.
11444b9ad928SBarry Smith 
1145ad4df100SBarry Smith    Logically Collective on PC
11464b9ad928SBarry Smith 
11474b9ad928SBarry Smith    Input Parameters:
11484b9ad928SBarry Smith +  mg - the multigrid context
11494b9ad928SBarry Smith -  n - the number of smoothing steps
11504b9ad928SBarry Smith 
11514b9ad928SBarry Smith    Options Database Key:
11524b9ad928SBarry Smith .  -pc_mg_smoothup <n> - Sets number of post-smoothing steps
11534b9ad928SBarry Smith 
11544b9ad928SBarry Smith    Level: advanced
11554b9ad928SBarry Smith 
11564b9ad928SBarry Smith    Note: this does not set a value on the coarsest grid, since we assume that
1157a8c7a070SBarry Smith     there is no separate smooth up on the coarsest grid.
11584b9ad928SBarry Smith 
11594b9ad928SBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid
11604b9ad928SBarry Smith 
116197177400SBarry Smith .seealso: PCMGSetNumberSmoothDown()
11624b9ad928SBarry Smith @*/
11637087cfbeSBarry Smith PetscErrorCode  PCMGSetNumberSmoothUp(PC pc,PetscInt n)
11644b9ad928SBarry Smith {
1165f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
1166f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
11676849ba73SBarry Smith   PetscErrorCode ierr;
116879416396SBarry Smith   PetscInt       i,levels;
11694b9ad928SBarry Smith 
11704b9ad928SBarry Smith   PetscFunctionBegin;
11710700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1172ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1173c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
1174f3fbd535SBarry Smith   levels = mglevels[0]->levels;
11754b9ad928SBarry Smith 
11764b9ad928SBarry Smith   for (i=1; i<levels; i++) {
11774b9ad928SBarry Smith     /* make sure smoother up and down are different */
11780298fd71SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,NULL);CHKERRQ(ierr);
1179f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
11802fa5cd67SKarl Rupp 
118131567311SBarry Smith     mg->default_smoothu = n;
11824b9ad928SBarry Smith   }
11834b9ad928SBarry Smith   PetscFunctionReturn(0);
11844b9ad928SBarry Smith }
11854b9ad928SBarry Smith 
11864b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/
11874b9ad928SBarry Smith 
11883b09bd56SBarry Smith /*MC
1189ccb205f8SBarry Smith    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
11903b09bd56SBarry Smith     information about the coarser grid matrices and restriction/interpolation operators.
11913b09bd56SBarry Smith 
11923b09bd56SBarry Smith    Options Database Keys:
11933b09bd56SBarry Smith +  -pc_mg_levels <nlevels> - number of levels including finest
1194b4519db0SPatrick Sanan .  -pc_mg_cycles <v,w> -
119579416396SBarry Smith .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
11963b09bd56SBarry Smith .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
11978c1c2452SJed Brown .  -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default
11983b09bd56SBarry Smith .  -pc_mg_log - log information about time spent on each level of the solver
11998cf6037eSBarry Smith .  -pc_mg_galerkin - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R'
12008cf6037eSBarry Smith .  -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1)
12018cf6037eSBarry Smith .  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1202e3c5b3baSBarry Smith                         to the Socket viewer for reading from MATLAB.
12038cf6037eSBarry Smith -  -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices
12048cf6037eSBarry Smith                         to the binary output file called binaryoutput
12053b09bd56SBarry Smith 
120624c3aa18SBarry 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
12073b09bd56SBarry Smith 
12088cf6037eSBarry Smith        When run with a single level the smoother options are used on that level NOT the coarse grid solver options
12098cf6037eSBarry Smith 
12103b09bd56SBarry Smith    Level: intermediate
12113b09bd56SBarry Smith 
12128f87f92bSBarry Smith    Concepts: multigrid/multilevel
12133b09bd56SBarry Smith 
12148cf6037eSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE
12150d353602SBarry Smith            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
121697177400SBarry Smith            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
121797177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
12180d353602SBarry Smith            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
12193b09bd56SBarry Smith M*/
12203b09bd56SBarry Smith 
12214b9ad928SBarry Smith #undef __FUNCT__
12224b9ad928SBarry Smith #define __FUNCT__ "PCCreate_MG"
12238cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc)
12244b9ad928SBarry Smith {
1225f3fbd535SBarry Smith   PC_MG          *mg;
1226f3fbd535SBarry Smith   PetscErrorCode ierr;
1227f3fbd535SBarry Smith 
12284b9ad928SBarry Smith   PetscFunctionBegin;
1229b00a9115SJed Brown   ierr        = PetscNewLog(pc,&mg);CHKERRQ(ierr);
1230f3fbd535SBarry Smith   pc->data    = (void*)mg;
1231f3fbd535SBarry Smith   mg->nlevels = -1;
123210eca3edSLisandro Dalcin   mg->am      = PC_MG_MULTIPLICATIVE;
1233f3fbd535SBarry Smith 
123437a44384SMark Adams   pc->useAmat = PETSC_TRUE;
123537a44384SMark Adams 
12364b9ad928SBarry Smith   pc->ops->apply          = PCApply_MG;
12374b9ad928SBarry Smith   pc->ops->setup          = PCSetUp_MG;
1238a06653b4SBarry Smith   pc->ops->reset          = PCReset_MG;
12394b9ad928SBarry Smith   pc->ops->destroy        = PCDestroy_MG;
12404b9ad928SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_MG;
12414b9ad928SBarry Smith   pc->ops->view           = PCView_MG;
1242fb15c04eSBarry Smith 
1243fb15c04eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",PCMGSetGalerkin_MG);CHKERRQ(ierr);
12444b9ad928SBarry Smith   PetscFunctionReturn(0);
12454b9ad928SBarry Smith }
1246