xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision f69a0ea3504314d029ee423e0de2950ece2dab86)
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   PetscScalar    zero = 0.0;
174b9ad928SBarry Smith 
184b9ad928SBarry Smith   PetscFunctionBegin;
194b9ad928SBarry Smith   if (converged) *converged = PETSC_FALSE;
204b9ad928SBarry Smith 
214b9ad928SBarry Smith   if (mg->eventsolve) {ierr = PetscLogEventBegin(mg->eventsolve,0,0,0,0);CHKERRQ(ierr);}
2223ce1328SBarry Smith   ierr = KSPSolve(mg->smoothd,mg->b,mg->x);CHKERRQ(ierr);
234b9ad928SBarry Smith   if (mg->eventsolve) {ierr = PetscLogEventEnd(mg->eventsolve,0,0,0,0);CHKERRQ(ierr);}
244b9ad928SBarry Smith   if (mg->level) {  /* not the coarsest grid */
254b9ad928SBarry Smith     ierr = (*mg->residual)(mg->A,mg->b,mg->x,mg->r);CHKERRQ(ierr);
264b9ad928SBarry Smith 
274b9ad928SBarry Smith     /* if on finest level and have convergence criteria set */
284b9ad928SBarry Smith     if (mg->level == mg->levels-1 && mg->ttol) {
294b9ad928SBarry Smith       PetscReal rnorm;
304b9ad928SBarry Smith       ierr = VecNorm(mg->r,NORM_2,&rnorm);CHKERRQ(ierr);
314b9ad928SBarry Smith       if (rnorm <= mg->ttol) {
324b9ad928SBarry Smith         *converged = PETSC_TRUE;
3370441072SBarry Smith         if (rnorm < mg->abstol) {
349dcbbd2bSBarry Smith           ierr = PetscLogInfo((0,"PCMGMCycle_Private:Linear solver has converged. Residual norm %g is less than absolute tolerance %g\n",rnorm,mg->abstol));CHKERRQ(ierr);
354b9ad928SBarry Smith         } else {
369dcbbd2bSBarry Smith           ierr = PetscLogInfo((0,"PCMGMCycle_Private:Linear solver has converged. Residual norm %g is less than relative tolerance times initial residual norm %g\n",rnorm,mg->ttol));CHKERRQ(ierr);
374b9ad928SBarry Smith         }
384b9ad928SBarry Smith         PetscFunctionReturn(0);
394b9ad928SBarry Smith       }
404b9ad928SBarry Smith     }
414b9ad928SBarry Smith 
424b9ad928SBarry Smith     mgc = *(mglevels - 1);
434b9ad928SBarry Smith     ierr = MatRestrict(mg->restrct,mg->r,mgc->b);CHKERRQ(ierr);
442dcb1b2aSMatthew Knepley     ierr = VecSet(mgc->x,zero);CHKERRQ(ierr);
454b9ad928SBarry Smith     while (cycles--) {
469dcbbd2bSBarry Smith       ierr = PCMGMCycle_Private(mglevels-1,converged);CHKERRQ(ierr);
474b9ad928SBarry Smith     }
484b9ad928SBarry Smith     ierr = MatInterpolateAdd(mg->interpolate,mgc->x,mg->x,mg->x);CHKERRQ(ierr);
494b9ad928SBarry Smith     if (mg->eventsolve) {ierr = PetscLogEventBegin(mg->eventsolve,0,0,0,0);CHKERRQ(ierr);}
5023ce1328SBarry Smith     ierr = KSPSolve(mg->smoothu,mg->b,mg->x);CHKERRQ(ierr);
514b9ad928SBarry Smith     if (mg->eventsolve) {ierr = PetscLogEventEnd(mg->eventsolve,0,0,0,0);CHKERRQ(ierr);}
524b9ad928SBarry Smith   }
534b9ad928SBarry Smith   PetscFunctionReturn(0);
544b9ad928SBarry Smith }
554b9ad928SBarry Smith 
564b9ad928SBarry Smith /*
579dcbbd2bSBarry Smith        PCMGCreate_Private - Creates a PC_MG structure for use with the
584b9ad928SBarry Smith                multigrid code. Level 0 is the coarsest. (But the
594b9ad928SBarry Smith                finest level is stored first in the array).
604b9ad928SBarry Smith 
614b9ad928SBarry Smith */
624b9ad928SBarry Smith #undef __FUNCT__
639dcbbd2bSBarry Smith #define __FUNCT__ "PCMGCreate_Private"
649dcbbd2bSBarry Smith static PetscErrorCode PCMGCreate_Private(MPI_Comm comm,PetscInt levels,PC pc,MPI_Comm *comms,PC_MG ***result)
654b9ad928SBarry Smith {
669dcbbd2bSBarry Smith   PC_MG          **mg;
676849ba73SBarry Smith   PetscErrorCode ierr;
6879416396SBarry Smith   PetscInt       i;
6979416396SBarry Smith   PetscMPIInt    size;
70*f69a0ea3SMatthew Knepley   const char     *prefix;
714b9ad928SBarry Smith   PC             ipc;
724b9ad928SBarry Smith 
734b9ad928SBarry Smith   PetscFunctionBegin;
749dcbbd2bSBarry Smith   ierr = PetscMalloc(levels*sizeof(PC_MG*),&mg);CHKERRQ(ierr);
759dcbbd2bSBarry Smith   ierr = PetscLogObjectMemory(pc,levels*(sizeof(PC_MG*)+sizeof(PC_MG)));CHKERRQ(ierr);
764b9ad928SBarry Smith 
774b9ad928SBarry Smith   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
784b9ad928SBarry Smith 
794b9ad928SBarry Smith   for (i=0; i<levels; i++) {
809dcbbd2bSBarry Smith     ierr = PetscNew(PC_MG,&mg[i]);CHKERRQ(ierr);
814b9ad928SBarry Smith     mg[i]->level           = i;
824b9ad928SBarry Smith     mg[i]->levels          = levels;
834b9ad928SBarry Smith     mg[i]->cycles          = 1;
84c2be2410SBarry Smith     mg[i]->galerkin        = PETSC_FALSE;
85c2be2410SBarry Smith     mg[i]->galerkinused    = PETSC_FALSE;
8679416396SBarry Smith     mg[i]->default_smoothu = 1;
8779416396SBarry Smith     mg[i]->default_smoothd = 1;
884b9ad928SBarry Smith 
894b9ad928SBarry Smith     if (comms) comm = comms[i];
904b9ad928SBarry Smith     ierr = KSPCreate(comm,&mg[i]->smoothd);CHKERRQ(ierr);
9179416396SBarry Smith     ierr = KSPSetTolerances(mg[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg[i]->default_smoothd);CHKERRQ(ierr);
924b9ad928SBarry Smith     ierr = KSPSetOptionsPrefix(mg[i]->smoothd,prefix);CHKERRQ(ierr);
934b9ad928SBarry Smith 
944b9ad928SBarry Smith     /* do special stuff for coarse grid */
954b9ad928SBarry Smith     if (!i && levels > 1) {
964b9ad928SBarry Smith       ierr = KSPAppendOptionsPrefix(mg[0]->smoothd,"mg_coarse_");CHKERRQ(ierr);
974b9ad928SBarry Smith 
984b9ad928SBarry Smith       /* coarse solve is (redundant) LU by default */
994b9ad928SBarry Smith       ierr = KSPSetType(mg[0]->smoothd,KSPPREONLY);CHKERRQ(ierr);
1004b9ad928SBarry Smith       ierr = KSPGetPC(mg[0]->smoothd,&ipc);CHKERRQ(ierr);
1014b9ad928SBarry Smith       ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1024b9ad928SBarry Smith       if (size > 1) {
1034b9ad928SBarry Smith         ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr);
1044b9ad928SBarry Smith         ierr = PCRedundantGetPC(ipc,&ipc);CHKERRQ(ierr);
1054b9ad928SBarry Smith       }
1064b9ad928SBarry Smith       ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr);
1074b9ad928SBarry Smith 
1084b9ad928SBarry Smith     } else {
1092e70eadcSBarry Smith       char tprefix[128];
11013f74950SBarry Smith       sprintf(tprefix,"mg_levels_%d_",(int)i);
1112e70eadcSBarry Smith       ierr = KSPAppendOptionsPrefix(mg[i]->smoothd,tprefix);CHKERRQ(ierr);
1124b9ad928SBarry Smith     }
11352e6d16bSBarry Smith     ierr = PetscLogObjectParent(pc,mg[i]->smoothd);CHKERRQ(ierr);
1144b9ad928SBarry Smith     mg[i]->smoothu         = mg[i]->smoothd;
1154b9ad928SBarry Smith     mg[i]->rtol = 0.0;
11670441072SBarry Smith     mg[i]->abstol = 0.0;
1174b9ad928SBarry Smith     mg[i]->dtol = 0.0;
1184b9ad928SBarry Smith     mg[i]->ttol = 0.0;
1194b9ad928SBarry Smith     mg[i]->eventsetup = 0;
1204b9ad928SBarry Smith     mg[i]->eventsolve = 0;
1214b9ad928SBarry Smith   }
1224b9ad928SBarry Smith   *result = mg;
1234b9ad928SBarry Smith   PetscFunctionReturn(0);
1244b9ad928SBarry Smith }
1254b9ad928SBarry Smith 
1264b9ad928SBarry Smith #undef __FUNCT__
1274b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_MG"
1286849ba73SBarry Smith static PetscErrorCode PCDestroy_MG(PC pc)
1294b9ad928SBarry Smith {
1309dcbbd2bSBarry Smith   PC_MG          **mg = (PC_MG**)pc->data;
1316849ba73SBarry Smith   PetscErrorCode ierr;
13279416396SBarry Smith   PetscInt       i,n = mg[0]->levels;
1334b9ad928SBarry Smith 
1344b9ad928SBarry Smith   PetscFunctionBegin;
135c2be2410SBarry Smith   if (mg[0]->galerkinused) {
136c2be2410SBarry Smith     Mat B;
137c2be2410SBarry Smith     for (i=0; i<n-1; i++) {
138c2be2410SBarry Smith       ierr = KSPGetOperators(mg[i]->smoothd,0,&B,0);CHKERRQ(ierr);
139c2be2410SBarry Smith       ierr = MatDestroy(B);CHKERRQ(ierr);
140c2be2410SBarry Smith     }
141c2be2410SBarry Smith   }
142c2be2410SBarry Smith 
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;
1804b9ad928SBarry Smith   PetscScalar    zero = 0.0;
1816849ba73SBarry Smith   PetscErrorCode ierr;
18279416396SBarry Smith   PetscInt       levels = mg[0]->levels;
1834b9ad928SBarry Smith 
1844b9ad928SBarry Smith   PetscFunctionBegin;
1854b9ad928SBarry Smith   mg[levels-1]->b = b;
1864b9ad928SBarry Smith   mg[levels-1]->x = x;
1879dcbbd2bSBarry Smith   if (!mg[levels-1]->r && mg[0]->am == PC_MG_ADDITIVE) {
1880a6bb862SBarry Smith     Vec tvec;
1890a6bb862SBarry Smith     ierr = VecDuplicate(mg[levels-1]->b,&tvec);CHKERRQ(ierr);
19097177400SBarry Smith     ierr = PCMGSetR(pc,levels-1,tvec);CHKERRQ(ierr);
1910a6bb862SBarry Smith     ierr = VecDestroy(tvec);CHKERRQ(ierr);
1920a6bb862SBarry Smith   }
1939dcbbd2bSBarry Smith   if (mg[0]->am == PC_MG_MULTIPLICATIVE) {
1942dcb1b2aSMatthew Knepley     ierr = VecSet(x,zero);CHKERRQ(ierr);
1959dcbbd2bSBarry Smith     ierr = PCMGMCycle_Private(mg+levels-1,PETSC_NULL);CHKERRQ(ierr);
1964b9ad928SBarry Smith   }
1979dcbbd2bSBarry Smith   else if (mg[0]->am == PC_MG_ADDITIVE) {
1989dcbbd2bSBarry Smith     ierr = PCMGACycle_Private(mg);CHKERRQ(ierr);
1994b9ad928SBarry Smith   }
2009dcbbd2bSBarry Smith   else if (mg[0]->am == PC_MG_KASKADE) {
2019dcbbd2bSBarry Smith     ierr = PCMGKCycle_Private(mg);CHKERRQ(ierr);
2024b9ad928SBarry Smith   }
2034b9ad928SBarry Smith   else {
2049dcbbd2bSBarry Smith     ierr = PCMGFCycle_Private(mg);CHKERRQ(ierr);
2054b9ad928SBarry Smith   }
2064b9ad928SBarry Smith   PetscFunctionReturn(0);
2074b9ad928SBarry Smith }
2084b9ad928SBarry Smith 
2094b9ad928SBarry Smith #undef __FUNCT__
2104b9ad928SBarry Smith #define __FUNCT__ "PCApplyRichardson_MG"
21179416396SBarry Smith static PetscErrorCode PCApplyRichardson_MG(PC pc,Vec b,Vec x,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its)
2124b9ad928SBarry Smith {
2139dcbbd2bSBarry Smith   PC_MG          **mg = (PC_MG**)pc->data;
214dfbe8321SBarry Smith   PetscErrorCode ierr;
21579416396SBarry Smith   PetscInt       levels = mg[0]->levels;
2164b9ad928SBarry Smith   PetscTruth     converged = PETSC_FALSE;
2174b9ad928SBarry Smith 
2184b9ad928SBarry Smith   PetscFunctionBegin;
2194b9ad928SBarry Smith   mg[levels-1]->b    = b;
2204b9ad928SBarry Smith   mg[levels-1]->x    = x;
2214b9ad928SBarry Smith 
2224b9ad928SBarry Smith   mg[levels-1]->rtol = rtol;
22370441072SBarry Smith   mg[levels-1]->abstol = abstol;
2244b9ad928SBarry Smith   mg[levels-1]->dtol = dtol;
2254b9ad928SBarry Smith   if (rtol) {
2264b9ad928SBarry Smith     /* compute initial residual norm for relative convergence test */
2274b9ad928SBarry Smith     PetscReal rnorm;
2284b9ad928SBarry Smith     ierr               = (*mg[levels-1]->residual)(mg[levels-1]->A,b,x,w);CHKERRQ(ierr);
2294b9ad928SBarry Smith     ierr               = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr);
23070441072SBarry Smith     mg[levels-1]->ttol = PetscMax(rtol*rnorm,abstol);
23170441072SBarry Smith   } else if (abstol) {
23270441072SBarry Smith     mg[levels-1]->ttol = abstol;
2334b9ad928SBarry Smith   } else {
2344b9ad928SBarry Smith     mg[levels-1]->ttol = 0.0;
2354b9ad928SBarry Smith   }
2364b9ad928SBarry Smith 
2374b9ad928SBarry Smith   while (its-- && !converged) {
2389dcbbd2bSBarry Smith     ierr = PCMGMCycle_Private(mg+levels-1,&converged);CHKERRQ(ierr);
2394b9ad928SBarry Smith   }
2404b9ad928SBarry Smith   PetscFunctionReturn(0);
2414b9ad928SBarry Smith }
2424b9ad928SBarry Smith 
2434b9ad928SBarry Smith #undef __FUNCT__
2444b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_MG"
2456ca4d86aSHong Zhang PetscErrorCode PCSetFromOptions_MG(PC pc)
2464b9ad928SBarry Smith {
247dfbe8321SBarry Smith   PetscErrorCode ierr;
2489dcbbd2bSBarry Smith   PetscInt       m,levels = 1;
2494b9ad928SBarry Smith   PetscTruth     flg;
2509dcbbd2bSBarry Smith   PC_MG          **mg = (PC_MG**)pc->data;
2519dcbbd2bSBarry Smith   PCMGType       mgtype = mg[0]->am;;
2524b9ad928SBarry Smith 
2534b9ad928SBarry Smith   PetscFunctionBegin;
2544b9ad928SBarry Smith 
2554b9ad928SBarry Smith   ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr);
2564b9ad928SBarry Smith     if (!pc->data) {
2579dcbbd2bSBarry Smith       ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr);
25897177400SBarry Smith       ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr);
2594b9ad928SBarry Smith     }
2609dcbbd2bSBarry Smith     ierr = PetscOptionsInt("-pc_mg_cycles","1 for V cycle, 2 for W-cycle","PCMGSetCycles",1,&m,&flg);CHKERRQ(ierr);
2614b9ad928SBarry Smith     if (flg) {
26297177400SBarry Smith       ierr = PCMGSetCycles(pc,m);CHKERRQ(ierr);
2634b9ad928SBarry Smith     }
2649dcbbd2bSBarry Smith     ierr = PetscOptionsName("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",&flg);CHKERRQ(ierr);
265c2be2410SBarry Smith     if (flg) {
26697177400SBarry Smith       ierr = PCMGSetGalerkin(pc);CHKERRQ(ierr);
267c2be2410SBarry Smith     }
2689dcbbd2bSBarry Smith     ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr);
2694b9ad928SBarry Smith     if (flg) {
27097177400SBarry Smith       ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr);
2714b9ad928SBarry Smith     }
2729dcbbd2bSBarry Smith     ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr);
2734b9ad928SBarry Smith     if (flg) {
27497177400SBarry Smith       ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr);
2754b9ad928SBarry Smith     }
2769dcbbd2bSBarry Smith     ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
2779dcbbd2bSBarry Smith     if (flg) {ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);}
2784b9ad928SBarry Smith     ierr = PetscOptionsName("-pc_mg_log","Log times for each multigrid level","None",&flg);CHKERRQ(ierr);
2794b9ad928SBarry Smith     if (flg) {
2804f5ab15aSBarry Smith       PetscInt i;
2814b9ad928SBarry Smith       char     eventname[128];
2824b9ad928SBarry Smith       if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
2834b9ad928SBarry Smith       levels = mg[0]->levels;
2844b9ad928SBarry Smith       for (i=0; i<levels; i++) {
28577431f27SBarry Smith         sprintf(eventname,"MSetup Level %d",(int)i);
2864b9ad928SBarry Smith         ierr = PetscLogEventRegister(&mg[i]->eventsetup,eventname,pc->cookie);CHKERRQ(ierr);
28725a62d46SBarry Smith         sprintf(eventname,"MGSolve Level %d to 0",(int)i);
2884b9ad928SBarry Smith         ierr = PetscLogEventRegister(&mg[i]->eventsolve,eventname,pc->cookie);CHKERRQ(ierr);
2894b9ad928SBarry Smith       }
2904b9ad928SBarry Smith     }
2914b9ad928SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
2924b9ad928SBarry Smith   PetscFunctionReturn(0);
2934b9ad928SBarry Smith }
2944b9ad928SBarry Smith 
2959dcbbd2bSBarry Smith const char *PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0};
2969dcbbd2bSBarry Smith 
2974b9ad928SBarry Smith #undef __FUNCT__
2984b9ad928SBarry Smith #define __FUNCT__ "PCView_MG"
2996849ba73SBarry Smith static PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
3004b9ad928SBarry Smith {
3019dcbbd2bSBarry Smith   PC_MG          **mg = (PC_MG**)pc->data;
302dfbe8321SBarry Smith   PetscErrorCode ierr;
30379416396SBarry Smith   PetscInt       levels = mg[0]->levels,i;
30432077d6dSBarry Smith   PetscTruth     iascii;
3054b9ad928SBarry Smith 
3064b9ad928SBarry Smith   PetscFunctionBegin;
30732077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
30832077d6dSBarry Smith   if (iascii) {
30977431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  MG: type is %s, levels=%D cycles=%D, pre-smooths=%D, post-smooths=%D\n",
3109dcbbd2bSBarry Smith                       PCMGTypes[mg[0]->am],levels,mg[0]->cycles,mg[0]->default_smoothd,mg[0]->default_smoothu);CHKERRQ(ierr);
311c2be2410SBarry Smith     if (mg[0]->galerkin) {
312c2be2410SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
313c2be2410SBarry Smith     }
3144b9ad928SBarry Smith     for (i=0; i<levels; i++) {
315b03c7568SBarry Smith       if (!i) {
316b03c7568SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Coarse gride solver -- level %D -------------------------------\n",i);CHKERRQ(ierr);
317b03c7568SBarry Smith       } else {
31877431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
319b03c7568SBarry Smith       }
3204b9ad928SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
3214b9ad928SBarry Smith       ierr = KSPView(mg[i]->smoothd,viewer);CHKERRQ(ierr);
3224b9ad928SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
323b03c7568SBarry Smith       if (i && mg[i]->smoothd == mg[i]->smoothu) {
3244b9ad928SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
325b03c7568SBarry Smith       } else if (i){
32677431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
3274b9ad928SBarry Smith         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
3284b9ad928SBarry Smith         ierr = KSPView(mg[i]->smoothu,viewer);CHKERRQ(ierr);
3294b9ad928SBarry Smith         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
3304b9ad928SBarry Smith       }
3314b9ad928SBarry Smith     }
3324b9ad928SBarry Smith   } else {
33379a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCMG",((PetscObject)viewer)->type_name);
3344b9ad928SBarry Smith   }
3354b9ad928SBarry Smith   PetscFunctionReturn(0);
3364b9ad928SBarry Smith }
3374b9ad928SBarry Smith 
3384b9ad928SBarry Smith /*
3394b9ad928SBarry Smith     Calls setup for the KSP on each level
3404b9ad928SBarry Smith */
3414b9ad928SBarry Smith #undef __FUNCT__
3424b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_MG"
3436849ba73SBarry Smith static PetscErrorCode PCSetUp_MG(PC pc)
3444b9ad928SBarry Smith {
3459dcbbd2bSBarry Smith   PC_MG          **mg = (PC_MG**)pc->data;
346dfbe8321SBarry Smith   PetscErrorCode ierr;
34779416396SBarry Smith   PetscInt       i,n = mg[0]->levels;
3484b9ad928SBarry Smith   PC             cpc;
34968eff7e6SBarry Smith   PetscTruth     preonly,lu,redundant,cholesky,monitor = PETSC_FALSE,dump;
3504b9ad928SBarry Smith   PetscViewer    ascii;
3514b9ad928SBarry Smith   MPI_Comm       comm;
352c2be2410SBarry Smith   Mat            dA,dB;
353c2be2410SBarry Smith   MatStructure   uflag;
3540a6bb862SBarry Smith   Vec            tvec;
3554b9ad928SBarry Smith 
3564b9ad928SBarry Smith   PetscFunctionBegin;
357958c9bccSBarry Smith   if (!pc->setupcalled) {
3584b9ad928SBarry Smith     ierr = PetscOptionsHasName(0,"-pc_mg_monitor",&monitor);CHKERRQ(ierr);
3594b9ad928SBarry Smith 
360b03c7568SBarry Smith     for (i=0; i<n; i++) {
3614b9ad928SBarry Smith       if (monitor) {
3624b9ad928SBarry Smith         ierr = PetscObjectGetComm((PetscObject)mg[i]->smoothd,&comm);CHKERRQ(ierr);
3634b9ad928SBarry Smith         ierr = PetscViewerASCIIOpen(comm,"stdout",&ascii);CHKERRQ(ierr);
3644b9ad928SBarry Smith         ierr = PetscViewerASCIISetTab(ascii,n-i);CHKERRQ(ierr);
3656849ba73SBarry Smith         ierr = KSPSetMonitor(mg[i]->smoothd,KSPDefaultMonitor,ascii,(PetscErrorCode(*)(void*))PetscViewerDestroy);CHKERRQ(ierr);
3664b9ad928SBarry Smith       }
3674b9ad928SBarry Smith       ierr = KSPSetFromOptions(mg[i]->smoothd);CHKERRQ(ierr);
3684b9ad928SBarry Smith     }
369b03c7568SBarry Smith     for (i=1; i<n; i++) {
3704b9ad928SBarry Smith       if (mg[i]->smoothu && mg[i]->smoothu != mg[i]->smoothd) {
3714b9ad928SBarry Smith         if (monitor) {
3724b9ad928SBarry Smith           ierr = PetscObjectGetComm((PetscObject)mg[i]->smoothu,&comm);CHKERRQ(ierr);
3734b9ad928SBarry Smith           ierr = PetscViewerASCIIOpen(comm,"stdout",&ascii);CHKERRQ(ierr);
3744b9ad928SBarry Smith           ierr = PetscViewerASCIISetTab(ascii,n-i);CHKERRQ(ierr);
3756849ba73SBarry Smith           ierr = KSPSetMonitor(mg[i]->smoothu,KSPDefaultMonitor,ascii,(PetscErrorCode(*)(void*))PetscViewerDestroy);CHKERRQ(ierr);
3764b9ad928SBarry Smith         }
3774b9ad928SBarry Smith         ierr = KSPSetFromOptions(mg[i]->smoothu);CHKERRQ(ierr);
3784b9ad928SBarry Smith       }
3794b9ad928SBarry Smith     }
380fccaa45eSBarry Smith     for (i=1; i<n; i++) {
381fccaa45eSBarry Smith       if (mg[i]->restrct && !mg[i]->interpolate) {
38297177400SBarry Smith         ierr = PCMGSetInterpolate(pc,i,mg[i]->restrct);CHKERRQ(ierr);
383fccaa45eSBarry Smith       }
384fccaa45eSBarry Smith       if (!mg[i]->restrct && mg[i]->interpolate) {
38597177400SBarry Smith         ierr = PCMGSetRestriction(pc,i,mg[i]->interpolate);CHKERRQ(ierr);
386fccaa45eSBarry Smith       }
387fccaa45eSBarry Smith #if defined(PETSC_USE_DEBUG)
388fccaa45eSBarry Smith       if (!mg[i]->restrct || !mg[i]->interpolate) {
389fccaa45eSBarry Smith         SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Need to set restriction or interpolation on level %d",(int)i);
390fccaa45eSBarry Smith       }
391fccaa45eSBarry Smith #endif
392fccaa45eSBarry Smith     }
3930a6bb862SBarry Smith     for (i=0; i<n-1; i++) {
3946ca4d86aSHong Zhang       if (!mg[i]->r && i) {
3950a6bb862SBarry Smith         ierr = VecDuplicate(mg[i]->b,&tvec);CHKERRQ(ierr);
39697177400SBarry Smith         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
3970a6bb862SBarry Smith         ierr = VecDestroy(tvec);CHKERRQ(ierr);
3980a6bb862SBarry Smith       }
3990a6bb862SBarry Smith       if (!mg[i]->x) {
4000a6bb862SBarry Smith         ierr = VecDuplicate(mg[i]->b,&tvec);CHKERRQ(ierr);
40197177400SBarry Smith         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
4020a6bb862SBarry Smith         ierr = VecDestroy(tvec);CHKERRQ(ierr);
4030a6bb862SBarry Smith       }
4040a6bb862SBarry Smith     }
4054b9ad928SBarry Smith   }
4064b9ad928SBarry Smith 
407c2be2410SBarry Smith   /* If user did not provide fine grid operators, use those from PC */
408c2be2410SBarry Smith   /* BUG BUG BUG This will work ONLY the first time called: hence if the user changes
409c2be2410SBarry Smith      the PC matrices between solves PCMG will continue to use first set provided */
410c2be2410SBarry Smith   ierr = KSPGetOperators(mg[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr);
411c2be2410SBarry Smith   if (!dA  && !dB) {
41297177400SBarry Smith     ierr = PetscLogInfo((pc,"PCSetUp_MG: Using outer operators to define finest grid operator \n  because PCMGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n"));
413c2be2410SBarry Smith     ierr = KSPSetOperators(mg[n-1]->smoothd,pc->mat,pc->pmat,uflag);CHKERRQ(ierr);
414c2be2410SBarry Smith   }
415c2be2410SBarry Smith 
416c2be2410SBarry Smith   if (mg[0]->galerkin) {
4175e8befabSKris Buschelman     Mat B;
418c2be2410SBarry Smith     mg[0]->galerkinused = PETSC_TRUE;
419c2be2410SBarry Smith     /* currently only handle case where mat and pmat are the same on coarser levels */
420c2be2410SBarry Smith     ierr = KSPGetOperators(mg[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr);
421c2be2410SBarry Smith     if (!pc->setupcalled) {
422c2be2410SBarry Smith       for (i=n-2; i>-1; i--) {
423c2be2410SBarry Smith         ierr = MatPtAP(dB,mg[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
424c2be2410SBarry Smith         ierr = KSPSetOperators(mg[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
425c2be2410SBarry Smith         dB   = B;
426c2be2410SBarry Smith       }
427c2be2410SBarry Smith     } else {
428c2be2410SBarry Smith       for (i=n-2; i>-1; i--) {
429c2be2410SBarry Smith         ierr = KSPGetOperators(mg[i]->smoothd,0,&B,0);CHKERRQ(ierr);
430c2be2410SBarry Smith         ierr = MatPtAP(dB,mg[i]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
431c2be2410SBarry Smith         ierr = KSPSetOperators(mg[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
432c2be2410SBarry Smith         dB   = B;
433c2be2410SBarry Smith       }
434c2be2410SBarry Smith     }
435c2be2410SBarry Smith   }
436c2be2410SBarry Smith 
4374b9ad928SBarry Smith   for (i=1; i<n; i++) {
438b03c7568SBarry Smith     if (mg[i]->smoothu == mg[i]->smoothd) {
439b03c7568SBarry Smith       /* if doing only down then initial guess is zero */
4404b9ad928SBarry Smith       ierr = KSPSetInitialGuessNonzero(mg[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
441b03c7568SBarry Smith     }
4424b9ad928SBarry Smith     if (mg[i]->eventsetup) {ierr = PetscLogEventBegin(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);}
4434b9ad928SBarry Smith     ierr = KSPSetUp(mg[i]->smoothd);CHKERRQ(ierr);
4444b9ad928SBarry Smith     if (mg[i]->eventsetup) {ierr = PetscLogEventEnd(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);}
4454b9ad928SBarry Smith   }
446b03c7568SBarry Smith   for (i=1; i<n; i++) {
4474b9ad928SBarry Smith     if (mg[i]->smoothu && mg[i]->smoothu != mg[i]->smoothd) {
448c2be2410SBarry Smith       PC           uppc,downpc;
44997f1f81fSBarry Smith       Mat          downmat,downpmat,upmat,uppmat;
45097f1f81fSBarry Smith       MatStructure matflag;
45197f1f81fSBarry Smith 
45297f1f81fSBarry Smith       /* check if operators have been set for up, if not use down operators to set them */
45397f1f81fSBarry Smith       ierr = KSPGetPC(mg[i]->smoothu,&uppc);CHKERRQ(ierr);
45497f1f81fSBarry Smith       ierr = PCGetOperators(uppc,&upmat,&uppmat,PETSC_NULL);CHKERRQ(ierr);
45597f1f81fSBarry Smith       if (!upmat) {
45697f1f81fSBarry Smith         ierr = KSPGetPC(mg[i]->smoothd,&downpc);CHKERRQ(ierr);
45797f1f81fSBarry Smith         ierr = PCGetOperators(downpc,&downmat,&downpmat,&matflag);CHKERRQ(ierr);
45897f1f81fSBarry Smith         ierr = KSPSetOperators(mg[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr);
45997f1f81fSBarry Smith       }
46097f1f81fSBarry Smith 
4614b9ad928SBarry Smith       ierr = KSPSetInitialGuessNonzero(mg[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
4624b9ad928SBarry Smith       if (mg[i]->eventsetup) {ierr = PetscLogEventBegin(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);}
4634b9ad928SBarry Smith       ierr = KSPSetUp(mg[i]->smoothu);CHKERRQ(ierr);
4644b9ad928SBarry Smith       if (mg[i]->eventsetup) {ierr = PetscLogEventEnd(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);}
4654b9ad928SBarry Smith     }
4664b9ad928SBarry Smith   }
4674b9ad928SBarry Smith 
4684b9ad928SBarry Smith   /*
4694b9ad928SBarry Smith       If coarse solver is not direct method then DO NOT USE preonly
4704b9ad928SBarry Smith   */
4714b9ad928SBarry Smith   ierr = PetscTypeCompare((PetscObject)mg[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr);
4724b9ad928SBarry Smith   if (preonly) {
4734b9ad928SBarry Smith     ierr = KSPGetPC(mg[0]->smoothd,&cpc);CHKERRQ(ierr);
4744b9ad928SBarry Smith     ierr = PetscTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr);
4754b9ad928SBarry Smith     ierr = PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr);
47668eff7e6SBarry Smith     ierr = PetscTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr);
47768eff7e6SBarry Smith     if (!lu && !redundant && !cholesky) {
4784b9ad928SBarry Smith       ierr = KSPSetType(mg[0]->smoothd,KSPGMRES);CHKERRQ(ierr);
4794b9ad928SBarry Smith     }
4804b9ad928SBarry Smith   }
4814b9ad928SBarry Smith 
482958c9bccSBarry Smith   if (!pc->setupcalled) {
4834b9ad928SBarry Smith     if (monitor) {
4844b9ad928SBarry Smith       ierr = PetscObjectGetComm((PetscObject)mg[0]->smoothd,&comm);CHKERRQ(ierr);
4854b9ad928SBarry Smith       ierr = PetscViewerASCIIOpen(comm,"stdout",&ascii);CHKERRQ(ierr);
4864b9ad928SBarry Smith       ierr = PetscViewerASCIISetTab(ascii,n);CHKERRQ(ierr);
4876849ba73SBarry Smith       ierr = KSPSetMonitor(mg[0]->smoothd,KSPDefaultMonitor,ascii,(PetscErrorCode(*)(void*))PetscViewerDestroy);CHKERRQ(ierr);
4884b9ad928SBarry Smith     }
4894b9ad928SBarry Smith     ierr = KSPSetFromOptions(mg[0]->smoothd);CHKERRQ(ierr);
4904b9ad928SBarry Smith   }
4914b9ad928SBarry Smith 
4924b9ad928SBarry Smith   if (mg[0]->eventsetup) {ierr = PetscLogEventBegin(mg[0]->eventsetup,0,0,0,0);CHKERRQ(ierr);}
4934b9ad928SBarry Smith   ierr = KSPSetUp(mg[0]->smoothd);CHKERRQ(ierr);
4944b9ad928SBarry Smith   if (mg[0]->eventsetup) {ierr = PetscLogEventEnd(mg[0]->eventsetup,0,0,0,0);CHKERRQ(ierr);}
4954b9ad928SBarry Smith 
4964cf18111SSatish Balay #if defined(PETSC_USE_SOCKET_VIEWER)
4974b9ad928SBarry Smith   /*
4984b9ad928SBarry Smith      Dump the interpolation/restriction matrices to matlab plus the
4994b9ad928SBarry Smith    Jacobian/stiffness on each level. This allows Matlab users to
5004b9ad928SBarry Smith    easily check if the Galerkin condition A_c = R A_f R^T is satisfied */
5014b9ad928SBarry Smith   ierr = PetscOptionsHasName(pc->prefix,"-pc_mg_dump_matlab",&dump);CHKERRQ(ierr);
5024b9ad928SBarry Smith   if (dump) {
5034b9ad928SBarry Smith     for (i=1; i<n; i++) {
5044b9ad928SBarry Smith       ierr = MatView(mg[i]->restrct,PETSC_VIEWER_SOCKET_(pc->comm));CHKERRQ(ierr);
5054b9ad928SBarry Smith     }
5064b9ad928SBarry Smith     for (i=0; i<n; i++) {
5074b9ad928SBarry Smith       ierr = KSPGetPC(mg[i]->smoothd,&pc);CHKERRQ(ierr);
5084b9ad928SBarry Smith       ierr = MatView(pc->mat,PETSC_VIEWER_SOCKET_(pc->comm));CHKERRQ(ierr);
5094b9ad928SBarry Smith     }
5104b9ad928SBarry Smith   }
511c45a1595SBarry Smith #endif
512c45a1595SBarry Smith 
513c2be2410SBarry Smith   ierr = PetscOptionsHasName(pc->prefix,"-pc_mg_dump_binary",&dump);CHKERRQ(ierr);
514c2be2410SBarry Smith   if (dump) {
515c2be2410SBarry Smith     for (i=1; i<n; i++) {
516c2be2410SBarry Smith       ierr = MatView(mg[i]->restrct,PETSC_VIEWER_BINARY_(pc->comm));CHKERRQ(ierr);
517c2be2410SBarry Smith     }
518c2be2410SBarry Smith     for (i=0; i<n; i++) {
519c2be2410SBarry Smith       ierr = KSPGetPC(mg[i]->smoothd,&pc);CHKERRQ(ierr);
520c2be2410SBarry Smith       ierr = MatView(pc->mat,PETSC_VIEWER_BINARY_(pc->comm));CHKERRQ(ierr);
521c2be2410SBarry Smith     }
522c2be2410SBarry Smith   }
5234b9ad928SBarry Smith   PetscFunctionReturn(0);
5244b9ad928SBarry Smith }
5254b9ad928SBarry Smith 
5264b9ad928SBarry Smith /* -------------------------------------------------------------------------------------*/
5274b9ad928SBarry Smith 
5284b9ad928SBarry Smith #undef __FUNCT__
5299dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetLevels"
5304b9ad928SBarry Smith /*@C
53197177400SBarry Smith    PCMGSetLevels - Sets the number of levels to use with MG.
5324b9ad928SBarry Smith    Must be called before any other MG routine.
5334b9ad928SBarry Smith 
5344b9ad928SBarry Smith    Collective on PC
5354b9ad928SBarry Smith 
5364b9ad928SBarry Smith    Input Parameters:
5374b9ad928SBarry Smith +  pc - the preconditioner context
5384b9ad928SBarry Smith .  levels - the number of levels
5394b9ad928SBarry Smith -  comms - optional communicators for each level; this is to allow solving the coarser problems
5404b9ad928SBarry Smith            on smaller sets of processors. Use PETSC_NULL_OBJECT for default in Fortran
5414b9ad928SBarry Smith 
5424b9ad928SBarry Smith    Level: intermediate
5434b9ad928SBarry Smith 
5444b9ad928SBarry Smith    Notes:
5454b9ad928SBarry Smith      If the number of levels is one then the multigrid uses the -mg_levels prefix
5464b9ad928SBarry Smith   for setting the level options rather than the -mg_coarse prefix.
5474b9ad928SBarry Smith 
5484b9ad928SBarry Smith .keywords: MG, set, levels, multigrid
5494b9ad928SBarry Smith 
55097177400SBarry Smith .seealso: PCMGSetType(), PCMGGetLevels()
5514b9ad928SBarry Smith @*/
55297177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
5534b9ad928SBarry Smith {
554dfbe8321SBarry Smith   PetscErrorCode ierr;
5559dcbbd2bSBarry Smith   PC_MG          **mg;
5564b9ad928SBarry Smith 
5574b9ad928SBarry Smith   PetscFunctionBegin;
5584482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
5594b9ad928SBarry Smith 
5604b9ad928SBarry Smith   if (pc->data) {
5611302d50aSBarry Smith     SETERRQ(PETSC_ERR_ORDER,"Number levels already set for MG\n\
56297177400SBarry Smith     make sure that you call PCMGSetLevels() before KSPSetFromOptions()");
5634b9ad928SBarry Smith   }
5649dcbbd2bSBarry Smith   ierr                     = PCMGCreate_Private(pc->comm,levels,pc,comms,&mg);CHKERRQ(ierr);
5659dcbbd2bSBarry Smith   mg[0]->am                = PC_MG_MULTIPLICATIVE;
5664b9ad928SBarry Smith   pc->data                 = (void*)mg;
5674b9ad928SBarry Smith   pc->ops->applyrichardson = PCApplyRichardson_MG;
5684b9ad928SBarry Smith   PetscFunctionReturn(0);
5694b9ad928SBarry Smith }
5704b9ad928SBarry Smith 
5714b9ad928SBarry Smith #undef __FUNCT__
5729dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetLevels"
5734b9ad928SBarry Smith /*@
57497177400SBarry Smith    PCMGGetLevels - Gets the number of levels to use with MG.
5754b9ad928SBarry Smith 
5764b9ad928SBarry Smith    Not Collective
5774b9ad928SBarry Smith 
5784b9ad928SBarry Smith    Input Parameter:
5794b9ad928SBarry Smith .  pc - the preconditioner context
5804b9ad928SBarry Smith 
5814b9ad928SBarry Smith    Output parameter:
5824b9ad928SBarry Smith .  levels - the number of levels
5834b9ad928SBarry Smith 
5844b9ad928SBarry Smith    Level: advanced
5854b9ad928SBarry Smith 
5864b9ad928SBarry Smith .keywords: MG, get, levels, multigrid
5874b9ad928SBarry Smith 
58897177400SBarry Smith .seealso: PCMGSetLevels()
5894b9ad928SBarry Smith @*/
59097177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetLevels(PC pc,PetscInt *levels)
5914b9ad928SBarry Smith {
5929dcbbd2bSBarry Smith   PC_MG  **mg;
5934b9ad928SBarry Smith 
5944b9ad928SBarry Smith   PetscFunctionBegin;
5954482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
5964482741eSBarry Smith   PetscValidIntPointer(levels,2);
5974b9ad928SBarry Smith 
5989dcbbd2bSBarry Smith   mg      = (PC_MG**)pc->data;
5994b9ad928SBarry Smith   *levels = mg[0]->levels;
6004b9ad928SBarry Smith   PetscFunctionReturn(0);
6014b9ad928SBarry Smith }
6024b9ad928SBarry Smith 
6034b9ad928SBarry Smith #undef __FUNCT__
6049dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetType"
6054b9ad928SBarry Smith /*@
60697177400SBarry Smith    PCMGSetType - Determines the form of multigrid to use:
6074b9ad928SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
6084b9ad928SBarry Smith 
6094b9ad928SBarry Smith    Collective on PC
6104b9ad928SBarry Smith 
6114b9ad928SBarry Smith    Input Parameters:
6124b9ad928SBarry Smith +  pc - the preconditioner context
6139dcbbd2bSBarry Smith -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
6149dcbbd2bSBarry Smith    PC_MG_FULL, PC_MG_KASKADE
6154b9ad928SBarry Smith 
6164b9ad928SBarry Smith    Options Database Key:
6174b9ad928SBarry Smith .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
6184b9ad928SBarry Smith    additive, full, kaskade
6194b9ad928SBarry Smith 
6204b9ad928SBarry Smith    Level: advanced
6214b9ad928SBarry Smith 
6224b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
6234b9ad928SBarry Smith 
62497177400SBarry Smith .seealso: PCMGSetLevels()
6254b9ad928SBarry Smith @*/
6269dcbbd2bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetType(PC pc,PCMGType form)
6274b9ad928SBarry Smith {
6289dcbbd2bSBarry Smith   PC_MG **mg;
6294b9ad928SBarry Smith 
6304b9ad928SBarry Smith   PetscFunctionBegin;
6314482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
6329dcbbd2bSBarry Smith   mg = (PC_MG**)pc->data;
6334b9ad928SBarry Smith 
6344b9ad928SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
6354b9ad928SBarry Smith   mg[0]->am = form;
6369dcbbd2bSBarry Smith   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
6374b9ad928SBarry Smith   else pc->ops->applyrichardson = 0;
6384b9ad928SBarry Smith   PetscFunctionReturn(0);
6394b9ad928SBarry Smith }
6404b9ad928SBarry Smith 
6414b9ad928SBarry Smith #undef __FUNCT__
6429dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetCycles"
6434b9ad928SBarry Smith /*@
64497177400SBarry Smith    PCMGSetCycles - Sets the type cycles to use.  Use PCMGSetCyclesOnLevel() for more
6454b9ad928SBarry Smith    complicated cycling.
6464b9ad928SBarry Smith 
6474b9ad928SBarry Smith    Collective on PC
6484b9ad928SBarry Smith 
6494b9ad928SBarry Smith    Input Parameters:
650c2be2410SBarry Smith +  pc - the multigrid context
6514b9ad928SBarry Smith -  n - the number of cycles
6524b9ad928SBarry Smith 
6534b9ad928SBarry Smith    Options Database Key:
6544b9ad928SBarry Smith $  -pc_mg_cycles n - 1 denotes a V-cycle; 2 denotes a W-cycle.
6554b9ad928SBarry Smith 
6564b9ad928SBarry Smith    Level: advanced
6574b9ad928SBarry Smith 
6584b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
6594b9ad928SBarry Smith 
66097177400SBarry Smith .seealso: PCMGSetCyclesOnLevel()
6614b9ad928SBarry Smith @*/
66297177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetCycles(PC pc,PetscInt n)
6634b9ad928SBarry Smith {
6649dcbbd2bSBarry Smith   PC_MG    **mg;
66579416396SBarry Smith   PetscInt i,levels;
6664b9ad928SBarry Smith 
6674b9ad928SBarry Smith   PetscFunctionBegin;
6684482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
6699dcbbd2bSBarry Smith   mg     = (PC_MG**)pc->data;
6704b9ad928SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
6714b9ad928SBarry Smith   levels = mg[0]->levels;
6724b9ad928SBarry Smith 
6734b9ad928SBarry Smith   for (i=0; i<levels; i++) {
6744b9ad928SBarry Smith     mg[i]->cycles  = n;
6754b9ad928SBarry Smith   }
6764b9ad928SBarry Smith   PetscFunctionReturn(0);
6774b9ad928SBarry Smith }
6784b9ad928SBarry Smith 
6794b9ad928SBarry Smith #undef __FUNCT__
6809dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetGalerkin"
681c2be2410SBarry Smith /*@
68297177400SBarry Smith    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
683c2be2410SBarry Smith       finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t
684c2be2410SBarry Smith 
685c2be2410SBarry Smith    Collective on PC
686c2be2410SBarry Smith 
687c2be2410SBarry Smith    Input Parameters:
688c2be2410SBarry Smith +  pc - the multigrid context
689c2be2410SBarry Smith -  n - the number of cycles
690c2be2410SBarry Smith 
691c2be2410SBarry Smith    Options Database Key:
692c2be2410SBarry Smith $  -pc_mg_galerkin
693c2be2410SBarry Smith 
694c2be2410SBarry Smith    Level: intermediate
695c2be2410SBarry Smith 
696c2be2410SBarry Smith .keywords: MG, set, Galerkin
697c2be2410SBarry Smith 
698c2be2410SBarry Smith @*/
69997177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetGalerkin(PC pc)
700c2be2410SBarry Smith {
7019dcbbd2bSBarry Smith   PC_MG    **mg;
702c2be2410SBarry Smith   PetscInt i,levels;
703c2be2410SBarry Smith 
704c2be2410SBarry Smith   PetscFunctionBegin;
705c2be2410SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
7069dcbbd2bSBarry Smith   mg     = (PC_MG**)pc->data;
707c2be2410SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
708c2be2410SBarry Smith   levels = mg[0]->levels;
709c2be2410SBarry Smith 
710c2be2410SBarry Smith   for (i=0; i<levels; i++) {
711c2be2410SBarry Smith     mg[i]->galerkin = PETSC_TRUE;
712c2be2410SBarry Smith   }
713c2be2410SBarry Smith   PetscFunctionReturn(0);
714c2be2410SBarry Smith }
715c2be2410SBarry Smith 
716c2be2410SBarry Smith #undef __FUNCT__
7179dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothDown"
7184b9ad928SBarry Smith /*@
71997177400SBarry Smith    PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
72097177400SBarry Smith    use on all levels. Use PCMGGetSmootherDown() to set different
7214b9ad928SBarry Smith    pre-smoothing steps on different levels.
7224b9ad928SBarry Smith 
7234b9ad928SBarry Smith    Collective on PC
7244b9ad928SBarry Smith 
7254b9ad928SBarry Smith    Input Parameters:
7264b9ad928SBarry Smith +  mg - the multigrid context
7274b9ad928SBarry Smith -  n - the number of smoothing steps
7284b9ad928SBarry Smith 
7294b9ad928SBarry Smith    Options Database Key:
7304b9ad928SBarry Smith .  -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps
7314b9ad928SBarry Smith 
7324b9ad928SBarry Smith    Level: advanced
7334b9ad928SBarry Smith 
7344b9ad928SBarry Smith .keywords: MG, smooth, down, pre-smoothing, steps, multigrid
7354b9ad928SBarry Smith 
73697177400SBarry Smith .seealso: PCMGSetNumberSmoothUp()
7374b9ad928SBarry Smith @*/
73897177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothDown(PC pc,PetscInt n)
7394b9ad928SBarry Smith {
7409dcbbd2bSBarry Smith   PC_MG          **mg;
7416849ba73SBarry Smith   PetscErrorCode ierr;
74279416396SBarry Smith   PetscInt       i,levels;
7434b9ad928SBarry Smith 
7444b9ad928SBarry Smith   PetscFunctionBegin;
7454482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
7469dcbbd2bSBarry Smith   mg     = (PC_MG**)pc->data;
7474b9ad928SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
7484b9ad928SBarry Smith   levels = mg[0]->levels;
7494b9ad928SBarry Smith 
7504b9ad928SBarry Smith   for (i=0; i<levels; i++) {
7514b9ad928SBarry Smith     /* make sure smoother up and down are different */
75297177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
7534b9ad928SBarry Smith     ierr = KSPSetTolerances(mg[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
7544b9ad928SBarry Smith     mg[i]->default_smoothd = n;
7554b9ad928SBarry Smith   }
7564b9ad928SBarry Smith   PetscFunctionReturn(0);
7574b9ad928SBarry Smith }
7584b9ad928SBarry Smith 
7594b9ad928SBarry Smith #undef __FUNCT__
7609dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothUp"
7614b9ad928SBarry Smith /*@
76297177400SBarry Smith    PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
76397177400SBarry Smith    on all levels. Use PCMGGetSmootherUp() to set different numbers of
7644b9ad928SBarry Smith    post-smoothing steps on different levels.
7654b9ad928SBarry Smith 
7664b9ad928SBarry Smith    Collective on PC
7674b9ad928SBarry Smith 
7684b9ad928SBarry Smith    Input Parameters:
7694b9ad928SBarry Smith +  mg - the multigrid context
7704b9ad928SBarry Smith -  n - the number of smoothing steps
7714b9ad928SBarry Smith 
7724b9ad928SBarry Smith    Options Database Key:
7734b9ad928SBarry Smith .  -pc_mg_smoothup <n> - Sets number of post-smoothing steps
7744b9ad928SBarry Smith 
7754b9ad928SBarry Smith    Level: advanced
7764b9ad928SBarry Smith 
7774b9ad928SBarry Smith    Note: this does not set a value on the coarsest grid, since we assume that
778b03c7568SBarry Smith     there is no seperate smooth up on the coarsest grid.
7794b9ad928SBarry Smith 
7804b9ad928SBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid
7814b9ad928SBarry Smith 
78297177400SBarry Smith .seealso: PCMGSetNumberSmoothDown()
7834b9ad928SBarry Smith @*/
78497177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothUp(PC pc,PetscInt n)
7854b9ad928SBarry Smith {
7869dcbbd2bSBarry Smith   PC_MG          **mg;
7876849ba73SBarry Smith   PetscErrorCode ierr;
78879416396SBarry Smith   PetscInt       i,levels;
7894b9ad928SBarry Smith 
7904b9ad928SBarry Smith   PetscFunctionBegin;
7914482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
7929dcbbd2bSBarry Smith   mg     = (PC_MG**)pc->data;
7934b9ad928SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
7944b9ad928SBarry Smith   levels = mg[0]->levels;
7954b9ad928SBarry Smith 
7964b9ad928SBarry Smith   for (i=1; i<levels; i++) {
7974b9ad928SBarry Smith     /* make sure smoother up and down are different */
79897177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
7994b9ad928SBarry Smith     ierr = KSPSetTolerances(mg[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
8004b9ad928SBarry Smith     mg[i]->default_smoothu = n;
8014b9ad928SBarry Smith   }
8024b9ad928SBarry Smith   PetscFunctionReturn(0);
8034b9ad928SBarry Smith }
8044b9ad928SBarry Smith 
8054b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/
8064b9ad928SBarry Smith 
8073b09bd56SBarry Smith /*MC
8083b09bd56SBarry Smith    PCMG - Use geometric multigrid preconditioning. This preconditioner requires you provide additional
8093b09bd56SBarry Smith     information about the coarser grid matrices and restriction/interpolation operators.
8103b09bd56SBarry Smith 
8113b09bd56SBarry Smith    Options Database Keys:
8123b09bd56SBarry Smith +  -pc_mg_levels <nlevels> - number of levels including finest
8133b09bd56SBarry Smith .  -pc_mg_cycles 1 or 2 - for V or W-cycle
81479416396SBarry Smith .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
8153b09bd56SBarry Smith .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
8163b09bd56SBarry Smith .  -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default
8173b09bd56SBarry Smith .  -pc_mg_log - log information about time spent on each level of the solver
8183b09bd56SBarry Smith .  -pc_mg_monitor - print information on the multigrid convergence
81968eff7e6SBarry Smith .  -pc_mg_galerkin - use Galerkin process to compute coarser operators
8203b09bd56SBarry Smith -  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
8213b09bd56SBarry Smith                         to the Socket viewer for reading from Matlab.
8223b09bd56SBarry Smith 
8233b09bd56SBarry Smith    Notes:
8243b09bd56SBarry Smith 
8253b09bd56SBarry Smith    Level: intermediate
8263b09bd56SBarry Smith 
8273b09bd56SBarry Smith    Concepts: multigrid
8283b09bd56SBarry Smith 
8293b09bd56SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
8309dcbbd2bSBarry Smith            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycles(), PCMGSetNumberSmoothDown(),
83197177400SBarry Smith            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
83297177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
83397177400SBarry Smith            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
8343b09bd56SBarry Smith M*/
8353b09bd56SBarry Smith 
8364b9ad928SBarry Smith EXTERN_C_BEGIN
8374b9ad928SBarry Smith #undef __FUNCT__
8384b9ad928SBarry Smith #define __FUNCT__ "PCCreate_MG"
839dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_MG(PC pc)
8404b9ad928SBarry Smith {
8414b9ad928SBarry Smith   PetscFunctionBegin;
8424b9ad928SBarry Smith   pc->ops->apply          = PCApply_MG;
8434b9ad928SBarry Smith   pc->ops->setup          = PCSetUp_MG;
8444b9ad928SBarry Smith   pc->ops->destroy        = PCDestroy_MG;
8454b9ad928SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_MG;
8464b9ad928SBarry Smith   pc->ops->view           = PCView_MG;
8474b9ad928SBarry Smith 
8484b9ad928SBarry Smith   pc->data                = (void*)0;
8494b9ad928SBarry Smith   PetscFunctionReturn(0);
8504b9ad928SBarry Smith }
8514b9ad928SBarry Smith EXTERN_C_END
852