xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision aea2a34e36487be060e5e9669bb39cba3c3f8e3b)
1dba47a55SKris Buschelman #define PETSCKSP_DLL
2dba47a55SKris Buschelman 
34b9ad928SBarry Smith /*
44b9ad928SBarry Smith     Defines the multigrid preconditioner interface.
54b9ad928SBarry Smith */
64b9ad928SBarry 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"
119dcbbd2bSBarry Smith PetscErrorCode PCMGMCycle_Private(PC_MG **mglevels,PetscTruth *converged)
124b9ad928SBarry Smith {
139dcbbd2bSBarry Smith   PC_MG          *mg = *mglevels,*mgc;
146849ba73SBarry Smith   PetscErrorCode ierr;
150d353602SBarry Smith   PetscInt       cycles = (PetscInt) mg->cycles;
164b9ad928SBarry Smith 
174b9ad928SBarry Smith   PetscFunctionBegin;
184b9ad928SBarry Smith   if (converged) *converged = PETSC_FALSE;
194b9ad928SBarry Smith 
2032cf1786SBarry Smith   if (mg->eventsmoothsolve) {ierr = PetscLogEventBegin(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
210d353602SBarry Smith   ierr = KSPSolve(mg->smoothd,mg->b,mg->x);CHKERRQ(ierr);  /* pre-smooth */
2232cf1786SBarry Smith   if (mg->eventsmoothsolve) {ierr = PetscLogEventEnd(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
234b9ad928SBarry Smith   if (mg->level) {  /* not the coarsest grid */
2432cf1786SBarry Smith     if (mg->eventresidual) {ierr = PetscLogEventBegin(mg->eventresidual,0,0,0,0);CHKERRQ(ierr);}
254b9ad928SBarry Smith     ierr = (*mg->residual)(mg->A,mg->b,mg->x,mg->r);CHKERRQ(ierr);
2632cf1786SBarry Smith     if (mg->eventresidual) {ierr = PetscLogEventEnd(mg->eventresidual,0,0,0,0);CHKERRQ(ierr);}
274b9ad928SBarry Smith 
284b9ad928SBarry Smith     /* if on finest level and have convergence criteria set */
294b9ad928SBarry Smith     if (mg->level == mg->levels-1 && mg->ttol) {
304b9ad928SBarry Smith       PetscReal rnorm;
314b9ad928SBarry Smith       ierr = VecNorm(mg->r,NORM_2,&rnorm);CHKERRQ(ierr);
324b9ad928SBarry Smith       if (rnorm <= mg->ttol) {
334b9ad928SBarry Smith         *converged = PETSC_TRUE;
3470441072SBarry Smith         if (rnorm < mg->abstol) {
35ae15b995SBarry Smith           ierr = PetscInfo2(0,"Linear solver has converged. Residual norm %G is less than absolute tolerance %G\n",rnorm,mg->abstol);CHKERRQ(ierr);
364b9ad928SBarry Smith         } else {
37ae15b995SBarry Smith           ierr = PetscInfo2(0,"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--) {
499dcbbd2bSBarry Smith       ierr = PCMGMCycle_Private(mglevels-1,converged);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);
809dcbbd2bSBarry Smith   ierr = PetscLogObjectMemory(pc,levels*(sizeof(PC_MG*)+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++) {
859dcbbd2bSBarry Smith     ierr = PetscNew(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);
9679416396SBarry Smith     ierr = KSPSetTolerances(mg[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg[i]->default_smoothd);CHKERRQ(ierr);
974b9ad928SBarry Smith     ierr = KSPSetOptionsPrefix(mg[i]->smoothd,prefix);CHKERRQ(ierr);
984b9ad928SBarry Smith 
994b9ad928SBarry Smith     /* do special stuff for coarse grid */
1004b9ad928SBarry Smith     if (!i && levels > 1) {
1014b9ad928SBarry Smith       ierr = KSPAppendOptionsPrefix(mg[0]->smoothd,"mg_coarse_");CHKERRQ(ierr);
1024b9ad928SBarry Smith 
1034b9ad928SBarry Smith       /* coarse solve is (redundant) LU by default */
1044b9ad928SBarry Smith       ierr = KSPSetType(mg[0]->smoothd,KSPPREONLY);CHKERRQ(ierr);
1054b9ad928SBarry Smith       ierr = KSPGetPC(mg[0]->smoothd,&ipc);CHKERRQ(ierr);
1064b9ad928SBarry Smith       ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1074b9ad928SBarry Smith       if (size > 1) {
1084b9ad928SBarry Smith         ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr);
1090d7810c8SBarry Smith       } else {
1104b9ad928SBarry Smith         ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr);
1110d7810c8SBarry Smith       }
1124b9ad928SBarry Smith 
1134b9ad928SBarry Smith     } else {
1142e70eadcSBarry Smith       char tprefix[128];
11513f74950SBarry Smith       sprintf(tprefix,"mg_levels_%d_",(int)i);
1162e70eadcSBarry Smith       ierr = KSPAppendOptionsPrefix(mg[i]->smoothd,tprefix);CHKERRQ(ierr);
1174b9ad928SBarry Smith     }
11852e6d16bSBarry Smith     ierr = PetscLogObjectParent(pc,mg[i]->smoothd);CHKERRQ(ierr);
1194b9ad928SBarry Smith     mg[i]->smoothu             = mg[i]->smoothd;
1204b9ad928SBarry Smith     mg[i]->rtol                = 0.0;
12170441072SBarry Smith     mg[i]->abstol              = 0.0;
1224b9ad928SBarry Smith     mg[i]->dtol                = 0.0;
1234b9ad928SBarry Smith     mg[i]->ttol                = 0.0;
12432cf1786SBarry Smith     mg[i]->eventsmoothsetup    = 0;
12532cf1786SBarry Smith     mg[i]->eventsmoothsolve    = 0;
12632cf1786SBarry Smith     mg[i]->eventresidual       = 0;
12732cf1786SBarry Smith     mg[i]->eventinterprestrict = 0;
1288cc2d5dfSBarry Smith     mg[i]->cyclesperpcapply    = 1;
1294b9ad928SBarry Smith   }
1304b9ad928SBarry Smith   *result = mg;
1314b9ad928SBarry Smith   PetscFunctionReturn(0);
1324b9ad928SBarry Smith }
1334b9ad928SBarry Smith 
1344b9ad928SBarry Smith #undef __FUNCT__
1354b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_MG"
1366849ba73SBarry Smith static PetscErrorCode PCDestroy_MG(PC pc)
1374b9ad928SBarry Smith {
1389dcbbd2bSBarry Smith   PC_MG          **mg = (PC_MG**)pc->data;
1396849ba73SBarry Smith   PetscErrorCode ierr;
14079416396SBarry Smith   PetscInt       i,n = mg[0]->levels;
1414b9ad928SBarry Smith 
1424b9ad928SBarry Smith   PetscFunctionBegin;
143fccaa45eSBarry Smith   for (i=0; i<n-1; i++) {
144630ba2eeSBarry Smith     if (mg[i+1]->r) {ierr = VecDestroy(mg[i+1]->r);CHKERRQ(ierr);}
145fccaa45eSBarry Smith     if (mg[i]->b) {ierr = VecDestroy(mg[i]->b);CHKERRQ(ierr);}
146fccaa45eSBarry Smith     if (mg[i]->x) {ierr = VecDestroy(mg[i]->x);CHKERRQ(ierr);}
147fccaa45eSBarry Smith     if (mg[i+1]->restrct) {ierr = MatDestroy(mg[i+1]->restrct);CHKERRQ(ierr);}
148fccaa45eSBarry Smith     if (mg[i+1]->interpolate) {ierr = MatDestroy(mg[i+1]->interpolate);CHKERRQ(ierr);}
149fccaa45eSBarry Smith   }
150fccaa45eSBarry Smith 
1514b9ad928SBarry Smith   for (i=0; i<n; i++) {
1524b9ad928SBarry Smith     if (mg[i]->smoothd != mg[i]->smoothu) {
1534b9ad928SBarry Smith       ierr = KSPDestroy(mg[i]->smoothd);CHKERRQ(ierr);
1544b9ad928SBarry Smith     }
1554b9ad928SBarry Smith     ierr = KSPDestroy(mg[i]->smoothu);CHKERRQ(ierr);
1564b9ad928SBarry Smith     ierr = PetscFree(mg[i]);CHKERRQ(ierr);
1574b9ad928SBarry Smith   }
1584b9ad928SBarry Smith   ierr = PetscFree(mg);CHKERRQ(ierr);
1594b9ad928SBarry Smith   PetscFunctionReturn(0);
1604b9ad928SBarry Smith }
1614b9ad928SBarry Smith 
1624b9ad928SBarry Smith 
1634b9ad928SBarry Smith 
1649dcbbd2bSBarry Smith EXTERN PetscErrorCode PCMGACycle_Private(PC_MG**);
1659dcbbd2bSBarry Smith EXTERN PetscErrorCode PCMGFCycle_Private(PC_MG**);
1669dcbbd2bSBarry Smith EXTERN PetscErrorCode PCMGKCycle_Private(PC_MG**);
1674b9ad928SBarry Smith 
1684b9ad928SBarry Smith /*
1694b9ad928SBarry Smith    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
1704b9ad928SBarry Smith              or full cycle of multigrid.
1714b9ad928SBarry Smith 
1724b9ad928SBarry Smith   Note:
1739dcbbd2bSBarry Smith   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
1744b9ad928SBarry Smith */
1754b9ad928SBarry Smith #undef __FUNCT__
1764b9ad928SBarry Smith #define __FUNCT__ "PCApply_MG"
1776849ba73SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
1784b9ad928SBarry Smith {
1799dcbbd2bSBarry Smith   PC_MG          **mg = (PC_MG**)pc->data;
1806849ba73SBarry Smith   PetscErrorCode ierr;
1818cc2d5dfSBarry Smith   PetscInt       levels = mg[0]->levels,i;
1824b9ad928SBarry Smith 
1834b9ad928SBarry Smith   PetscFunctionBegin;
1844b9ad928SBarry Smith   mg[levels-1]->b = b;
1854b9ad928SBarry Smith   mg[levels-1]->x = x;
186ef70c39aSMatthew Knepley   if (!mg[levels-1]->r && mg[0]->am != PC_MG_ADDITIVE && levels > 1) {
1870a6bb862SBarry Smith     Vec tvec;
1880a6bb862SBarry Smith     ierr = VecDuplicate(mg[levels-1]->b,&tvec);CHKERRQ(ierr);
18997177400SBarry Smith     ierr = PCMGSetR(pc,levels-1,tvec);CHKERRQ(ierr);
1900a6bb862SBarry Smith     ierr = VecDestroy(tvec);CHKERRQ(ierr);
1910a6bb862SBarry Smith   }
1929dcbbd2bSBarry Smith   if (mg[0]->am == PC_MG_MULTIPLICATIVE) {
193efb30889SBarry Smith     ierr = VecSet(x,0.0);CHKERRQ(ierr);
1948cc2d5dfSBarry Smith     for (i=0; i<mg[0]->cyclesperpcapply; i++) {
1959dcbbd2bSBarry Smith       ierr = PCMGMCycle_Private(mg+levels-1,PETSC_NULL);CHKERRQ(ierr);
1964b9ad928SBarry Smith     }
1978cc2d5dfSBarry Smith   }
1989dcbbd2bSBarry Smith   else if (mg[0]->am == PC_MG_ADDITIVE) {
1999dcbbd2bSBarry Smith     ierr = PCMGACycle_Private(mg);CHKERRQ(ierr);
2004b9ad928SBarry Smith   }
2019dcbbd2bSBarry Smith   else if (mg[0]->am == PC_MG_KASKADE) {
2029dcbbd2bSBarry Smith     ierr = PCMGKCycle_Private(mg);CHKERRQ(ierr);
2034b9ad928SBarry Smith   }
2044b9ad928SBarry Smith   else {
2059dcbbd2bSBarry Smith     ierr = PCMGFCycle_Private(mg);CHKERRQ(ierr);
2064b9ad928SBarry Smith   }
2074b9ad928SBarry Smith   PetscFunctionReturn(0);
2084b9ad928SBarry Smith }
2094b9ad928SBarry Smith 
2104b9ad928SBarry Smith #undef __FUNCT__
2114b9ad928SBarry Smith #define __FUNCT__ "PCApplyRichardson_MG"
21279416396SBarry Smith static PetscErrorCode PCApplyRichardson_MG(PC pc,Vec b,Vec x,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its)
2134b9ad928SBarry Smith {
2149dcbbd2bSBarry Smith   PC_MG          **mg = (PC_MG**)pc->data;
215dfbe8321SBarry Smith   PetscErrorCode ierr;
21679416396SBarry Smith   PetscInt       levels = mg[0]->levels;
2174b9ad928SBarry Smith   PetscTruth     converged = PETSC_FALSE;
2184b9ad928SBarry Smith 
2194b9ad928SBarry Smith   PetscFunctionBegin;
2204b9ad928SBarry Smith   mg[levels-1]->b    = b;
2214b9ad928SBarry Smith   mg[levels-1]->x    = x;
2224b9ad928SBarry Smith 
2234b9ad928SBarry Smith   mg[levels-1]->rtol = rtol;
22470441072SBarry Smith   mg[levels-1]->abstol = abstol;
2254b9ad928SBarry Smith   mg[levels-1]->dtol = dtol;
2264b9ad928SBarry Smith   if (rtol) {
2274b9ad928SBarry Smith     /* compute initial residual norm for relative convergence test */
2284b9ad928SBarry Smith     PetscReal rnorm;
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);
23170441072SBarry Smith     mg[levels-1]->ttol = PetscMax(rtol*rnorm,abstol);
23270441072SBarry Smith   } else if (abstol) {
23370441072SBarry Smith     mg[levels-1]->ttol = abstol;
2344b9ad928SBarry Smith   } else {
2354b9ad928SBarry Smith     mg[levels-1]->ttol = 0.0;
2364b9ad928SBarry Smith   }
2374b9ad928SBarry Smith 
2384b9ad928SBarry Smith   while (its-- && !converged) {
2399dcbbd2bSBarry Smith     ierr = PCMGMCycle_Private(mg+levels-1,&converged);CHKERRQ(ierr);
2404b9ad928SBarry Smith   }
2414b9ad928SBarry Smith   PetscFunctionReturn(0);
2424b9ad928SBarry Smith }
2434b9ad928SBarry Smith 
2444b9ad928SBarry Smith #undef __FUNCT__
2454b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_MG"
2466ca4d86aSHong Zhang PetscErrorCode PCSetFromOptions_MG(PC pc)
2474b9ad928SBarry Smith {
248dfbe8321SBarry Smith   PetscErrorCode ierr;
2498cc2d5dfSBarry Smith   PetscInt       m,levels = 1,cycles;
2504b9ad928SBarry Smith   PetscTruth     flg;
2519dcbbd2bSBarry Smith   PC_MG          **mg = (PC_MG**)pc->data;
25290da95b0SMatthew Knepley   PCMGType       mgtype = PC_MG_ADDITIVE;
2531aa34eeaSBarry Smith   PCMGCycleType  mgctype;
2544b9ad928SBarry Smith 
2554b9ad928SBarry Smith   PetscFunctionBegin;
2564b9ad928SBarry Smith   ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr);
2574b9ad928SBarry Smith     if (!pc->data) {
2589dcbbd2bSBarry Smith       ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr);
25997177400SBarry Smith       ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr);
260cf502942SBarry Smith       mg = (PC_MG**)pc->data;
2614b9ad928SBarry Smith     }
262f2070a76SMatthew Knepley     mgctype = (PCMGCycleType) mg[0]->cycles;
2630d353602SBarry Smith     ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr);
2644b9ad928SBarry Smith     if (flg) {
2650d353602SBarry Smith       ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr);
2660d353602SBarry Smith     };
2679dcbbd2bSBarry Smith     ierr = PetscOptionsName("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",&flg);CHKERRQ(ierr);
268c2be2410SBarry Smith     if (flg) {
26997177400SBarry Smith       ierr = PCMGSetGalerkin(pc);CHKERRQ(ierr);
270c2be2410SBarry Smith     }
2719dcbbd2bSBarry Smith     ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr);
2724b9ad928SBarry Smith     if (flg) {
27397177400SBarry Smith       ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr);
2744b9ad928SBarry Smith     }
2759dcbbd2bSBarry Smith     ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr);
2764b9ad928SBarry Smith     if (flg) {
27797177400SBarry Smith       ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr);
2784b9ad928SBarry Smith     }
2799dcbbd2bSBarry Smith     ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
2808cc2d5dfSBarry Smith     if (flg) {
2818cc2d5dfSBarry Smith       ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
2828cc2d5dfSBarry Smith     }
2838cc2d5dfSBarry Smith     if (mg[0]->am == PC_MG_MULTIPLICATIVE) {
2848cc2d5dfSBarry Smith       ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGSetLevels",mg[0]->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr);
2858cc2d5dfSBarry Smith       if (flg) {
2868cc2d5dfSBarry Smith 	ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr);
2878cc2d5dfSBarry Smith       }
2888cc2d5dfSBarry Smith     }
2894b9ad928SBarry Smith     ierr = PetscOptionsName("-pc_mg_log","Log times for each multigrid level","None",&flg);CHKERRQ(ierr);
2904b9ad928SBarry Smith     if (flg) {
2914f5ab15aSBarry Smith       PetscInt i;
2924b9ad928SBarry Smith       char     eventname[128];
2934b9ad928SBarry Smith       if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
2944b9ad928SBarry Smith       levels = mg[0]->levels;
2954b9ad928SBarry Smith       for (i=0; i<levels; i++) {
29632cf1786SBarry Smith         sprintf(eventname,"MGSetup Level %d",(int)i);
29732cf1786SBarry Smith         ierr = PetscLogEventRegister(&mg[i]->eventsmoothsetup,eventname,pc->cookie);CHKERRQ(ierr);
29832cf1786SBarry Smith         sprintf(eventname,"MGSmooth Level %d",(int)i);
29932cf1786SBarry Smith         ierr = PetscLogEventRegister(&mg[i]->eventsmoothsolve,eventname,pc->cookie);CHKERRQ(ierr);
30032cf1786SBarry Smith         if (i) {
30132cf1786SBarry Smith           sprintf(eventname,"MGResid Level %d",(int)i);
30232cf1786SBarry Smith           ierr = PetscLogEventRegister(&mg[i]->eventresidual,eventname,pc->cookie);CHKERRQ(ierr);
30332cf1786SBarry Smith           sprintf(eventname,"MGInterp Level %d",(int)i);
30432cf1786SBarry Smith           ierr = PetscLogEventRegister(&mg[i]->eventinterprestrict,eventname,pc->cookie);CHKERRQ(ierr);
30532cf1786SBarry Smith         }
3064b9ad928SBarry Smith       }
3074b9ad928SBarry Smith     }
3084b9ad928SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
3094b9ad928SBarry Smith   PetscFunctionReturn(0);
3104b9ad928SBarry Smith }
3114b9ad928SBarry Smith 
3129dcbbd2bSBarry Smith const char *PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0};
3130d353602SBarry Smith const char *PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0};
3149dcbbd2bSBarry Smith 
3154b9ad928SBarry Smith #undef __FUNCT__
3164b9ad928SBarry Smith #define __FUNCT__ "PCView_MG"
3176849ba73SBarry Smith static PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
3184b9ad928SBarry Smith {
3199dcbbd2bSBarry Smith   PC_MG          **mg = (PC_MG**)pc->data;
320dfbe8321SBarry Smith   PetscErrorCode ierr;
32179416396SBarry Smith   PetscInt       levels = mg[0]->levels,i;
32232077d6dSBarry Smith   PetscTruth     iascii;
3234b9ad928SBarry Smith 
3244b9ad928SBarry Smith   PetscFunctionBegin;
32532077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
32632077d6dSBarry Smith   if (iascii) {
3270d353602SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  MG: type is %s, levels=%D cycles=%s, pre-smooths=%D, post-smooths=%D\n",
3280d353602SBarry Smith 				  PCMGTypes[mg[0]->am],levels,(mg[0]->cycles == PC_MG_CYCLE_V) ? "v" : "w",
3290d353602SBarry Smith                                   mg[0]->default_smoothd,mg[0]->default_smoothu);CHKERRQ(ierr);
330c2be2410SBarry Smith     if (mg[0]->galerkin) {
331c2be2410SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
332c2be2410SBarry Smith     }
3334b9ad928SBarry Smith     for (i=0; i<levels; i++) {
334b03c7568SBarry Smith       if (!i) {
335b03c7568SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Coarse gride solver -- level %D -------------------------------\n",i);CHKERRQ(ierr);
336b03c7568SBarry Smith       } else {
33777431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
338b03c7568SBarry Smith       }
3394b9ad928SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
3404b9ad928SBarry Smith       ierr = KSPView(mg[i]->smoothd,viewer);CHKERRQ(ierr);
3414b9ad928SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
342b03c7568SBarry Smith       if (i && mg[i]->smoothd == mg[i]->smoothu) {
3434b9ad928SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
344b03c7568SBarry Smith       } else if (i){
34577431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
3464b9ad928SBarry Smith         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
3474b9ad928SBarry Smith         ierr = KSPView(mg[i]->smoothu,viewer);CHKERRQ(ierr);
3484b9ad928SBarry Smith         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
3494b9ad928SBarry Smith       }
3504b9ad928SBarry Smith     }
3514b9ad928SBarry Smith   } else {
35279a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCMG",((PetscObject)viewer)->type_name);
3534b9ad928SBarry Smith   }
3544b9ad928SBarry Smith   PetscFunctionReturn(0);
3554b9ad928SBarry Smith }
3564b9ad928SBarry Smith 
3574b9ad928SBarry Smith /*
3584b9ad928SBarry Smith     Calls setup for the KSP on each level
3594b9ad928SBarry Smith */
3604b9ad928SBarry Smith #undef __FUNCT__
3614b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_MG"
3626849ba73SBarry Smith static PetscErrorCode PCSetUp_MG(PC pc)
3634b9ad928SBarry Smith {
3649dcbbd2bSBarry Smith   PC_MG                   **mg = (PC_MG**)pc->data;
365dfbe8321SBarry Smith   PetscErrorCode          ierr;
36679416396SBarry Smith   PetscInt                i,n = mg[0]->levels;
36777122347SBarry Smith   PC                      cpc,mpc;
368906ed7ccSBarry Smith   PetscTruth              preonly,lu,redundant,cholesky,monitor = PETSC_FALSE,dump,opsset;
36923d894e5SBarry Smith   PetscViewerASCIIMonitor ascii;
37023d894e5SBarry Smith   PetscViewer             viewer = PETSC_NULL;
3714b9ad928SBarry Smith   MPI_Comm                comm;
372c2be2410SBarry Smith   Mat                     dA,dB;
373c2be2410SBarry Smith   MatStructure            uflag;
3740a6bb862SBarry Smith   Vec                     tvec;
3754b9ad928SBarry Smith 
3764b9ad928SBarry Smith   PetscFunctionBegin;
377852665d3SBarry Smith 
37843fb2f97SBarry Smith   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
37943fb2f97SBarry Smith   /* so use those from global PC */
38043fb2f97SBarry Smith   /* Is this what we always want? What if user wants to keep old one? */
381906ed7ccSBarry Smith   ierr = KSPGetOperatorsSet(mg[n-1]->smoothd,PETSC_NULL,&opsset);CHKERRQ(ierr);
38243fb2f97SBarry Smith   ierr = KSPGetPC(mg[0]->smoothd,&cpc);CHKERRQ(ierr);
38377122347SBarry Smith   ierr = KSPGetPC(mg[n-1]->smoothd,&mpc);CHKERRQ(ierr);
38477122347SBarry Smith   if (!opsset || ((cpc->setupcalled == 1) && (mpc->setupcalled == 2))) {
385852665d3SBarry 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);
3861cfe3bddSBarry Smith     ierr = KSPSetOperators(mg[n-1]->smoothd,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr);
387852665d3SBarry Smith   }
388852665d3SBarry Smith 
389852665d3SBarry Smith   if (mg[0]->galerkin) {
390852665d3SBarry Smith     Mat B;
391852665d3SBarry Smith     mg[0]->galerkinused = PETSC_TRUE;
392852665d3SBarry Smith     /* currently only handle case where mat and pmat are the same on coarser levels */
393852665d3SBarry Smith     ierr = KSPGetOperators(mg[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr);
394852665d3SBarry Smith     if (!pc->setupcalled) {
395852665d3SBarry Smith       for (i=n-2; i>-1; i--) {
396852665d3SBarry Smith         ierr = MatPtAP(dB,mg[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
397852665d3SBarry Smith         ierr = KSPSetOperators(mg[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
398906ed7ccSBarry Smith 	if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
399852665d3SBarry Smith         dB   = B;
400852665d3SBarry Smith       }
401906ed7ccSBarry Smith       ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);
402852665d3SBarry Smith     } else {
403852665d3SBarry Smith       for (i=n-2; i>-1; i--) {
404906ed7ccSBarry Smith         ierr = KSPGetOperators(mg[i]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
405852665d3SBarry Smith         ierr = MatPtAP(dB,mg[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
406852665d3SBarry Smith         ierr = KSPSetOperators(mg[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
407852665d3SBarry Smith         dB   = B;
408852665d3SBarry Smith       }
409852665d3SBarry Smith     }
410852665d3SBarry Smith   }
411852665d3SBarry Smith 
412958c9bccSBarry Smith   if (!pc->setupcalled) {
4134b9ad928SBarry Smith     ierr = PetscOptionsHasName(0,"-pc_mg_monitor",&monitor);CHKERRQ(ierr);
4144b9ad928SBarry Smith 
415b03c7568SBarry Smith     for (i=0; i<n; i++) {
4164b9ad928SBarry Smith       if (monitor) {
4174b9ad928SBarry Smith         ierr = PetscObjectGetComm((PetscObject)mg[i]->smoothd,&comm);CHKERRQ(ierr);
41823d894e5SBarry Smith         ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr);
41923d894e5SBarry Smith         ierr = KSPMonitorSet(mg[i]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
4204b9ad928SBarry Smith       }
4214b9ad928SBarry Smith       ierr = KSPSetFromOptions(mg[i]->smoothd);CHKERRQ(ierr);
4224b9ad928SBarry Smith     }
423b03c7568SBarry Smith     for (i=1; i<n; i++) {
424a98fc643SBarry Smith       if (mg[i]->smoothu && (mg[i]->smoothu != mg[i]->smoothd)) {
4254b9ad928SBarry Smith         if (monitor) {
4264b9ad928SBarry Smith           ierr = PetscObjectGetComm((PetscObject)mg[i]->smoothu,&comm);CHKERRQ(ierr);
42723d894e5SBarry Smith           ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr);
42823d894e5SBarry Smith           ierr = KSPMonitorSet(mg[i]->smoothu,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
4294b9ad928SBarry Smith         }
4304b9ad928SBarry Smith         ierr = KSPSetFromOptions(mg[i]->smoothu);CHKERRQ(ierr);
4314b9ad928SBarry Smith       }
4324b9ad928SBarry Smith     }
433fccaa45eSBarry Smith     for (i=1; i<n; i++) {
4340cace4b0SBarry Smith       if (!mg[i]->residual) {
4350cace4b0SBarry Smith         Mat mat;
4360cace4b0SBarry Smith         ierr = KSPGetOperators(mg[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr);
4370cace4b0SBarry Smith         ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr);
4380cace4b0SBarry Smith       }
439fccaa45eSBarry Smith       if (mg[i]->restrct && !mg[i]->interpolate) {
440*aea2a34eSBarry Smith         ierr = PCMGSetInterpolation(pc,i,mg[i]->restrct);CHKERRQ(ierr);
441fccaa45eSBarry Smith       }
442fccaa45eSBarry Smith       if (!mg[i]->restrct && mg[i]->interpolate) {
44397177400SBarry Smith         ierr = PCMGSetRestriction(pc,i,mg[i]->interpolate);CHKERRQ(ierr);
444fccaa45eSBarry Smith       }
445fccaa45eSBarry Smith #if defined(PETSC_USE_DEBUG)
446fccaa45eSBarry Smith       if (!mg[i]->restrct || !mg[i]->interpolate) {
447fccaa45eSBarry Smith         SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Need to set restriction or interpolation on level %d",(int)i);
448fccaa45eSBarry Smith       }
449fccaa45eSBarry Smith #endif
450fccaa45eSBarry Smith     }
4510a6bb862SBarry Smith     for (i=0; i<n-1; i++) {
45237b0e6c0SBarry Smith       if (!mg[i]->b) {
453906ed7ccSBarry Smith         Vec *vec;
454906ed7ccSBarry Smith         ierr = KSPGetVecs(mg[i]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
455906ed7ccSBarry Smith         ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr);
456906ed7ccSBarry Smith         ierr = PetscFree(vec);CHKERRQ(ierr);
45737b0e6c0SBarry Smith       }
4586ca4d86aSHong Zhang       if (!mg[i]->r && i) {
4590a6bb862SBarry Smith         ierr = VecDuplicate(mg[i]->b,&tvec);CHKERRQ(ierr);
46097177400SBarry Smith         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
4610a6bb862SBarry Smith         ierr = VecDestroy(tvec);CHKERRQ(ierr);
4620a6bb862SBarry Smith       }
4630a6bb862SBarry Smith       if (!mg[i]->x) {
4640a6bb862SBarry Smith         ierr = VecDuplicate(mg[i]->b,&tvec);CHKERRQ(ierr);
46597177400SBarry Smith         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
4660a6bb862SBarry Smith         ierr = VecDestroy(tvec);CHKERRQ(ierr);
4670a6bb862SBarry Smith       }
4680a6bb862SBarry Smith     }
4694b9ad928SBarry Smith   }
4704b9ad928SBarry Smith 
471c2be2410SBarry Smith 
4724b9ad928SBarry Smith   for (i=1; i<n; i++) {
473b03c7568SBarry Smith     if (mg[i]->smoothu == mg[i]->smoothd) {
474b03c7568SBarry Smith       /* if doing only down then initial guess is zero */
4754b9ad928SBarry Smith       ierr = KSPSetInitialGuessNonzero(mg[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
476b03c7568SBarry Smith     }
47732cf1786SBarry Smith     if (mg[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mg[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
4784b9ad928SBarry Smith     ierr = KSPSetUp(mg[i]->smoothd);CHKERRQ(ierr);
47932cf1786SBarry Smith     if (mg[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mg[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
4804b9ad928SBarry Smith   }
481b03c7568SBarry Smith   for (i=1; i<n; i++) {
4824b9ad928SBarry Smith     if (mg[i]->smoothu && mg[i]->smoothu != mg[i]->smoothd) {
483906ed7ccSBarry Smith       Mat          downmat,downpmat;
48497f1f81fSBarry Smith       MatStructure matflag;
485906ed7ccSBarry Smith       PetscTruth   opsset;
48697f1f81fSBarry Smith 
48797f1f81fSBarry Smith       /* check if operators have been set for up, if not use down operators to set them */
488906ed7ccSBarry Smith       ierr = KSPGetOperatorsSet(mg[i]->smoothu,&opsset,PETSC_NULL);CHKERRQ(ierr);
489906ed7ccSBarry Smith       if (!opsset) {
490906ed7ccSBarry Smith         ierr = KSPGetOperators(mg[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr);
49197f1f81fSBarry Smith         ierr = KSPSetOperators(mg[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr);
49297f1f81fSBarry Smith       }
49397f1f81fSBarry Smith 
4944b9ad928SBarry Smith       ierr = KSPSetInitialGuessNonzero(mg[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
49532cf1786SBarry Smith       if (mg[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mg[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
4964b9ad928SBarry Smith       ierr = KSPSetUp(mg[i]->smoothu);CHKERRQ(ierr);
49732cf1786SBarry Smith       if (mg[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mg[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
4984b9ad928SBarry Smith     }
4994b9ad928SBarry Smith   }
5004b9ad928SBarry Smith 
5014b9ad928SBarry Smith   /*
5024b9ad928SBarry Smith       If coarse solver is not direct method then DO NOT USE preonly
5034b9ad928SBarry Smith   */
5044b9ad928SBarry Smith   ierr = PetscTypeCompare((PetscObject)mg[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr);
5054b9ad928SBarry Smith   if (preonly) {
5064b9ad928SBarry Smith     ierr = PetscTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr);
5074b9ad928SBarry Smith     ierr = PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr);
50868eff7e6SBarry Smith     ierr = PetscTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr);
50968eff7e6SBarry Smith     if (!lu && !redundant && !cholesky) {
5104b9ad928SBarry Smith       ierr = KSPSetType(mg[0]->smoothd,KSPGMRES);CHKERRQ(ierr);
5114b9ad928SBarry Smith     }
5124b9ad928SBarry Smith   }
5134b9ad928SBarry Smith 
514958c9bccSBarry Smith   if (!pc->setupcalled) {
5154b9ad928SBarry Smith     if (monitor) {
5164b9ad928SBarry Smith       ierr = PetscObjectGetComm((PetscObject)mg[0]->smoothd,&comm);CHKERRQ(ierr);
51723d894e5SBarry Smith       ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n,&ascii);CHKERRQ(ierr);
51823d894e5SBarry Smith       ierr = KSPMonitorSet(mg[0]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
5194b9ad928SBarry Smith     }
5204b9ad928SBarry Smith     ierr = KSPSetFromOptions(mg[0]->smoothd);CHKERRQ(ierr);
5214b9ad928SBarry Smith   }
5224b9ad928SBarry Smith 
52332cf1786SBarry Smith   if (mg[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mg[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
5244b9ad928SBarry Smith   ierr = KSPSetUp(mg[0]->smoothd);CHKERRQ(ierr);
52532cf1786SBarry Smith   if (mg[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mg[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
5264b9ad928SBarry Smith 
5274b9ad928SBarry Smith   /*
5286805f65bSBarry Smith      Dump the interpolation/restriction matrices plus the
5294b9ad928SBarry Smith    Jacobian/stiffness on each level. This allows Matlab users to
5306805f65bSBarry Smith    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
5316805f65bSBarry Smith 
5326805f65bSBarry Smith    Only support one or the other at the same time.
5336805f65bSBarry Smith   */
5346805f65bSBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
5354b9ad928SBarry Smith   ierr = PetscOptionsHasName(pc->prefix,"-pc_mg_dump_matlab",&dump);CHKERRQ(ierr);
5364b9ad928SBarry Smith   if (dump) {
5376805f65bSBarry Smith     viewer = PETSC_VIEWER_SOCKET_(pc->comm);
5384b9ad928SBarry Smith   }
539c45a1595SBarry Smith #endif
540c2be2410SBarry Smith   ierr = PetscOptionsHasName(pc->prefix,"-pc_mg_dump_binary",&dump);CHKERRQ(ierr);
541c2be2410SBarry Smith   if (dump) {
5426805f65bSBarry Smith     viewer = PETSC_VIEWER_BINARY_(pc->comm);
5436805f65bSBarry Smith   }
5446805f65bSBarry Smith 
5456805f65bSBarry Smith   if (viewer) {
546c2be2410SBarry Smith     for (i=1; i<n; i++) {
5476805f65bSBarry Smith       ierr = MatView(mg[i]->restrct,viewer);CHKERRQ(ierr);
548c2be2410SBarry Smith     }
549c2be2410SBarry Smith     for (i=0; i<n; i++) {
550c2be2410SBarry Smith       ierr = KSPGetPC(mg[i]->smoothd,&pc);CHKERRQ(ierr);
5516805f65bSBarry Smith       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
552c2be2410SBarry Smith     }
553c2be2410SBarry Smith   }
5544b9ad928SBarry Smith   PetscFunctionReturn(0);
5554b9ad928SBarry Smith }
5564b9ad928SBarry Smith 
5574b9ad928SBarry Smith /* -------------------------------------------------------------------------------------*/
5584b9ad928SBarry Smith 
5594b9ad928SBarry Smith #undef __FUNCT__
5609dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetLevels"
5614b9ad928SBarry Smith /*@C
56297177400SBarry Smith    PCMGSetLevels - Sets the number of levels to use with MG.
5634b9ad928SBarry Smith    Must be called before any other MG routine.
5644b9ad928SBarry Smith 
5654b9ad928SBarry Smith    Collective on PC
5664b9ad928SBarry Smith 
5674b9ad928SBarry Smith    Input Parameters:
5684b9ad928SBarry Smith +  pc - the preconditioner context
5694b9ad928SBarry Smith .  levels - the number of levels
5704b9ad928SBarry Smith -  comms - optional communicators for each level; this is to allow solving the coarser problems
5714b9ad928SBarry Smith            on smaller sets of processors. Use PETSC_NULL_OBJECT for default in Fortran
5724b9ad928SBarry Smith 
5734b9ad928SBarry Smith    Level: intermediate
5744b9ad928SBarry Smith 
5754b9ad928SBarry Smith    Notes:
5764b9ad928SBarry Smith      If the number of levels is one then the multigrid uses the -mg_levels prefix
5774b9ad928SBarry Smith   for setting the level options rather than the -mg_coarse prefix.
5784b9ad928SBarry Smith 
5794b9ad928SBarry Smith .keywords: MG, set, levels, multigrid
5804b9ad928SBarry Smith 
58197177400SBarry Smith .seealso: PCMGSetType(), PCMGGetLevels()
5824b9ad928SBarry Smith @*/
58397177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
5844b9ad928SBarry Smith {
585dfbe8321SBarry Smith   PetscErrorCode ierr;
586ada7143aSSatish Balay   PC_MG          **mg=0;
5874b9ad928SBarry Smith 
5884b9ad928SBarry Smith   PetscFunctionBegin;
5894482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
5904b9ad928SBarry Smith 
5914b9ad928SBarry Smith   if (pc->data) {
5921302d50aSBarry Smith     SETERRQ(PETSC_ERR_ORDER,"Number levels already set for MG\n\
59397177400SBarry Smith     make sure that you call PCMGSetLevels() before KSPSetFromOptions()");
5944b9ad928SBarry Smith   }
5959dcbbd2bSBarry Smith   ierr                     = PCMGCreate_Private(pc->comm,levels,pc,comms,&mg);CHKERRQ(ierr);
5969dcbbd2bSBarry Smith   mg[0]->am                = PC_MG_MULTIPLICATIVE;
5974b9ad928SBarry Smith   pc->data                 = (void*)mg;
5984b9ad928SBarry Smith   pc->ops->applyrichardson = PCApplyRichardson_MG;
5994b9ad928SBarry Smith   PetscFunctionReturn(0);
6004b9ad928SBarry Smith }
6014b9ad928SBarry Smith 
6024b9ad928SBarry Smith #undef __FUNCT__
6039dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetLevels"
6044b9ad928SBarry Smith /*@
60597177400SBarry Smith    PCMGGetLevels - Gets the number of levels to use with MG.
6064b9ad928SBarry Smith 
6074b9ad928SBarry Smith    Not Collective
6084b9ad928SBarry Smith 
6094b9ad928SBarry Smith    Input Parameter:
6104b9ad928SBarry Smith .  pc - the preconditioner context
6114b9ad928SBarry Smith 
6124b9ad928SBarry Smith    Output parameter:
6134b9ad928SBarry Smith .  levels - the number of levels
6144b9ad928SBarry Smith 
6154b9ad928SBarry Smith    Level: advanced
6164b9ad928SBarry Smith 
6174b9ad928SBarry Smith .keywords: MG, get, levels, multigrid
6184b9ad928SBarry Smith 
61997177400SBarry Smith .seealso: PCMGSetLevels()
6204b9ad928SBarry Smith @*/
62197177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetLevels(PC pc,PetscInt *levels)
6224b9ad928SBarry Smith {
6239dcbbd2bSBarry Smith   PC_MG  **mg;
6244b9ad928SBarry Smith 
6254b9ad928SBarry Smith   PetscFunctionBegin;
6264482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
6274482741eSBarry Smith   PetscValidIntPointer(levels,2);
6284b9ad928SBarry Smith 
6299dcbbd2bSBarry Smith   mg      = (PC_MG**)pc->data;
6304b9ad928SBarry Smith   *levels = mg[0]->levels;
6314b9ad928SBarry Smith   PetscFunctionReturn(0);
6324b9ad928SBarry Smith }
6334b9ad928SBarry Smith 
6344b9ad928SBarry Smith #undef __FUNCT__
6359dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetType"
6364b9ad928SBarry Smith /*@
63797177400SBarry Smith    PCMGSetType - Determines the form of multigrid to use:
6384b9ad928SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
6394b9ad928SBarry Smith 
6404b9ad928SBarry Smith    Collective on PC
6414b9ad928SBarry Smith 
6424b9ad928SBarry Smith    Input Parameters:
6434b9ad928SBarry Smith +  pc - the preconditioner context
6449dcbbd2bSBarry Smith -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
6459dcbbd2bSBarry Smith    PC_MG_FULL, PC_MG_KASKADE
6464b9ad928SBarry Smith 
6474b9ad928SBarry Smith    Options Database Key:
6484b9ad928SBarry Smith .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
6494b9ad928SBarry Smith    additive, full, kaskade
6504b9ad928SBarry Smith 
6514b9ad928SBarry Smith    Level: advanced
6524b9ad928SBarry Smith 
6534b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
6544b9ad928SBarry Smith 
65597177400SBarry Smith .seealso: PCMGSetLevels()
6564b9ad928SBarry Smith @*/
6579dcbbd2bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetType(PC pc,PCMGType form)
6584b9ad928SBarry Smith {
6599dcbbd2bSBarry Smith   PC_MG **mg;
6604b9ad928SBarry Smith 
6614b9ad928SBarry Smith   PetscFunctionBegin;
6624482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
6639dcbbd2bSBarry Smith   mg = (PC_MG**)pc->data;
6644b9ad928SBarry Smith 
6654b9ad928SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
6664b9ad928SBarry Smith   mg[0]->am = form;
6679dcbbd2bSBarry Smith   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
6684b9ad928SBarry Smith   else pc->ops->applyrichardson = 0;
6694b9ad928SBarry Smith   PetscFunctionReturn(0);
6704b9ad928SBarry Smith }
6714b9ad928SBarry Smith 
6724b9ad928SBarry Smith #undef __FUNCT__
6730d353602SBarry Smith #define __FUNCT__ "PCMGSetCycleType"
6744b9ad928SBarry Smith /*@
6750d353602SBarry Smith    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
6764b9ad928SBarry Smith    complicated cycling.
6774b9ad928SBarry Smith 
6784b9ad928SBarry Smith    Collective on PC
6794b9ad928SBarry Smith 
6804b9ad928SBarry Smith    Input Parameters:
681c2be2410SBarry Smith +  pc - the multigrid context
6820d353602SBarry Smith -  PC_MG_CYCLE_V or PC_MG_CYCLE_W
6834b9ad928SBarry Smith 
6844b9ad928SBarry Smith    Options Database Key:
6850d353602SBarry Smith $  -pc_mg_cycle_type v or w
6864b9ad928SBarry Smith 
6874b9ad928SBarry Smith    Level: advanced
6884b9ad928SBarry Smith 
6894b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
6904b9ad928SBarry Smith 
6910d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel()
6924b9ad928SBarry Smith @*/
6930d353602SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetCycleType(PC pc,PCMGCycleType n)
6944b9ad928SBarry Smith {
6959dcbbd2bSBarry Smith   PC_MG    **mg;
69679416396SBarry Smith   PetscInt i,levels;
6974b9ad928SBarry Smith 
6984b9ad928SBarry Smith   PetscFunctionBegin;
6994482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
7009dcbbd2bSBarry Smith   mg     = (PC_MG**)pc->data;
7014b9ad928SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
7024b9ad928SBarry Smith   levels = mg[0]->levels;
7034b9ad928SBarry Smith 
7044b9ad928SBarry Smith   for (i=0; i<levels; i++) {
7054b9ad928SBarry Smith     mg[i]->cycles  = n;
7064b9ad928SBarry Smith   }
7074b9ad928SBarry Smith   PetscFunctionReturn(0);
7084b9ad928SBarry Smith }
7094b9ad928SBarry Smith 
7104b9ad928SBarry Smith #undef __FUNCT__
7118cc2d5dfSBarry Smith #define __FUNCT__ "PCMGMultiplicativeSetCycles"
7128cc2d5dfSBarry Smith /*@
7138cc2d5dfSBarry Smith    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
7148cc2d5dfSBarry Smith          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
7158cc2d5dfSBarry Smith 
7168cc2d5dfSBarry Smith    Collective on PC
7178cc2d5dfSBarry Smith 
7188cc2d5dfSBarry Smith    Input Parameters:
7198cc2d5dfSBarry Smith +  pc - the multigrid context
7208cc2d5dfSBarry Smith -  n - number of cycles (default is 1)
7218cc2d5dfSBarry Smith 
7228cc2d5dfSBarry Smith    Options Database Key:
7238cc2d5dfSBarry Smith $  -pc_mg_multiplicative_cycles n
7248cc2d5dfSBarry Smith 
7258cc2d5dfSBarry Smith    Level: advanced
7268cc2d5dfSBarry Smith 
7278cc2d5dfSBarry Smith    Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
7288cc2d5dfSBarry Smith 
7298cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
7308cc2d5dfSBarry Smith 
7318cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
7328cc2d5dfSBarry Smith @*/
7338cc2d5dfSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
7348cc2d5dfSBarry Smith {
7358cc2d5dfSBarry Smith   PC_MG    **mg;
7368cc2d5dfSBarry Smith   PetscInt i,levels;
7378cc2d5dfSBarry Smith 
7388cc2d5dfSBarry Smith   PetscFunctionBegin;
7398cc2d5dfSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
7408cc2d5dfSBarry Smith   mg     = (PC_MG**)pc->data;
7418cc2d5dfSBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
7428cc2d5dfSBarry Smith   levels = mg[0]->levels;
7438cc2d5dfSBarry Smith 
7448cc2d5dfSBarry Smith   for (i=0; i<levels; i++) {
7458cc2d5dfSBarry Smith     mg[i]->cyclesperpcapply  = n;
7468cc2d5dfSBarry Smith   }
7478cc2d5dfSBarry Smith   PetscFunctionReturn(0);
7488cc2d5dfSBarry Smith }
7498cc2d5dfSBarry Smith 
7508cc2d5dfSBarry Smith #undef __FUNCT__
7519dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetGalerkin"
752c2be2410SBarry Smith /*@
75397177400SBarry Smith    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
754c2be2410SBarry Smith       finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t
755c2be2410SBarry Smith 
756c2be2410SBarry Smith    Collective on PC
757c2be2410SBarry Smith 
758c2be2410SBarry Smith    Input Parameters:
7593fc8bf9cSBarry Smith .  pc - the multigrid context
760c2be2410SBarry Smith 
761c2be2410SBarry Smith    Options Database Key:
762c2be2410SBarry Smith $  -pc_mg_galerkin
763c2be2410SBarry Smith 
764c2be2410SBarry Smith    Level: intermediate
765c2be2410SBarry Smith 
766c2be2410SBarry Smith .keywords: MG, set, Galerkin
767c2be2410SBarry Smith 
7683fc8bf9cSBarry Smith .seealso: PCMGGetGalerkin()
7693fc8bf9cSBarry Smith 
770c2be2410SBarry Smith @*/
77197177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetGalerkin(PC pc)
772c2be2410SBarry Smith {
7739dcbbd2bSBarry Smith   PC_MG    **mg;
774c2be2410SBarry Smith   PetscInt i,levels;
775c2be2410SBarry Smith 
776c2be2410SBarry Smith   PetscFunctionBegin;
777c2be2410SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
7789dcbbd2bSBarry Smith   mg     = (PC_MG**)pc->data;
779c2be2410SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
780c2be2410SBarry Smith   levels = mg[0]->levels;
781c2be2410SBarry Smith 
782c2be2410SBarry Smith   for (i=0; i<levels; i++) {
783c2be2410SBarry Smith     mg[i]->galerkin = PETSC_TRUE;
784c2be2410SBarry Smith   }
785c2be2410SBarry Smith   PetscFunctionReturn(0);
786c2be2410SBarry Smith }
787c2be2410SBarry Smith 
788c2be2410SBarry Smith #undef __FUNCT__
7893fc8bf9cSBarry Smith #define __FUNCT__ "PCMGGetGalerkin"
7903fc8bf9cSBarry Smith /*@
7913fc8bf9cSBarry Smith    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
7923fc8bf9cSBarry Smith       A_i-1 = r_i * A_i * r_i^t
7933fc8bf9cSBarry Smith 
7943fc8bf9cSBarry Smith    Not Collective
7953fc8bf9cSBarry Smith 
7963fc8bf9cSBarry Smith    Input Parameter:
7973fc8bf9cSBarry Smith .  pc - the multigrid context
7983fc8bf9cSBarry Smith 
7993fc8bf9cSBarry Smith    Output Parameter:
8003fc8bf9cSBarry Smith .  gelerkin - PETSC_TRUE or PETSC_FALSE
8013fc8bf9cSBarry Smith 
8023fc8bf9cSBarry Smith    Options Database Key:
8033fc8bf9cSBarry Smith $  -pc_mg_galerkin
8043fc8bf9cSBarry Smith 
8053fc8bf9cSBarry Smith    Level: intermediate
8063fc8bf9cSBarry Smith 
8073fc8bf9cSBarry Smith .keywords: MG, set, Galerkin
8083fc8bf9cSBarry Smith 
8093fc8bf9cSBarry Smith .seealso: PCMGSetGalerkin()
8103fc8bf9cSBarry Smith 
8113fc8bf9cSBarry Smith @*/
8123fc8bf9cSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetGalerkin(PC pc,PetscTruth *galerkin)
8133fc8bf9cSBarry Smith {
8143fc8bf9cSBarry Smith   PC_MG    **mg;
8153fc8bf9cSBarry Smith 
8163fc8bf9cSBarry Smith   PetscFunctionBegin;
8173fc8bf9cSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
8183fc8bf9cSBarry Smith   mg     = (PC_MG**)pc->data;
8193fc8bf9cSBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
8203fc8bf9cSBarry Smith   *galerkin = mg[0]->galerkin;
8213fc8bf9cSBarry Smith   PetscFunctionReturn(0);
8223fc8bf9cSBarry Smith }
8233fc8bf9cSBarry Smith 
8243fc8bf9cSBarry Smith #undef __FUNCT__
8259dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothDown"
8264b9ad928SBarry Smith /*@
82797177400SBarry Smith    PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
82897177400SBarry Smith    use on all levels. Use PCMGGetSmootherDown() to set different
8294b9ad928SBarry Smith    pre-smoothing steps on different levels.
8304b9ad928SBarry Smith 
8314b9ad928SBarry Smith    Collective on PC
8324b9ad928SBarry Smith 
8334b9ad928SBarry Smith    Input Parameters:
8344b9ad928SBarry Smith +  mg - the multigrid context
8354b9ad928SBarry Smith -  n - the number of smoothing steps
8364b9ad928SBarry Smith 
8374b9ad928SBarry Smith    Options Database Key:
8384b9ad928SBarry Smith .  -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps
8394b9ad928SBarry Smith 
8404b9ad928SBarry Smith    Level: advanced
8414b9ad928SBarry Smith 
8424b9ad928SBarry Smith .keywords: MG, smooth, down, pre-smoothing, steps, multigrid
8434b9ad928SBarry Smith 
84497177400SBarry Smith .seealso: PCMGSetNumberSmoothUp()
8454b9ad928SBarry Smith @*/
84697177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothDown(PC pc,PetscInt n)
8474b9ad928SBarry Smith {
8489dcbbd2bSBarry Smith   PC_MG          **mg;
8496849ba73SBarry Smith   PetscErrorCode ierr;
85079416396SBarry Smith   PetscInt       i,levels;
8514b9ad928SBarry Smith 
8524b9ad928SBarry Smith   PetscFunctionBegin;
8534482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
8549dcbbd2bSBarry Smith   mg     = (PC_MG**)pc->data;
8554b9ad928SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
8564b9ad928SBarry Smith   levels = mg[0]->levels;
8574b9ad928SBarry Smith 
858b05257ddSBarry Smith   for (i=1; i<levels; i++) {
8594b9ad928SBarry Smith     /* make sure smoother up and down are different */
86097177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
8614b9ad928SBarry Smith     ierr = KSPSetTolerances(mg[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
8624b9ad928SBarry Smith     mg[i]->default_smoothd = n;
8634b9ad928SBarry Smith   }
8644b9ad928SBarry Smith   PetscFunctionReturn(0);
8654b9ad928SBarry Smith }
8664b9ad928SBarry Smith 
8674b9ad928SBarry Smith #undef __FUNCT__
8689dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothUp"
8694b9ad928SBarry Smith /*@
87097177400SBarry Smith    PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
87197177400SBarry Smith    on all levels. Use PCMGGetSmootherUp() to set different numbers of
8724b9ad928SBarry Smith    post-smoothing steps on different levels.
8734b9ad928SBarry Smith 
8744b9ad928SBarry Smith    Collective on PC
8754b9ad928SBarry Smith 
8764b9ad928SBarry Smith    Input Parameters:
8774b9ad928SBarry Smith +  mg - the multigrid context
8784b9ad928SBarry Smith -  n - the number of smoothing steps
8794b9ad928SBarry Smith 
8804b9ad928SBarry Smith    Options Database Key:
8814b9ad928SBarry Smith .  -pc_mg_smoothup <n> - Sets number of post-smoothing steps
8824b9ad928SBarry Smith 
8834b9ad928SBarry Smith    Level: advanced
8844b9ad928SBarry Smith 
8854b9ad928SBarry Smith    Note: this does not set a value on the coarsest grid, since we assume that
886a8c7a070SBarry Smith     there is no separate smooth up on the coarsest grid.
8874b9ad928SBarry Smith 
8884b9ad928SBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid
8894b9ad928SBarry Smith 
89097177400SBarry Smith .seealso: PCMGSetNumberSmoothDown()
8914b9ad928SBarry Smith @*/
89297177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothUp(PC pc,PetscInt n)
8934b9ad928SBarry Smith {
8949dcbbd2bSBarry Smith   PC_MG          **mg;
8956849ba73SBarry Smith   PetscErrorCode ierr;
89679416396SBarry Smith   PetscInt       i,levels;
8974b9ad928SBarry Smith 
8984b9ad928SBarry Smith   PetscFunctionBegin;
8994482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
9009dcbbd2bSBarry Smith   mg     = (PC_MG**)pc->data;
9014b9ad928SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
9024b9ad928SBarry Smith   levels = mg[0]->levels;
9034b9ad928SBarry Smith 
9044b9ad928SBarry Smith   for (i=1; i<levels; i++) {
9054b9ad928SBarry Smith     /* make sure smoother up and down are different */
90697177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
9074b9ad928SBarry Smith     ierr = KSPSetTolerances(mg[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
9084b9ad928SBarry Smith     mg[i]->default_smoothu = n;
9094b9ad928SBarry Smith   }
9104b9ad928SBarry Smith   PetscFunctionReturn(0);
9114b9ad928SBarry Smith }
9124b9ad928SBarry Smith 
9134b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/
9144b9ad928SBarry Smith 
9153b09bd56SBarry Smith /*MC
916ccb205f8SBarry Smith    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
9173b09bd56SBarry Smith     information about the coarser grid matrices and restriction/interpolation operators.
9183b09bd56SBarry Smith 
9193b09bd56SBarry Smith    Options Database Keys:
9203b09bd56SBarry Smith +  -pc_mg_levels <nlevels> - number of levels including finest
9210d353602SBarry Smith .  -pc_mg_cycles v or w
92279416396SBarry Smith .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
9233b09bd56SBarry Smith .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
9243b09bd56SBarry Smith .  -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default
9253b09bd56SBarry Smith .  -pc_mg_log - log information about time spent on each level of the solver
9263b09bd56SBarry Smith .  -pc_mg_monitor - print information on the multigrid convergence
92768eff7e6SBarry Smith .  -pc_mg_galerkin - use Galerkin process to compute coarser operators
9283b09bd56SBarry Smith -  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
9293b09bd56SBarry Smith                         to the Socket viewer for reading from Matlab.
9303b09bd56SBarry Smith 
9313b09bd56SBarry Smith    Notes:
9323b09bd56SBarry Smith 
9333b09bd56SBarry Smith    Level: intermediate
9343b09bd56SBarry Smith 
9353b09bd56SBarry Smith    Concepts: multigrid
9363b09bd56SBarry Smith 
9373b09bd56SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
9380d353602SBarry Smith            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
93997177400SBarry Smith            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
94097177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
9410d353602SBarry Smith            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
9423b09bd56SBarry Smith M*/
9433b09bd56SBarry Smith 
9444b9ad928SBarry Smith EXTERN_C_BEGIN
9454b9ad928SBarry Smith #undef __FUNCT__
9464b9ad928SBarry Smith #define __FUNCT__ "PCCreate_MG"
947dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_MG(PC pc)
9484b9ad928SBarry Smith {
9494b9ad928SBarry Smith   PetscFunctionBegin;
9504b9ad928SBarry Smith   pc->ops->apply          = PCApply_MG;
9514b9ad928SBarry Smith   pc->ops->setup          = PCSetUp_MG;
9524b9ad928SBarry Smith   pc->ops->destroy        = PCDestroy_MG;
9534b9ad928SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_MG;
9544b9ad928SBarry Smith   pc->ops->view           = PCView_MG;
9554b9ad928SBarry Smith 
9564b9ad928SBarry Smith   pc->data                = (void*)0;
9574b9ad928SBarry Smith   PetscFunctionReturn(0);
9584b9ad928SBarry Smith }
9594b9ad928SBarry Smith EXTERN_C_END
960