1dba47a55SKris Buschelman #define PETSCKSP_DLL 2dba47a55SKris Buschelman 34b9ad928SBarry Smith /* 44b9ad928SBarry Smith Defines the multigrid preconditioner interface. 54b9ad928SBarry Smith */ 64b9ad928SBarry Smith #include "src/ksp/pc/impls/mg/mgimpl.h" /*I "petscmg.h" I*/ 74b9ad928SBarry Smith 84b9ad928SBarry Smith 94b9ad928SBarry Smith #undef __FUNCT__ 109dcbbd2bSBarry Smith #define __FUNCT__ "PCMGMCycle_Private" 119dcbbd2bSBarry Smith PetscErrorCode PCMGMCycle_Private(PC_MG **mglevels,PetscTruth *converged) 124b9ad928SBarry Smith { 139dcbbd2bSBarry Smith PC_MG *mg = *mglevels,*mgc; 146849ba73SBarry Smith PetscErrorCode ierr; 150d353602SBarry Smith PetscInt cycles = (PetscInt) mg->cycles; 164b9ad928SBarry Smith 174b9ad928SBarry Smith PetscFunctionBegin; 184b9ad928SBarry Smith if (converged) *converged = PETSC_FALSE; 194b9ad928SBarry Smith 204b9ad928SBarry Smith if (mg->eventsolve) {ierr = PetscLogEventBegin(mg->eventsolve,0,0,0,0);CHKERRQ(ierr);} 210d353602SBarry Smith ierr = KSPSolve(mg->smoothd,mg->b,mg->x);CHKERRQ(ierr); /* pre-smooth */ 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);} 490d353602SBarry Smith ierr = KSPSolve(mg->smoothu,mg->b,mg->x);CHKERRQ(ierr); /* post smooth */ 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; 820d353602SBarry Smith mg[i]->cycles = PC_MG_CYCLE_V; 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; 120*8cc2d5dfSBarry Smith mg[i]->cyclesperpcapply = 1; 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; 135fccaa45eSBarry Smith for (i=0; i<n-1; i++) { 136630ba2eeSBarry Smith if (mg[i+1]->r) {ierr = VecDestroy(mg[i+1]->r);CHKERRQ(ierr);} 137fccaa45eSBarry Smith if (mg[i]->b) {ierr = VecDestroy(mg[i]->b);CHKERRQ(ierr);} 138fccaa45eSBarry Smith if (mg[i]->x) {ierr = VecDestroy(mg[i]->x);CHKERRQ(ierr);} 139fccaa45eSBarry Smith if (mg[i+1]->restrct) {ierr = MatDestroy(mg[i+1]->restrct);CHKERRQ(ierr);} 140fccaa45eSBarry Smith if (mg[i+1]->interpolate) {ierr = MatDestroy(mg[i+1]->interpolate);CHKERRQ(ierr);} 141fccaa45eSBarry Smith } 142fccaa45eSBarry Smith 1434b9ad928SBarry Smith for (i=0; i<n; i++) { 1444b9ad928SBarry Smith if (mg[i]->smoothd != mg[i]->smoothu) { 1454b9ad928SBarry Smith ierr = KSPDestroy(mg[i]->smoothd);CHKERRQ(ierr); 1464b9ad928SBarry Smith } 1474b9ad928SBarry Smith ierr = KSPDestroy(mg[i]->smoothu);CHKERRQ(ierr); 1484b9ad928SBarry Smith ierr = PetscFree(mg[i]);CHKERRQ(ierr); 1494b9ad928SBarry Smith } 1504b9ad928SBarry Smith ierr = PetscFree(mg);CHKERRQ(ierr); 1514b9ad928SBarry Smith PetscFunctionReturn(0); 1524b9ad928SBarry Smith } 1534b9ad928SBarry Smith 1544b9ad928SBarry Smith 1554b9ad928SBarry Smith 1569dcbbd2bSBarry Smith EXTERN PetscErrorCode PCMGACycle_Private(PC_MG**); 1579dcbbd2bSBarry Smith EXTERN PetscErrorCode PCMGFCycle_Private(PC_MG**); 1589dcbbd2bSBarry Smith EXTERN PetscErrorCode PCMGKCycle_Private(PC_MG**); 1594b9ad928SBarry Smith 1604b9ad928SBarry Smith /* 1614b9ad928SBarry Smith PCApply_MG - Runs either an additive, multiplicative, Kaskadic 1624b9ad928SBarry Smith or full cycle of multigrid. 1634b9ad928SBarry Smith 1644b9ad928SBarry Smith Note: 1659dcbbd2bSBarry Smith A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle(). 1664b9ad928SBarry Smith */ 1674b9ad928SBarry Smith #undef __FUNCT__ 1684b9ad928SBarry Smith #define __FUNCT__ "PCApply_MG" 1696849ba73SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x) 1704b9ad928SBarry Smith { 1719dcbbd2bSBarry Smith PC_MG **mg = (PC_MG**)pc->data; 1726849ba73SBarry Smith PetscErrorCode ierr; 173*8cc2d5dfSBarry Smith PetscInt levels = mg[0]->levels,i; 1744b9ad928SBarry Smith 1754b9ad928SBarry Smith PetscFunctionBegin; 1764b9ad928SBarry Smith mg[levels-1]->b = b; 1774b9ad928SBarry Smith mg[levels-1]->x = x; 178ef70c39aSMatthew Knepley if (!mg[levels-1]->r && mg[0]->am != PC_MG_ADDITIVE && levels > 1) { 1790a6bb862SBarry Smith Vec tvec; 1800a6bb862SBarry Smith ierr = VecDuplicate(mg[levels-1]->b,&tvec);CHKERRQ(ierr); 18197177400SBarry Smith ierr = PCMGSetR(pc,levels-1,tvec);CHKERRQ(ierr); 1820a6bb862SBarry Smith ierr = VecDestroy(tvec);CHKERRQ(ierr); 1830a6bb862SBarry Smith } 1849dcbbd2bSBarry Smith if (mg[0]->am == PC_MG_MULTIPLICATIVE) { 185efb30889SBarry Smith ierr = VecSet(x,0.0);CHKERRQ(ierr); 186*8cc2d5dfSBarry Smith for (i=0; i<mg[0]->cyclesperpcapply; i++) { 1879dcbbd2bSBarry Smith ierr = PCMGMCycle_Private(mg+levels-1,PETSC_NULL);CHKERRQ(ierr); 1884b9ad928SBarry Smith } 189*8cc2d5dfSBarry Smith } 1909dcbbd2bSBarry Smith else if (mg[0]->am == PC_MG_ADDITIVE) { 1919dcbbd2bSBarry Smith ierr = PCMGACycle_Private(mg);CHKERRQ(ierr); 1924b9ad928SBarry Smith } 1939dcbbd2bSBarry Smith else if (mg[0]->am == PC_MG_KASKADE) { 1949dcbbd2bSBarry Smith ierr = PCMGKCycle_Private(mg);CHKERRQ(ierr); 1954b9ad928SBarry Smith } 1964b9ad928SBarry Smith else { 1979dcbbd2bSBarry Smith ierr = PCMGFCycle_Private(mg);CHKERRQ(ierr); 1984b9ad928SBarry Smith } 1994b9ad928SBarry Smith PetscFunctionReturn(0); 2004b9ad928SBarry Smith } 2014b9ad928SBarry Smith 2024b9ad928SBarry Smith #undef __FUNCT__ 2034b9ad928SBarry Smith #define __FUNCT__ "PCApplyRichardson_MG" 20479416396SBarry Smith static PetscErrorCode PCApplyRichardson_MG(PC pc,Vec b,Vec x,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its) 2054b9ad928SBarry Smith { 2069dcbbd2bSBarry Smith PC_MG **mg = (PC_MG**)pc->data; 207dfbe8321SBarry Smith PetscErrorCode ierr; 20879416396SBarry Smith PetscInt levels = mg[0]->levels; 2094b9ad928SBarry Smith PetscTruth converged = PETSC_FALSE; 2104b9ad928SBarry Smith 2114b9ad928SBarry Smith PetscFunctionBegin; 2124b9ad928SBarry Smith mg[levels-1]->b = b; 2134b9ad928SBarry Smith mg[levels-1]->x = x; 2144b9ad928SBarry Smith 2154b9ad928SBarry Smith mg[levels-1]->rtol = rtol; 21670441072SBarry Smith mg[levels-1]->abstol = abstol; 2174b9ad928SBarry Smith mg[levels-1]->dtol = dtol; 2184b9ad928SBarry Smith if (rtol) { 2194b9ad928SBarry Smith /* compute initial residual norm for relative convergence test */ 2204b9ad928SBarry Smith PetscReal rnorm; 2214b9ad928SBarry Smith ierr = (*mg[levels-1]->residual)(mg[levels-1]->A,b,x,w);CHKERRQ(ierr); 2224b9ad928SBarry Smith ierr = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr); 22370441072SBarry Smith mg[levels-1]->ttol = PetscMax(rtol*rnorm,abstol); 22470441072SBarry Smith } else if (abstol) { 22570441072SBarry Smith mg[levels-1]->ttol = abstol; 2264b9ad928SBarry Smith } else { 2274b9ad928SBarry Smith mg[levels-1]->ttol = 0.0; 2284b9ad928SBarry Smith } 2294b9ad928SBarry Smith 2304b9ad928SBarry Smith while (its-- && !converged) { 2319dcbbd2bSBarry Smith ierr = PCMGMCycle_Private(mg+levels-1,&converged);CHKERRQ(ierr); 2324b9ad928SBarry Smith } 2334b9ad928SBarry Smith PetscFunctionReturn(0); 2344b9ad928SBarry Smith } 2354b9ad928SBarry Smith 2364b9ad928SBarry Smith #undef __FUNCT__ 2374b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_MG" 2386ca4d86aSHong Zhang PetscErrorCode PCSetFromOptions_MG(PC pc) 2394b9ad928SBarry Smith { 240dfbe8321SBarry Smith PetscErrorCode ierr; 241*8cc2d5dfSBarry Smith PetscInt m,levels = 1,cycles; 2424b9ad928SBarry Smith PetscTruth flg; 2439dcbbd2bSBarry Smith PC_MG **mg = (PC_MG**)pc->data; 244cf502942SBarry Smith PCMGType mgtype; 2450d353602SBarry Smith PCMGCycleType mgctype; 2464b9ad928SBarry Smith 2474b9ad928SBarry Smith PetscFunctionBegin; 2484b9ad928SBarry Smith ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr); 2494b9ad928SBarry Smith if (!pc->data) { 2509dcbbd2bSBarry Smith ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr); 25197177400SBarry Smith ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr); 252cf502942SBarry Smith mg = (PC_MG**)pc->data; 2534b9ad928SBarry Smith } 254*8cc2d5dfSBarry Smith mgctype = mg[0]->cycles; 2550d353602SBarry Smith ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr); 2564b9ad928SBarry Smith if (flg) { 2570d353602SBarry Smith ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr); 2580d353602SBarry Smith }; 2599dcbbd2bSBarry Smith ierr = PetscOptionsName("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",&flg);CHKERRQ(ierr); 260c2be2410SBarry Smith if (flg) { 26197177400SBarry Smith ierr = PCMGSetGalerkin(pc);CHKERRQ(ierr); 262c2be2410SBarry Smith } 2639dcbbd2bSBarry Smith ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr); 2644b9ad928SBarry Smith if (flg) { 26597177400SBarry Smith ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr); 2664b9ad928SBarry Smith } 2679dcbbd2bSBarry Smith ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr); 2684b9ad928SBarry Smith if (flg) { 26997177400SBarry Smith ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr); 2704b9ad928SBarry Smith } 2719dcbbd2bSBarry Smith ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr); 272*8cc2d5dfSBarry Smith if (flg) { 273*8cc2d5dfSBarry Smith ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr); 274*8cc2d5dfSBarry Smith } 275*8cc2d5dfSBarry Smith if (mg[0]->am == PC_MG_MULTIPLICATIVE) { 276*8cc2d5dfSBarry Smith ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGSetLevels",mg[0]->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr); 277*8cc2d5dfSBarry Smith if (flg) { 278*8cc2d5dfSBarry Smith ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr); 279*8cc2d5dfSBarry Smith } 280*8cc2d5dfSBarry Smith } 2814b9ad928SBarry Smith ierr = PetscOptionsName("-pc_mg_log","Log times for each multigrid level","None",&flg);CHKERRQ(ierr); 2824b9ad928SBarry Smith if (flg) { 2834f5ab15aSBarry Smith PetscInt i; 2844b9ad928SBarry Smith char eventname[128]; 2854b9ad928SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 2864b9ad928SBarry Smith levels = mg[0]->levels; 2874b9ad928SBarry Smith for (i=0; i<levels; i++) { 28877431f27SBarry Smith sprintf(eventname,"MSetup Level %d",(int)i); 2894b9ad928SBarry Smith ierr = PetscLogEventRegister(&mg[i]->eventsetup,eventname,pc->cookie);CHKERRQ(ierr); 29025a62d46SBarry Smith sprintf(eventname,"MGSolve Level %d to 0",(int)i); 2914b9ad928SBarry Smith ierr = PetscLogEventRegister(&mg[i]->eventsolve,eventname,pc->cookie);CHKERRQ(ierr); 2924b9ad928SBarry Smith } 2934b9ad928SBarry Smith } 2944b9ad928SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 2954b9ad928SBarry Smith PetscFunctionReturn(0); 2964b9ad928SBarry Smith } 2974b9ad928SBarry Smith 2989dcbbd2bSBarry Smith const char *PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0}; 2990d353602SBarry Smith const char *PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0}; 3009dcbbd2bSBarry Smith 3014b9ad928SBarry Smith #undef __FUNCT__ 3024b9ad928SBarry Smith #define __FUNCT__ "PCView_MG" 3036849ba73SBarry Smith static PetscErrorCode PCView_MG(PC pc,PetscViewer viewer) 3044b9ad928SBarry Smith { 3059dcbbd2bSBarry Smith PC_MG **mg = (PC_MG**)pc->data; 306dfbe8321SBarry Smith PetscErrorCode ierr; 30779416396SBarry Smith PetscInt levels = mg[0]->levels,i; 30832077d6dSBarry Smith PetscTruth iascii; 3094b9ad928SBarry Smith 3104b9ad928SBarry Smith PetscFunctionBegin; 31132077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 31232077d6dSBarry Smith if (iascii) { 3130d353602SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," MG: type is %s, levels=%D cycles=%s, pre-smooths=%D, post-smooths=%D\n", 3140d353602SBarry Smith PCMGTypes[mg[0]->am],levels,(mg[0]->cycles == PC_MG_CYCLE_V) ? "v" : "w", 3150d353602SBarry Smith mg[0]->default_smoothd,mg[0]->default_smoothu);CHKERRQ(ierr); 316c2be2410SBarry Smith if (mg[0]->galerkin) { 317c2be2410SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 318c2be2410SBarry Smith } 3194b9ad928SBarry Smith for (i=0; i<levels; i++) { 320b03c7568SBarry Smith if (!i) { 321b03c7568SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Coarse gride solver -- level %D -------------------------------\n",i);CHKERRQ(ierr); 322b03c7568SBarry Smith } else { 32377431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 324b03c7568SBarry Smith } 3254b9ad928SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 3264b9ad928SBarry Smith ierr = KSPView(mg[i]->smoothd,viewer);CHKERRQ(ierr); 3274b9ad928SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 328b03c7568SBarry Smith if (i && mg[i]->smoothd == mg[i]->smoothu) { 3294b9ad928SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 330b03c7568SBarry Smith } else if (i){ 33177431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 3324b9ad928SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 3334b9ad928SBarry Smith ierr = KSPView(mg[i]->smoothu,viewer);CHKERRQ(ierr); 3344b9ad928SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 3354b9ad928SBarry Smith } 3364b9ad928SBarry Smith } 3374b9ad928SBarry Smith } else { 33879a5c55eSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCMG",((PetscObject)viewer)->type_name); 3394b9ad928SBarry Smith } 3404b9ad928SBarry Smith PetscFunctionReturn(0); 3414b9ad928SBarry Smith } 3424b9ad928SBarry Smith 3434b9ad928SBarry Smith /* 3444b9ad928SBarry Smith Calls setup for the KSP on each level 3454b9ad928SBarry Smith */ 3464b9ad928SBarry Smith #undef __FUNCT__ 3474b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_MG" 3486849ba73SBarry Smith static PetscErrorCode PCSetUp_MG(PC pc) 3494b9ad928SBarry Smith { 3509dcbbd2bSBarry Smith PC_MG **mg = (PC_MG**)pc->data; 351dfbe8321SBarry Smith PetscErrorCode ierr; 35279416396SBarry Smith PetscInt i,n = mg[0]->levels; 3534b9ad928SBarry Smith PC cpc; 354906ed7ccSBarry Smith PetscTruth preonly,lu,redundant,cholesky,monitor = PETSC_FALSE,dump,opsset; 35523d894e5SBarry Smith PetscViewerASCIIMonitor ascii; 35623d894e5SBarry Smith PetscViewer viewer = PETSC_NULL; 3574b9ad928SBarry Smith MPI_Comm comm; 358c2be2410SBarry Smith Mat dA,dB; 359c2be2410SBarry Smith MatStructure uflag; 3600a6bb862SBarry Smith Vec tvec; 3614b9ad928SBarry Smith 3624b9ad928SBarry Smith PetscFunctionBegin; 363852665d3SBarry Smith 36443fb2f97SBarry Smith /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */ 36543fb2f97SBarry Smith /* so use those from global PC */ 36643fb2f97SBarry Smith /* Is this what we always want? What if user wants to keep old one? */ 367906ed7ccSBarry Smith ierr = KSPGetOperatorsSet(mg[n-1]->smoothd,PETSC_NULL,&opsset);CHKERRQ(ierr); 36843fb2f97SBarry Smith ierr = KSPGetPC(mg[0]->smoothd,&cpc);CHKERRQ(ierr); 36943fb2f97SBarry Smith if (!opsset || cpc->setupcalled == 2) { 370852665d3SBarry 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); 3711cfe3bddSBarry Smith ierr = KSPSetOperators(mg[n-1]->smoothd,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr); 372852665d3SBarry Smith } 373852665d3SBarry Smith 374852665d3SBarry Smith if (mg[0]->galerkin) { 375852665d3SBarry Smith Mat B; 376852665d3SBarry Smith mg[0]->galerkinused = PETSC_TRUE; 377852665d3SBarry Smith /* currently only handle case where mat and pmat are the same on coarser levels */ 378852665d3SBarry Smith ierr = KSPGetOperators(mg[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr); 379852665d3SBarry Smith if (!pc->setupcalled) { 380852665d3SBarry Smith for (i=n-2; i>-1; i--) { 381852665d3SBarry Smith ierr = MatPtAP(dB,mg[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 382852665d3SBarry Smith ierr = KSPSetOperators(mg[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 383906ed7ccSBarry Smith if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);} 384852665d3SBarry Smith dB = B; 385852665d3SBarry Smith } 386906ed7ccSBarry Smith ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr); 387852665d3SBarry Smith } else { 388852665d3SBarry Smith for (i=n-2; i>-1; i--) { 389906ed7ccSBarry Smith ierr = KSPGetOperators(mg[i]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr); 390852665d3SBarry Smith ierr = MatPtAP(dB,mg[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 391852665d3SBarry Smith ierr = KSPSetOperators(mg[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 392852665d3SBarry Smith dB = B; 393852665d3SBarry Smith } 394852665d3SBarry Smith } 395852665d3SBarry Smith } 396852665d3SBarry Smith 397958c9bccSBarry Smith if (!pc->setupcalled) { 3984b9ad928SBarry Smith ierr = PetscOptionsHasName(0,"-pc_mg_monitor",&monitor);CHKERRQ(ierr); 3994b9ad928SBarry Smith 400b03c7568SBarry Smith for (i=0; i<n; i++) { 4014b9ad928SBarry Smith if (monitor) { 4024b9ad928SBarry Smith ierr = PetscObjectGetComm((PetscObject)mg[i]->smoothd,&comm);CHKERRQ(ierr); 40323d894e5SBarry Smith ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr); 40423d894e5SBarry Smith ierr = KSPMonitorSet(mg[i]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 4054b9ad928SBarry Smith } 4064b9ad928SBarry Smith ierr = KSPSetFromOptions(mg[i]->smoothd);CHKERRQ(ierr); 4074b9ad928SBarry Smith } 408b03c7568SBarry Smith for (i=1; i<n; i++) { 409a98fc643SBarry Smith if (mg[i]->smoothu && (mg[i]->smoothu != mg[i]->smoothd)) { 4104b9ad928SBarry Smith if (monitor) { 4114b9ad928SBarry Smith ierr = PetscObjectGetComm((PetscObject)mg[i]->smoothu,&comm);CHKERRQ(ierr); 41223d894e5SBarry Smith ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr); 41323d894e5SBarry Smith ierr = KSPMonitorSet(mg[i]->smoothu,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 4144b9ad928SBarry Smith } 4154b9ad928SBarry Smith ierr = KSPSetFromOptions(mg[i]->smoothu);CHKERRQ(ierr); 4164b9ad928SBarry Smith } 4174b9ad928SBarry Smith } 418fccaa45eSBarry Smith for (i=1; i<n; i++) { 4190cace4b0SBarry Smith if (!mg[i]->residual) { 4200cace4b0SBarry Smith Mat mat; 4210cace4b0SBarry Smith ierr = KSPGetOperators(mg[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr); 4220cace4b0SBarry Smith ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr); 4230cace4b0SBarry Smith } 424fccaa45eSBarry Smith if (mg[i]->restrct && !mg[i]->interpolate) { 42597177400SBarry Smith ierr = PCMGSetInterpolate(pc,i,mg[i]->restrct);CHKERRQ(ierr); 426fccaa45eSBarry Smith } 427fccaa45eSBarry Smith if (!mg[i]->restrct && mg[i]->interpolate) { 42897177400SBarry Smith ierr = PCMGSetRestriction(pc,i,mg[i]->interpolate);CHKERRQ(ierr); 429fccaa45eSBarry Smith } 430fccaa45eSBarry Smith #if defined(PETSC_USE_DEBUG) 431fccaa45eSBarry Smith if (!mg[i]->restrct || !mg[i]->interpolate) { 432fccaa45eSBarry Smith SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Need to set restriction or interpolation on level %d",(int)i); 433fccaa45eSBarry Smith } 434fccaa45eSBarry Smith #endif 435fccaa45eSBarry Smith } 4360a6bb862SBarry Smith for (i=0; i<n-1; i++) { 43737b0e6c0SBarry Smith if (!mg[i]->b) { 438906ed7ccSBarry Smith Vec *vec; 439906ed7ccSBarry Smith ierr = KSPGetVecs(mg[i]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr); 440906ed7ccSBarry Smith ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr); 441906ed7ccSBarry Smith ierr = PetscFree(vec);CHKERRQ(ierr); 44237b0e6c0SBarry Smith } 4436ca4d86aSHong Zhang if (!mg[i]->r && i) { 4440a6bb862SBarry Smith ierr = VecDuplicate(mg[i]->b,&tvec);CHKERRQ(ierr); 44597177400SBarry Smith ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr); 4460a6bb862SBarry Smith ierr = VecDestroy(tvec);CHKERRQ(ierr); 4470a6bb862SBarry Smith } 4480a6bb862SBarry Smith if (!mg[i]->x) { 4490a6bb862SBarry Smith ierr = VecDuplicate(mg[i]->b,&tvec);CHKERRQ(ierr); 45097177400SBarry Smith ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr); 4510a6bb862SBarry Smith ierr = VecDestroy(tvec);CHKERRQ(ierr); 4520a6bb862SBarry Smith } 4530a6bb862SBarry Smith } 4544b9ad928SBarry Smith } 4554b9ad928SBarry Smith 456c2be2410SBarry Smith 4574b9ad928SBarry Smith for (i=1; i<n; i++) { 458b03c7568SBarry Smith if (mg[i]->smoothu == mg[i]->smoothd) { 459b03c7568SBarry Smith /* if doing only down then initial guess is zero */ 4604b9ad928SBarry Smith ierr = KSPSetInitialGuessNonzero(mg[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 461b03c7568SBarry Smith } 4624b9ad928SBarry Smith if (mg[i]->eventsetup) {ierr = PetscLogEventBegin(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 4634b9ad928SBarry Smith ierr = KSPSetUp(mg[i]->smoothd);CHKERRQ(ierr); 4644b9ad928SBarry Smith if (mg[i]->eventsetup) {ierr = PetscLogEventEnd(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 4654b9ad928SBarry Smith } 466b03c7568SBarry Smith for (i=1; i<n; i++) { 4674b9ad928SBarry Smith if (mg[i]->smoothu && mg[i]->smoothu != mg[i]->smoothd) { 468906ed7ccSBarry Smith Mat downmat,downpmat; 46997f1f81fSBarry Smith MatStructure matflag; 470906ed7ccSBarry Smith PetscTruth opsset; 47197f1f81fSBarry Smith 47297f1f81fSBarry Smith /* check if operators have been set for up, if not use down operators to set them */ 473906ed7ccSBarry Smith ierr = KSPGetOperatorsSet(mg[i]->smoothu,&opsset,PETSC_NULL);CHKERRQ(ierr); 474906ed7ccSBarry Smith if (!opsset) { 475906ed7ccSBarry Smith ierr = KSPGetOperators(mg[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr); 47697f1f81fSBarry Smith ierr = KSPSetOperators(mg[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr); 47797f1f81fSBarry Smith } 47897f1f81fSBarry Smith 4794b9ad928SBarry Smith ierr = KSPSetInitialGuessNonzero(mg[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr); 4804b9ad928SBarry Smith if (mg[i]->eventsetup) {ierr = PetscLogEventBegin(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 4814b9ad928SBarry Smith ierr = KSPSetUp(mg[i]->smoothu);CHKERRQ(ierr); 4824b9ad928SBarry Smith if (mg[i]->eventsetup) {ierr = PetscLogEventEnd(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 4834b9ad928SBarry Smith } 4844b9ad928SBarry Smith } 4854b9ad928SBarry Smith 4864b9ad928SBarry Smith /* 4874b9ad928SBarry Smith If coarse solver is not direct method then DO NOT USE preonly 4884b9ad928SBarry Smith */ 4894b9ad928SBarry Smith ierr = PetscTypeCompare((PetscObject)mg[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr); 4904b9ad928SBarry Smith if (preonly) { 4914b9ad928SBarry Smith ierr = PetscTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr); 4924b9ad928SBarry Smith ierr = PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr); 49368eff7e6SBarry Smith ierr = PetscTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr); 49468eff7e6SBarry Smith if (!lu && !redundant && !cholesky) { 4954b9ad928SBarry Smith ierr = KSPSetType(mg[0]->smoothd,KSPGMRES);CHKERRQ(ierr); 4964b9ad928SBarry Smith } 4974b9ad928SBarry Smith } 4984b9ad928SBarry Smith 499958c9bccSBarry Smith if (!pc->setupcalled) { 5004b9ad928SBarry Smith if (monitor) { 5014b9ad928SBarry Smith ierr = PetscObjectGetComm((PetscObject)mg[0]->smoothd,&comm);CHKERRQ(ierr); 50223d894e5SBarry Smith ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n,&ascii);CHKERRQ(ierr); 50323d894e5SBarry Smith ierr = KSPMonitorSet(mg[0]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 5044b9ad928SBarry Smith } 5054b9ad928SBarry Smith ierr = KSPSetFromOptions(mg[0]->smoothd);CHKERRQ(ierr); 5064b9ad928SBarry Smith } 5074b9ad928SBarry Smith 5084b9ad928SBarry Smith if (mg[0]->eventsetup) {ierr = PetscLogEventBegin(mg[0]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 5094b9ad928SBarry Smith ierr = KSPSetUp(mg[0]->smoothd);CHKERRQ(ierr); 5104b9ad928SBarry Smith if (mg[0]->eventsetup) {ierr = PetscLogEventEnd(mg[0]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 5114b9ad928SBarry Smith 5124b9ad928SBarry Smith /* 5136805f65bSBarry Smith Dump the interpolation/restriction matrices plus the 5144b9ad928SBarry Smith Jacobian/stiffness on each level. This allows Matlab users to 5156805f65bSBarry Smith easily check if the Galerkin condition A_c = R A_f R^T is satisfied. 5166805f65bSBarry Smith 5176805f65bSBarry Smith Only support one or the other at the same time. 5186805f65bSBarry Smith */ 5196805f65bSBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER) 5204b9ad928SBarry Smith ierr = PetscOptionsHasName(pc->prefix,"-pc_mg_dump_matlab",&dump);CHKERRQ(ierr); 5214b9ad928SBarry Smith if (dump) { 5226805f65bSBarry Smith viewer = PETSC_VIEWER_SOCKET_(pc->comm); 5234b9ad928SBarry Smith } 524c45a1595SBarry Smith #endif 525c2be2410SBarry Smith ierr = PetscOptionsHasName(pc->prefix,"-pc_mg_dump_binary",&dump);CHKERRQ(ierr); 526c2be2410SBarry Smith if (dump) { 5276805f65bSBarry Smith viewer = PETSC_VIEWER_BINARY_(pc->comm); 5286805f65bSBarry Smith } 5296805f65bSBarry Smith 5306805f65bSBarry Smith if (viewer) { 531c2be2410SBarry Smith for (i=1; i<n; i++) { 5326805f65bSBarry Smith ierr = MatView(mg[i]->restrct,viewer);CHKERRQ(ierr); 533c2be2410SBarry Smith } 534c2be2410SBarry Smith for (i=0; i<n; i++) { 535c2be2410SBarry Smith ierr = KSPGetPC(mg[i]->smoothd,&pc);CHKERRQ(ierr); 5366805f65bSBarry Smith ierr = MatView(pc->mat,viewer);CHKERRQ(ierr); 537c2be2410SBarry Smith } 538c2be2410SBarry Smith } 5394b9ad928SBarry Smith PetscFunctionReturn(0); 5404b9ad928SBarry Smith } 5414b9ad928SBarry Smith 5424b9ad928SBarry Smith /* -------------------------------------------------------------------------------------*/ 5434b9ad928SBarry Smith 5444b9ad928SBarry Smith #undef __FUNCT__ 5459dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetLevels" 5464b9ad928SBarry Smith /*@C 54797177400SBarry Smith PCMGSetLevels - Sets the number of levels to use with MG. 5484b9ad928SBarry Smith Must be called before any other MG routine. 5494b9ad928SBarry Smith 5504b9ad928SBarry Smith Collective on PC 5514b9ad928SBarry Smith 5524b9ad928SBarry Smith Input Parameters: 5534b9ad928SBarry Smith + pc - the preconditioner context 5544b9ad928SBarry Smith . levels - the number of levels 5554b9ad928SBarry Smith - comms - optional communicators for each level; this is to allow solving the coarser problems 5564b9ad928SBarry Smith on smaller sets of processors. Use PETSC_NULL_OBJECT for default in Fortran 5574b9ad928SBarry Smith 5584b9ad928SBarry Smith Level: intermediate 5594b9ad928SBarry Smith 5604b9ad928SBarry Smith Notes: 5614b9ad928SBarry Smith If the number of levels is one then the multigrid uses the -mg_levels prefix 5624b9ad928SBarry Smith for setting the level options rather than the -mg_coarse prefix. 5634b9ad928SBarry Smith 5644b9ad928SBarry Smith .keywords: MG, set, levels, multigrid 5654b9ad928SBarry Smith 56697177400SBarry Smith .seealso: PCMGSetType(), PCMGGetLevels() 5674b9ad928SBarry Smith @*/ 56897177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms) 5694b9ad928SBarry Smith { 570dfbe8321SBarry Smith PetscErrorCode ierr; 571ada7143aSSatish Balay PC_MG **mg=0; 5724b9ad928SBarry Smith 5734b9ad928SBarry Smith PetscFunctionBegin; 5744482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 5754b9ad928SBarry Smith 5764b9ad928SBarry Smith if (pc->data) { 5771302d50aSBarry Smith SETERRQ(PETSC_ERR_ORDER,"Number levels already set for MG\n\ 57897177400SBarry Smith make sure that you call PCMGSetLevels() before KSPSetFromOptions()"); 5794b9ad928SBarry Smith } 5809dcbbd2bSBarry Smith ierr = PCMGCreate_Private(pc->comm,levels,pc,comms,&mg);CHKERRQ(ierr); 5819dcbbd2bSBarry Smith mg[0]->am = PC_MG_MULTIPLICATIVE; 5824b9ad928SBarry Smith pc->data = (void*)mg; 5834b9ad928SBarry Smith pc->ops->applyrichardson = PCApplyRichardson_MG; 5844b9ad928SBarry Smith PetscFunctionReturn(0); 5854b9ad928SBarry Smith } 5864b9ad928SBarry Smith 5874b9ad928SBarry Smith #undef __FUNCT__ 5889dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetLevels" 5894b9ad928SBarry Smith /*@ 59097177400SBarry Smith PCMGGetLevels - Gets the number of levels to use with MG. 5914b9ad928SBarry Smith 5924b9ad928SBarry Smith Not Collective 5934b9ad928SBarry Smith 5944b9ad928SBarry Smith Input Parameter: 5954b9ad928SBarry Smith . pc - the preconditioner context 5964b9ad928SBarry Smith 5974b9ad928SBarry Smith Output parameter: 5984b9ad928SBarry Smith . levels - the number of levels 5994b9ad928SBarry Smith 6004b9ad928SBarry Smith Level: advanced 6014b9ad928SBarry Smith 6024b9ad928SBarry Smith .keywords: MG, get, levels, multigrid 6034b9ad928SBarry Smith 60497177400SBarry Smith .seealso: PCMGSetLevels() 6054b9ad928SBarry Smith @*/ 60697177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetLevels(PC pc,PetscInt *levels) 6074b9ad928SBarry Smith { 6089dcbbd2bSBarry Smith PC_MG **mg; 6094b9ad928SBarry Smith 6104b9ad928SBarry Smith PetscFunctionBegin; 6114482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 6124482741eSBarry Smith PetscValidIntPointer(levels,2); 6134b9ad928SBarry Smith 6149dcbbd2bSBarry Smith mg = (PC_MG**)pc->data; 6154b9ad928SBarry Smith *levels = mg[0]->levels; 6164b9ad928SBarry Smith PetscFunctionReturn(0); 6174b9ad928SBarry Smith } 6184b9ad928SBarry Smith 6194b9ad928SBarry Smith #undef __FUNCT__ 6209dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetType" 6214b9ad928SBarry Smith /*@ 62297177400SBarry Smith PCMGSetType - Determines the form of multigrid to use: 6234b9ad928SBarry Smith multiplicative, additive, full, or the Kaskade algorithm. 6244b9ad928SBarry Smith 6254b9ad928SBarry Smith Collective on PC 6264b9ad928SBarry Smith 6274b9ad928SBarry Smith Input Parameters: 6284b9ad928SBarry Smith + pc - the preconditioner context 6299dcbbd2bSBarry Smith - form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE, 6309dcbbd2bSBarry Smith PC_MG_FULL, PC_MG_KASKADE 6314b9ad928SBarry Smith 6324b9ad928SBarry Smith Options Database Key: 6334b9ad928SBarry Smith . -pc_mg_type <form> - Sets <form>, one of multiplicative, 6344b9ad928SBarry Smith additive, full, kaskade 6354b9ad928SBarry Smith 6364b9ad928SBarry Smith Level: advanced 6374b9ad928SBarry Smith 6384b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid 6394b9ad928SBarry Smith 64097177400SBarry Smith .seealso: PCMGSetLevels() 6414b9ad928SBarry Smith @*/ 6429dcbbd2bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetType(PC pc,PCMGType form) 6434b9ad928SBarry Smith { 6449dcbbd2bSBarry Smith PC_MG **mg; 6454b9ad928SBarry Smith 6464b9ad928SBarry Smith PetscFunctionBegin; 6474482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 6489dcbbd2bSBarry Smith mg = (PC_MG**)pc->data; 6494b9ad928SBarry Smith 6504b9ad928SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 6514b9ad928SBarry Smith mg[0]->am = form; 6529dcbbd2bSBarry Smith if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG; 6534b9ad928SBarry Smith else pc->ops->applyrichardson = 0; 6544b9ad928SBarry Smith PetscFunctionReturn(0); 6554b9ad928SBarry Smith } 6564b9ad928SBarry Smith 6574b9ad928SBarry Smith #undef __FUNCT__ 6580d353602SBarry Smith #define __FUNCT__ "PCMGSetCycleType" 6594b9ad928SBarry Smith /*@ 6600d353602SBarry Smith PCMGSetCycleType - Sets the type cycles to use. Use PCMGSetCycleTypeOnLevel() for more 6614b9ad928SBarry Smith complicated cycling. 6624b9ad928SBarry Smith 6634b9ad928SBarry Smith Collective on PC 6644b9ad928SBarry Smith 6654b9ad928SBarry Smith Input Parameters: 666c2be2410SBarry Smith + pc - the multigrid context 6670d353602SBarry Smith - PC_MG_CYCLE_V or PC_MG_CYCLE_W 6684b9ad928SBarry Smith 6694b9ad928SBarry Smith Options Database Key: 6700d353602SBarry Smith $ -pc_mg_cycle_type v or w 6714b9ad928SBarry Smith 6724b9ad928SBarry Smith Level: advanced 6734b9ad928SBarry Smith 6744b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 6754b9ad928SBarry Smith 6760d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel() 6774b9ad928SBarry Smith @*/ 6780d353602SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetCycleType(PC pc,PCMGCycleType n) 6794b9ad928SBarry Smith { 6809dcbbd2bSBarry Smith PC_MG **mg; 68179416396SBarry Smith PetscInt i,levels; 6824b9ad928SBarry Smith 6834b9ad928SBarry Smith PetscFunctionBegin; 6844482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 6859dcbbd2bSBarry Smith mg = (PC_MG**)pc->data; 6864b9ad928SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 6874b9ad928SBarry Smith levels = mg[0]->levels; 6884b9ad928SBarry Smith 6894b9ad928SBarry Smith for (i=0; i<levels; i++) { 6904b9ad928SBarry Smith mg[i]->cycles = n; 6914b9ad928SBarry Smith } 6924b9ad928SBarry Smith PetscFunctionReturn(0); 6934b9ad928SBarry Smith } 6944b9ad928SBarry Smith 6954b9ad928SBarry Smith #undef __FUNCT__ 696*8cc2d5dfSBarry Smith #define __FUNCT__ "PCMGMultiplicativeSetCycles" 697*8cc2d5dfSBarry Smith /*@ 698*8cc2d5dfSBarry Smith PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step 699*8cc2d5dfSBarry Smith of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used 700*8cc2d5dfSBarry Smith 701*8cc2d5dfSBarry Smith Collective on PC 702*8cc2d5dfSBarry Smith 703*8cc2d5dfSBarry Smith Input Parameters: 704*8cc2d5dfSBarry Smith + pc - the multigrid context 705*8cc2d5dfSBarry Smith - n - number of cycles (default is 1) 706*8cc2d5dfSBarry Smith 707*8cc2d5dfSBarry Smith Options Database Key: 708*8cc2d5dfSBarry Smith $ -pc_mg_multiplicative_cycles n 709*8cc2d5dfSBarry Smith 710*8cc2d5dfSBarry Smith Level: advanced 711*8cc2d5dfSBarry Smith 712*8cc2d5dfSBarry Smith Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType() 713*8cc2d5dfSBarry Smith 714*8cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 715*8cc2d5dfSBarry Smith 716*8cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType() 717*8cc2d5dfSBarry Smith @*/ 718*8cc2d5dfSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGMultiplicativeSetCycles(PC pc,PetscInt n) 719*8cc2d5dfSBarry Smith { 720*8cc2d5dfSBarry Smith PC_MG **mg; 721*8cc2d5dfSBarry Smith PetscInt i,levels; 722*8cc2d5dfSBarry Smith 723*8cc2d5dfSBarry Smith PetscFunctionBegin; 724*8cc2d5dfSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 725*8cc2d5dfSBarry Smith mg = (PC_MG**)pc->data; 726*8cc2d5dfSBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 727*8cc2d5dfSBarry Smith levels = mg[0]->levels; 728*8cc2d5dfSBarry Smith 729*8cc2d5dfSBarry Smith for (i=0; i<levels; i++) { 730*8cc2d5dfSBarry Smith mg[i]->cyclesperpcapply = n; 731*8cc2d5dfSBarry Smith } 732*8cc2d5dfSBarry Smith PetscFunctionReturn(0); 733*8cc2d5dfSBarry Smith } 734*8cc2d5dfSBarry Smith 735*8cc2d5dfSBarry Smith #undef __FUNCT__ 7369dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetGalerkin" 737c2be2410SBarry Smith /*@ 73897177400SBarry Smith PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the 739c2be2410SBarry Smith finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t 740c2be2410SBarry Smith 741c2be2410SBarry Smith Collective on PC 742c2be2410SBarry Smith 743c2be2410SBarry Smith Input Parameters: 7443fc8bf9cSBarry Smith . pc - the multigrid context 745c2be2410SBarry Smith 746c2be2410SBarry Smith Options Database Key: 747c2be2410SBarry Smith $ -pc_mg_galerkin 748c2be2410SBarry Smith 749c2be2410SBarry Smith Level: intermediate 750c2be2410SBarry Smith 751c2be2410SBarry Smith .keywords: MG, set, Galerkin 752c2be2410SBarry Smith 7533fc8bf9cSBarry Smith .seealso: PCMGGetGalerkin() 7543fc8bf9cSBarry Smith 755c2be2410SBarry Smith @*/ 75697177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetGalerkin(PC pc) 757c2be2410SBarry Smith { 7589dcbbd2bSBarry Smith PC_MG **mg; 759c2be2410SBarry Smith PetscInt i,levels; 760c2be2410SBarry Smith 761c2be2410SBarry Smith PetscFunctionBegin; 762c2be2410SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 7639dcbbd2bSBarry Smith mg = (PC_MG**)pc->data; 764c2be2410SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 765c2be2410SBarry Smith levels = mg[0]->levels; 766c2be2410SBarry Smith 767c2be2410SBarry Smith for (i=0; i<levels; i++) { 768c2be2410SBarry Smith mg[i]->galerkin = PETSC_TRUE; 769c2be2410SBarry Smith } 770c2be2410SBarry Smith PetscFunctionReturn(0); 771c2be2410SBarry Smith } 772c2be2410SBarry Smith 773c2be2410SBarry Smith #undef __FUNCT__ 7743fc8bf9cSBarry Smith #define __FUNCT__ "PCMGGetGalerkin" 7753fc8bf9cSBarry Smith /*@ 7763fc8bf9cSBarry Smith PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. 7773fc8bf9cSBarry Smith A_i-1 = r_i * A_i * r_i^t 7783fc8bf9cSBarry Smith 7793fc8bf9cSBarry Smith Not Collective 7803fc8bf9cSBarry Smith 7813fc8bf9cSBarry Smith Input Parameter: 7823fc8bf9cSBarry Smith . pc - the multigrid context 7833fc8bf9cSBarry Smith 7843fc8bf9cSBarry Smith Output Parameter: 7853fc8bf9cSBarry Smith . gelerkin - PETSC_TRUE or PETSC_FALSE 7863fc8bf9cSBarry Smith 7873fc8bf9cSBarry Smith Options Database Key: 7883fc8bf9cSBarry Smith $ -pc_mg_galerkin 7893fc8bf9cSBarry Smith 7903fc8bf9cSBarry Smith Level: intermediate 7913fc8bf9cSBarry Smith 7923fc8bf9cSBarry Smith .keywords: MG, set, Galerkin 7933fc8bf9cSBarry Smith 7943fc8bf9cSBarry Smith .seealso: PCMGSetGalerkin() 7953fc8bf9cSBarry Smith 7963fc8bf9cSBarry Smith @*/ 7973fc8bf9cSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetGalerkin(PC pc,PetscTruth *galerkin) 7983fc8bf9cSBarry Smith { 7993fc8bf9cSBarry Smith PC_MG **mg; 8003fc8bf9cSBarry Smith 8013fc8bf9cSBarry Smith PetscFunctionBegin; 8023fc8bf9cSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 8033fc8bf9cSBarry Smith mg = (PC_MG**)pc->data; 8043fc8bf9cSBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 8053fc8bf9cSBarry Smith *galerkin = mg[0]->galerkin; 8063fc8bf9cSBarry Smith PetscFunctionReturn(0); 8073fc8bf9cSBarry Smith } 8083fc8bf9cSBarry Smith 8093fc8bf9cSBarry Smith #undef __FUNCT__ 8109dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothDown" 8114b9ad928SBarry Smith /*@ 81297177400SBarry Smith PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to 81397177400SBarry Smith use on all levels. Use PCMGGetSmootherDown() to set different 8144b9ad928SBarry Smith pre-smoothing steps on different levels. 8154b9ad928SBarry Smith 8164b9ad928SBarry Smith Collective on PC 8174b9ad928SBarry Smith 8184b9ad928SBarry Smith Input Parameters: 8194b9ad928SBarry Smith + mg - the multigrid context 8204b9ad928SBarry Smith - n - the number of smoothing steps 8214b9ad928SBarry Smith 8224b9ad928SBarry Smith Options Database Key: 8234b9ad928SBarry Smith . -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps 8244b9ad928SBarry Smith 8254b9ad928SBarry Smith Level: advanced 8264b9ad928SBarry Smith 8274b9ad928SBarry Smith .keywords: MG, smooth, down, pre-smoothing, steps, multigrid 8284b9ad928SBarry Smith 82997177400SBarry Smith .seealso: PCMGSetNumberSmoothUp() 8304b9ad928SBarry Smith @*/ 83197177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothDown(PC pc,PetscInt n) 8324b9ad928SBarry Smith { 8339dcbbd2bSBarry Smith PC_MG **mg; 8346849ba73SBarry Smith PetscErrorCode ierr; 83579416396SBarry Smith PetscInt i,levels; 8364b9ad928SBarry Smith 8374b9ad928SBarry Smith PetscFunctionBegin; 8384482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 8399dcbbd2bSBarry Smith mg = (PC_MG**)pc->data; 8404b9ad928SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 8414b9ad928SBarry Smith levels = mg[0]->levels; 8424b9ad928SBarry Smith 843b05257ddSBarry Smith for (i=1; i<levels; i++) { 8444b9ad928SBarry Smith /* make sure smoother up and down are different */ 84597177400SBarry Smith ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 8464b9ad928SBarry Smith ierr = KSPSetTolerances(mg[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 8474b9ad928SBarry Smith mg[i]->default_smoothd = n; 8484b9ad928SBarry Smith } 8494b9ad928SBarry Smith PetscFunctionReturn(0); 8504b9ad928SBarry Smith } 8514b9ad928SBarry Smith 8524b9ad928SBarry Smith #undef __FUNCT__ 8539dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothUp" 8544b9ad928SBarry Smith /*@ 85597177400SBarry Smith PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use 85697177400SBarry Smith on all levels. Use PCMGGetSmootherUp() to set different numbers of 8574b9ad928SBarry Smith post-smoothing steps on different levels. 8584b9ad928SBarry Smith 8594b9ad928SBarry Smith Collective on PC 8604b9ad928SBarry Smith 8614b9ad928SBarry Smith Input Parameters: 8624b9ad928SBarry Smith + mg - the multigrid context 8634b9ad928SBarry Smith - n - the number of smoothing steps 8644b9ad928SBarry Smith 8654b9ad928SBarry Smith Options Database Key: 8664b9ad928SBarry Smith . -pc_mg_smoothup <n> - Sets number of post-smoothing steps 8674b9ad928SBarry Smith 8684b9ad928SBarry Smith Level: advanced 8694b9ad928SBarry Smith 8704b9ad928SBarry Smith Note: this does not set a value on the coarsest grid, since we assume that 871a8c7a070SBarry Smith there is no separate smooth up on the coarsest grid. 8724b9ad928SBarry Smith 8734b9ad928SBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid 8744b9ad928SBarry Smith 87597177400SBarry Smith .seealso: PCMGSetNumberSmoothDown() 8764b9ad928SBarry Smith @*/ 87797177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothUp(PC pc,PetscInt n) 8784b9ad928SBarry Smith { 8799dcbbd2bSBarry Smith PC_MG **mg; 8806849ba73SBarry Smith PetscErrorCode ierr; 88179416396SBarry Smith PetscInt i,levels; 8824b9ad928SBarry Smith 8834b9ad928SBarry Smith PetscFunctionBegin; 8844482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 8859dcbbd2bSBarry Smith mg = (PC_MG**)pc->data; 8864b9ad928SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 8874b9ad928SBarry Smith levels = mg[0]->levels; 8884b9ad928SBarry Smith 8894b9ad928SBarry Smith for (i=1; i<levels; i++) { 8904b9ad928SBarry Smith /* make sure smoother up and down are different */ 89197177400SBarry Smith ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 8924b9ad928SBarry Smith ierr = KSPSetTolerances(mg[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 8934b9ad928SBarry Smith mg[i]->default_smoothu = n; 8944b9ad928SBarry Smith } 8954b9ad928SBarry Smith PetscFunctionReturn(0); 8964b9ad928SBarry Smith } 8974b9ad928SBarry Smith 8984b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/ 8994b9ad928SBarry Smith 9003b09bd56SBarry Smith /*MC 901ccb205f8SBarry Smith PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional 9023b09bd56SBarry Smith information about the coarser grid matrices and restriction/interpolation operators. 9033b09bd56SBarry Smith 9043b09bd56SBarry Smith Options Database Keys: 9053b09bd56SBarry Smith + -pc_mg_levels <nlevels> - number of levels including finest 9060d353602SBarry Smith . -pc_mg_cycles v or w 90779416396SBarry Smith . -pc_mg_smoothup <n> - number of smoothing steps after interpolation 9083b09bd56SBarry Smith . -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator 9093b09bd56SBarry Smith . -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default 9103b09bd56SBarry Smith . -pc_mg_log - log information about time spent on each level of the solver 9113b09bd56SBarry Smith . -pc_mg_monitor - print information on the multigrid convergence 91268eff7e6SBarry Smith . -pc_mg_galerkin - use Galerkin process to compute coarser operators 9133b09bd56SBarry Smith - -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 9143b09bd56SBarry Smith to the Socket viewer for reading from Matlab. 9153b09bd56SBarry Smith 9163b09bd56SBarry Smith Notes: 9173b09bd56SBarry Smith 9183b09bd56SBarry Smith Level: intermediate 9193b09bd56SBarry Smith 9203b09bd56SBarry Smith Concepts: multigrid 9213b09bd56SBarry Smith 9223b09bd56SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 9230d353602SBarry Smith PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(), 92497177400SBarry Smith PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 92597177400SBarry Smith PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 9260d353602SBarry Smith PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 9273b09bd56SBarry Smith M*/ 9283b09bd56SBarry Smith 9294b9ad928SBarry Smith EXTERN_C_BEGIN 9304b9ad928SBarry Smith #undef __FUNCT__ 9314b9ad928SBarry Smith #define __FUNCT__ "PCCreate_MG" 932dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_MG(PC pc) 9334b9ad928SBarry Smith { 9344b9ad928SBarry Smith PetscFunctionBegin; 9354b9ad928SBarry Smith pc->ops->apply = PCApply_MG; 9364b9ad928SBarry Smith pc->ops->setup = PCSetUp_MG; 9374b9ad928SBarry Smith pc->ops->destroy = PCDestroy_MG; 9384b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_MG; 9394b9ad928SBarry Smith pc->ops->view = PCView_MG; 9404b9ad928SBarry Smith 9414b9ad928SBarry Smith pc->data = (void*)0; 9424b9ad928SBarry Smith PetscFunctionReturn(0); 9434b9ad928SBarry Smith } 9444b9ad928SBarry Smith EXTERN_C_END 945