xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision 906ed7cc33fecbafab01746eec64dcdcc8a4842f)
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;
1579416396SBarry Smith   PetscInt       cycles = mg->cycles;
164b9ad928SBarry Smith 
174b9ad928SBarry Smith   PetscFunctionBegin;
184b9ad928SBarry Smith   if (converged) *converged = PETSC_FALSE;
194b9ad928SBarry Smith 
204b9ad928SBarry Smith   if (mg->eventsolve) {ierr = PetscLogEventBegin(mg->eventsolve,0,0,0,0);CHKERRQ(ierr);}
2123ce1328SBarry Smith   ierr = KSPSolve(mg->smoothd,mg->b,mg->x);CHKERRQ(ierr);
224b9ad928SBarry Smith   if (mg->eventsolve) {ierr = PetscLogEventEnd(mg->eventsolve,0,0,0,0);CHKERRQ(ierr);}
234b9ad928SBarry Smith   if (mg->level) {  /* not the coarsest grid */
244b9ad928SBarry Smith     ierr = (*mg->residual)(mg->A,mg->b,mg->x,mg->r);CHKERRQ(ierr);
254b9ad928SBarry Smith 
264b9ad928SBarry Smith     /* if on finest level and have convergence criteria set */
274b9ad928SBarry Smith     if (mg->level == mg->levels-1 && mg->ttol) {
284b9ad928SBarry Smith       PetscReal rnorm;
294b9ad928SBarry Smith       ierr = VecNorm(mg->r,NORM_2,&rnorm);CHKERRQ(ierr);
304b9ad928SBarry Smith       if (rnorm <= mg->ttol) {
314b9ad928SBarry Smith         *converged = PETSC_TRUE;
3270441072SBarry Smith         if (rnorm < mg->abstol) {
33ae15b995SBarry Smith           ierr = PetscInfo2(0,"Linear solver has converged. Residual norm %G is less than absolute tolerance %G\n",rnorm,mg->abstol);CHKERRQ(ierr);
344b9ad928SBarry Smith         } else {
35ae15b995SBarry 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);
364b9ad928SBarry Smith         }
374b9ad928SBarry Smith         PetscFunctionReturn(0);
384b9ad928SBarry Smith       }
394b9ad928SBarry Smith     }
404b9ad928SBarry Smith 
414b9ad928SBarry Smith     mgc = *(mglevels - 1);
424b9ad928SBarry Smith     ierr = MatRestrict(mg->restrct,mg->r,mgc->b);CHKERRQ(ierr);
43efb30889SBarry Smith     ierr = VecSet(mgc->x,0.0);CHKERRQ(ierr);
444b9ad928SBarry Smith     while (cycles--) {
459dcbbd2bSBarry Smith       ierr = PCMGMCycle_Private(mglevels-1,converged);CHKERRQ(ierr);
464b9ad928SBarry Smith     }
474b9ad928SBarry Smith     ierr = MatInterpolateAdd(mg->interpolate,mgc->x,mg->x,mg->x);CHKERRQ(ierr);
484b9ad928SBarry Smith     if (mg->eventsolve) {ierr = PetscLogEventBegin(mg->eventsolve,0,0,0,0);CHKERRQ(ierr);}
4923ce1328SBarry Smith     ierr = KSPSolve(mg->smoothu,mg->b,mg->x);CHKERRQ(ierr);
504b9ad928SBarry Smith     if (mg->eventsolve) {ierr = PetscLogEventEnd(mg->eventsolve,0,0,0,0);CHKERRQ(ierr);}
514b9ad928SBarry Smith   }
524b9ad928SBarry Smith   PetscFunctionReturn(0);
534b9ad928SBarry Smith }
544b9ad928SBarry Smith 
554b9ad928SBarry Smith /*
569dcbbd2bSBarry Smith        PCMGCreate_Private - Creates a PC_MG structure for use with the
574b9ad928SBarry Smith                multigrid code. Level 0 is the coarsest. (But the
584b9ad928SBarry Smith                finest level is stored first in the array).
594b9ad928SBarry Smith 
604b9ad928SBarry Smith */
614b9ad928SBarry Smith #undef __FUNCT__
629dcbbd2bSBarry Smith #define __FUNCT__ "PCMGCreate_Private"
639dcbbd2bSBarry Smith static PetscErrorCode PCMGCreate_Private(MPI_Comm comm,PetscInt levels,PC pc,MPI_Comm *comms,PC_MG ***result)
644b9ad928SBarry Smith {
659dcbbd2bSBarry Smith   PC_MG          **mg;
666849ba73SBarry Smith   PetscErrorCode ierr;
6779416396SBarry Smith   PetscInt       i;
6879416396SBarry Smith   PetscMPIInt    size;
69f69a0ea3SMatthew Knepley   const char     *prefix;
704b9ad928SBarry Smith   PC             ipc;
714b9ad928SBarry Smith 
724b9ad928SBarry Smith   PetscFunctionBegin;
739dcbbd2bSBarry Smith   ierr = PetscMalloc(levels*sizeof(PC_MG*),&mg);CHKERRQ(ierr);
749dcbbd2bSBarry Smith   ierr = PetscLogObjectMemory(pc,levels*(sizeof(PC_MG*)+sizeof(PC_MG)));CHKERRQ(ierr);
754b9ad928SBarry Smith 
764b9ad928SBarry Smith   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
774b9ad928SBarry Smith 
784b9ad928SBarry Smith   for (i=0; i<levels; i++) {
799dcbbd2bSBarry Smith     ierr = PetscNew(PC_MG,&mg[i]);CHKERRQ(ierr);
804b9ad928SBarry Smith     mg[i]->level           = i;
814b9ad928SBarry Smith     mg[i]->levels          = levels;
824b9ad928SBarry Smith     mg[i]->cycles          = 1;
83c2be2410SBarry Smith     mg[i]->galerkin        = PETSC_FALSE;
84c2be2410SBarry Smith     mg[i]->galerkinused    = PETSC_FALSE;
8579416396SBarry Smith     mg[i]->default_smoothu = 1;
8679416396SBarry Smith     mg[i]->default_smoothd = 1;
874b9ad928SBarry Smith 
884b9ad928SBarry Smith     if (comms) comm = comms[i];
894b9ad928SBarry Smith     ierr = KSPCreate(comm,&mg[i]->smoothd);CHKERRQ(ierr);
9079416396SBarry Smith     ierr = KSPSetTolerances(mg[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg[i]->default_smoothd);CHKERRQ(ierr);
914b9ad928SBarry Smith     ierr = KSPSetOptionsPrefix(mg[i]->smoothd,prefix);CHKERRQ(ierr);
924b9ad928SBarry Smith 
934b9ad928SBarry Smith     /* do special stuff for coarse grid */
944b9ad928SBarry Smith     if (!i && levels > 1) {
954b9ad928SBarry Smith       ierr = KSPAppendOptionsPrefix(mg[0]->smoothd,"mg_coarse_");CHKERRQ(ierr);
964b9ad928SBarry Smith 
974b9ad928SBarry Smith       /* coarse solve is (redundant) LU by default */
984b9ad928SBarry Smith       ierr = KSPSetType(mg[0]->smoothd,KSPPREONLY);CHKERRQ(ierr);
994b9ad928SBarry Smith       ierr = KSPGetPC(mg[0]->smoothd,&ipc);CHKERRQ(ierr);
1004b9ad928SBarry Smith       ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1014b9ad928SBarry Smith       if (size > 1) {
1024b9ad928SBarry Smith         ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr);
1034b9ad928SBarry Smith         ierr = PCRedundantGetPC(ipc,&ipc);CHKERRQ(ierr);
1044b9ad928SBarry Smith       }
1054b9ad928SBarry Smith       ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr);
1064b9ad928SBarry Smith 
1074b9ad928SBarry Smith     } else {
1082e70eadcSBarry Smith       char tprefix[128];
10913f74950SBarry Smith       sprintf(tprefix,"mg_levels_%d_",(int)i);
1102e70eadcSBarry Smith       ierr = KSPAppendOptionsPrefix(mg[i]->smoothd,tprefix);CHKERRQ(ierr);
1114b9ad928SBarry Smith     }
11252e6d16bSBarry Smith     ierr = PetscLogObjectParent(pc,mg[i]->smoothd);CHKERRQ(ierr);
1134b9ad928SBarry Smith     mg[i]->smoothu         = mg[i]->smoothd;
1144b9ad928SBarry Smith     mg[i]->rtol = 0.0;
11570441072SBarry Smith     mg[i]->abstol = 0.0;
1164b9ad928SBarry Smith     mg[i]->dtol = 0.0;
1174b9ad928SBarry Smith     mg[i]->ttol = 0.0;
1184b9ad928SBarry Smith     mg[i]->eventsetup = 0;
1194b9ad928SBarry Smith     mg[i]->eventsolve = 0;
1204b9ad928SBarry Smith   }
1214b9ad928SBarry Smith   *result = mg;
1224b9ad928SBarry Smith   PetscFunctionReturn(0);
1234b9ad928SBarry Smith }
1244b9ad928SBarry Smith 
1254b9ad928SBarry Smith #undef __FUNCT__
1264b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_MG"
1276849ba73SBarry Smith static PetscErrorCode PCDestroy_MG(PC pc)
1284b9ad928SBarry Smith {
1299dcbbd2bSBarry Smith   PC_MG          **mg = (PC_MG**)pc->data;
1306849ba73SBarry Smith   PetscErrorCode ierr;
13179416396SBarry Smith   PetscInt       i,n = mg[0]->levels;
1324b9ad928SBarry Smith 
1334b9ad928SBarry Smith   PetscFunctionBegin;
134fccaa45eSBarry Smith   for (i=0; i<n-1; i++) {
135630ba2eeSBarry Smith     if (mg[i+1]->r) {ierr = VecDestroy(mg[i+1]->r);CHKERRQ(ierr);}
136fccaa45eSBarry Smith     if (mg[i]->b) {ierr = VecDestroy(mg[i]->b);CHKERRQ(ierr);}
137fccaa45eSBarry Smith     if (mg[i]->x) {ierr = VecDestroy(mg[i]->x);CHKERRQ(ierr);}
138fccaa45eSBarry Smith     if (mg[i+1]->restrct) {ierr = MatDestroy(mg[i+1]->restrct);CHKERRQ(ierr);}
139fccaa45eSBarry Smith     if (mg[i+1]->interpolate) {ierr = MatDestroy(mg[i+1]->interpolate);CHKERRQ(ierr);}
140fccaa45eSBarry Smith   }
141fccaa45eSBarry Smith 
1424b9ad928SBarry Smith   for (i=0; i<n; i++) {
1434b9ad928SBarry Smith     if (mg[i]->smoothd != mg[i]->smoothu) {
1444b9ad928SBarry Smith       ierr = KSPDestroy(mg[i]->smoothd);CHKERRQ(ierr);
1454b9ad928SBarry Smith     }
1464b9ad928SBarry Smith     ierr = KSPDestroy(mg[i]->smoothu);CHKERRQ(ierr);
1474b9ad928SBarry Smith     ierr = PetscFree(mg[i]);CHKERRQ(ierr);
1484b9ad928SBarry Smith   }
1494b9ad928SBarry Smith   ierr = PetscFree(mg);CHKERRQ(ierr);
1504b9ad928SBarry Smith   PetscFunctionReturn(0);
1514b9ad928SBarry Smith }
1524b9ad928SBarry Smith 
1534b9ad928SBarry Smith 
1544b9ad928SBarry Smith 
1559dcbbd2bSBarry Smith EXTERN PetscErrorCode PCMGACycle_Private(PC_MG**);
1569dcbbd2bSBarry Smith EXTERN PetscErrorCode PCMGFCycle_Private(PC_MG**);
1579dcbbd2bSBarry Smith EXTERN PetscErrorCode PCMGKCycle_Private(PC_MG**);
1584b9ad928SBarry Smith 
1594b9ad928SBarry Smith /*
1604b9ad928SBarry Smith    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
1614b9ad928SBarry Smith              or full cycle of multigrid.
1624b9ad928SBarry Smith 
1634b9ad928SBarry Smith   Note:
1649dcbbd2bSBarry Smith   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
1654b9ad928SBarry Smith */
1664b9ad928SBarry Smith #undef __FUNCT__
1674b9ad928SBarry Smith #define __FUNCT__ "PCApply_MG"
1686849ba73SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
1694b9ad928SBarry Smith {
1709dcbbd2bSBarry Smith   PC_MG          **mg = (PC_MG**)pc->data;
1716849ba73SBarry Smith   PetscErrorCode ierr;
17279416396SBarry Smith   PetscInt       levels = mg[0]->levels;
1734b9ad928SBarry Smith 
1744b9ad928SBarry Smith   PetscFunctionBegin;
1754b9ad928SBarry Smith   mg[levels-1]->b = b;
1764b9ad928SBarry Smith   mg[levels-1]->x = x;
177ef70c39aSMatthew Knepley   if (!mg[levels-1]->r && mg[0]->am != PC_MG_ADDITIVE && levels > 1) {
1780a6bb862SBarry Smith     Vec tvec;
1790a6bb862SBarry Smith     ierr = VecDuplicate(mg[levels-1]->b,&tvec);CHKERRQ(ierr);
18097177400SBarry Smith     ierr = PCMGSetR(pc,levels-1,tvec);CHKERRQ(ierr);
1810a6bb862SBarry Smith     ierr = VecDestroy(tvec);CHKERRQ(ierr);
1820a6bb862SBarry Smith   }
1839dcbbd2bSBarry Smith   if (mg[0]->am == PC_MG_MULTIPLICATIVE) {
184efb30889SBarry Smith     ierr = VecSet(x,0.0);CHKERRQ(ierr);
1859dcbbd2bSBarry Smith     ierr = PCMGMCycle_Private(mg+levels-1,PETSC_NULL);CHKERRQ(ierr);
1864b9ad928SBarry Smith   }
1879dcbbd2bSBarry Smith   else if (mg[0]->am == PC_MG_ADDITIVE) {
1889dcbbd2bSBarry Smith     ierr = PCMGACycle_Private(mg);CHKERRQ(ierr);
1894b9ad928SBarry Smith   }
1909dcbbd2bSBarry Smith   else if (mg[0]->am == PC_MG_KASKADE) {
1919dcbbd2bSBarry Smith     ierr = PCMGKCycle_Private(mg);CHKERRQ(ierr);
1924b9ad928SBarry Smith   }
1934b9ad928SBarry Smith   else {
1949dcbbd2bSBarry Smith     ierr = PCMGFCycle_Private(mg);CHKERRQ(ierr);
1954b9ad928SBarry Smith   }
1964b9ad928SBarry Smith   PetscFunctionReturn(0);
1974b9ad928SBarry Smith }
1984b9ad928SBarry Smith 
1994b9ad928SBarry Smith #undef __FUNCT__
2004b9ad928SBarry Smith #define __FUNCT__ "PCApplyRichardson_MG"
20179416396SBarry Smith static PetscErrorCode PCApplyRichardson_MG(PC pc,Vec b,Vec x,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its)
2024b9ad928SBarry Smith {
2039dcbbd2bSBarry Smith   PC_MG          **mg = (PC_MG**)pc->data;
204dfbe8321SBarry Smith   PetscErrorCode ierr;
20579416396SBarry Smith   PetscInt       levels = mg[0]->levels;
2064b9ad928SBarry Smith   PetscTruth     converged = PETSC_FALSE;
2074b9ad928SBarry Smith 
2084b9ad928SBarry Smith   PetscFunctionBegin;
2094b9ad928SBarry Smith   mg[levels-1]->b    = b;
2104b9ad928SBarry Smith   mg[levels-1]->x    = x;
2114b9ad928SBarry Smith 
2124b9ad928SBarry Smith   mg[levels-1]->rtol = rtol;
21370441072SBarry Smith   mg[levels-1]->abstol = abstol;
2144b9ad928SBarry Smith   mg[levels-1]->dtol = dtol;
2154b9ad928SBarry Smith   if (rtol) {
2164b9ad928SBarry Smith     /* compute initial residual norm for relative convergence test */
2174b9ad928SBarry Smith     PetscReal rnorm;
2184b9ad928SBarry Smith     ierr               = (*mg[levels-1]->residual)(mg[levels-1]->A,b,x,w);CHKERRQ(ierr);
2194b9ad928SBarry Smith     ierr               = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr);
22070441072SBarry Smith     mg[levels-1]->ttol = PetscMax(rtol*rnorm,abstol);
22170441072SBarry Smith   } else if (abstol) {
22270441072SBarry Smith     mg[levels-1]->ttol = abstol;
2234b9ad928SBarry Smith   } else {
2244b9ad928SBarry Smith     mg[levels-1]->ttol = 0.0;
2254b9ad928SBarry Smith   }
2264b9ad928SBarry Smith 
2274b9ad928SBarry Smith   while (its-- && !converged) {
2289dcbbd2bSBarry Smith     ierr = PCMGMCycle_Private(mg+levels-1,&converged);CHKERRQ(ierr);
2294b9ad928SBarry Smith   }
2304b9ad928SBarry Smith   PetscFunctionReturn(0);
2314b9ad928SBarry Smith }
2324b9ad928SBarry Smith 
2334b9ad928SBarry Smith #undef __FUNCT__
2344b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_MG"
2356ca4d86aSHong Zhang PetscErrorCode PCSetFromOptions_MG(PC pc)
2364b9ad928SBarry Smith {
237dfbe8321SBarry Smith   PetscErrorCode ierr;
2389dcbbd2bSBarry Smith   PetscInt       m,levels = 1;
2394b9ad928SBarry Smith   PetscTruth     flg;
2409dcbbd2bSBarry Smith   PC_MG          **mg = (PC_MG**)pc->data;
241cf502942SBarry Smith   PCMGType       mgtype;
2424b9ad928SBarry Smith 
2434b9ad928SBarry Smith   PetscFunctionBegin;
2444b9ad928SBarry Smith 
2454b9ad928SBarry Smith   ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr);
2464b9ad928SBarry Smith     if (!pc->data) {
2479dcbbd2bSBarry Smith       ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr);
24897177400SBarry Smith       ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr);
249cf502942SBarry Smith       mg = (PC_MG**)pc->data;
2504b9ad928SBarry Smith     }
251cf502942SBarry Smith     mgtype = mg[0]->am;
2529dcbbd2bSBarry Smith     ierr = PetscOptionsInt("-pc_mg_cycles","1 for V cycle, 2 for W-cycle","PCMGSetCycles",1,&m,&flg);CHKERRQ(ierr);
2534b9ad928SBarry Smith     if (flg) {
25497177400SBarry Smith       ierr = PCMGSetCycles(pc,m);CHKERRQ(ierr);
2554b9ad928SBarry Smith     }
2569dcbbd2bSBarry Smith     ierr = PetscOptionsName("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",&flg);CHKERRQ(ierr);
257c2be2410SBarry Smith     if (flg) {
25897177400SBarry Smith       ierr = PCMGSetGalerkin(pc);CHKERRQ(ierr);
259c2be2410SBarry Smith     }
2609dcbbd2bSBarry Smith     ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr);
2614b9ad928SBarry Smith     if (flg) {
26297177400SBarry Smith       ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr);
2634b9ad928SBarry Smith     }
2649dcbbd2bSBarry Smith     ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr);
2654b9ad928SBarry Smith     if (flg) {
26697177400SBarry Smith       ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr);
2674b9ad928SBarry Smith     }
2689dcbbd2bSBarry Smith     ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
2699dcbbd2bSBarry Smith     if (flg) {ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);}
2704b9ad928SBarry Smith     ierr = PetscOptionsName("-pc_mg_log","Log times for each multigrid level","None",&flg);CHKERRQ(ierr);
2714b9ad928SBarry Smith     if (flg) {
2724f5ab15aSBarry Smith       PetscInt i;
2734b9ad928SBarry Smith       char     eventname[128];
2744b9ad928SBarry Smith       if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
2754b9ad928SBarry Smith       levels = mg[0]->levels;
2764b9ad928SBarry Smith       for (i=0; i<levels; i++) {
27777431f27SBarry Smith         sprintf(eventname,"MSetup Level %d",(int)i);
2784b9ad928SBarry Smith         ierr = PetscLogEventRegister(&mg[i]->eventsetup,eventname,pc->cookie);CHKERRQ(ierr);
27925a62d46SBarry Smith         sprintf(eventname,"MGSolve Level %d to 0",(int)i);
2804b9ad928SBarry Smith         ierr = PetscLogEventRegister(&mg[i]->eventsolve,eventname,pc->cookie);CHKERRQ(ierr);
2814b9ad928SBarry Smith       }
2824b9ad928SBarry Smith     }
2834b9ad928SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
2844b9ad928SBarry Smith   PetscFunctionReturn(0);
2854b9ad928SBarry Smith }
2864b9ad928SBarry Smith 
2879dcbbd2bSBarry Smith const char *PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0};
2889dcbbd2bSBarry Smith 
2894b9ad928SBarry Smith #undef __FUNCT__
2904b9ad928SBarry Smith #define __FUNCT__ "PCView_MG"
2916849ba73SBarry Smith static PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
2924b9ad928SBarry Smith {
2939dcbbd2bSBarry Smith   PC_MG          **mg = (PC_MG**)pc->data;
294dfbe8321SBarry Smith   PetscErrorCode ierr;
29579416396SBarry Smith   PetscInt       levels = mg[0]->levels,i;
29632077d6dSBarry Smith   PetscTruth     iascii;
2974b9ad928SBarry Smith 
2984b9ad928SBarry Smith   PetscFunctionBegin;
29932077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
30032077d6dSBarry Smith   if (iascii) {
30177431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  MG: type is %s, levels=%D cycles=%D, pre-smooths=%D, post-smooths=%D\n",
3029dcbbd2bSBarry Smith                       PCMGTypes[mg[0]->am],levels,mg[0]->cycles,mg[0]->default_smoothd,mg[0]->default_smoothu);CHKERRQ(ierr);
303c2be2410SBarry Smith     if (mg[0]->galerkin) {
304c2be2410SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
305c2be2410SBarry Smith     }
3064b9ad928SBarry Smith     for (i=0; i<levels; i++) {
307b03c7568SBarry Smith       if (!i) {
308b03c7568SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Coarse gride solver -- level %D -------------------------------\n",i);CHKERRQ(ierr);
309b03c7568SBarry Smith       } else {
31077431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
311b03c7568SBarry Smith       }
3124b9ad928SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
3134b9ad928SBarry Smith       ierr = KSPView(mg[i]->smoothd,viewer);CHKERRQ(ierr);
3144b9ad928SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
315b03c7568SBarry Smith       if (i && mg[i]->smoothd == mg[i]->smoothu) {
3164b9ad928SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
317b03c7568SBarry Smith       } else if (i){
31877431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
3194b9ad928SBarry Smith         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
3204b9ad928SBarry Smith         ierr = KSPView(mg[i]->smoothu,viewer);CHKERRQ(ierr);
3214b9ad928SBarry Smith         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
3224b9ad928SBarry Smith       }
3234b9ad928SBarry Smith     }
3244b9ad928SBarry Smith   } else {
32579a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCMG",((PetscObject)viewer)->type_name);
3264b9ad928SBarry Smith   }
3274b9ad928SBarry Smith   PetscFunctionReturn(0);
3284b9ad928SBarry Smith }
3294b9ad928SBarry Smith 
3304b9ad928SBarry Smith /*
3314b9ad928SBarry Smith     Calls setup for the KSP on each level
3324b9ad928SBarry Smith */
3334b9ad928SBarry Smith #undef __FUNCT__
3344b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_MG"
3356849ba73SBarry Smith static PetscErrorCode PCSetUp_MG(PC pc)
3364b9ad928SBarry Smith {
3379dcbbd2bSBarry Smith   PC_MG          **mg = (PC_MG**)pc->data;
338dfbe8321SBarry Smith   PetscErrorCode ierr;
33979416396SBarry Smith   PetscInt       i,n = mg[0]->levels;
3404b9ad928SBarry Smith   PC             cpc;
341*906ed7ccSBarry Smith   PetscTruth     preonly,lu,redundant,cholesky,monitor = PETSC_FALSE,dump,opsset;
3424b9ad928SBarry Smith   PetscViewer    ascii;
3434b9ad928SBarry Smith   MPI_Comm       comm;
344c2be2410SBarry Smith   Mat            dA,dB;
345c2be2410SBarry Smith   MatStructure   uflag;
3460a6bb862SBarry Smith   Vec            tvec;
3474b9ad928SBarry Smith 
3484b9ad928SBarry Smith   PetscFunctionBegin;
349852665d3SBarry Smith 
350852665d3SBarry Smith   /* If user did not provide fine grid operators, use those from PC */
351852665d3SBarry Smith   /* BUG BUG BUG This will work ONLY the first time called: hence if the user changes
352852665d3SBarry Smith      the PC matrices between solves PCMG will continue to use first set provided */
353*906ed7ccSBarry Smith   ierr = KSPGetOperatorsSet(mg[n-1]->smoothd,PETSC_NULL,&opsset);CHKERRQ(ierr);
354*906ed7ccSBarry Smith   if (!opsset) {
355852665d3SBarry 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);
356852665d3SBarry Smith     ierr = KSPSetOperators(mg[n-1]->smoothd,pc->mat,pc->pmat,uflag);CHKERRQ(ierr);
357852665d3SBarry Smith   }
358852665d3SBarry Smith 
359852665d3SBarry Smith   if (mg[0]->galerkin) {
360852665d3SBarry Smith     Mat B;
361852665d3SBarry Smith     mg[0]->galerkinused = PETSC_TRUE;
362852665d3SBarry Smith     /* currently only handle case where mat and pmat are the same on coarser levels */
363852665d3SBarry Smith     ierr = KSPGetOperators(mg[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr);
364852665d3SBarry Smith     if (!pc->setupcalled) {
365852665d3SBarry Smith       for (i=n-2; i>-1; i--) {
366852665d3SBarry Smith         ierr = MatPtAP(dB,mg[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
367852665d3SBarry Smith         ierr = KSPSetOperators(mg[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
368*906ed7ccSBarry Smith 	if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
369852665d3SBarry Smith         dB   = B;
370852665d3SBarry Smith       }
371*906ed7ccSBarry Smith       ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);
372852665d3SBarry Smith     } else {
373852665d3SBarry Smith       for (i=n-2; i>-1; i--) {
374*906ed7ccSBarry Smith         ierr = KSPGetOperators(mg[i]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
375852665d3SBarry Smith         ierr = MatPtAP(dB,mg[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
376852665d3SBarry Smith         ierr = KSPSetOperators(mg[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
377852665d3SBarry Smith         dB   = B;
378852665d3SBarry Smith       }
379852665d3SBarry Smith     }
380852665d3SBarry Smith   }
381852665d3SBarry Smith 
382958c9bccSBarry Smith   if (!pc->setupcalled) {
3834b9ad928SBarry Smith     ierr = PetscOptionsHasName(0,"-pc_mg_monitor",&monitor);CHKERRQ(ierr);
3844b9ad928SBarry Smith 
385b03c7568SBarry Smith     for (i=0; i<n; i++) {
3864b9ad928SBarry Smith       if (monitor) {
3874b9ad928SBarry Smith         ierr = PetscObjectGetComm((PetscObject)mg[i]->smoothd,&comm);CHKERRQ(ierr);
3884b9ad928SBarry Smith         ierr = PetscViewerASCIIOpen(comm,"stdout",&ascii);CHKERRQ(ierr);
3894b9ad928SBarry Smith         ierr = PetscViewerASCIISetTab(ascii,n-i);CHKERRQ(ierr);
3906849ba73SBarry Smith         ierr = KSPSetMonitor(mg[i]->smoothd,KSPDefaultMonitor,ascii,(PetscErrorCode(*)(void*))PetscViewerDestroy);CHKERRQ(ierr);
3914b9ad928SBarry Smith       }
3924b9ad928SBarry Smith       ierr = KSPSetFromOptions(mg[i]->smoothd);CHKERRQ(ierr);
3934b9ad928SBarry Smith     }
394b03c7568SBarry Smith     for (i=1; i<n; i++) {
395a98fc643SBarry Smith       if (mg[i]->smoothu && (mg[i]->smoothu != mg[i]->smoothd)) {
3964b9ad928SBarry Smith         if (monitor) {
3974b9ad928SBarry Smith           ierr = PetscObjectGetComm((PetscObject)mg[i]->smoothu,&comm);CHKERRQ(ierr);
3984b9ad928SBarry Smith           ierr = PetscViewerASCIIOpen(comm,"stdout",&ascii);CHKERRQ(ierr);
3994b9ad928SBarry Smith           ierr = PetscViewerASCIISetTab(ascii,n-i);CHKERRQ(ierr);
4006849ba73SBarry Smith           ierr = KSPSetMonitor(mg[i]->smoothu,KSPDefaultMonitor,ascii,(PetscErrorCode(*)(void*))PetscViewerDestroy);CHKERRQ(ierr);
4014b9ad928SBarry Smith         }
4024b9ad928SBarry Smith         ierr = KSPSetFromOptions(mg[i]->smoothu);CHKERRQ(ierr);
4034b9ad928SBarry Smith       }
4044b9ad928SBarry Smith     }
405fccaa45eSBarry Smith     for (i=1; i<n; i++) {
4060cace4b0SBarry Smith       if (!mg[i]->residual) {
4070cace4b0SBarry Smith         Mat mat;
4080cace4b0SBarry Smith         ierr = KSPGetOperators(mg[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr);
4090cace4b0SBarry Smith         ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr);
4100cace4b0SBarry Smith       }
411fccaa45eSBarry Smith       if (mg[i]->restrct && !mg[i]->interpolate) {
41297177400SBarry Smith         ierr = PCMGSetInterpolate(pc,i,mg[i]->restrct);CHKERRQ(ierr);
413fccaa45eSBarry Smith       }
414fccaa45eSBarry Smith       if (!mg[i]->restrct && mg[i]->interpolate) {
41597177400SBarry Smith         ierr = PCMGSetRestriction(pc,i,mg[i]->interpolate);CHKERRQ(ierr);
416fccaa45eSBarry Smith       }
417fccaa45eSBarry Smith #if defined(PETSC_USE_DEBUG)
418fccaa45eSBarry Smith       if (!mg[i]->restrct || !mg[i]->interpolate) {
419fccaa45eSBarry Smith         SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Need to set restriction or interpolation on level %d",(int)i);
420fccaa45eSBarry Smith       }
421fccaa45eSBarry Smith #endif
422fccaa45eSBarry Smith     }
4230a6bb862SBarry Smith     for (i=0; i<n-1; i++) {
42437b0e6c0SBarry Smith       if (!mg[i]->b) {
425*906ed7ccSBarry Smith         Vec *vec;
426*906ed7ccSBarry Smith         ierr = KSPGetVecs(mg[i]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
427*906ed7ccSBarry Smith         ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr);
428*906ed7ccSBarry Smith         ierr = PetscFree(vec);CHKERRQ(ierr);
42937b0e6c0SBarry Smith       }
4306ca4d86aSHong Zhang       if (!mg[i]->r && i) {
4310a6bb862SBarry Smith         ierr = VecDuplicate(mg[i]->b,&tvec);CHKERRQ(ierr);
43297177400SBarry Smith         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
4330a6bb862SBarry Smith         ierr = VecDestroy(tvec);CHKERRQ(ierr);
4340a6bb862SBarry Smith       }
4350a6bb862SBarry Smith       if (!mg[i]->x) {
4360a6bb862SBarry Smith         ierr = VecDuplicate(mg[i]->b,&tvec);CHKERRQ(ierr);
43797177400SBarry Smith         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
4380a6bb862SBarry Smith         ierr = VecDestroy(tvec);CHKERRQ(ierr);
4390a6bb862SBarry Smith       }
4400a6bb862SBarry Smith     }
4414b9ad928SBarry Smith   }
4424b9ad928SBarry Smith 
443c2be2410SBarry Smith 
4444b9ad928SBarry Smith   for (i=1; i<n; i++) {
445b03c7568SBarry Smith     if (mg[i]->smoothu == mg[i]->smoothd) {
446b03c7568SBarry Smith       /* if doing only down then initial guess is zero */
4474b9ad928SBarry Smith       ierr = KSPSetInitialGuessNonzero(mg[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
448b03c7568SBarry Smith     }
4494b9ad928SBarry Smith     if (mg[i]->eventsetup) {ierr = PetscLogEventBegin(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);}
4504b9ad928SBarry Smith     ierr = KSPSetUp(mg[i]->smoothd);CHKERRQ(ierr);
4514b9ad928SBarry Smith     if (mg[i]->eventsetup) {ierr = PetscLogEventEnd(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);}
4524b9ad928SBarry Smith   }
453b03c7568SBarry Smith   for (i=1; i<n; i++) {
4544b9ad928SBarry Smith     if (mg[i]->smoothu && mg[i]->smoothu != mg[i]->smoothd) {
455*906ed7ccSBarry Smith       Mat          downmat,downpmat;
45697f1f81fSBarry Smith       MatStructure matflag;
457*906ed7ccSBarry Smith       PetscTruth   opsset;
45897f1f81fSBarry Smith 
45997f1f81fSBarry Smith       /* check if operators have been set for up, if not use down operators to set them */
460*906ed7ccSBarry Smith       ierr = KSPGetOperatorsSet(mg[i]->smoothu,&opsset,PETSC_NULL);CHKERRQ(ierr);
461*906ed7ccSBarry Smith       if (!opsset) {
462*906ed7ccSBarry Smith         ierr = KSPGetOperators(mg[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr);
46397f1f81fSBarry Smith         ierr = KSPSetOperators(mg[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr);
46497f1f81fSBarry Smith       }
46597f1f81fSBarry Smith 
4664b9ad928SBarry Smith       ierr = KSPSetInitialGuessNonzero(mg[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
4674b9ad928SBarry Smith       if (mg[i]->eventsetup) {ierr = PetscLogEventBegin(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);}
4684b9ad928SBarry Smith       ierr = KSPSetUp(mg[i]->smoothu);CHKERRQ(ierr);
4694b9ad928SBarry Smith       if (mg[i]->eventsetup) {ierr = PetscLogEventEnd(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);}
4704b9ad928SBarry Smith     }
4714b9ad928SBarry Smith   }
4724b9ad928SBarry Smith 
4734b9ad928SBarry Smith   /*
4744b9ad928SBarry Smith       If coarse solver is not direct method then DO NOT USE preonly
4754b9ad928SBarry Smith   */
4764b9ad928SBarry Smith   ierr = PetscTypeCompare((PetscObject)mg[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr);
4774b9ad928SBarry Smith   if (preonly) {
4784b9ad928SBarry Smith     ierr = KSPGetPC(mg[0]->smoothd,&cpc);CHKERRQ(ierr);
4794b9ad928SBarry Smith     ierr = PetscTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr);
4804b9ad928SBarry Smith     ierr = PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr);
48168eff7e6SBarry Smith     ierr = PetscTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr);
48268eff7e6SBarry Smith     if (!lu && !redundant && !cholesky) {
4834b9ad928SBarry Smith       ierr = KSPSetType(mg[0]->smoothd,KSPGMRES);CHKERRQ(ierr);
4844b9ad928SBarry Smith     }
4854b9ad928SBarry Smith   }
4864b9ad928SBarry Smith 
487958c9bccSBarry Smith   if (!pc->setupcalled) {
4884b9ad928SBarry Smith     if (monitor) {
4894b9ad928SBarry Smith       ierr = PetscObjectGetComm((PetscObject)mg[0]->smoothd,&comm);CHKERRQ(ierr);
4904b9ad928SBarry Smith       ierr = PetscViewerASCIIOpen(comm,"stdout",&ascii);CHKERRQ(ierr);
4914b9ad928SBarry Smith       ierr = PetscViewerASCIISetTab(ascii,n);CHKERRQ(ierr);
4926849ba73SBarry Smith       ierr = KSPSetMonitor(mg[0]->smoothd,KSPDefaultMonitor,ascii,(PetscErrorCode(*)(void*))PetscViewerDestroy);CHKERRQ(ierr);
4934b9ad928SBarry Smith     }
4944b9ad928SBarry Smith     ierr = KSPSetFromOptions(mg[0]->smoothd);CHKERRQ(ierr);
4954b9ad928SBarry Smith   }
4964b9ad928SBarry Smith 
4974b9ad928SBarry Smith   if (mg[0]->eventsetup) {ierr = PetscLogEventBegin(mg[0]->eventsetup,0,0,0,0);CHKERRQ(ierr);}
4984b9ad928SBarry Smith   ierr = KSPSetUp(mg[0]->smoothd);CHKERRQ(ierr);
4994b9ad928SBarry Smith   if (mg[0]->eventsetup) {ierr = PetscLogEventEnd(mg[0]->eventsetup,0,0,0,0);CHKERRQ(ierr);}
5004b9ad928SBarry Smith 
5014cf18111SSatish Balay #if defined(PETSC_USE_SOCKET_VIEWER)
5024b9ad928SBarry Smith   /*
5034b9ad928SBarry Smith      Dump the interpolation/restriction matrices to matlab plus the
5044b9ad928SBarry Smith    Jacobian/stiffness on each level. This allows Matlab users to
5054b9ad928SBarry Smith    easily check if the Galerkin condition A_c = R A_f R^T is satisfied */
5064b9ad928SBarry Smith   ierr = PetscOptionsHasName(pc->prefix,"-pc_mg_dump_matlab",&dump);CHKERRQ(ierr);
5074b9ad928SBarry Smith   if (dump) {
5084b9ad928SBarry Smith     for (i=1; i<n; i++) {
5094b9ad928SBarry Smith       ierr = MatView(mg[i]->restrct,PETSC_VIEWER_SOCKET_(pc->comm));CHKERRQ(ierr);
5104b9ad928SBarry Smith     }
5114b9ad928SBarry Smith     for (i=0; i<n; i++) {
5124b9ad928SBarry Smith       ierr = KSPGetPC(mg[i]->smoothd,&pc);CHKERRQ(ierr);
5134b9ad928SBarry Smith       ierr = MatView(pc->mat,PETSC_VIEWER_SOCKET_(pc->comm));CHKERRQ(ierr);
5144b9ad928SBarry Smith     }
5154b9ad928SBarry Smith   }
516c45a1595SBarry Smith #endif
517c45a1595SBarry Smith 
518c2be2410SBarry Smith   ierr = PetscOptionsHasName(pc->prefix,"-pc_mg_dump_binary",&dump);CHKERRQ(ierr);
519c2be2410SBarry Smith   if (dump) {
520c2be2410SBarry Smith     for (i=1; i<n; i++) {
521c2be2410SBarry Smith       ierr = MatView(mg[i]->restrct,PETSC_VIEWER_BINARY_(pc->comm));CHKERRQ(ierr);
522c2be2410SBarry Smith     }
523c2be2410SBarry Smith     for (i=0; i<n; i++) {
524c2be2410SBarry Smith       ierr = KSPGetPC(mg[i]->smoothd,&pc);CHKERRQ(ierr);
525c2be2410SBarry Smith       ierr = MatView(pc->mat,PETSC_VIEWER_BINARY_(pc->comm));CHKERRQ(ierr);
526c2be2410SBarry Smith     }
527c2be2410SBarry Smith   }
5284b9ad928SBarry Smith   PetscFunctionReturn(0);
5294b9ad928SBarry Smith }
5304b9ad928SBarry Smith 
5314b9ad928SBarry Smith /* -------------------------------------------------------------------------------------*/
5324b9ad928SBarry Smith 
5334b9ad928SBarry Smith #undef __FUNCT__
5349dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetLevels"
5354b9ad928SBarry Smith /*@C
53697177400SBarry Smith    PCMGSetLevels - Sets the number of levels to use with MG.
5374b9ad928SBarry Smith    Must be called before any other MG routine.
5384b9ad928SBarry Smith 
5394b9ad928SBarry Smith    Collective on PC
5404b9ad928SBarry Smith 
5414b9ad928SBarry Smith    Input Parameters:
5424b9ad928SBarry Smith +  pc - the preconditioner context
5434b9ad928SBarry Smith .  levels - the number of levels
5444b9ad928SBarry Smith -  comms - optional communicators for each level; this is to allow solving the coarser problems
5454b9ad928SBarry Smith            on smaller sets of processors. Use PETSC_NULL_OBJECT for default in Fortran
5464b9ad928SBarry Smith 
5474b9ad928SBarry Smith    Level: intermediate
5484b9ad928SBarry Smith 
5494b9ad928SBarry Smith    Notes:
5504b9ad928SBarry Smith      If the number of levels is one then the multigrid uses the -mg_levels prefix
5514b9ad928SBarry Smith   for setting the level options rather than the -mg_coarse prefix.
5524b9ad928SBarry Smith 
5534b9ad928SBarry Smith .keywords: MG, set, levels, multigrid
5544b9ad928SBarry Smith 
55597177400SBarry Smith .seealso: PCMGSetType(), PCMGGetLevels()
5564b9ad928SBarry Smith @*/
55797177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
5584b9ad928SBarry Smith {
559dfbe8321SBarry Smith   PetscErrorCode ierr;
560ada7143aSSatish Balay   PC_MG          **mg=0;
5614b9ad928SBarry Smith 
5624b9ad928SBarry Smith   PetscFunctionBegin;
5634482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
5644b9ad928SBarry Smith 
5654b9ad928SBarry Smith   if (pc->data) {
5661302d50aSBarry Smith     SETERRQ(PETSC_ERR_ORDER,"Number levels already set for MG\n\
56797177400SBarry Smith     make sure that you call PCMGSetLevels() before KSPSetFromOptions()");
5684b9ad928SBarry Smith   }
5699dcbbd2bSBarry Smith   ierr                     = PCMGCreate_Private(pc->comm,levels,pc,comms,&mg);CHKERRQ(ierr);
5709dcbbd2bSBarry Smith   mg[0]->am                = PC_MG_MULTIPLICATIVE;
5714b9ad928SBarry Smith   pc->data                 = (void*)mg;
5724b9ad928SBarry Smith   pc->ops->applyrichardson = PCApplyRichardson_MG;
5734b9ad928SBarry Smith   PetscFunctionReturn(0);
5744b9ad928SBarry Smith }
5754b9ad928SBarry Smith 
5764b9ad928SBarry Smith #undef __FUNCT__
5779dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetLevels"
5784b9ad928SBarry Smith /*@
57997177400SBarry Smith    PCMGGetLevels - Gets the number of levels to use with MG.
5804b9ad928SBarry Smith 
5814b9ad928SBarry Smith    Not Collective
5824b9ad928SBarry Smith 
5834b9ad928SBarry Smith    Input Parameter:
5844b9ad928SBarry Smith .  pc - the preconditioner context
5854b9ad928SBarry Smith 
5864b9ad928SBarry Smith    Output parameter:
5874b9ad928SBarry Smith .  levels - the number of levels
5884b9ad928SBarry Smith 
5894b9ad928SBarry Smith    Level: advanced
5904b9ad928SBarry Smith 
5914b9ad928SBarry Smith .keywords: MG, get, levels, multigrid
5924b9ad928SBarry Smith 
59397177400SBarry Smith .seealso: PCMGSetLevels()
5944b9ad928SBarry Smith @*/
59597177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetLevels(PC pc,PetscInt *levels)
5964b9ad928SBarry Smith {
5979dcbbd2bSBarry Smith   PC_MG  **mg;
5984b9ad928SBarry Smith 
5994b9ad928SBarry Smith   PetscFunctionBegin;
6004482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
6014482741eSBarry Smith   PetscValidIntPointer(levels,2);
6024b9ad928SBarry Smith 
6039dcbbd2bSBarry Smith   mg      = (PC_MG**)pc->data;
6044b9ad928SBarry Smith   *levels = mg[0]->levels;
6054b9ad928SBarry Smith   PetscFunctionReturn(0);
6064b9ad928SBarry Smith }
6074b9ad928SBarry Smith 
6084b9ad928SBarry Smith #undef __FUNCT__
6099dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetType"
6104b9ad928SBarry Smith /*@
61197177400SBarry Smith    PCMGSetType - Determines the form of multigrid to use:
6124b9ad928SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
6134b9ad928SBarry Smith 
6144b9ad928SBarry Smith    Collective on PC
6154b9ad928SBarry Smith 
6164b9ad928SBarry Smith    Input Parameters:
6174b9ad928SBarry Smith +  pc - the preconditioner context
6189dcbbd2bSBarry Smith -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
6199dcbbd2bSBarry Smith    PC_MG_FULL, PC_MG_KASKADE
6204b9ad928SBarry Smith 
6214b9ad928SBarry Smith    Options Database Key:
6224b9ad928SBarry Smith .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
6234b9ad928SBarry Smith    additive, full, kaskade
6244b9ad928SBarry Smith 
6254b9ad928SBarry Smith    Level: advanced
6264b9ad928SBarry Smith 
6274b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
6284b9ad928SBarry Smith 
62997177400SBarry Smith .seealso: PCMGSetLevels()
6304b9ad928SBarry Smith @*/
6319dcbbd2bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetType(PC pc,PCMGType form)
6324b9ad928SBarry Smith {
6339dcbbd2bSBarry Smith   PC_MG **mg;
6344b9ad928SBarry Smith 
6354b9ad928SBarry Smith   PetscFunctionBegin;
6364482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
6379dcbbd2bSBarry Smith   mg = (PC_MG**)pc->data;
6384b9ad928SBarry Smith 
6394b9ad928SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
6404b9ad928SBarry Smith   mg[0]->am = form;
6419dcbbd2bSBarry Smith   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
6424b9ad928SBarry Smith   else pc->ops->applyrichardson = 0;
6434b9ad928SBarry Smith   PetscFunctionReturn(0);
6444b9ad928SBarry Smith }
6454b9ad928SBarry Smith 
6464b9ad928SBarry Smith #undef __FUNCT__
6479dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetCycles"
6484b9ad928SBarry Smith /*@
64997177400SBarry Smith    PCMGSetCycles - Sets the type cycles to use.  Use PCMGSetCyclesOnLevel() for more
6504b9ad928SBarry Smith    complicated cycling.
6514b9ad928SBarry Smith 
6524b9ad928SBarry Smith    Collective on PC
6534b9ad928SBarry Smith 
6544b9ad928SBarry Smith    Input Parameters:
655c2be2410SBarry Smith +  pc - the multigrid context
6564b9ad928SBarry Smith -  n - the number of cycles
6574b9ad928SBarry Smith 
6584b9ad928SBarry Smith    Options Database Key:
6594b9ad928SBarry Smith $  -pc_mg_cycles n - 1 denotes a V-cycle; 2 denotes a W-cycle.
6604b9ad928SBarry Smith 
6614b9ad928SBarry Smith    Level: advanced
6624b9ad928SBarry Smith 
6634b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
6644b9ad928SBarry Smith 
66597177400SBarry Smith .seealso: PCMGSetCyclesOnLevel()
6664b9ad928SBarry Smith @*/
66797177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetCycles(PC pc,PetscInt n)
6684b9ad928SBarry Smith {
6699dcbbd2bSBarry Smith   PC_MG    **mg;
67079416396SBarry Smith   PetscInt i,levels;
6714b9ad928SBarry Smith 
6724b9ad928SBarry Smith   PetscFunctionBegin;
6734482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
6749dcbbd2bSBarry Smith   mg     = (PC_MG**)pc->data;
6754b9ad928SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
6764b9ad928SBarry Smith   levels = mg[0]->levels;
6774b9ad928SBarry Smith 
6784b9ad928SBarry Smith   for (i=0; i<levels; i++) {
6794b9ad928SBarry Smith     mg[i]->cycles  = n;
6804b9ad928SBarry Smith   }
6814b9ad928SBarry Smith   PetscFunctionReturn(0);
6824b9ad928SBarry Smith }
6834b9ad928SBarry Smith 
6844b9ad928SBarry Smith #undef __FUNCT__
6859dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetGalerkin"
686c2be2410SBarry Smith /*@
68797177400SBarry Smith    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
688c2be2410SBarry Smith       finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t
689c2be2410SBarry Smith 
690c2be2410SBarry Smith    Collective on PC
691c2be2410SBarry Smith 
692c2be2410SBarry Smith    Input Parameters:
6933fc8bf9cSBarry Smith .  pc - the multigrid context
694c2be2410SBarry Smith 
695c2be2410SBarry Smith    Options Database Key:
696c2be2410SBarry Smith $  -pc_mg_galerkin
697c2be2410SBarry Smith 
698c2be2410SBarry Smith    Level: intermediate
699c2be2410SBarry Smith 
700c2be2410SBarry Smith .keywords: MG, set, Galerkin
701c2be2410SBarry Smith 
7023fc8bf9cSBarry Smith .seealso: PCMGGetGalerkin()
7033fc8bf9cSBarry Smith 
704c2be2410SBarry Smith @*/
70597177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetGalerkin(PC pc)
706c2be2410SBarry Smith {
7079dcbbd2bSBarry Smith   PC_MG    **mg;
708c2be2410SBarry Smith   PetscInt i,levels;
709c2be2410SBarry Smith 
710c2be2410SBarry Smith   PetscFunctionBegin;
711c2be2410SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
7129dcbbd2bSBarry Smith   mg     = (PC_MG**)pc->data;
713c2be2410SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
714c2be2410SBarry Smith   levels = mg[0]->levels;
715c2be2410SBarry Smith 
716c2be2410SBarry Smith   for (i=0; i<levels; i++) {
717c2be2410SBarry Smith     mg[i]->galerkin = PETSC_TRUE;
718c2be2410SBarry Smith   }
719c2be2410SBarry Smith   PetscFunctionReturn(0);
720c2be2410SBarry Smith }
721c2be2410SBarry Smith 
722c2be2410SBarry Smith #undef __FUNCT__
7233fc8bf9cSBarry Smith #define __FUNCT__ "PCMGGetGalerkin"
7243fc8bf9cSBarry Smith /*@
7253fc8bf9cSBarry Smith    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
7263fc8bf9cSBarry Smith       A_i-1 = r_i * A_i * r_i^t
7273fc8bf9cSBarry Smith 
7283fc8bf9cSBarry Smith    Not Collective
7293fc8bf9cSBarry Smith 
7303fc8bf9cSBarry Smith    Input Parameter:
7313fc8bf9cSBarry Smith .  pc - the multigrid context
7323fc8bf9cSBarry Smith 
7333fc8bf9cSBarry Smith    Output Parameter:
7343fc8bf9cSBarry Smith .  gelerkin - PETSC_TRUE or PETSC_FALSE
7353fc8bf9cSBarry Smith 
7363fc8bf9cSBarry Smith    Options Database Key:
7373fc8bf9cSBarry Smith $  -pc_mg_galerkin
7383fc8bf9cSBarry Smith 
7393fc8bf9cSBarry Smith    Level: intermediate
7403fc8bf9cSBarry Smith 
7413fc8bf9cSBarry Smith .keywords: MG, set, Galerkin
7423fc8bf9cSBarry Smith 
7433fc8bf9cSBarry Smith .seealso: PCMGSetGalerkin()
7443fc8bf9cSBarry Smith 
7453fc8bf9cSBarry Smith @*/
7463fc8bf9cSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetGalerkin(PC pc,PetscTruth *galerkin)
7473fc8bf9cSBarry Smith {
7483fc8bf9cSBarry Smith   PC_MG    **mg;
7493fc8bf9cSBarry Smith 
7503fc8bf9cSBarry Smith   PetscFunctionBegin;
7513fc8bf9cSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
7523fc8bf9cSBarry Smith   mg     = (PC_MG**)pc->data;
7533fc8bf9cSBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
7543fc8bf9cSBarry Smith   *galerkin = mg[0]->galerkin;
7553fc8bf9cSBarry Smith   PetscFunctionReturn(0);
7563fc8bf9cSBarry Smith }
7573fc8bf9cSBarry Smith 
7583fc8bf9cSBarry Smith #undef __FUNCT__
7599dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothDown"
7604b9ad928SBarry Smith /*@
76197177400SBarry Smith    PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
76297177400SBarry Smith    use on all levels. Use PCMGGetSmootherDown() to set different
7634b9ad928SBarry Smith    pre-smoothing steps on different levels.
7644b9ad928SBarry Smith 
7654b9ad928SBarry Smith    Collective on PC
7664b9ad928SBarry Smith 
7674b9ad928SBarry Smith    Input Parameters:
7684b9ad928SBarry Smith +  mg - the multigrid context
7694b9ad928SBarry Smith -  n - the number of smoothing steps
7704b9ad928SBarry Smith 
7714b9ad928SBarry Smith    Options Database Key:
7724b9ad928SBarry Smith .  -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps
7734b9ad928SBarry Smith 
7744b9ad928SBarry Smith    Level: advanced
7754b9ad928SBarry Smith 
7764b9ad928SBarry Smith .keywords: MG, smooth, down, pre-smoothing, steps, multigrid
7774b9ad928SBarry Smith 
77897177400SBarry Smith .seealso: PCMGSetNumberSmoothUp()
7794b9ad928SBarry Smith @*/
78097177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothDown(PC pc,PetscInt n)
7814b9ad928SBarry Smith {
7829dcbbd2bSBarry Smith   PC_MG          **mg;
7836849ba73SBarry Smith   PetscErrorCode ierr;
78479416396SBarry Smith   PetscInt       i,levels;
7854b9ad928SBarry Smith 
7864b9ad928SBarry Smith   PetscFunctionBegin;
7874482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
7889dcbbd2bSBarry Smith   mg     = (PC_MG**)pc->data;
7894b9ad928SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
7904b9ad928SBarry Smith   levels = mg[0]->levels;
7914b9ad928SBarry Smith 
792b05257ddSBarry Smith   for (i=1; i<levels; i++) {
7934b9ad928SBarry Smith     /* make sure smoother up and down are different */
79497177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
7954b9ad928SBarry Smith     ierr = KSPSetTolerances(mg[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
7964b9ad928SBarry Smith     mg[i]->default_smoothd = n;
7974b9ad928SBarry Smith   }
7984b9ad928SBarry Smith   PetscFunctionReturn(0);
7994b9ad928SBarry Smith }
8004b9ad928SBarry Smith 
8014b9ad928SBarry Smith #undef __FUNCT__
8029dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothUp"
8034b9ad928SBarry Smith /*@
80497177400SBarry Smith    PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
80597177400SBarry Smith    on all levels. Use PCMGGetSmootherUp() to set different numbers of
8064b9ad928SBarry Smith    post-smoothing steps on different levels.
8074b9ad928SBarry Smith 
8084b9ad928SBarry Smith    Collective on PC
8094b9ad928SBarry Smith 
8104b9ad928SBarry Smith    Input Parameters:
8114b9ad928SBarry Smith +  mg - the multigrid context
8124b9ad928SBarry Smith -  n - the number of smoothing steps
8134b9ad928SBarry Smith 
8144b9ad928SBarry Smith    Options Database Key:
8154b9ad928SBarry Smith .  -pc_mg_smoothup <n> - Sets number of post-smoothing steps
8164b9ad928SBarry Smith 
8174b9ad928SBarry Smith    Level: advanced
8184b9ad928SBarry Smith 
8194b9ad928SBarry Smith    Note: this does not set a value on the coarsest grid, since we assume that
820a8c7a070SBarry Smith     there is no separate smooth up on the coarsest grid.
8214b9ad928SBarry Smith 
8224b9ad928SBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid
8234b9ad928SBarry Smith 
82497177400SBarry Smith .seealso: PCMGSetNumberSmoothDown()
8254b9ad928SBarry Smith @*/
82697177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothUp(PC pc,PetscInt n)
8274b9ad928SBarry Smith {
8289dcbbd2bSBarry Smith   PC_MG          **mg;
8296849ba73SBarry Smith   PetscErrorCode ierr;
83079416396SBarry Smith   PetscInt       i,levels;
8314b9ad928SBarry Smith 
8324b9ad928SBarry Smith   PetscFunctionBegin;
8334482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
8349dcbbd2bSBarry Smith   mg     = (PC_MG**)pc->data;
8354b9ad928SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
8364b9ad928SBarry Smith   levels = mg[0]->levels;
8374b9ad928SBarry Smith 
8384b9ad928SBarry Smith   for (i=1; i<levels; i++) {
8394b9ad928SBarry Smith     /* make sure smoother up and down are different */
84097177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
8414b9ad928SBarry Smith     ierr = KSPSetTolerances(mg[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
8424b9ad928SBarry Smith     mg[i]->default_smoothu = n;
8434b9ad928SBarry Smith   }
8444b9ad928SBarry Smith   PetscFunctionReturn(0);
8454b9ad928SBarry Smith }
8464b9ad928SBarry Smith 
8474b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/
8484b9ad928SBarry Smith 
8493b09bd56SBarry Smith /*MC
8503b09bd56SBarry Smith    PCMG - Use geometric multigrid preconditioning. This preconditioner requires you provide additional
8513b09bd56SBarry Smith     information about the coarser grid matrices and restriction/interpolation operators.
8523b09bd56SBarry Smith 
8533b09bd56SBarry Smith    Options Database Keys:
8543b09bd56SBarry Smith +  -pc_mg_levels <nlevels> - number of levels including finest
8553b09bd56SBarry Smith .  -pc_mg_cycles 1 or 2 - for V or W-cycle
85679416396SBarry Smith .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
8573b09bd56SBarry Smith .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
8583b09bd56SBarry Smith .  -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default
8593b09bd56SBarry Smith .  -pc_mg_log - log information about time spent on each level of the solver
8603b09bd56SBarry Smith .  -pc_mg_monitor - print information on the multigrid convergence
86168eff7e6SBarry Smith .  -pc_mg_galerkin - use Galerkin process to compute coarser operators
8623b09bd56SBarry Smith -  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
8633b09bd56SBarry Smith                         to the Socket viewer for reading from Matlab.
8643b09bd56SBarry Smith 
8653b09bd56SBarry Smith    Notes:
8663b09bd56SBarry Smith 
8673b09bd56SBarry Smith    Level: intermediate
8683b09bd56SBarry Smith 
8693b09bd56SBarry Smith    Concepts: multigrid
8703b09bd56SBarry Smith 
8713b09bd56SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
8729dcbbd2bSBarry Smith            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycles(), PCMGSetNumberSmoothDown(),
87397177400SBarry Smith            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
87497177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
87597177400SBarry Smith            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
8763b09bd56SBarry Smith M*/
8773b09bd56SBarry Smith 
8784b9ad928SBarry Smith EXTERN_C_BEGIN
8794b9ad928SBarry Smith #undef __FUNCT__
8804b9ad928SBarry Smith #define __FUNCT__ "PCCreate_MG"
881dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_MG(PC pc)
8824b9ad928SBarry Smith {
8834b9ad928SBarry Smith   PetscFunctionBegin;
8844b9ad928SBarry Smith   pc->ops->apply          = PCApply_MG;
8854b9ad928SBarry Smith   pc->ops->setup          = PCSetUp_MG;
8864b9ad928SBarry Smith   pc->ops->destroy        = PCDestroy_MG;
8874b9ad928SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_MG;
8884b9ad928SBarry Smith   pc->ops->view           = PCView_MG;
8894b9ad928SBarry Smith 
8904b9ad928SBarry Smith   pc->data                = (void*)0;
8914b9ad928SBarry Smith   PetscFunctionReturn(0);
8924b9ad928SBarry Smith }
8934b9ad928SBarry Smith EXTERN_C_END
894