xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision 7319c6545cdf078a8212d19c3cdefe116c1a5fdd)
1dba47a55SKris Buschelman #define PETSCKSP_DLL
2dba47a55SKris Buschelman 
34b9ad928SBarry Smith /*
44b9ad928SBarry Smith     Defines the multigrid preconditioner interface.
54b9ad928SBarry Smith */
67c4f633dSBarry Smith #include "../src/ksp/pc/impls/mg/mgimpl.h"                    /*I "petscmg.h" I*/
74b9ad928SBarry Smith 
84b9ad928SBarry Smith 
94b9ad928SBarry Smith #undef __FUNCT__
109dcbbd2bSBarry Smith #define __FUNCT__ "PCMGMCycle_Private"
114d0a8057SBarry Smith PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG **mglevels,PCRichardsonConvergedReason *reason)
124b9ad928SBarry Smith {
139dcbbd2bSBarry Smith   PC_MG          *mg = *mglevels,*mgc;
146849ba73SBarry Smith   PetscErrorCode ierr;
15e0ef7d32SBarry Smith   PetscInt       cycles = (mg->level == 1) ? 1 : (PetscInt) mg->cycles;
164b9ad928SBarry Smith 
174b9ad928SBarry Smith   PetscFunctionBegin;
184b9ad928SBarry Smith 
1932cf1786SBarry Smith   if (mg->eventsmoothsolve) {ierr = PetscLogEventBegin(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
200d353602SBarry Smith   ierr = KSPSolve(mg->smoothd,mg->b,mg->x);CHKERRQ(ierr);  /* pre-smooth */
2132cf1786SBarry Smith   if (mg->eventsmoothsolve) {ierr = PetscLogEventEnd(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
224b9ad928SBarry Smith   if (mg->level) {  /* not the coarsest grid */
2332cf1786SBarry Smith     if (mg->eventresidual) {ierr = PetscLogEventBegin(mg->eventresidual,0,0,0,0);CHKERRQ(ierr);}
244b9ad928SBarry Smith     ierr = (*mg->residual)(mg->A,mg->b,mg->x,mg->r);CHKERRQ(ierr);
2532cf1786SBarry Smith     if (mg->eventresidual) {ierr = PetscLogEventEnd(mg->eventresidual,0,0,0,0);CHKERRQ(ierr);}
264b9ad928SBarry Smith 
274b9ad928SBarry Smith     /* if on finest level and have convergence criteria set */
284d0a8057SBarry Smith     if (mg->level == mg->levels-1 && mg->ttol && reason) {
294b9ad928SBarry Smith       PetscReal rnorm;
304b9ad928SBarry Smith       ierr = VecNorm(mg->r,NORM_2,&rnorm);CHKERRQ(ierr);
314b9ad928SBarry Smith       if (rnorm <= mg->ttol) {
3270441072SBarry Smith         if (rnorm < mg->abstol) {
334d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_ATOL;
341e2582c4SBarry Smith           ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %G is less than absolute tolerance %G\n",rnorm,mg->abstol);CHKERRQ(ierr);
354b9ad928SBarry Smith         } else {
364d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_RTOL;
371e2582c4SBarry Smith           ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %G is less than relative tolerance times initial residual norm %G\n",rnorm,mg->ttol);CHKERRQ(ierr);
384b9ad928SBarry Smith         }
394b9ad928SBarry Smith         PetscFunctionReturn(0);
404b9ad928SBarry Smith       }
414b9ad928SBarry Smith     }
424b9ad928SBarry Smith 
434b9ad928SBarry Smith     mgc = *(mglevels - 1);
4432cf1786SBarry Smith     if (mg->eventinterprestrict) {ierr = PetscLogEventBegin(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
454b9ad928SBarry Smith     ierr = MatRestrict(mg->restrct,mg->r,mgc->b);CHKERRQ(ierr);
4632cf1786SBarry Smith     if (mg->eventinterprestrict) {ierr = PetscLogEventEnd(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
47efb30889SBarry Smith     ierr = VecSet(mgc->x,0.0);CHKERRQ(ierr);
484b9ad928SBarry Smith     while (cycles--) {
494d0a8057SBarry Smith       ierr = PCMGMCycle_Private(pc,mglevels-1,reason);CHKERRQ(ierr);
504b9ad928SBarry Smith     }
5132cf1786SBarry Smith     if (mg->eventinterprestrict) {ierr = PetscLogEventBegin(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
524b9ad928SBarry Smith     ierr = MatInterpolateAdd(mg->interpolate,mgc->x,mg->x,mg->x);CHKERRQ(ierr);
5332cf1786SBarry Smith     if (mg->eventinterprestrict) {ierr = PetscLogEventEnd(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
5432cf1786SBarry Smith     if (mg->eventsmoothsolve) {ierr = PetscLogEventBegin(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
550d353602SBarry Smith     ierr = KSPSolve(mg->smoothu,mg->b,mg->x);CHKERRQ(ierr);    /* post smooth */
5632cf1786SBarry Smith     if (mg->eventsmoothsolve) {ierr = PetscLogEventEnd(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
574b9ad928SBarry Smith   }
584b9ad928SBarry Smith   PetscFunctionReturn(0);
594b9ad928SBarry Smith }
604b9ad928SBarry Smith 
614b9ad928SBarry Smith /*
629dcbbd2bSBarry Smith        PCMGCreate_Private - Creates a PC_MG structure for use with the
634b9ad928SBarry Smith                multigrid code. Level 0 is the coarsest. (But the
644b9ad928SBarry Smith                finest level is stored first in the array).
654b9ad928SBarry Smith 
664b9ad928SBarry Smith */
674b9ad928SBarry Smith #undef __FUNCT__
689dcbbd2bSBarry Smith #define __FUNCT__ "PCMGCreate_Private"
699dcbbd2bSBarry Smith static PetscErrorCode PCMGCreate_Private(MPI_Comm comm,PetscInt levels,PC pc,MPI_Comm *comms,PC_MG ***result)
704b9ad928SBarry Smith {
719dcbbd2bSBarry Smith   PC_MG          **mg;
726849ba73SBarry Smith   PetscErrorCode ierr;
7379416396SBarry Smith   PetscInt       i;
7479416396SBarry Smith   PetscMPIInt    size;
75f69a0ea3SMatthew Knepley   const char     *prefix;
764b9ad928SBarry Smith   PC             ipc;
774b9ad928SBarry Smith 
784b9ad928SBarry Smith   PetscFunctionBegin;
799dcbbd2bSBarry Smith   ierr = PetscMalloc(levels*sizeof(PC_MG*),&mg);CHKERRQ(ierr);
8038f2d2fdSLisandro Dalcin   ierr = PetscLogObjectMemory(pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr);
814b9ad928SBarry Smith 
824b9ad928SBarry Smith   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
834b9ad928SBarry Smith 
844b9ad928SBarry Smith   for (i=0; i<levels; i++) {
8538f2d2fdSLisandro Dalcin     ierr = PetscNewLog(pc,PC_MG,&mg[i]);CHKERRQ(ierr);
864b9ad928SBarry Smith     mg[i]->level           = i;
874b9ad928SBarry Smith     mg[i]->levels          = levels;
880d353602SBarry Smith     mg[i]->cycles          = PC_MG_CYCLE_V;
89c2be2410SBarry Smith     mg[i]->galerkin        = PETSC_FALSE;
90c2be2410SBarry Smith     mg[i]->galerkinused    = PETSC_FALSE;
9179416396SBarry Smith     mg[i]->default_smoothu = 1;
9279416396SBarry Smith     mg[i]->default_smoothd = 1;
934b9ad928SBarry Smith 
944b9ad928SBarry Smith     if (comms) comm = comms[i];
954b9ad928SBarry Smith     ierr = KSPCreate(comm,&mg[i]->smoothd);CHKERRQ(ierr);
961cee3971SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)mg[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr);
9779416396SBarry Smith     ierr = KSPSetTolerances(mg[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg[i]->default_smoothd);CHKERRQ(ierr);
984b9ad928SBarry Smith     ierr = KSPSetOptionsPrefix(mg[i]->smoothd,prefix);CHKERRQ(ierr);
994b9ad928SBarry Smith 
1004b9ad928SBarry Smith     /* do special stuff for coarse grid */
1014b9ad928SBarry Smith     if (!i && levels > 1) {
1024b9ad928SBarry Smith       ierr = KSPAppendOptionsPrefix(mg[0]->smoothd,"mg_coarse_");CHKERRQ(ierr);
1034b9ad928SBarry Smith 
1044b9ad928SBarry Smith       /* coarse solve is (redundant) LU by default */
1054b9ad928SBarry Smith       ierr = KSPSetType(mg[0]->smoothd,KSPPREONLY);CHKERRQ(ierr);
1064b9ad928SBarry Smith       ierr = KSPGetPC(mg[0]->smoothd,&ipc);CHKERRQ(ierr);
1074b9ad928SBarry Smith       ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1084b9ad928SBarry Smith       if (size > 1) {
1094b9ad928SBarry Smith         ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr);
1100d7810c8SBarry Smith       } else {
1114b9ad928SBarry Smith         ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr);
1120d7810c8SBarry Smith       }
1134b9ad928SBarry Smith 
1144b9ad928SBarry Smith     } else {
1152e70eadcSBarry Smith       char tprefix[128];
11613f74950SBarry Smith       sprintf(tprefix,"mg_levels_%d_",(int)i);
1172e70eadcSBarry Smith       ierr = KSPAppendOptionsPrefix(mg[i]->smoothd,tprefix);CHKERRQ(ierr);
1184b9ad928SBarry Smith     }
11952e6d16bSBarry Smith     ierr = PetscLogObjectParent(pc,mg[i]->smoothd);CHKERRQ(ierr);
1204b9ad928SBarry Smith     mg[i]->smoothu             = mg[i]->smoothd;
1214b9ad928SBarry Smith     mg[i]->rtol                = 0.0;
12270441072SBarry Smith     mg[i]->abstol              = 0.0;
1234b9ad928SBarry Smith     mg[i]->dtol                = 0.0;
1244b9ad928SBarry Smith     mg[i]->ttol                = 0.0;
12532cf1786SBarry Smith     mg[i]->eventsmoothsetup    = 0;
12632cf1786SBarry Smith     mg[i]->eventsmoothsolve    = 0;
12732cf1786SBarry Smith     mg[i]->eventresidual       = 0;
12832cf1786SBarry Smith     mg[i]->eventinterprestrict = 0;
1298cc2d5dfSBarry Smith     mg[i]->cyclesperpcapply    = 1;
1304b9ad928SBarry Smith   }
1314b9ad928SBarry Smith   *result = mg;
1324b9ad928SBarry Smith   PetscFunctionReturn(0);
1334b9ad928SBarry Smith }
1344b9ad928SBarry Smith 
1354b9ad928SBarry Smith #undef __FUNCT__
1364b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_MG"
1377233f9f0SBarry Smith PetscErrorCode PCDestroy_MG(PC pc)
1384b9ad928SBarry Smith {
1399dcbbd2bSBarry Smith   PC_MG          **mg = (PC_MG**)pc->data;
1406849ba73SBarry Smith   PetscErrorCode ierr;
14138f2d2fdSLisandro Dalcin   PetscInt       i,n;
1424b9ad928SBarry Smith 
1434b9ad928SBarry Smith   PetscFunctionBegin;
14438f2d2fdSLisandro Dalcin   if (!mg) PetscFunctionReturn(0);
14514ca34e6SBarry Smith 
14638f2d2fdSLisandro Dalcin   n = mg[0]->levels;
147fccaa45eSBarry Smith   for (i=0; i<n-1; i++) {
148630ba2eeSBarry Smith     if (mg[i+1]->r) {ierr = VecDestroy(mg[i+1]->r);CHKERRQ(ierr);}
149fccaa45eSBarry Smith     if (mg[i]->b) {ierr = VecDestroy(mg[i]->b);CHKERRQ(ierr);}
150fccaa45eSBarry Smith     if (mg[i]->x) {ierr = VecDestroy(mg[i]->x);CHKERRQ(ierr);}
151fccaa45eSBarry Smith     if (mg[i+1]->restrct) {ierr = MatDestroy(mg[i+1]->restrct);CHKERRQ(ierr);}
152fccaa45eSBarry Smith     if (mg[i+1]->interpolate) {ierr = MatDestroy(mg[i+1]->interpolate);CHKERRQ(ierr);}
153fccaa45eSBarry Smith   }
154fccaa45eSBarry Smith 
1554b9ad928SBarry Smith   for (i=0; i<n; i++) {
1564b9ad928SBarry Smith     if (mg[i]->smoothd != mg[i]->smoothu) {
1574b9ad928SBarry Smith       ierr = KSPDestroy(mg[i]->smoothd);CHKERRQ(ierr);
1584b9ad928SBarry Smith     }
1594b9ad928SBarry Smith     ierr = KSPDestroy(mg[i]->smoothu);CHKERRQ(ierr);
1604b9ad928SBarry Smith     ierr = PetscFree(mg[i]);CHKERRQ(ierr);
1614b9ad928SBarry Smith   }
1624b9ad928SBarry Smith   ierr = PetscFree(mg);CHKERRQ(ierr);
1634b9ad928SBarry Smith   PetscFunctionReturn(0);
1644b9ad928SBarry Smith }
1654b9ad928SBarry Smith 
1664b9ad928SBarry Smith 
1674b9ad928SBarry Smith 
1689dcbbd2bSBarry Smith EXTERN PetscErrorCode PCMGACycle_Private(PC_MG**);
1691e2582c4SBarry Smith EXTERN PetscErrorCode PCMGFCycle_Private(PC,PC_MG**);
1709dcbbd2bSBarry Smith EXTERN PetscErrorCode PCMGKCycle_Private(PC_MG**);
1714b9ad928SBarry Smith 
1724b9ad928SBarry Smith /*
1734b9ad928SBarry Smith    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
1744b9ad928SBarry Smith              or full cycle of multigrid.
1754b9ad928SBarry Smith 
1764b9ad928SBarry Smith   Note:
1779dcbbd2bSBarry Smith   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
1784b9ad928SBarry Smith */
1794b9ad928SBarry Smith #undef __FUNCT__
1804b9ad928SBarry Smith #define __FUNCT__ "PCApply_MG"
1816849ba73SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
1824b9ad928SBarry Smith {
1839dcbbd2bSBarry Smith   PC_MG          **mg = (PC_MG**)pc->data;
1846849ba73SBarry Smith   PetscErrorCode ierr;
1858cc2d5dfSBarry Smith   PetscInt       levels = mg[0]->levels,i;
1864b9ad928SBarry Smith 
1874b9ad928SBarry Smith   PetscFunctionBegin;
1884b9ad928SBarry Smith   mg[levels-1]->b = b;
1894b9ad928SBarry Smith   mg[levels-1]->x = x;
1909dcbbd2bSBarry Smith   if (mg[0]->am == PC_MG_MULTIPLICATIVE) {
191efb30889SBarry Smith     ierr = VecSet(x,0.0);CHKERRQ(ierr);
1928cc2d5dfSBarry Smith     for (i=0; i<mg[0]->cyclesperpcapply; i++) {
1931e2582c4SBarry Smith       ierr = PCMGMCycle_Private(pc,mg+levels-1,PETSC_NULL);CHKERRQ(ierr);
1944b9ad928SBarry Smith     }
1958cc2d5dfSBarry Smith   }
1969dcbbd2bSBarry Smith   else if (mg[0]->am == PC_MG_ADDITIVE) {
1979dcbbd2bSBarry Smith     ierr = PCMGACycle_Private(mg);CHKERRQ(ierr);
1984b9ad928SBarry Smith   }
1999dcbbd2bSBarry Smith   else if (mg[0]->am == PC_MG_KASKADE) {
2009dcbbd2bSBarry Smith     ierr = PCMGKCycle_Private(mg);CHKERRQ(ierr);
2014b9ad928SBarry Smith   }
2024b9ad928SBarry Smith   else {
2031e2582c4SBarry Smith     ierr = PCMGFCycle_Private(pc,mg);CHKERRQ(ierr);
2044b9ad928SBarry Smith   }
2054b9ad928SBarry Smith   PetscFunctionReturn(0);
2064b9ad928SBarry Smith }
2074b9ad928SBarry Smith 
2084b9ad928SBarry Smith #undef __FUNCT__
2094b9ad928SBarry Smith #define __FUNCT__ "PCApplyRichardson_MG"
210*7319c654SBarry Smith static PetscErrorCode PCApplyRichardson_MG(PC pc,Vec b,Vec x,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscTruth zeroguess,PetscInt *outits,PCRichardsonConvergedReason *reason)
2114b9ad928SBarry Smith {
2129dcbbd2bSBarry Smith   PC_MG          **mg = (PC_MG**)pc->data;
213dfbe8321SBarry Smith   PetscErrorCode ierr;
2144d0a8057SBarry Smith   PetscInt       levels = mg[0]->levels,i;
2154b9ad928SBarry Smith 
2164b9ad928SBarry Smith   PetscFunctionBegin;
2174b9ad928SBarry Smith   mg[levels-1]->b    = b;
2184b9ad928SBarry Smith   mg[levels-1]->x    = x;
2194b9ad928SBarry Smith 
2204b9ad928SBarry Smith   mg[levels-1]->rtol = rtol;
22170441072SBarry Smith   mg[levels-1]->abstol = abstol;
2224b9ad928SBarry Smith   mg[levels-1]->dtol = dtol;
2234b9ad928SBarry Smith   if (rtol) {
2244b9ad928SBarry Smith     /* compute initial residual norm for relative convergence test */
2254b9ad928SBarry Smith     PetscReal rnorm;
226*7319c654SBarry Smith     if (zeroguess) {
227*7319c654SBarry Smith       ierr               = VecNorm(b,NORM_2,&rnorm);CHKERRQ(ierr);
228*7319c654SBarry Smith     } else {
2294b9ad928SBarry Smith       ierr               = (*mg[levels-1]->residual)(mg[levels-1]->A,b,x,w);CHKERRQ(ierr);
2304b9ad928SBarry Smith       ierr               = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr);
231*7319c654SBarry Smith     }
23270441072SBarry Smith     mg[levels-1]->ttol = PetscMax(rtol*rnorm,abstol);
23370441072SBarry Smith   } else if (abstol) {
23470441072SBarry Smith     mg[levels-1]->ttol = abstol;
2354b9ad928SBarry Smith   } else {
2364b9ad928SBarry Smith     mg[levels-1]->ttol = 0.0;
2374b9ad928SBarry Smith   }
2384b9ad928SBarry Smith 
2394d0a8057SBarry Smith   /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't
2404d0a8057SBarry Smith      stop prematurely do to small residual */
2414d0a8057SBarry Smith   for (i=1; i<levels; i++) {
2424d0a8057SBarry Smith     ierr = KSPSetTolerances(mg[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
2434d0a8057SBarry Smith     if (mg[i]->smoothu != mg[i]->smoothd) {
2444d0a8057SBarry Smith       ierr = KSPSetTolerances(mg[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
2454b9ad928SBarry Smith     }
2464d0a8057SBarry Smith   }
2474d0a8057SBarry Smith 
2484d0a8057SBarry Smith   *reason = (PCRichardsonConvergedReason)0;
2494d0a8057SBarry Smith   for (i=0; i<its; i++) {
2504d0a8057SBarry Smith     ierr = PCMGMCycle_Private(pc,mg+levels-1,reason);CHKERRQ(ierr);
2514d0a8057SBarry Smith     if (*reason) break;
2524d0a8057SBarry Smith   }
2534d0a8057SBarry Smith   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
2544d0a8057SBarry Smith   *outits = i;
2554b9ad928SBarry Smith   PetscFunctionReturn(0);
2564b9ad928SBarry Smith }
2574b9ad928SBarry Smith 
2584b9ad928SBarry Smith #undef __FUNCT__
2594b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_MG"
2606ca4d86aSHong Zhang PetscErrorCode PCSetFromOptions_MG(PC pc)
2614b9ad928SBarry Smith {
262dfbe8321SBarry Smith   PetscErrorCode ierr;
2638cc2d5dfSBarry Smith   PetscInt       m,levels = 1,cycles;
2644b9ad928SBarry Smith   PetscTruth     flg;
2659dcbbd2bSBarry Smith   PC_MG          **mg = (PC_MG**)pc->data;
26690da95b0SMatthew Knepley   PCMGType       mgtype = PC_MG_ADDITIVE;
2671aa34eeaSBarry Smith   PCMGCycleType  mgctype;
2684b9ad928SBarry Smith 
2694b9ad928SBarry Smith   PetscFunctionBegin;
2704b9ad928SBarry Smith   ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr);
2714b9ad928SBarry Smith     if (!pc->data) {
2729dcbbd2bSBarry Smith       ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr);
27397177400SBarry Smith       ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr);
274cf502942SBarry Smith       mg = (PC_MG**)pc->data;
2754b9ad928SBarry Smith     }
276f2070a76SMatthew Knepley     mgctype = (PCMGCycleType) mg[0]->cycles;
2770d353602SBarry Smith     ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr);
2784b9ad928SBarry Smith     if (flg) {
2790d353602SBarry Smith       ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr);
2800d353602SBarry Smith     };
28190d69ab7SBarry Smith     flg  = PETSC_FALSE;
28290d69ab7SBarry Smith     ierr = PetscOptionsTruth("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
283c2be2410SBarry Smith     if (flg) {
28497177400SBarry Smith       ierr = PCMGSetGalerkin(pc);CHKERRQ(ierr);
285c2be2410SBarry Smith     }
2869dcbbd2bSBarry Smith     ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr);
2874b9ad928SBarry Smith     if (flg) {
28897177400SBarry Smith       ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr);
2894b9ad928SBarry Smith     }
2909dcbbd2bSBarry Smith     ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr);
2914b9ad928SBarry Smith     if (flg) {
29297177400SBarry Smith       ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr);
2934b9ad928SBarry Smith     }
2949dcbbd2bSBarry Smith     ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
2958cc2d5dfSBarry Smith     if (flg) {
2968cc2d5dfSBarry Smith       ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
2978cc2d5dfSBarry Smith     }
2988cc2d5dfSBarry Smith     if (mg[0]->am == PC_MG_MULTIPLICATIVE) {
2998cc2d5dfSBarry Smith       ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGSetLevels",mg[0]->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr);
3008cc2d5dfSBarry Smith       if (flg) {
3018cc2d5dfSBarry Smith 	ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr);
3028cc2d5dfSBarry Smith       }
3038cc2d5dfSBarry Smith     }
30490d69ab7SBarry Smith     flg  = PETSC_FALSE;
30590d69ab7SBarry Smith     ierr = PetscOptionsTruth("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
3064b9ad928SBarry Smith     if (flg) {
3074f5ab15aSBarry Smith       PetscInt i;
3084b9ad928SBarry Smith       char     eventname[128];
3094b9ad928SBarry Smith       if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
3104b9ad928SBarry Smith       levels = mg[0]->levels;
3114b9ad928SBarry Smith       for (i=0; i<levels; i++) {
31232cf1786SBarry Smith         sprintf(eventname,"MGSetup Level %d",(int)i);
3138cbcd9ccSBarry Smith         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->cookie,&mg[i]->eventsmoothsetup);CHKERRQ(ierr);
31432cf1786SBarry Smith         sprintf(eventname,"MGSmooth Level %d",(int)i);
3158cbcd9ccSBarry Smith         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->cookie,&mg[i]->eventsmoothsolve);CHKERRQ(ierr);
31632cf1786SBarry Smith         if (i) {
31732cf1786SBarry Smith           sprintf(eventname,"MGResid Level %d",(int)i);
318a3bc4eb9SBarry Smith           ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->cookie,&mg[i]->eventresidual);CHKERRQ(ierr);
31932cf1786SBarry Smith           sprintf(eventname,"MGInterp Level %d",(int)i);
3208cbcd9ccSBarry Smith           ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->cookie,&mg[i]->eventinterprestrict);CHKERRQ(ierr);
32132cf1786SBarry Smith         }
3224b9ad928SBarry Smith       }
3234b9ad928SBarry Smith     }
3244b9ad928SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
3254b9ad928SBarry Smith   PetscFunctionReturn(0);
3264b9ad928SBarry Smith }
3274b9ad928SBarry Smith 
3289dcbbd2bSBarry Smith const char *PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0};
3290d353602SBarry Smith const char *PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0};
3309dcbbd2bSBarry Smith 
3314b9ad928SBarry Smith #undef __FUNCT__
3324b9ad928SBarry Smith #define __FUNCT__ "PCView_MG"
3337233f9f0SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
3344b9ad928SBarry Smith {
3359dcbbd2bSBarry Smith   PC_MG          **mg = (PC_MG**)pc->data;
336dfbe8321SBarry Smith   PetscErrorCode ierr;
33779416396SBarry Smith   PetscInt       levels = mg[0]->levels,i;
33832077d6dSBarry Smith   PetscTruth     iascii;
3394b9ad928SBarry Smith 
3404b9ad928SBarry Smith   PetscFunctionBegin;
34132077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
34232077d6dSBarry Smith   if (iascii) {
343a8e6722dSMatthew Knepley     ierr = PetscViewerASCIIPrintf(viewer,"  MG: type is %s, levels=%D cycles=%s\n", PCMGTypes[mg[0]->am],levels,(mg[0]->cycles == PC_MG_CYCLE_V) ? "v" : "w");CHKERRQ(ierr);
344d792f543SMatthew Knepley     if (mg[0]->am == PC_MG_MULTIPLICATIVE) {
345d792f543SMatthew Knepley       ierr = PetscViewerASCIIPrintf(viewer,"    Cycles per PCApply=%d\n",mg[0]->cyclesperpcapply);CHKERRQ(ierr);
346d792f543SMatthew Knepley     }
347c2be2410SBarry Smith     if (mg[0]->galerkin) {
348c2be2410SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
349c2be2410SBarry Smith     }
3504b9ad928SBarry Smith     for (i=0; i<levels; i++) {
351b03c7568SBarry Smith       if (!i) {
3525b3b9d15SJed Brown         ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %D presmooths=%D postsmooths=%D -----\n",i,mg[0]->default_smoothd,mg[0]->default_smoothu);CHKERRQ(ierr);
353b03c7568SBarry Smith       } else {
354a8e6722dSMatthew Knepley         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D smooths=%D --------------------\n",i,mg[i]->default_smoothd);CHKERRQ(ierr);
355b03c7568SBarry Smith       }
3564b9ad928SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
3574b9ad928SBarry Smith       ierr = KSPView(mg[i]->smoothd,viewer);CHKERRQ(ierr);
3584b9ad928SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
359b03c7568SBarry Smith       if (i && mg[i]->smoothd == mg[i]->smoothu) {
3604b9ad928SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
361b03c7568SBarry Smith       } else if (i){
362a8e6722dSMatthew Knepley         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D smooths=%D --------------------\n",i,mg[i]->default_smoothu);CHKERRQ(ierr);
3634b9ad928SBarry Smith         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
3644b9ad928SBarry Smith         ierr = KSPView(mg[i]->smoothu,viewer);CHKERRQ(ierr);
3654b9ad928SBarry Smith         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
3664b9ad928SBarry Smith       }
3674b9ad928SBarry Smith     }
3684b9ad928SBarry Smith   } else {
36979a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCMG",((PetscObject)viewer)->type_name);
3704b9ad928SBarry Smith   }
3714b9ad928SBarry Smith   PetscFunctionReturn(0);
3724b9ad928SBarry Smith }
3734b9ad928SBarry Smith 
3744b9ad928SBarry Smith /*
3754b9ad928SBarry Smith     Calls setup for the KSP on each level
3764b9ad928SBarry Smith */
3774b9ad928SBarry Smith #undef __FUNCT__
3784b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_MG"
3797233f9f0SBarry Smith PetscErrorCode PCSetUp_MG(PC pc)
3804b9ad928SBarry Smith {
3819dcbbd2bSBarry Smith   PC_MG                   **mg = (PC_MG**)pc->data;
382dfbe8321SBarry Smith   PetscErrorCode          ierr;
38379416396SBarry Smith   PetscInt                i,n = mg[0]->levels;
38477122347SBarry Smith   PC                      cpc,mpc;
38590d69ab7SBarry Smith   PetscTruth              preonly,lu,redundant,cholesky,monitor = PETSC_FALSE,dump = PETSC_FALSE,opsset;
38623d894e5SBarry Smith   PetscViewerASCIIMonitor ascii;
38723d894e5SBarry Smith   PetscViewer             viewer = PETSC_NULL;
3884b9ad928SBarry Smith   MPI_Comm                comm;
389c2be2410SBarry Smith   Mat                     dA,dB;
390c2be2410SBarry Smith   MatStructure            uflag;
3910a6bb862SBarry Smith   Vec                     tvec;
3924b9ad928SBarry Smith 
3934b9ad928SBarry Smith   PetscFunctionBegin;
394852665d3SBarry Smith 
39543fb2f97SBarry Smith   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
39643fb2f97SBarry Smith   /* so use those from global PC */
39743fb2f97SBarry Smith   /* Is this what we always want? What if user wants to keep old one? */
398906ed7ccSBarry Smith   ierr = KSPGetOperatorsSet(mg[n-1]->smoothd,PETSC_NULL,&opsset);CHKERRQ(ierr);
39943fb2f97SBarry Smith   ierr = KSPGetPC(mg[0]->smoothd,&cpc);CHKERRQ(ierr);
40077122347SBarry Smith   ierr = KSPGetPC(mg[n-1]->smoothd,&mpc);CHKERRQ(ierr);
40177122347SBarry Smith   if (!opsset || ((cpc->setupcalled == 1) && (mpc->setupcalled == 2))) {
402852665d3SBarry 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);
4031cfe3bddSBarry Smith     ierr = KSPSetOperators(mg[n-1]->smoothd,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr);
404852665d3SBarry Smith   }
405852665d3SBarry Smith 
406852665d3SBarry Smith   if (mg[0]->galerkin) {
407852665d3SBarry Smith     Mat B;
408852665d3SBarry Smith     mg[0]->galerkinused = PETSC_TRUE;
409852665d3SBarry Smith     /* currently only handle case where mat and pmat are the same on coarser levels */
410852665d3SBarry Smith     ierr = KSPGetOperators(mg[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr);
411852665d3SBarry Smith     if (!pc->setupcalled) {
412852665d3SBarry Smith       for (i=n-2; i>-1; i--) {
413852665d3SBarry Smith         ierr = MatPtAP(dB,mg[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
414852665d3SBarry Smith         ierr = KSPSetOperators(mg[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
415906ed7ccSBarry Smith 	if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
416852665d3SBarry Smith         dB   = B;
417852665d3SBarry Smith       }
418906ed7ccSBarry Smith       ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);
419852665d3SBarry Smith     } else {
420852665d3SBarry Smith       for (i=n-2; i>-1; i--) {
421906ed7ccSBarry Smith         ierr = KSPGetOperators(mg[i]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
422852665d3SBarry Smith         ierr = MatPtAP(dB,mg[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
423852665d3SBarry Smith         ierr = KSPSetOperators(mg[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
424852665d3SBarry Smith         dB   = B;
425852665d3SBarry Smith       }
426852665d3SBarry Smith     }
427852665d3SBarry Smith   }
428852665d3SBarry Smith 
429958c9bccSBarry Smith   if (!pc->setupcalled) {
43090d69ab7SBarry Smith     ierr = PetscOptionsGetTruth(0,"-pc_mg_monitor",&monitor,PETSC_NULL);CHKERRQ(ierr);
4314b9ad928SBarry Smith 
432b03c7568SBarry Smith     for (i=0; i<n; i++) {
4334b9ad928SBarry Smith       if (monitor) {
4344b9ad928SBarry Smith         ierr = PetscObjectGetComm((PetscObject)mg[i]->smoothd,&comm);CHKERRQ(ierr);
43523d894e5SBarry Smith         ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr);
43623d894e5SBarry Smith         ierr = KSPMonitorSet(mg[i]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
4374b9ad928SBarry Smith       }
4384b9ad928SBarry Smith       ierr = KSPSetFromOptions(mg[i]->smoothd);CHKERRQ(ierr);
4394b9ad928SBarry Smith     }
440b03c7568SBarry Smith     for (i=1; i<n; i++) {
441a98fc643SBarry Smith       if (mg[i]->smoothu && (mg[i]->smoothu != mg[i]->smoothd)) {
4424b9ad928SBarry Smith         if (monitor) {
4434b9ad928SBarry Smith           ierr = PetscObjectGetComm((PetscObject)mg[i]->smoothu,&comm);CHKERRQ(ierr);
44423d894e5SBarry Smith           ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr);
44523d894e5SBarry Smith           ierr = KSPMonitorSet(mg[i]->smoothu,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
4464b9ad928SBarry Smith         }
4474b9ad928SBarry Smith         ierr = KSPSetFromOptions(mg[i]->smoothu);CHKERRQ(ierr);
4484b9ad928SBarry Smith       }
4494b9ad928SBarry Smith     }
450fccaa45eSBarry Smith     for (i=1; i<n; i++) {
4510cace4b0SBarry Smith       if (!mg[i]->residual) {
4520cace4b0SBarry Smith         Mat mat;
4530cace4b0SBarry Smith         ierr = KSPGetOperators(mg[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr);
4540cace4b0SBarry Smith         ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr);
4550cace4b0SBarry Smith       }
456fccaa45eSBarry Smith       if (mg[i]->restrct && !mg[i]->interpolate) {
457aea2a34eSBarry Smith         ierr = PCMGSetInterpolation(pc,i,mg[i]->restrct);CHKERRQ(ierr);
458fccaa45eSBarry Smith       }
459fccaa45eSBarry Smith       if (!mg[i]->restrct && mg[i]->interpolate) {
46097177400SBarry Smith         ierr = PCMGSetRestriction(pc,i,mg[i]->interpolate);CHKERRQ(ierr);
461fccaa45eSBarry Smith       }
462fccaa45eSBarry Smith #if defined(PETSC_USE_DEBUG)
463fccaa45eSBarry Smith       if (!mg[i]->restrct || !mg[i]->interpolate) {
464fccaa45eSBarry Smith         SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Need to set restriction or interpolation on level %d",(int)i);
465fccaa45eSBarry Smith       }
466fccaa45eSBarry Smith #endif
467fccaa45eSBarry Smith     }
4680a6bb862SBarry Smith     for (i=0; i<n-1; i++) {
46937b0e6c0SBarry Smith       if (!mg[i]->b) {
470906ed7ccSBarry Smith         Vec *vec;
471906ed7ccSBarry Smith         ierr = KSPGetVecs(mg[i]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
472906ed7ccSBarry Smith         ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr);
4734fdc791cSSatish Balay         ierr = VecDestroy(*vec);CHKERRQ(ierr);
474906ed7ccSBarry Smith         ierr = PetscFree(vec);CHKERRQ(ierr);
47537b0e6c0SBarry Smith       }
4766ca4d86aSHong Zhang       if (!mg[i]->r && i) {
4770a6bb862SBarry Smith         ierr = VecDuplicate(mg[i]->b,&tvec);CHKERRQ(ierr);
47897177400SBarry Smith         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
4790a6bb862SBarry Smith         ierr = VecDestroy(tvec);CHKERRQ(ierr);
4800a6bb862SBarry Smith       }
4810a6bb862SBarry Smith       if (!mg[i]->x) {
4820a6bb862SBarry Smith         ierr = VecDuplicate(mg[i]->b,&tvec);CHKERRQ(ierr);
48397177400SBarry Smith         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
4840a6bb862SBarry Smith         ierr = VecDestroy(tvec);CHKERRQ(ierr);
4850a6bb862SBarry Smith       }
4860a6bb862SBarry Smith     }
487dfef70bdSHong Zhang     if (n != 1 && !mg[n-1]->r) {
488dfef70bdSHong Zhang       /* PCMGSetR() on the finest level if user did not supply it */
4890f77fa87SBarry Smith       Vec *vec;
4900f77fa87SBarry Smith       ierr = KSPGetVecs(mg[n-1]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
4910f77fa87SBarry Smith       ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr);
492797d2ccbSBarry Smith       ierr = VecDestroy(*vec);CHKERRQ(ierr);
4930f77fa87SBarry Smith       ierr = PetscFree(vec);CHKERRQ(ierr);
4940f77fa87SBarry Smith     }
4954b9ad928SBarry Smith   }
4964b9ad928SBarry Smith 
497c2be2410SBarry Smith 
4984b9ad928SBarry Smith   for (i=1; i<n; i++) {
499b03c7568SBarry Smith     if (mg[i]->smoothu == mg[i]->smoothd) {
500b03c7568SBarry Smith       /* if doing only down then initial guess is zero */
5014b9ad928SBarry Smith       ierr = KSPSetInitialGuessNonzero(mg[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
502b03c7568SBarry Smith     }
50332cf1786SBarry Smith     if (mg[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mg[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
5044b9ad928SBarry Smith     ierr = KSPSetUp(mg[i]->smoothd);CHKERRQ(ierr);
50532cf1786SBarry Smith     if (mg[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mg[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
5064b9ad928SBarry Smith   }
507b03c7568SBarry Smith   for (i=1; i<n; i++) {
5084b9ad928SBarry Smith     if (mg[i]->smoothu && mg[i]->smoothu != mg[i]->smoothd) {
509906ed7ccSBarry Smith       Mat          downmat,downpmat;
51097f1f81fSBarry Smith       MatStructure matflag;
511906ed7ccSBarry Smith       PetscTruth   opsset;
51297f1f81fSBarry Smith 
51397f1f81fSBarry Smith       /* check if operators have been set for up, if not use down operators to set them */
514906ed7ccSBarry Smith       ierr = KSPGetOperatorsSet(mg[i]->smoothu,&opsset,PETSC_NULL);CHKERRQ(ierr);
515906ed7ccSBarry Smith       if (!opsset) {
516906ed7ccSBarry Smith         ierr = KSPGetOperators(mg[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr);
51797f1f81fSBarry Smith         ierr = KSPSetOperators(mg[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr);
51897f1f81fSBarry Smith       }
51997f1f81fSBarry Smith 
5204b9ad928SBarry Smith       ierr = KSPSetInitialGuessNonzero(mg[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
52132cf1786SBarry Smith       if (mg[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mg[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
5224b9ad928SBarry Smith       ierr = KSPSetUp(mg[i]->smoothu);CHKERRQ(ierr);
52332cf1786SBarry Smith       if (mg[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mg[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
5244b9ad928SBarry Smith     }
5254b9ad928SBarry Smith   }
5264b9ad928SBarry Smith 
5274b9ad928SBarry Smith   /*
5284b9ad928SBarry Smith       If coarse solver is not direct method then DO NOT USE preonly
5294b9ad928SBarry Smith   */
5304b9ad928SBarry Smith   ierr = PetscTypeCompare((PetscObject)mg[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr);
5314b9ad928SBarry Smith   if (preonly) {
5324b9ad928SBarry Smith     ierr = PetscTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr);
5334b9ad928SBarry Smith     ierr = PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr);
53468eff7e6SBarry Smith     ierr = PetscTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr);
53568eff7e6SBarry Smith     if (!lu && !redundant && !cholesky) {
5364b9ad928SBarry Smith       ierr = KSPSetType(mg[0]->smoothd,KSPGMRES);CHKERRQ(ierr);
5374b9ad928SBarry Smith     }
5384b9ad928SBarry Smith   }
5394b9ad928SBarry Smith 
540958c9bccSBarry Smith   if (!pc->setupcalled) {
5414b9ad928SBarry Smith     if (monitor) {
5424b9ad928SBarry Smith       ierr = PetscObjectGetComm((PetscObject)mg[0]->smoothd,&comm);CHKERRQ(ierr);
54323d894e5SBarry Smith       ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n,&ascii);CHKERRQ(ierr);
54423d894e5SBarry Smith       ierr = KSPMonitorSet(mg[0]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
5454b9ad928SBarry Smith     }
5464b9ad928SBarry Smith     ierr = KSPSetFromOptions(mg[0]->smoothd);CHKERRQ(ierr);
5474b9ad928SBarry Smith   }
5484b9ad928SBarry Smith 
54932cf1786SBarry Smith   if (mg[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mg[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
5504b9ad928SBarry Smith   ierr = KSPSetUp(mg[0]->smoothd);CHKERRQ(ierr);
55132cf1786SBarry Smith   if (mg[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mg[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
5524b9ad928SBarry Smith 
5534b9ad928SBarry Smith   /*
5546805f65bSBarry Smith      Dump the interpolation/restriction matrices plus the
5554b9ad928SBarry Smith    Jacobian/stiffness on each level. This allows Matlab users to
5566805f65bSBarry Smith    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
5576805f65bSBarry Smith 
5586805f65bSBarry Smith    Only support one or the other at the same time.
5596805f65bSBarry Smith   */
5606805f65bSBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
56190d69ab7SBarry Smith   ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,PETSC_NULL);CHKERRQ(ierr);
5624b9ad928SBarry Smith   if (dump) {
5637adad957SLisandro Dalcin     viewer = PETSC_VIEWER_SOCKET_(((PetscObject)pc)->comm);
5644b9ad928SBarry Smith   }
56590d69ab7SBarry Smith   dump = PETSC_FALSE;
566c45a1595SBarry Smith #endif
56790d69ab7SBarry Smith   ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,PETSC_NULL);CHKERRQ(ierr);
568c2be2410SBarry Smith   if (dump) {
5697adad957SLisandro Dalcin     viewer = PETSC_VIEWER_BINARY_(((PetscObject)pc)->comm);
5706805f65bSBarry Smith   }
5716805f65bSBarry Smith 
5726805f65bSBarry Smith   if (viewer) {
573c2be2410SBarry Smith     for (i=1; i<n; i++) {
5746805f65bSBarry Smith       ierr = MatView(mg[i]->restrct,viewer);CHKERRQ(ierr);
575c2be2410SBarry Smith     }
576c2be2410SBarry Smith     for (i=0; i<n; i++) {
577c2be2410SBarry Smith       ierr = KSPGetPC(mg[i]->smoothd,&pc);CHKERRQ(ierr);
5786805f65bSBarry Smith       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
579c2be2410SBarry Smith     }
580c2be2410SBarry Smith   }
5814b9ad928SBarry Smith   PetscFunctionReturn(0);
5824b9ad928SBarry Smith }
5834b9ad928SBarry Smith 
5844b9ad928SBarry Smith /* -------------------------------------------------------------------------------------*/
5854b9ad928SBarry Smith 
5864b9ad928SBarry Smith #undef __FUNCT__
5879dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetLevels"
5884b9ad928SBarry Smith /*@C
58997177400SBarry Smith    PCMGSetLevels - Sets the number of levels to use with MG.
5904b9ad928SBarry Smith    Must be called before any other MG routine.
5914b9ad928SBarry Smith 
5924b9ad928SBarry Smith    Collective on PC
5934b9ad928SBarry Smith 
5944b9ad928SBarry Smith    Input Parameters:
5954b9ad928SBarry Smith +  pc - the preconditioner context
5964b9ad928SBarry Smith .  levels - the number of levels
5974b9ad928SBarry Smith -  comms - optional communicators for each level; this is to allow solving the coarser problems
5984b9ad928SBarry Smith            on smaller sets of processors. Use PETSC_NULL_OBJECT for default in Fortran
5994b9ad928SBarry Smith 
6004b9ad928SBarry Smith    Level: intermediate
6014b9ad928SBarry Smith 
6024b9ad928SBarry Smith    Notes:
6034b9ad928SBarry Smith      If the number of levels is one then the multigrid uses the -mg_levels prefix
6044b9ad928SBarry Smith   for setting the level options rather than the -mg_coarse prefix.
6054b9ad928SBarry Smith 
6064b9ad928SBarry Smith .keywords: MG, set, levels, multigrid
6074b9ad928SBarry Smith 
60897177400SBarry Smith .seealso: PCMGSetType(), PCMGGetLevels()
6094b9ad928SBarry Smith @*/
61097177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
6114b9ad928SBarry Smith {
612dfbe8321SBarry Smith   PetscErrorCode ierr;
613ada7143aSSatish Balay   PC_MG          **mg=0;
6144b9ad928SBarry Smith 
6154b9ad928SBarry Smith   PetscFunctionBegin;
6164482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
6174b9ad928SBarry Smith 
6184b9ad928SBarry Smith   if (pc->data) {
6191302d50aSBarry Smith     SETERRQ(PETSC_ERR_ORDER,"Number levels already set for MG\n\
62097177400SBarry Smith     make sure that you call PCMGSetLevels() before KSPSetFromOptions()");
6214b9ad928SBarry Smith   }
6227adad957SLisandro Dalcin   ierr                     = PCMGCreate_Private(((PetscObject)pc)->comm,levels,pc,comms,&mg);CHKERRQ(ierr);
6239dcbbd2bSBarry Smith   mg[0]->am                = PC_MG_MULTIPLICATIVE;
6244b9ad928SBarry Smith   pc->data                 = (void*)mg;
6254b9ad928SBarry Smith   pc->ops->applyrichardson = PCApplyRichardson_MG;
6264b9ad928SBarry Smith   PetscFunctionReturn(0);
6274b9ad928SBarry Smith }
6284b9ad928SBarry Smith 
6294b9ad928SBarry Smith #undef __FUNCT__
6309dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetLevels"
6314b9ad928SBarry Smith /*@
63297177400SBarry Smith    PCMGGetLevels - Gets the number of levels to use with MG.
6334b9ad928SBarry Smith 
6344b9ad928SBarry Smith    Not Collective
6354b9ad928SBarry Smith 
6364b9ad928SBarry Smith    Input Parameter:
6374b9ad928SBarry Smith .  pc - the preconditioner context
6384b9ad928SBarry Smith 
6394b9ad928SBarry Smith    Output parameter:
6404b9ad928SBarry Smith .  levels - the number of levels
6414b9ad928SBarry Smith 
6424b9ad928SBarry Smith    Level: advanced
6434b9ad928SBarry Smith 
6444b9ad928SBarry Smith .keywords: MG, get, levels, multigrid
6454b9ad928SBarry Smith 
64697177400SBarry Smith .seealso: PCMGSetLevels()
6474b9ad928SBarry Smith @*/
64897177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetLevels(PC pc,PetscInt *levels)
6494b9ad928SBarry Smith {
6509dcbbd2bSBarry Smith   PC_MG  **mg;
6514b9ad928SBarry Smith 
6524b9ad928SBarry Smith   PetscFunctionBegin;
6534482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
6544482741eSBarry Smith   PetscValidIntPointer(levels,2);
6554b9ad928SBarry Smith 
6569dcbbd2bSBarry Smith   mg      = (PC_MG**)pc->data;
6574b9ad928SBarry Smith   *levels = mg[0]->levels;
6584b9ad928SBarry Smith   PetscFunctionReturn(0);
6594b9ad928SBarry Smith }
6604b9ad928SBarry Smith 
6614b9ad928SBarry Smith #undef __FUNCT__
6629dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetType"
6634b9ad928SBarry Smith /*@
66497177400SBarry Smith    PCMGSetType - Determines the form of multigrid to use:
6654b9ad928SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
6664b9ad928SBarry Smith 
6674b9ad928SBarry Smith    Collective on PC
6684b9ad928SBarry Smith 
6694b9ad928SBarry Smith    Input Parameters:
6704b9ad928SBarry Smith +  pc - the preconditioner context
6719dcbbd2bSBarry Smith -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
6729dcbbd2bSBarry Smith    PC_MG_FULL, PC_MG_KASKADE
6734b9ad928SBarry Smith 
6744b9ad928SBarry Smith    Options Database Key:
6754b9ad928SBarry Smith .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
6764b9ad928SBarry Smith    additive, full, kaskade
6774b9ad928SBarry Smith 
6784b9ad928SBarry Smith    Level: advanced
6794b9ad928SBarry Smith 
6804b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
6814b9ad928SBarry Smith 
68297177400SBarry Smith .seealso: PCMGSetLevels()
6834b9ad928SBarry Smith @*/
6849dcbbd2bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetType(PC pc,PCMGType form)
6854b9ad928SBarry Smith {
6869dcbbd2bSBarry Smith   PC_MG **mg;
6874b9ad928SBarry Smith 
6884b9ad928SBarry Smith   PetscFunctionBegin;
6894482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
6909dcbbd2bSBarry Smith   mg = (PC_MG**)pc->data;
6914b9ad928SBarry Smith 
6924b9ad928SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
6934b9ad928SBarry Smith   mg[0]->am = form;
6949dcbbd2bSBarry Smith   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
6954b9ad928SBarry Smith   else pc->ops->applyrichardson = 0;
6964b9ad928SBarry Smith   PetscFunctionReturn(0);
6974b9ad928SBarry Smith }
6984b9ad928SBarry Smith 
6994b9ad928SBarry Smith #undef __FUNCT__
7000d353602SBarry Smith #define __FUNCT__ "PCMGSetCycleType"
7014b9ad928SBarry Smith /*@
7020d353602SBarry Smith    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
7034b9ad928SBarry Smith    complicated cycling.
7044b9ad928SBarry Smith 
7054b9ad928SBarry Smith    Collective on PC
7064b9ad928SBarry Smith 
7074b9ad928SBarry Smith    Input Parameters:
708c2be2410SBarry Smith +  pc - the multigrid context
7090d353602SBarry Smith -  PC_MG_CYCLE_V or PC_MG_CYCLE_W
7104b9ad928SBarry Smith 
7114b9ad928SBarry Smith    Options Database Key:
7120d353602SBarry Smith $  -pc_mg_cycle_type v or w
7134b9ad928SBarry Smith 
7144b9ad928SBarry Smith    Level: advanced
7154b9ad928SBarry Smith 
7164b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
7174b9ad928SBarry Smith 
7180d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel()
7194b9ad928SBarry Smith @*/
7200d353602SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetCycleType(PC pc,PCMGCycleType n)
7214b9ad928SBarry Smith {
7229dcbbd2bSBarry Smith   PC_MG    **mg;
72379416396SBarry Smith   PetscInt i,levels;
7244b9ad928SBarry Smith 
7254b9ad928SBarry Smith   PetscFunctionBegin;
7264482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
7279dcbbd2bSBarry Smith   mg     = (PC_MG**)pc->data;
7284b9ad928SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
7294b9ad928SBarry Smith   levels = mg[0]->levels;
7304b9ad928SBarry Smith 
7314b9ad928SBarry Smith   for (i=0; i<levels; i++) {
7324b9ad928SBarry Smith     mg[i]->cycles  = n;
7334b9ad928SBarry Smith   }
7344b9ad928SBarry Smith   PetscFunctionReturn(0);
7354b9ad928SBarry Smith }
7364b9ad928SBarry Smith 
7374b9ad928SBarry Smith #undef __FUNCT__
7388cc2d5dfSBarry Smith #define __FUNCT__ "PCMGMultiplicativeSetCycles"
7398cc2d5dfSBarry Smith /*@
7408cc2d5dfSBarry Smith    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
7418cc2d5dfSBarry Smith          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
7428cc2d5dfSBarry Smith 
7438cc2d5dfSBarry Smith    Collective on PC
7448cc2d5dfSBarry Smith 
7458cc2d5dfSBarry Smith    Input Parameters:
7468cc2d5dfSBarry Smith +  pc - the multigrid context
7478cc2d5dfSBarry Smith -  n - number of cycles (default is 1)
7488cc2d5dfSBarry Smith 
7498cc2d5dfSBarry Smith    Options Database Key:
7508cc2d5dfSBarry Smith $  -pc_mg_multiplicative_cycles n
7518cc2d5dfSBarry Smith 
7528cc2d5dfSBarry Smith    Level: advanced
7538cc2d5dfSBarry Smith 
7548cc2d5dfSBarry Smith    Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
7558cc2d5dfSBarry Smith 
7568cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
7578cc2d5dfSBarry Smith 
7588cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
7598cc2d5dfSBarry Smith @*/
7608cc2d5dfSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
7618cc2d5dfSBarry Smith {
7628cc2d5dfSBarry Smith   PC_MG    **mg;
7638cc2d5dfSBarry Smith   PetscInt i,levels;
7648cc2d5dfSBarry Smith 
7658cc2d5dfSBarry Smith   PetscFunctionBegin;
7668cc2d5dfSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
7678cc2d5dfSBarry Smith   mg     = (PC_MG**)pc->data;
7688cc2d5dfSBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
7698cc2d5dfSBarry Smith   levels = mg[0]->levels;
7708cc2d5dfSBarry Smith 
7718cc2d5dfSBarry Smith   for (i=0; i<levels; i++) {
7728cc2d5dfSBarry Smith     mg[i]->cyclesperpcapply  = n;
7738cc2d5dfSBarry Smith   }
7748cc2d5dfSBarry Smith   PetscFunctionReturn(0);
7758cc2d5dfSBarry Smith }
7768cc2d5dfSBarry Smith 
7778cc2d5dfSBarry Smith #undef __FUNCT__
7789dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetGalerkin"
779c2be2410SBarry Smith /*@
78097177400SBarry Smith    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
781c2be2410SBarry Smith       finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t
782c2be2410SBarry Smith 
783c2be2410SBarry Smith    Collective on PC
784c2be2410SBarry Smith 
785c2be2410SBarry Smith    Input Parameters:
7863fc8bf9cSBarry Smith .  pc - the multigrid context
787c2be2410SBarry Smith 
788c2be2410SBarry Smith    Options Database Key:
789c2be2410SBarry Smith $  -pc_mg_galerkin
790c2be2410SBarry Smith 
791c2be2410SBarry Smith    Level: intermediate
792c2be2410SBarry Smith 
793c2be2410SBarry Smith .keywords: MG, set, Galerkin
794c2be2410SBarry Smith 
7953fc8bf9cSBarry Smith .seealso: PCMGGetGalerkin()
7963fc8bf9cSBarry Smith 
797c2be2410SBarry Smith @*/
79897177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetGalerkin(PC pc)
799c2be2410SBarry Smith {
8009dcbbd2bSBarry Smith   PC_MG    **mg;
801c2be2410SBarry Smith   PetscInt i,levels;
802c2be2410SBarry Smith 
803c2be2410SBarry Smith   PetscFunctionBegin;
804c2be2410SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
8059dcbbd2bSBarry Smith   mg     = (PC_MG**)pc->data;
806c2be2410SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
807c2be2410SBarry Smith   levels = mg[0]->levels;
808c2be2410SBarry Smith 
809c2be2410SBarry Smith   for (i=0; i<levels; i++) {
810c2be2410SBarry Smith     mg[i]->galerkin = PETSC_TRUE;
811c2be2410SBarry Smith   }
812c2be2410SBarry Smith   PetscFunctionReturn(0);
813c2be2410SBarry Smith }
814c2be2410SBarry Smith 
815c2be2410SBarry Smith #undef __FUNCT__
8163fc8bf9cSBarry Smith #define __FUNCT__ "PCMGGetGalerkin"
8173fc8bf9cSBarry Smith /*@
8183fc8bf9cSBarry Smith    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
8193fc8bf9cSBarry Smith       A_i-1 = r_i * A_i * r_i^t
8203fc8bf9cSBarry Smith 
8213fc8bf9cSBarry Smith    Not Collective
8223fc8bf9cSBarry Smith 
8233fc8bf9cSBarry Smith    Input Parameter:
8243fc8bf9cSBarry Smith .  pc - the multigrid context
8253fc8bf9cSBarry Smith 
8263fc8bf9cSBarry Smith    Output Parameter:
8273fc8bf9cSBarry Smith .  gelerkin - PETSC_TRUE or PETSC_FALSE
8283fc8bf9cSBarry Smith 
8293fc8bf9cSBarry Smith    Options Database Key:
8303fc8bf9cSBarry Smith $  -pc_mg_galerkin
8313fc8bf9cSBarry Smith 
8323fc8bf9cSBarry Smith    Level: intermediate
8333fc8bf9cSBarry Smith 
8343fc8bf9cSBarry Smith .keywords: MG, set, Galerkin
8353fc8bf9cSBarry Smith 
8363fc8bf9cSBarry Smith .seealso: PCMGSetGalerkin()
8373fc8bf9cSBarry Smith 
8383fc8bf9cSBarry Smith @*/
8393fc8bf9cSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetGalerkin(PC pc,PetscTruth *galerkin)
8403fc8bf9cSBarry Smith {
8413fc8bf9cSBarry Smith   PC_MG    **mg;
8423fc8bf9cSBarry Smith 
8433fc8bf9cSBarry Smith   PetscFunctionBegin;
8443fc8bf9cSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
8453fc8bf9cSBarry Smith   mg     = (PC_MG**)pc->data;
8463fc8bf9cSBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
8473fc8bf9cSBarry Smith   *galerkin = mg[0]->galerkin;
8483fc8bf9cSBarry Smith   PetscFunctionReturn(0);
8493fc8bf9cSBarry Smith }
8503fc8bf9cSBarry Smith 
8513fc8bf9cSBarry Smith #undef __FUNCT__
8529dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothDown"
8534b9ad928SBarry Smith /*@
85497177400SBarry Smith    PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
85597177400SBarry Smith    use on all levels. Use PCMGGetSmootherDown() to set different
8564b9ad928SBarry Smith    pre-smoothing steps on different levels.
8574b9ad928SBarry Smith 
8584b9ad928SBarry Smith    Collective on PC
8594b9ad928SBarry Smith 
8604b9ad928SBarry Smith    Input Parameters:
8614b9ad928SBarry Smith +  mg - the multigrid context
8624b9ad928SBarry Smith -  n - the number of smoothing steps
8634b9ad928SBarry Smith 
8644b9ad928SBarry Smith    Options Database Key:
8654b9ad928SBarry Smith .  -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps
8664b9ad928SBarry Smith 
8674b9ad928SBarry Smith    Level: advanced
8684b9ad928SBarry Smith 
8694b9ad928SBarry Smith .keywords: MG, smooth, down, pre-smoothing, steps, multigrid
8704b9ad928SBarry Smith 
87197177400SBarry Smith .seealso: PCMGSetNumberSmoothUp()
8724b9ad928SBarry Smith @*/
87397177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothDown(PC pc,PetscInt n)
8744b9ad928SBarry Smith {
8759dcbbd2bSBarry Smith   PC_MG          **mg;
8766849ba73SBarry Smith   PetscErrorCode ierr;
87779416396SBarry Smith   PetscInt       i,levels;
8784b9ad928SBarry Smith 
8794b9ad928SBarry Smith   PetscFunctionBegin;
8804482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
8819dcbbd2bSBarry Smith   mg     = (PC_MG**)pc->data;
8824b9ad928SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
8834b9ad928SBarry Smith   levels = mg[0]->levels;
8844b9ad928SBarry Smith 
885b05257ddSBarry Smith   for (i=1; i<levels; i++) {
8864b9ad928SBarry Smith     /* make sure smoother up and down are different */
88797177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
8884b9ad928SBarry Smith     ierr = KSPSetTolerances(mg[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
8894b9ad928SBarry Smith     mg[i]->default_smoothd = n;
8904b9ad928SBarry Smith   }
8914b9ad928SBarry Smith   PetscFunctionReturn(0);
8924b9ad928SBarry Smith }
8934b9ad928SBarry Smith 
8944b9ad928SBarry Smith #undef __FUNCT__
8959dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothUp"
8964b9ad928SBarry Smith /*@
89797177400SBarry Smith    PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
89897177400SBarry Smith    on all levels. Use PCMGGetSmootherUp() to set different numbers of
8994b9ad928SBarry Smith    post-smoothing steps on different levels.
9004b9ad928SBarry Smith 
9014b9ad928SBarry Smith    Collective on PC
9024b9ad928SBarry Smith 
9034b9ad928SBarry Smith    Input Parameters:
9044b9ad928SBarry Smith +  mg - the multigrid context
9054b9ad928SBarry Smith -  n - the number of smoothing steps
9064b9ad928SBarry Smith 
9074b9ad928SBarry Smith    Options Database Key:
9084b9ad928SBarry Smith .  -pc_mg_smoothup <n> - Sets number of post-smoothing steps
9094b9ad928SBarry Smith 
9104b9ad928SBarry Smith    Level: advanced
9114b9ad928SBarry Smith 
9124b9ad928SBarry Smith    Note: this does not set a value on the coarsest grid, since we assume that
913a8c7a070SBarry Smith     there is no separate smooth up on the coarsest grid.
9144b9ad928SBarry Smith 
9154b9ad928SBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid
9164b9ad928SBarry Smith 
91797177400SBarry Smith .seealso: PCMGSetNumberSmoothDown()
9184b9ad928SBarry Smith @*/
91997177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothUp(PC pc,PetscInt n)
9204b9ad928SBarry Smith {
9219dcbbd2bSBarry Smith   PC_MG          **mg;
9226849ba73SBarry Smith   PetscErrorCode ierr;
92379416396SBarry Smith   PetscInt       i,levels;
9244b9ad928SBarry Smith 
9254b9ad928SBarry Smith   PetscFunctionBegin;
9264482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
9279dcbbd2bSBarry Smith   mg     = (PC_MG**)pc->data;
9284b9ad928SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
9294b9ad928SBarry Smith   levels = mg[0]->levels;
9304b9ad928SBarry Smith 
9314b9ad928SBarry Smith   for (i=1; i<levels; i++) {
9324b9ad928SBarry Smith     /* make sure smoother up and down are different */
93397177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
9344b9ad928SBarry Smith     ierr = KSPSetTolerances(mg[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
9354b9ad928SBarry Smith     mg[i]->default_smoothu = n;
9364b9ad928SBarry Smith   }
9374b9ad928SBarry Smith   PetscFunctionReturn(0);
9384b9ad928SBarry Smith }
9394b9ad928SBarry Smith 
9404b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/
9414b9ad928SBarry Smith 
9423b09bd56SBarry Smith /*MC
943ccb205f8SBarry Smith    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
9443b09bd56SBarry Smith     information about the coarser grid matrices and restriction/interpolation operators.
9453b09bd56SBarry Smith 
9463b09bd56SBarry Smith    Options Database Keys:
9473b09bd56SBarry Smith +  -pc_mg_levels <nlevels> - number of levels including finest
9480d353602SBarry Smith .  -pc_mg_cycles v or w
94979416396SBarry Smith .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
9503b09bd56SBarry Smith .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
9513b09bd56SBarry Smith .  -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default
9523b09bd56SBarry Smith .  -pc_mg_log - log information about time spent on each level of the solver
9533b09bd56SBarry Smith .  -pc_mg_monitor - print information on the multigrid convergence
95468eff7e6SBarry Smith .  -pc_mg_galerkin - use Galerkin process to compute coarser operators
9553b09bd56SBarry Smith -  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
9563b09bd56SBarry Smith                         to the Socket viewer for reading from Matlab.
9573b09bd56SBarry Smith 
9583b09bd56SBarry Smith    Notes:
9593b09bd56SBarry Smith 
9603b09bd56SBarry Smith    Level: intermediate
9613b09bd56SBarry Smith 
9628f87f92bSBarry Smith    Concepts: multigrid/multilevel
9633b09bd56SBarry Smith 
9643b09bd56SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
9650d353602SBarry Smith            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
96697177400SBarry Smith            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
96797177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
9680d353602SBarry Smith            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
9693b09bd56SBarry Smith M*/
9703b09bd56SBarry Smith 
9714b9ad928SBarry Smith EXTERN_C_BEGIN
9724b9ad928SBarry Smith #undef __FUNCT__
9734b9ad928SBarry Smith #define __FUNCT__ "PCCreate_MG"
974dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_MG(PC pc)
9754b9ad928SBarry Smith {
9764b9ad928SBarry Smith   PetscFunctionBegin;
9774b9ad928SBarry Smith   pc->ops->apply          = PCApply_MG;
9784b9ad928SBarry Smith   pc->ops->setup          = PCSetUp_MG;
9794b9ad928SBarry Smith   pc->ops->destroy        = PCDestroy_MG;
9804b9ad928SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_MG;
9814b9ad928SBarry Smith   pc->ops->view           = PCView_MG;
9824b9ad928SBarry Smith 
9834b9ad928SBarry Smith   pc->data                = (void*)0;
9844b9ad928SBarry Smith   PetscFunctionReturn(0);
9854b9ad928SBarry Smith }
9864b9ad928SBarry Smith EXTERN_C_END
987