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; 134c2be2410SBarry Smith if (mg[0]->galerkinused) { 135c2be2410SBarry Smith Mat B; 136c2be2410SBarry Smith for (i=0; i<n-1; i++) { 137c2be2410SBarry Smith ierr = KSPGetOperators(mg[i]->smoothd,0,&B,0);CHKERRQ(ierr); 138c2be2410SBarry Smith ierr = MatDestroy(B);CHKERRQ(ierr); 139c2be2410SBarry Smith } 140c2be2410SBarry Smith } 141c2be2410SBarry Smith 142fccaa45eSBarry Smith for (i=0; i<n-1; i++) { 143630ba2eeSBarry Smith if (mg[i+1]->r) {ierr = VecDestroy(mg[i+1]->r);CHKERRQ(ierr);} 144fccaa45eSBarry Smith if (mg[i]->b) {ierr = VecDestroy(mg[i]->b);CHKERRQ(ierr);} 145fccaa45eSBarry Smith if (mg[i]->x) {ierr = VecDestroy(mg[i]->x);CHKERRQ(ierr);} 146fccaa45eSBarry Smith if (mg[i+1]->restrct) {ierr = MatDestroy(mg[i+1]->restrct);CHKERRQ(ierr);} 147fccaa45eSBarry Smith if (mg[i+1]->interpolate) {ierr = MatDestroy(mg[i+1]->interpolate);CHKERRQ(ierr);} 148fccaa45eSBarry Smith } 149fccaa45eSBarry Smith 1504b9ad928SBarry Smith for (i=0; i<n; i++) { 1514b9ad928SBarry Smith if (mg[i]->smoothd != mg[i]->smoothu) { 1524b9ad928SBarry Smith ierr = KSPDestroy(mg[i]->smoothd);CHKERRQ(ierr); 1534b9ad928SBarry Smith } 1544b9ad928SBarry Smith ierr = KSPDestroy(mg[i]->smoothu);CHKERRQ(ierr); 1554b9ad928SBarry Smith ierr = PetscFree(mg[i]);CHKERRQ(ierr); 1564b9ad928SBarry Smith } 1574b9ad928SBarry Smith ierr = PetscFree(mg);CHKERRQ(ierr); 1584b9ad928SBarry Smith PetscFunctionReturn(0); 1594b9ad928SBarry Smith } 1604b9ad928SBarry Smith 1614b9ad928SBarry Smith 1624b9ad928SBarry Smith 1639dcbbd2bSBarry Smith EXTERN PetscErrorCode PCMGACycle_Private(PC_MG**); 1649dcbbd2bSBarry Smith EXTERN PetscErrorCode PCMGFCycle_Private(PC_MG**); 1659dcbbd2bSBarry Smith EXTERN PetscErrorCode PCMGKCycle_Private(PC_MG**); 1664b9ad928SBarry Smith 1674b9ad928SBarry Smith /* 1684b9ad928SBarry Smith PCApply_MG - Runs either an additive, multiplicative, Kaskadic 1694b9ad928SBarry Smith or full cycle of multigrid. 1704b9ad928SBarry Smith 1714b9ad928SBarry Smith Note: 1729dcbbd2bSBarry Smith A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle(). 1734b9ad928SBarry Smith */ 1744b9ad928SBarry Smith #undef __FUNCT__ 1754b9ad928SBarry Smith #define __FUNCT__ "PCApply_MG" 1766849ba73SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x) 1774b9ad928SBarry Smith { 1789dcbbd2bSBarry Smith PC_MG **mg = (PC_MG**)pc->data; 1796849ba73SBarry Smith PetscErrorCode ierr; 18079416396SBarry Smith PetscInt levels = mg[0]->levels; 1814b9ad928SBarry Smith 1824b9ad928SBarry Smith PetscFunctionBegin; 1834b9ad928SBarry Smith mg[levels-1]->b = b; 1844b9ad928SBarry Smith mg[levels-1]->x = x; 185ef70c39aSMatthew Knepley if (!mg[levels-1]->r && mg[0]->am != PC_MG_ADDITIVE && levels > 1) { 1860a6bb862SBarry Smith Vec tvec; 1870a6bb862SBarry Smith ierr = VecDuplicate(mg[levels-1]->b,&tvec);CHKERRQ(ierr); 18897177400SBarry Smith ierr = PCMGSetR(pc,levels-1,tvec);CHKERRQ(ierr); 1890a6bb862SBarry Smith ierr = VecDestroy(tvec);CHKERRQ(ierr); 1900a6bb862SBarry Smith } 1919dcbbd2bSBarry Smith if (mg[0]->am == PC_MG_MULTIPLICATIVE) { 192efb30889SBarry Smith ierr = VecSet(x,0.0);CHKERRQ(ierr); 1939dcbbd2bSBarry Smith ierr = PCMGMCycle_Private(mg+levels-1,PETSC_NULL);CHKERRQ(ierr); 1944b9ad928SBarry Smith } 1959dcbbd2bSBarry Smith else if (mg[0]->am == PC_MG_ADDITIVE) { 1969dcbbd2bSBarry Smith ierr = PCMGACycle_Private(mg);CHKERRQ(ierr); 1974b9ad928SBarry Smith } 1989dcbbd2bSBarry Smith else if (mg[0]->am == PC_MG_KASKADE) { 1999dcbbd2bSBarry Smith ierr = PCMGKCycle_Private(mg);CHKERRQ(ierr); 2004b9ad928SBarry Smith } 2014b9ad928SBarry Smith else { 2029dcbbd2bSBarry Smith ierr = PCMGFCycle_Private(mg);CHKERRQ(ierr); 2034b9ad928SBarry Smith } 2044b9ad928SBarry Smith PetscFunctionReturn(0); 2054b9ad928SBarry Smith } 2064b9ad928SBarry Smith 2074b9ad928SBarry Smith #undef __FUNCT__ 2084b9ad928SBarry Smith #define __FUNCT__ "PCApplyRichardson_MG" 20979416396SBarry Smith static PetscErrorCode PCApplyRichardson_MG(PC pc,Vec b,Vec x,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its) 2104b9ad928SBarry Smith { 2119dcbbd2bSBarry Smith PC_MG **mg = (PC_MG**)pc->data; 212dfbe8321SBarry Smith PetscErrorCode ierr; 21379416396SBarry Smith PetscInt levels = mg[0]->levels; 2144b9ad928SBarry Smith PetscTruth converged = PETSC_FALSE; 2154b9ad928SBarry Smith 2164b9ad928SBarry Smith PetscFunctionBegin; 2174b9ad928SBarry Smith mg[levels-1]->b = b; 2184b9ad928SBarry Smith mg[levels-1]->x = x; 2194b9ad928SBarry Smith 2204b9ad928SBarry Smith mg[levels-1]->rtol = rtol; 22170441072SBarry Smith mg[levels-1]->abstol = abstol; 2224b9ad928SBarry Smith mg[levels-1]->dtol = dtol; 2234b9ad928SBarry Smith if (rtol) { 2244b9ad928SBarry Smith /* compute initial residual norm for relative convergence test */ 2254b9ad928SBarry Smith PetscReal rnorm; 2264b9ad928SBarry Smith ierr = (*mg[levels-1]->residual)(mg[levels-1]->A,b,x,w);CHKERRQ(ierr); 2274b9ad928SBarry Smith ierr = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr); 22870441072SBarry Smith mg[levels-1]->ttol = PetscMax(rtol*rnorm,abstol); 22970441072SBarry Smith } else if (abstol) { 23070441072SBarry Smith mg[levels-1]->ttol = abstol; 2314b9ad928SBarry Smith } else { 2324b9ad928SBarry Smith mg[levels-1]->ttol = 0.0; 2334b9ad928SBarry Smith } 2344b9ad928SBarry Smith 2354b9ad928SBarry Smith while (its-- && !converged) { 2369dcbbd2bSBarry Smith ierr = PCMGMCycle_Private(mg+levels-1,&converged);CHKERRQ(ierr); 2374b9ad928SBarry Smith } 2384b9ad928SBarry Smith PetscFunctionReturn(0); 2394b9ad928SBarry Smith } 2404b9ad928SBarry Smith 2414b9ad928SBarry Smith #undef __FUNCT__ 2424b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_MG" 2436ca4d86aSHong Zhang PetscErrorCode PCSetFromOptions_MG(PC pc) 2444b9ad928SBarry Smith { 245dfbe8321SBarry Smith PetscErrorCode ierr; 2469dcbbd2bSBarry Smith PetscInt m,levels = 1; 2474b9ad928SBarry Smith PetscTruth flg; 2489dcbbd2bSBarry Smith PC_MG **mg = (PC_MG**)pc->data; 249*cf502942SBarry Smith PCMGType mgtype; 2504b9ad928SBarry Smith 2514b9ad928SBarry Smith PetscFunctionBegin; 2524b9ad928SBarry Smith 2534b9ad928SBarry Smith ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr); 2544b9ad928SBarry Smith if (!pc->data) { 2559dcbbd2bSBarry Smith ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr); 25697177400SBarry Smith ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr); 257*cf502942SBarry Smith mg = (PC_MG**)pc->data; 2584b9ad928SBarry Smith } 259*cf502942SBarry Smith mgtype = mg[0]->am; 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++) { 370a98fc643SBarry 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++) { 3810cace4b0SBarry Smith if (!mg[i]->residual) { 3820cace4b0SBarry Smith Mat mat; 3830cace4b0SBarry Smith ierr = KSPGetOperators(mg[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr); 3840cace4b0SBarry Smith ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr); 3850cace4b0SBarry Smith } 386fccaa45eSBarry Smith if (mg[i]->restrct && !mg[i]->interpolate) { 38797177400SBarry Smith ierr = PCMGSetInterpolate(pc,i,mg[i]->restrct);CHKERRQ(ierr); 388fccaa45eSBarry Smith } 389fccaa45eSBarry Smith if (!mg[i]->restrct && mg[i]->interpolate) { 39097177400SBarry Smith ierr = PCMGSetRestriction(pc,i,mg[i]->interpolate);CHKERRQ(ierr); 391fccaa45eSBarry Smith } 392fccaa45eSBarry Smith #if defined(PETSC_USE_DEBUG) 393fccaa45eSBarry Smith if (!mg[i]->restrct || !mg[i]->interpolate) { 394fccaa45eSBarry Smith SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Need to set restriction or interpolation on level %d",(int)i); 395fccaa45eSBarry Smith } 396fccaa45eSBarry Smith #endif 397fccaa45eSBarry Smith } 3980a6bb862SBarry Smith for (i=0; i<n-1; i++) { 39937b0e6c0SBarry Smith if (!mg[i]->b) { 40037b0e6c0SBarry Smith Mat mat; 4010cace4b0SBarry Smith Vec vec; 40237b0e6c0SBarry Smith ierr = KSPGetOperators(mg[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr); 4030cace4b0SBarry Smith ierr = MatGetVecs(mat,&vec,PETSC_NULL);CHKERRQ(ierr); 4040cace4b0SBarry Smith ierr = PCMGSetRhs(pc,i,vec);CHKERRQ(ierr); 40537b0e6c0SBarry Smith } 4066ca4d86aSHong Zhang if (!mg[i]->r && i) { 4070a6bb862SBarry Smith ierr = VecDuplicate(mg[i]->b,&tvec);CHKERRQ(ierr); 40897177400SBarry Smith ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr); 4090a6bb862SBarry Smith ierr = VecDestroy(tvec);CHKERRQ(ierr); 4100a6bb862SBarry Smith } 4110a6bb862SBarry Smith if (!mg[i]->x) { 4120a6bb862SBarry Smith ierr = VecDuplicate(mg[i]->b,&tvec);CHKERRQ(ierr); 41397177400SBarry Smith ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr); 4140a6bb862SBarry Smith ierr = VecDestroy(tvec);CHKERRQ(ierr); 4150a6bb862SBarry Smith } 4160a6bb862SBarry Smith } 4174b9ad928SBarry Smith } 4184b9ad928SBarry Smith 419c2be2410SBarry Smith /* If user did not provide fine grid operators, use those from PC */ 420c2be2410SBarry Smith /* BUG BUG BUG This will work ONLY the first time called: hence if the user changes 421c2be2410SBarry Smith the PC matrices between solves PCMG will continue to use first set provided */ 422c2be2410SBarry Smith ierr = KSPGetOperators(mg[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr); 423c2be2410SBarry Smith if (!dA && !dB) { 424ae15b995SBarry 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); 425c2be2410SBarry Smith ierr = KSPSetOperators(mg[n-1]->smoothd,pc->mat,pc->pmat,uflag);CHKERRQ(ierr); 426c2be2410SBarry Smith } 427c2be2410SBarry Smith 428c2be2410SBarry Smith if (mg[0]->galerkin) { 4295e8befabSKris Buschelman Mat B; 430c2be2410SBarry Smith mg[0]->galerkinused = PETSC_TRUE; 431c2be2410SBarry Smith /* currently only handle case where mat and pmat are the same on coarser levels */ 432c2be2410SBarry Smith ierr = KSPGetOperators(mg[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr); 433c2be2410SBarry Smith if (!pc->setupcalled) { 434c2be2410SBarry Smith for (i=n-2; i>-1; i--) { 435c2be2410SBarry Smith ierr = MatPtAP(dB,mg[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 436c2be2410SBarry Smith ierr = KSPSetOperators(mg[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 437c2be2410SBarry Smith dB = B; 438c2be2410SBarry Smith } 439c2be2410SBarry Smith } else { 440c2be2410SBarry Smith for (i=n-2; i>-1; i--) { 441c2be2410SBarry Smith ierr = KSPGetOperators(mg[i]->smoothd,0,&B,0);CHKERRQ(ierr); 442a98fc643SBarry Smith ierr = MatPtAP(dB,mg[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 443c2be2410SBarry Smith ierr = KSPSetOperators(mg[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 444c2be2410SBarry Smith dB = B; 445c2be2410SBarry Smith } 446c2be2410SBarry Smith } 447c2be2410SBarry Smith } 448c2be2410SBarry Smith 4494b9ad928SBarry Smith for (i=1; i<n; i++) { 450b03c7568SBarry Smith if (mg[i]->smoothu == mg[i]->smoothd) { 451b03c7568SBarry Smith /* if doing only down then initial guess is zero */ 4524b9ad928SBarry Smith ierr = KSPSetInitialGuessNonzero(mg[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 453b03c7568SBarry Smith } 4544b9ad928SBarry Smith if (mg[i]->eventsetup) {ierr = PetscLogEventBegin(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 4554b9ad928SBarry Smith ierr = KSPSetUp(mg[i]->smoothd);CHKERRQ(ierr); 4564b9ad928SBarry Smith if (mg[i]->eventsetup) {ierr = PetscLogEventEnd(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 4574b9ad928SBarry Smith } 458b03c7568SBarry Smith for (i=1; i<n; i++) { 4594b9ad928SBarry Smith if (mg[i]->smoothu && mg[i]->smoothu != mg[i]->smoothd) { 460c2be2410SBarry Smith PC uppc,downpc; 46197f1f81fSBarry Smith Mat downmat,downpmat,upmat,uppmat; 46297f1f81fSBarry Smith MatStructure matflag; 46397f1f81fSBarry Smith 46497f1f81fSBarry Smith /* check if operators have been set for up, if not use down operators to set them */ 46597f1f81fSBarry Smith ierr = KSPGetPC(mg[i]->smoothu,&uppc);CHKERRQ(ierr); 46697f1f81fSBarry Smith ierr = PCGetOperators(uppc,&upmat,&uppmat,PETSC_NULL);CHKERRQ(ierr); 46797f1f81fSBarry Smith if (!upmat) { 46897f1f81fSBarry Smith ierr = KSPGetPC(mg[i]->smoothd,&downpc);CHKERRQ(ierr); 46997f1f81fSBarry Smith ierr = PCGetOperators(downpc,&downmat,&downpmat,&matflag);CHKERRQ(ierr); 47097f1f81fSBarry Smith ierr = KSPSetOperators(mg[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr); 47197f1f81fSBarry Smith } 47297f1f81fSBarry Smith 4734b9ad928SBarry Smith ierr = KSPSetInitialGuessNonzero(mg[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr); 4744b9ad928SBarry Smith if (mg[i]->eventsetup) {ierr = PetscLogEventBegin(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 4754b9ad928SBarry Smith ierr = KSPSetUp(mg[i]->smoothu);CHKERRQ(ierr); 4764b9ad928SBarry Smith if (mg[i]->eventsetup) {ierr = PetscLogEventEnd(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 4774b9ad928SBarry Smith } 4784b9ad928SBarry Smith } 4794b9ad928SBarry Smith 4804b9ad928SBarry Smith /* 4814b9ad928SBarry Smith If coarse solver is not direct method then DO NOT USE preonly 4824b9ad928SBarry Smith */ 4834b9ad928SBarry Smith ierr = PetscTypeCompare((PetscObject)mg[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr); 4844b9ad928SBarry Smith if (preonly) { 4854b9ad928SBarry Smith ierr = KSPGetPC(mg[0]->smoothd,&cpc);CHKERRQ(ierr); 4864b9ad928SBarry Smith ierr = PetscTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr); 4874b9ad928SBarry Smith ierr = PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr); 48868eff7e6SBarry Smith ierr = PetscTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr); 48968eff7e6SBarry Smith if (!lu && !redundant && !cholesky) { 4904b9ad928SBarry Smith ierr = KSPSetType(mg[0]->smoothd,KSPGMRES);CHKERRQ(ierr); 4914b9ad928SBarry Smith } 4924b9ad928SBarry Smith } 4934b9ad928SBarry Smith 494958c9bccSBarry Smith if (!pc->setupcalled) { 4954b9ad928SBarry Smith if (monitor) { 4964b9ad928SBarry Smith ierr = PetscObjectGetComm((PetscObject)mg[0]->smoothd,&comm);CHKERRQ(ierr); 4974b9ad928SBarry Smith ierr = PetscViewerASCIIOpen(comm,"stdout",&ascii);CHKERRQ(ierr); 4984b9ad928SBarry Smith ierr = PetscViewerASCIISetTab(ascii,n);CHKERRQ(ierr); 4996849ba73SBarry Smith ierr = KSPSetMonitor(mg[0]->smoothd,KSPDefaultMonitor,ascii,(PetscErrorCode(*)(void*))PetscViewerDestroy);CHKERRQ(ierr); 5004b9ad928SBarry Smith } 5014b9ad928SBarry Smith ierr = KSPSetFromOptions(mg[0]->smoothd);CHKERRQ(ierr); 5024b9ad928SBarry Smith } 5034b9ad928SBarry Smith 5044b9ad928SBarry Smith if (mg[0]->eventsetup) {ierr = PetscLogEventBegin(mg[0]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 5054b9ad928SBarry Smith ierr = KSPSetUp(mg[0]->smoothd);CHKERRQ(ierr); 5064b9ad928SBarry Smith if (mg[0]->eventsetup) {ierr = PetscLogEventEnd(mg[0]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 5074b9ad928SBarry Smith 5084cf18111SSatish Balay #if defined(PETSC_USE_SOCKET_VIEWER) 5094b9ad928SBarry Smith /* 5104b9ad928SBarry Smith Dump the interpolation/restriction matrices to matlab plus the 5114b9ad928SBarry Smith Jacobian/stiffness on each level. This allows Matlab users to 5124b9ad928SBarry Smith easily check if the Galerkin condition A_c = R A_f R^T is satisfied */ 5134b9ad928SBarry Smith ierr = PetscOptionsHasName(pc->prefix,"-pc_mg_dump_matlab",&dump);CHKERRQ(ierr); 5144b9ad928SBarry Smith if (dump) { 5154b9ad928SBarry Smith for (i=1; i<n; i++) { 5164b9ad928SBarry Smith ierr = MatView(mg[i]->restrct,PETSC_VIEWER_SOCKET_(pc->comm));CHKERRQ(ierr); 5174b9ad928SBarry Smith } 5184b9ad928SBarry Smith for (i=0; i<n; i++) { 5194b9ad928SBarry Smith ierr = KSPGetPC(mg[i]->smoothd,&pc);CHKERRQ(ierr); 5204b9ad928SBarry Smith ierr = MatView(pc->mat,PETSC_VIEWER_SOCKET_(pc->comm));CHKERRQ(ierr); 5214b9ad928SBarry Smith } 5224b9ad928SBarry Smith } 523c45a1595SBarry Smith #endif 524c45a1595SBarry Smith 525c2be2410SBarry Smith ierr = PetscOptionsHasName(pc->prefix,"-pc_mg_dump_binary",&dump);CHKERRQ(ierr); 526c2be2410SBarry Smith if (dump) { 527c2be2410SBarry Smith for (i=1; i<n; i++) { 528c2be2410SBarry Smith ierr = MatView(mg[i]->restrct,PETSC_VIEWER_BINARY_(pc->comm));CHKERRQ(ierr); 529c2be2410SBarry Smith } 530c2be2410SBarry Smith for (i=0; i<n; i++) { 531c2be2410SBarry Smith ierr = KSPGetPC(mg[i]->smoothd,&pc);CHKERRQ(ierr); 532c2be2410SBarry Smith ierr = MatView(pc->mat,PETSC_VIEWER_BINARY_(pc->comm));CHKERRQ(ierr); 533c2be2410SBarry Smith } 534c2be2410SBarry Smith } 5354b9ad928SBarry Smith PetscFunctionReturn(0); 5364b9ad928SBarry Smith } 5374b9ad928SBarry Smith 5384b9ad928SBarry Smith /* -------------------------------------------------------------------------------------*/ 5394b9ad928SBarry Smith 5404b9ad928SBarry Smith #undef __FUNCT__ 5419dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetLevels" 5424b9ad928SBarry Smith /*@C 54397177400SBarry Smith PCMGSetLevels - Sets the number of levels to use with MG. 5444b9ad928SBarry Smith Must be called before any other MG routine. 5454b9ad928SBarry Smith 5464b9ad928SBarry Smith Collective on PC 5474b9ad928SBarry Smith 5484b9ad928SBarry Smith Input Parameters: 5494b9ad928SBarry Smith + pc - the preconditioner context 5504b9ad928SBarry Smith . levels - the number of levels 5514b9ad928SBarry Smith - comms - optional communicators for each level; this is to allow solving the coarser problems 5524b9ad928SBarry Smith on smaller sets of processors. Use PETSC_NULL_OBJECT for default in Fortran 5534b9ad928SBarry Smith 5544b9ad928SBarry Smith Level: intermediate 5554b9ad928SBarry Smith 5564b9ad928SBarry Smith Notes: 5574b9ad928SBarry Smith If the number of levels is one then the multigrid uses the -mg_levels prefix 5584b9ad928SBarry Smith for setting the level options rather than the -mg_coarse prefix. 5594b9ad928SBarry Smith 5604b9ad928SBarry Smith .keywords: MG, set, levels, multigrid 5614b9ad928SBarry Smith 56297177400SBarry Smith .seealso: PCMGSetType(), PCMGGetLevels() 5634b9ad928SBarry Smith @*/ 56497177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms) 5654b9ad928SBarry Smith { 566dfbe8321SBarry Smith PetscErrorCode ierr; 567ada7143aSSatish Balay PC_MG **mg=0; 5684b9ad928SBarry Smith 5694b9ad928SBarry Smith PetscFunctionBegin; 5704482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 5714b9ad928SBarry Smith 5724b9ad928SBarry Smith if (pc->data) { 5731302d50aSBarry Smith SETERRQ(PETSC_ERR_ORDER,"Number levels already set for MG\n\ 57497177400SBarry Smith make sure that you call PCMGSetLevels() before KSPSetFromOptions()"); 5754b9ad928SBarry Smith } 5769dcbbd2bSBarry Smith ierr = PCMGCreate_Private(pc->comm,levels,pc,comms,&mg);CHKERRQ(ierr); 5779dcbbd2bSBarry Smith mg[0]->am = PC_MG_MULTIPLICATIVE; 5784b9ad928SBarry Smith pc->data = (void*)mg; 5794b9ad928SBarry Smith pc->ops->applyrichardson = PCApplyRichardson_MG; 5804b9ad928SBarry Smith PetscFunctionReturn(0); 5814b9ad928SBarry Smith } 5824b9ad928SBarry Smith 5834b9ad928SBarry Smith #undef __FUNCT__ 5849dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetLevels" 5854b9ad928SBarry Smith /*@ 58697177400SBarry Smith PCMGGetLevels - Gets the number of levels to use with MG. 5874b9ad928SBarry Smith 5884b9ad928SBarry Smith Not Collective 5894b9ad928SBarry Smith 5904b9ad928SBarry Smith Input Parameter: 5914b9ad928SBarry Smith . pc - the preconditioner context 5924b9ad928SBarry Smith 5934b9ad928SBarry Smith Output parameter: 5944b9ad928SBarry Smith . levels - the number of levels 5954b9ad928SBarry Smith 5964b9ad928SBarry Smith Level: advanced 5974b9ad928SBarry Smith 5984b9ad928SBarry Smith .keywords: MG, get, levels, multigrid 5994b9ad928SBarry Smith 60097177400SBarry Smith .seealso: PCMGSetLevels() 6014b9ad928SBarry Smith @*/ 60297177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetLevels(PC pc,PetscInt *levels) 6034b9ad928SBarry Smith { 6049dcbbd2bSBarry Smith PC_MG **mg; 6054b9ad928SBarry Smith 6064b9ad928SBarry Smith PetscFunctionBegin; 6074482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 6084482741eSBarry Smith PetscValidIntPointer(levels,2); 6094b9ad928SBarry Smith 6109dcbbd2bSBarry Smith mg = (PC_MG**)pc->data; 6114b9ad928SBarry Smith *levels = mg[0]->levels; 6124b9ad928SBarry Smith PetscFunctionReturn(0); 6134b9ad928SBarry Smith } 6144b9ad928SBarry Smith 6154b9ad928SBarry Smith #undef __FUNCT__ 6169dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetType" 6174b9ad928SBarry Smith /*@ 61897177400SBarry Smith PCMGSetType - Determines the form of multigrid to use: 6194b9ad928SBarry Smith multiplicative, additive, full, or the Kaskade algorithm. 6204b9ad928SBarry Smith 6214b9ad928SBarry Smith Collective on PC 6224b9ad928SBarry Smith 6234b9ad928SBarry Smith Input Parameters: 6244b9ad928SBarry Smith + pc - the preconditioner context 6259dcbbd2bSBarry Smith - form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE, 6269dcbbd2bSBarry Smith PC_MG_FULL, PC_MG_KASKADE 6274b9ad928SBarry Smith 6284b9ad928SBarry Smith Options Database Key: 6294b9ad928SBarry Smith . -pc_mg_type <form> - Sets <form>, one of multiplicative, 6304b9ad928SBarry Smith additive, full, kaskade 6314b9ad928SBarry Smith 6324b9ad928SBarry Smith Level: advanced 6334b9ad928SBarry Smith 6344b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid 6354b9ad928SBarry Smith 63697177400SBarry Smith .seealso: PCMGSetLevels() 6374b9ad928SBarry Smith @*/ 6389dcbbd2bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetType(PC pc,PCMGType form) 6394b9ad928SBarry Smith { 6409dcbbd2bSBarry Smith PC_MG **mg; 6414b9ad928SBarry Smith 6424b9ad928SBarry Smith PetscFunctionBegin; 6434482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 6449dcbbd2bSBarry Smith mg = (PC_MG**)pc->data; 6454b9ad928SBarry Smith 6464b9ad928SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 6474b9ad928SBarry Smith mg[0]->am = form; 6489dcbbd2bSBarry Smith if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG; 6494b9ad928SBarry Smith else pc->ops->applyrichardson = 0; 6504b9ad928SBarry Smith PetscFunctionReturn(0); 6514b9ad928SBarry Smith } 6524b9ad928SBarry Smith 6534b9ad928SBarry Smith #undef __FUNCT__ 6549dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetCycles" 6554b9ad928SBarry Smith /*@ 65697177400SBarry Smith PCMGSetCycles - Sets the type cycles to use. Use PCMGSetCyclesOnLevel() for more 6574b9ad928SBarry Smith complicated cycling. 6584b9ad928SBarry Smith 6594b9ad928SBarry Smith Collective on PC 6604b9ad928SBarry Smith 6614b9ad928SBarry Smith Input Parameters: 662c2be2410SBarry Smith + pc - the multigrid context 6634b9ad928SBarry Smith - n - the number of cycles 6644b9ad928SBarry Smith 6654b9ad928SBarry Smith Options Database Key: 6664b9ad928SBarry Smith $ -pc_mg_cycles n - 1 denotes a V-cycle; 2 denotes a W-cycle. 6674b9ad928SBarry Smith 6684b9ad928SBarry Smith Level: advanced 6694b9ad928SBarry Smith 6704b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 6714b9ad928SBarry Smith 67297177400SBarry Smith .seealso: PCMGSetCyclesOnLevel() 6734b9ad928SBarry Smith @*/ 67497177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetCycles(PC pc,PetscInt n) 6754b9ad928SBarry Smith { 6769dcbbd2bSBarry Smith PC_MG **mg; 67779416396SBarry Smith PetscInt i,levels; 6784b9ad928SBarry Smith 6794b9ad928SBarry Smith PetscFunctionBegin; 6804482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 6819dcbbd2bSBarry Smith mg = (PC_MG**)pc->data; 6824b9ad928SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 6834b9ad928SBarry Smith levels = mg[0]->levels; 6844b9ad928SBarry Smith 6854b9ad928SBarry Smith for (i=0; i<levels; i++) { 6864b9ad928SBarry Smith mg[i]->cycles = n; 6874b9ad928SBarry Smith } 6884b9ad928SBarry Smith PetscFunctionReturn(0); 6894b9ad928SBarry Smith } 6904b9ad928SBarry Smith 6914b9ad928SBarry Smith #undef __FUNCT__ 6929dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetGalerkin" 693c2be2410SBarry Smith /*@ 69497177400SBarry Smith PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the 695c2be2410SBarry Smith finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t 696c2be2410SBarry Smith 697c2be2410SBarry Smith Collective on PC 698c2be2410SBarry Smith 699c2be2410SBarry Smith Input Parameters: 7003fc8bf9cSBarry Smith . pc - the multigrid context 701c2be2410SBarry Smith 702c2be2410SBarry Smith Options Database Key: 703c2be2410SBarry Smith $ -pc_mg_galerkin 704c2be2410SBarry Smith 705c2be2410SBarry Smith Level: intermediate 706c2be2410SBarry Smith 707c2be2410SBarry Smith .keywords: MG, set, Galerkin 708c2be2410SBarry Smith 7093fc8bf9cSBarry Smith .seealso: PCMGGetGalerkin() 7103fc8bf9cSBarry Smith 711c2be2410SBarry Smith @*/ 71297177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetGalerkin(PC pc) 713c2be2410SBarry Smith { 7149dcbbd2bSBarry Smith PC_MG **mg; 715c2be2410SBarry Smith PetscInt i,levels; 716c2be2410SBarry Smith 717c2be2410SBarry Smith PetscFunctionBegin; 718c2be2410SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 7199dcbbd2bSBarry Smith mg = (PC_MG**)pc->data; 720c2be2410SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 721c2be2410SBarry Smith levels = mg[0]->levels; 722c2be2410SBarry Smith 723c2be2410SBarry Smith for (i=0; i<levels; i++) { 724c2be2410SBarry Smith mg[i]->galerkin = PETSC_TRUE; 725c2be2410SBarry Smith } 726c2be2410SBarry Smith PetscFunctionReturn(0); 727c2be2410SBarry Smith } 728c2be2410SBarry Smith 729c2be2410SBarry Smith #undef __FUNCT__ 7303fc8bf9cSBarry Smith #define __FUNCT__ "PCMGGetGalerkin" 7313fc8bf9cSBarry Smith /*@ 7323fc8bf9cSBarry Smith PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. 7333fc8bf9cSBarry Smith A_i-1 = r_i * A_i * r_i^t 7343fc8bf9cSBarry Smith 7353fc8bf9cSBarry Smith Not Collective 7363fc8bf9cSBarry Smith 7373fc8bf9cSBarry Smith Input Parameter: 7383fc8bf9cSBarry Smith . pc - the multigrid context 7393fc8bf9cSBarry Smith 7403fc8bf9cSBarry Smith Output Parameter: 7413fc8bf9cSBarry Smith . gelerkin - PETSC_TRUE or PETSC_FALSE 7423fc8bf9cSBarry Smith 7433fc8bf9cSBarry Smith Options Database Key: 7443fc8bf9cSBarry Smith $ -pc_mg_galerkin 7453fc8bf9cSBarry Smith 7463fc8bf9cSBarry Smith Level: intermediate 7473fc8bf9cSBarry Smith 7483fc8bf9cSBarry Smith .keywords: MG, set, Galerkin 7493fc8bf9cSBarry Smith 7503fc8bf9cSBarry Smith .seealso: PCMGSetGalerkin() 7513fc8bf9cSBarry Smith 7523fc8bf9cSBarry Smith @*/ 7533fc8bf9cSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetGalerkin(PC pc,PetscTruth *galerkin) 7543fc8bf9cSBarry Smith { 7553fc8bf9cSBarry Smith PC_MG **mg; 7563fc8bf9cSBarry Smith 7573fc8bf9cSBarry Smith PetscFunctionBegin; 7583fc8bf9cSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 7593fc8bf9cSBarry Smith mg = (PC_MG**)pc->data; 7603fc8bf9cSBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 7613fc8bf9cSBarry Smith *galerkin = mg[0]->galerkin; 7623fc8bf9cSBarry Smith PetscFunctionReturn(0); 7633fc8bf9cSBarry Smith } 7643fc8bf9cSBarry Smith 7653fc8bf9cSBarry Smith #undef __FUNCT__ 7669dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothDown" 7674b9ad928SBarry Smith /*@ 76897177400SBarry Smith PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to 76997177400SBarry Smith use on all levels. Use PCMGGetSmootherDown() to set different 7704b9ad928SBarry Smith pre-smoothing steps on different levels. 7714b9ad928SBarry Smith 7724b9ad928SBarry Smith Collective on PC 7734b9ad928SBarry Smith 7744b9ad928SBarry Smith Input Parameters: 7754b9ad928SBarry Smith + mg - the multigrid context 7764b9ad928SBarry Smith - n - the number of smoothing steps 7774b9ad928SBarry Smith 7784b9ad928SBarry Smith Options Database Key: 7794b9ad928SBarry Smith . -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps 7804b9ad928SBarry Smith 7814b9ad928SBarry Smith Level: advanced 7824b9ad928SBarry Smith 7834b9ad928SBarry Smith .keywords: MG, smooth, down, pre-smoothing, steps, multigrid 7844b9ad928SBarry Smith 78597177400SBarry Smith .seealso: PCMGSetNumberSmoothUp() 7864b9ad928SBarry Smith @*/ 78797177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothDown(PC pc,PetscInt n) 7884b9ad928SBarry Smith { 7899dcbbd2bSBarry Smith PC_MG **mg; 7906849ba73SBarry Smith PetscErrorCode ierr; 79179416396SBarry Smith PetscInt i,levels; 7924b9ad928SBarry Smith 7934b9ad928SBarry Smith PetscFunctionBegin; 7944482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 7959dcbbd2bSBarry Smith mg = (PC_MG**)pc->data; 7964b9ad928SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 7974b9ad928SBarry Smith levels = mg[0]->levels; 7984b9ad928SBarry Smith 799b05257ddSBarry Smith for (i=1; i<levels; i++) { 8004b9ad928SBarry Smith /* make sure smoother up and down are different */ 80197177400SBarry Smith ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 8024b9ad928SBarry Smith ierr = KSPSetTolerances(mg[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 8034b9ad928SBarry Smith mg[i]->default_smoothd = n; 8044b9ad928SBarry Smith } 8054b9ad928SBarry Smith PetscFunctionReturn(0); 8064b9ad928SBarry Smith } 8074b9ad928SBarry Smith 8084b9ad928SBarry Smith #undef __FUNCT__ 8099dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothUp" 8104b9ad928SBarry Smith /*@ 81197177400SBarry Smith PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use 81297177400SBarry Smith on all levels. Use PCMGGetSmootherUp() to set different numbers of 8134b9ad928SBarry Smith post-smoothing steps on different levels. 8144b9ad928SBarry Smith 8154b9ad928SBarry Smith Collective on PC 8164b9ad928SBarry Smith 8174b9ad928SBarry Smith Input Parameters: 8184b9ad928SBarry Smith + mg - the multigrid context 8194b9ad928SBarry Smith - n - the number of smoothing steps 8204b9ad928SBarry Smith 8214b9ad928SBarry Smith Options Database Key: 8224b9ad928SBarry Smith . -pc_mg_smoothup <n> - Sets number of post-smoothing steps 8234b9ad928SBarry Smith 8244b9ad928SBarry Smith Level: advanced 8254b9ad928SBarry Smith 8264b9ad928SBarry Smith Note: this does not set a value on the coarsest grid, since we assume that 827a8c7a070SBarry Smith there is no separate smooth up on the coarsest grid. 8284b9ad928SBarry Smith 8294b9ad928SBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid 8304b9ad928SBarry Smith 83197177400SBarry Smith .seealso: PCMGSetNumberSmoothDown() 8324b9ad928SBarry Smith @*/ 83397177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothUp(PC pc,PetscInt n) 8344b9ad928SBarry Smith { 8359dcbbd2bSBarry Smith PC_MG **mg; 8366849ba73SBarry Smith PetscErrorCode ierr; 83779416396SBarry Smith PetscInt i,levels; 8384b9ad928SBarry Smith 8394b9ad928SBarry Smith PetscFunctionBegin; 8404482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 8419dcbbd2bSBarry Smith mg = (PC_MG**)pc->data; 8424b9ad928SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 8434b9ad928SBarry Smith levels = mg[0]->levels; 8444b9ad928SBarry Smith 8454b9ad928SBarry Smith for (i=1; i<levels; i++) { 8464b9ad928SBarry Smith /* make sure smoother up and down are different */ 84797177400SBarry Smith ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 8484b9ad928SBarry Smith ierr = KSPSetTolerances(mg[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 8494b9ad928SBarry Smith mg[i]->default_smoothu = n; 8504b9ad928SBarry Smith } 8514b9ad928SBarry Smith PetscFunctionReturn(0); 8524b9ad928SBarry Smith } 8534b9ad928SBarry Smith 8544b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/ 8554b9ad928SBarry Smith 8563b09bd56SBarry Smith /*MC 8573b09bd56SBarry Smith PCMG - Use geometric multigrid preconditioning. This preconditioner requires you provide additional 8583b09bd56SBarry Smith information about the coarser grid matrices and restriction/interpolation operators. 8593b09bd56SBarry Smith 8603b09bd56SBarry Smith Options Database Keys: 8613b09bd56SBarry Smith + -pc_mg_levels <nlevels> - number of levels including finest 8623b09bd56SBarry Smith . -pc_mg_cycles 1 or 2 - for V or W-cycle 86379416396SBarry Smith . -pc_mg_smoothup <n> - number of smoothing steps after interpolation 8643b09bd56SBarry Smith . -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator 8653b09bd56SBarry Smith . -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default 8663b09bd56SBarry Smith . -pc_mg_log - log information about time spent on each level of the solver 8673b09bd56SBarry Smith . -pc_mg_monitor - print information on the multigrid convergence 86868eff7e6SBarry Smith . -pc_mg_galerkin - use Galerkin process to compute coarser operators 8693b09bd56SBarry Smith - -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 8703b09bd56SBarry Smith to the Socket viewer for reading from Matlab. 8713b09bd56SBarry Smith 8723b09bd56SBarry Smith Notes: 8733b09bd56SBarry Smith 8743b09bd56SBarry Smith Level: intermediate 8753b09bd56SBarry Smith 8763b09bd56SBarry Smith Concepts: multigrid 8773b09bd56SBarry Smith 8783b09bd56SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 8799dcbbd2bSBarry Smith PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycles(), PCMGSetNumberSmoothDown(), 88097177400SBarry Smith PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 88197177400SBarry Smith PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 88297177400SBarry Smith PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 8833b09bd56SBarry Smith M*/ 8843b09bd56SBarry Smith 8854b9ad928SBarry Smith EXTERN_C_BEGIN 8864b9ad928SBarry Smith #undef __FUNCT__ 8874b9ad928SBarry Smith #define __FUNCT__ "PCCreate_MG" 888dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_MG(PC pc) 8894b9ad928SBarry Smith { 8904b9ad928SBarry Smith PetscFunctionBegin; 8914b9ad928SBarry Smith pc->ops->apply = PCApply_MG; 8924b9ad928SBarry Smith pc->ops->setup = PCSetUp_MG; 8934b9ad928SBarry Smith pc->ops->destroy = PCDestroy_MG; 8944b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_MG; 8954b9ad928SBarry Smith pc->ops->view = PCView_MG; 8964b9ad928SBarry Smith 8974b9ad928SBarry Smith pc->data = (void*)0; 8984b9ad928SBarry Smith PetscFunctionReturn(0); 8994b9ad928SBarry Smith } 9004b9ad928SBarry Smith EXTERN_C_END 901