1dba47a55SKris Buschelman 24b9ad928SBarry Smith /* 34b9ad928SBarry Smith Defines the multigrid preconditioner interface. 44b9ad928SBarry Smith */ 5c6db04a5SJed Brown #include <../src/ksp/pc/impls/mg/mgimpl.h> /*I "petscpcmg.h" I*/ 64b9ad928SBarry Smith 74b9ad928SBarry Smith 84b9ad928SBarry Smith #undef __FUNCT__ 99dcbbd2bSBarry Smith #define __FUNCT__ "PCMGMCycle_Private" 1031567311SBarry Smith PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG_Levels **mglevelsin,PCRichardsonConvergedReason *reason) 114b9ad928SBarry Smith { 1231567311SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1331567311SBarry Smith PC_MG_Levels *mgc,*mglevels = *mglevelsin; 146849ba73SBarry Smith PetscErrorCode ierr; 1531567311SBarry Smith PetscInt cycles = (mglevels->level == 1) ? 1 : (PetscInt) mglevels->cycles; 164b9ad928SBarry Smith 174b9ad928SBarry Smith PetscFunctionBegin; 1863e6d426SJed Brown if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 1931567311SBarry Smith ierr = KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x);CHKERRQ(ierr); /* pre-smooth */ 2063e6d426SJed Brown if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 2131567311SBarry Smith if (mglevels->level) { /* not the coarsest grid */ 2263e6d426SJed Brown if (mglevels->eventresidual) {ierr = PetscLogEventBegin(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);} 2331567311SBarry Smith ierr = (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);CHKERRQ(ierr); 2463e6d426SJed Brown if (mglevels->eventresidual) {ierr = PetscLogEventEnd(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);} 254b9ad928SBarry Smith 264b9ad928SBarry Smith /* if on finest level and have convergence criteria set */ 2731567311SBarry Smith if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) { 284b9ad928SBarry Smith PetscReal rnorm; 2931567311SBarry Smith ierr = VecNorm(mglevels->r,NORM_2,&rnorm);CHKERRQ(ierr); 304b9ad928SBarry Smith if (rnorm <= mg->ttol) { 3170441072SBarry Smith if (rnorm < mg->abstol) { 324d0a8057SBarry Smith *reason = PCRICHARDSON_CONVERGED_ATOL; 331e2582c4SBarry Smith ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %G is less than absolute tolerance %G\n",rnorm,mg->abstol);CHKERRQ(ierr); 344b9ad928SBarry Smith } else { 354d0a8057SBarry Smith *reason = PCRICHARDSON_CONVERGED_RTOL; 361e2582c4SBarry Smith ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %G is less than relative tolerance times initial residual norm %G\n",rnorm,mg->ttol);CHKERRQ(ierr); 374b9ad928SBarry Smith } 384b9ad928SBarry Smith PetscFunctionReturn(0); 394b9ad928SBarry Smith } 404b9ad928SBarry Smith } 414b9ad928SBarry Smith 4231567311SBarry Smith mgc = *(mglevelsin - 1); 4363e6d426SJed Brown if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 4431567311SBarry Smith ierr = MatRestrict(mglevels->restrct,mglevels->r,mgc->b);CHKERRQ(ierr); 4563e6d426SJed Brown if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 46efb30889SBarry Smith ierr = VecSet(mgc->x,0.0);CHKERRQ(ierr); 474b9ad928SBarry Smith while (cycles--) { 4831567311SBarry Smith ierr = PCMGMCycle_Private(pc,mglevelsin-1,reason);CHKERRQ(ierr); 494b9ad928SBarry Smith } 5063e6d426SJed Brown if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 5131567311SBarry Smith ierr = MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x);CHKERRQ(ierr); 5263e6d426SJed Brown if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 5363e6d426SJed Brown if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 5431567311SBarry Smith ierr = KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x);CHKERRQ(ierr); /* post smooth */ 5563e6d426SJed Brown if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 564b9ad928SBarry Smith } 574b9ad928SBarry Smith PetscFunctionReturn(0); 584b9ad928SBarry Smith } 594b9ad928SBarry Smith 604b9ad928SBarry Smith #undef __FUNCT__ 614b9ad928SBarry Smith #define __FUNCT__ "PCApplyRichardson_MG" 62ace3abfcSBarry Smith static PetscErrorCode PCApplyRichardson_MG(PC pc,Vec b,Vec x,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscBool zeroguess,PetscInt *outits,PCRichardsonConvergedReason *reason) 634b9ad928SBarry Smith { 64f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 65f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 66dfbe8321SBarry Smith PetscErrorCode ierr; 67f3fbd535SBarry Smith PetscInt levels = mglevels[0]->levels,i; 684b9ad928SBarry Smith 694b9ad928SBarry Smith PetscFunctionBegin; 70f3fbd535SBarry Smith mglevels[levels-1]->b = b; 71f3fbd535SBarry Smith mglevels[levels-1]->x = x; 724b9ad928SBarry Smith 7331567311SBarry Smith mg->rtol = rtol; 7431567311SBarry Smith mg->abstol = abstol; 7531567311SBarry Smith mg->dtol = dtol; 764b9ad928SBarry Smith if (rtol) { 774b9ad928SBarry Smith /* compute initial residual norm for relative convergence test */ 784b9ad928SBarry Smith PetscReal rnorm; 797319c654SBarry Smith if (zeroguess) { 807319c654SBarry Smith ierr = VecNorm(b,NORM_2,&rnorm);CHKERRQ(ierr); 817319c654SBarry Smith } else { 82f3fbd535SBarry Smith ierr = (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);CHKERRQ(ierr); 834b9ad928SBarry Smith ierr = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr); 847319c654SBarry Smith } 8531567311SBarry Smith mg->ttol = PetscMax(rtol*rnorm,abstol); 86*2fa5cd67SKarl Rupp } else if (abstol) mg->ttol = abstol; 87*2fa5cd67SKarl Rupp else mg->ttol = 0.0; 884b9ad928SBarry Smith 894d0a8057SBarry Smith /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't 90336babb1SJed Brown stop prematurely due to small residual */ 914d0a8057SBarry Smith for (i=1; i<levels; i++) { 92f3fbd535SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 93f3fbd535SBarry Smith if (mglevels[i]->smoothu != mglevels[i]->smoothd) { 94f3fbd535SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 954b9ad928SBarry Smith } 964d0a8057SBarry Smith } 974d0a8057SBarry Smith 984d0a8057SBarry Smith *reason = (PCRichardsonConvergedReason)0; 994d0a8057SBarry Smith for (i=0; i<its; i++) { 100f3fbd535SBarry Smith ierr = PCMGMCycle_Private(pc,mglevels+levels-1,reason);CHKERRQ(ierr); 1014d0a8057SBarry Smith if (*reason) break; 1024d0a8057SBarry Smith } 1034d0a8057SBarry Smith if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS; 1044d0a8057SBarry Smith *outits = i; 1054b9ad928SBarry Smith PetscFunctionReturn(0); 1064b9ad928SBarry Smith } 1074b9ad928SBarry Smith 1084b9ad928SBarry Smith #undef __FUNCT__ 1093aeaf226SBarry Smith #define __FUNCT__ "PCReset_MG" 1103aeaf226SBarry Smith PetscErrorCode PCReset_MG(PC pc) 1113aeaf226SBarry Smith { 1123aeaf226SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1133aeaf226SBarry Smith PC_MG_Levels **mglevels = mg->levels; 1143aeaf226SBarry Smith PetscErrorCode ierr; 1153aeaf226SBarry Smith PetscInt i,n; 1163aeaf226SBarry Smith 1173aeaf226SBarry Smith PetscFunctionBegin; 1183aeaf226SBarry Smith if (mglevels) { 1193aeaf226SBarry Smith n = mglevels[0]->levels; 1203aeaf226SBarry Smith for (i=0; i<n-1; i++) { 1213aeaf226SBarry Smith ierr = VecDestroy(&mglevels[i+1]->r);CHKERRQ(ierr); 1223aeaf226SBarry Smith ierr = VecDestroy(&mglevels[i]->b);CHKERRQ(ierr); 1233aeaf226SBarry Smith ierr = VecDestroy(&mglevels[i]->x);CHKERRQ(ierr); 1243aeaf226SBarry Smith ierr = MatDestroy(&mglevels[i+1]->restrct);CHKERRQ(ierr); 1253aeaf226SBarry Smith ierr = MatDestroy(&mglevels[i+1]->interpolate);CHKERRQ(ierr); 12673250ac0SBarry Smith ierr = VecDestroy(&mglevels[i+1]->rscale);CHKERRQ(ierr); 1273aeaf226SBarry Smith } 1283aeaf226SBarry Smith 1293aeaf226SBarry Smith for (i=0; i<n; i++) { 1303aeaf226SBarry Smith ierr = MatDestroy(&mglevels[i]->A);CHKERRQ(ierr); 1313aeaf226SBarry Smith if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 1323aeaf226SBarry Smith ierr = KSPReset(mglevels[i]->smoothd);CHKERRQ(ierr); 1333aeaf226SBarry Smith } 1343aeaf226SBarry Smith ierr = KSPReset(mglevels[i]->smoothu);CHKERRQ(ierr); 1353aeaf226SBarry Smith } 1363aeaf226SBarry Smith } 1373aeaf226SBarry Smith PetscFunctionReturn(0); 1383aeaf226SBarry Smith } 1393aeaf226SBarry Smith 1403aeaf226SBarry Smith #undef __FUNCT__ 1419dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetLevels" 1424b9ad928SBarry Smith /*@C 14397177400SBarry Smith PCMGSetLevels - Sets the number of levels to use with MG. 1444b9ad928SBarry Smith Must be called before any other MG routine. 1454b9ad928SBarry Smith 146ad4df100SBarry Smith Logically Collective on PC 1474b9ad928SBarry Smith 1484b9ad928SBarry Smith Input Parameters: 1494b9ad928SBarry Smith + pc - the preconditioner context 1504b9ad928SBarry Smith . levels - the number of levels 1514b9ad928SBarry Smith - comms - optional communicators for each level; this is to allow solving the coarser problems 1524b9ad928SBarry Smith on smaller sets of processors. Use PETSC_NULL_OBJECT for default in Fortran 1534b9ad928SBarry Smith 1544b9ad928SBarry Smith Level: intermediate 1554b9ad928SBarry Smith 1564b9ad928SBarry Smith Notes: 1574b9ad928SBarry Smith If the number of levels is one then the multigrid uses the -mg_levels prefix 1584b9ad928SBarry Smith for setting the level options rather than the -mg_coarse prefix. 1594b9ad928SBarry Smith 1604b9ad928SBarry Smith .keywords: MG, set, levels, multigrid 1614b9ad928SBarry Smith 16297177400SBarry Smith .seealso: PCMGSetType(), PCMGGetLevels() 1634b9ad928SBarry Smith @*/ 1647087cfbeSBarry Smith PetscErrorCode PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms) 1654b9ad928SBarry Smith { 166dfbe8321SBarry Smith PetscErrorCode ierr; 167f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 168f3fbd535SBarry Smith MPI_Comm comm = ((PetscObject)pc)->comm; 1693aeaf226SBarry Smith PC_MG_Levels **mglevels = mg->levels; 170f3fbd535SBarry Smith PetscInt i; 171f3fbd535SBarry Smith PetscMPIInt size; 172f3fbd535SBarry Smith const char *prefix; 173f3fbd535SBarry Smith PC ipc; 1743aeaf226SBarry Smith PetscInt n; 1754b9ad928SBarry Smith 1764b9ad928SBarry Smith PetscFunctionBegin; 1770700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 178548767e0SJed Brown PetscValidLogicalCollectiveInt(pc,levels,2); 179548767e0SJed Brown if (mg->nlevels == levels) PetscFunctionReturn(0); 1803aeaf226SBarry Smith if (mglevels) { 1813aeaf226SBarry Smith /* changing the number of levels so free up the previous stuff */ 1823aeaf226SBarry Smith ierr = PCReset_MG(pc);CHKERRQ(ierr); 1833aeaf226SBarry Smith n = mglevels[0]->levels; 1843aeaf226SBarry Smith for (i=0; i<n; i++) { 1853aeaf226SBarry Smith if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 1863aeaf226SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr); 1873aeaf226SBarry Smith } 1883aeaf226SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr); 1893aeaf226SBarry Smith ierr = PetscFree(mglevels[i]);CHKERRQ(ierr); 1903aeaf226SBarry Smith } 1913aeaf226SBarry Smith ierr = PetscFree(mg->levels);CHKERRQ(ierr); 1923aeaf226SBarry Smith } 193f3fbd535SBarry Smith 194f3fbd535SBarry Smith mg->nlevels = levels; 195f3fbd535SBarry Smith 196f3fbd535SBarry Smith ierr = PetscMalloc(levels*sizeof(PC_MG*),&mglevels);CHKERRQ(ierr); 197f3fbd535SBarry Smith ierr = PetscLogObjectMemory(pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr); 198f3fbd535SBarry Smith 199f3fbd535SBarry Smith ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 200f3fbd535SBarry Smith 201a9db87a2SMatthew G Knepley mg->stageApply = 0; 202f3fbd535SBarry Smith for (i=0; i<levels; i++) { 203f3fbd535SBarry Smith ierr = PetscNewLog(pc,PC_MG_Levels,&mglevels[i]);CHKERRQ(ierr); 204*2fa5cd67SKarl Rupp 205f3fbd535SBarry Smith mglevels[i]->level = i; 206f3fbd535SBarry Smith mglevels[i]->levels = levels; 207f3fbd535SBarry Smith mglevels[i]->cycles = PC_MG_CYCLE_V; 208336babb1SJed Brown mg->default_smoothu = 2; 209336babb1SJed Brown mg->default_smoothd = 2; 21063e6d426SJed Brown mglevels[i]->eventsmoothsetup = 0; 21163e6d426SJed Brown mglevels[i]->eventsmoothsolve = 0; 21263e6d426SJed Brown mglevels[i]->eventresidual = 0; 21363e6d426SJed Brown mglevels[i]->eventinterprestrict = 0; 214f3fbd535SBarry Smith 215f3fbd535SBarry Smith if (comms) comm = comms[i]; 216f3fbd535SBarry Smith ierr = KSPCreate(comm,&mglevels[i]->smoothd);CHKERRQ(ierr); 217336babb1SJed Brown ierr = KSPSetType(mglevels[i]->smoothd,KSPCHEBYSHEV);CHKERRQ(ierr); 218336babb1SJed Brown ierr = KSPSetConvergenceTest(mglevels[i]->smoothd,KSPSkipConverged,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 219336babb1SJed Brown ierr = KSPSetNormType(mglevels[i]->smoothd,KSP_NORM_NONE);CHKERRQ(ierr); 220336babb1SJed Brown ierr = KSPGetPC(mglevels[i]->smoothd,&ipc);CHKERRQ(ierr); 221336babb1SJed Brown ierr = PCSetType(ipc,PCSOR);CHKERRQ(ierr); 222f3fbd535SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr); 223336babb1SJed Brown ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, i ? mg->default_smoothd : 1);CHKERRQ(ierr); 224f3fbd535SBarry Smith ierr = KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);CHKERRQ(ierr); 225f3fbd535SBarry Smith 226f3fbd535SBarry Smith /* do special stuff for coarse grid */ 227f3fbd535SBarry Smith if (!i && levels > 1) { 228f3fbd535SBarry Smith ierr = KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");CHKERRQ(ierr); 229f3fbd535SBarry Smith 2309dbfc187SHong Zhang /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */ 231f3fbd535SBarry Smith ierr = KSPSetType(mglevels[0]->smoothd,KSPPREONLY);CHKERRQ(ierr); 232f3fbd535SBarry Smith ierr = KSPGetPC(mglevels[0]->smoothd,&ipc);CHKERRQ(ierr); 233f3fbd535SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 234f3fbd535SBarry Smith if (size > 1) { 2359dbfc187SHong Zhang KSP innerksp; 2369dbfc187SHong Zhang PC innerpc; 237f3fbd535SBarry Smith ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr); 2389dbfc187SHong Zhang ierr = PCRedundantGetKSP(ipc,&innerksp);CHKERRQ(ierr); 2399dbfc187SHong Zhang ierr = KSPGetPC(innerksp,&innerpc);CHKERRQ(ierr); 2409dbfc187SHong Zhang ierr = PCFactorSetShiftType(innerpc,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr); 241f3fbd535SBarry Smith } else { 242f3fbd535SBarry Smith ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr); 2439dbfc187SHong Zhang ierr = PCFactorSetShiftType(ipc,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr); 244f3fbd535SBarry Smith } 245f3fbd535SBarry Smith } else { 246f3fbd535SBarry Smith char tprefix[128]; 247f3fbd535SBarry Smith sprintf(tprefix,"mg_levels_%d_",(int)i); 248f3fbd535SBarry Smith ierr = KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);CHKERRQ(ierr); 249f3fbd535SBarry Smith } 250f3fbd535SBarry Smith ierr = PetscLogObjectParent(pc,mglevels[i]->smoothd);CHKERRQ(ierr); 251*2fa5cd67SKarl Rupp 252f3fbd535SBarry Smith mglevels[i]->smoothu = mglevels[i]->smoothd; 25331567311SBarry Smith mg->rtol = 0.0; 25431567311SBarry Smith mg->abstol = 0.0; 25531567311SBarry Smith mg->dtol = 0.0; 25631567311SBarry Smith mg->ttol = 0.0; 25731567311SBarry Smith mg->cyclesperpcapply = 1; 258f3fbd535SBarry Smith } 25931567311SBarry Smith mg->am = PC_MG_MULTIPLICATIVE; 260f3fbd535SBarry Smith mg->levels = mglevels; 2614b9ad928SBarry Smith pc->ops->applyrichardson = PCApplyRichardson_MG; 2624b9ad928SBarry Smith PetscFunctionReturn(0); 2634b9ad928SBarry Smith } 2644b9ad928SBarry Smith 265c07bf074SBarry Smith 266c07bf074SBarry Smith #undef __FUNCT__ 267c07bf074SBarry Smith #define __FUNCT__ "PCDestroy_MG" 268c07bf074SBarry Smith PetscErrorCode PCDestroy_MG(PC pc) 269c07bf074SBarry Smith { 270c07bf074SBarry Smith PetscErrorCode ierr; 271a06653b4SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 272a06653b4SBarry Smith PC_MG_Levels **mglevels = mg->levels; 273a06653b4SBarry Smith PetscInt i,n; 274c07bf074SBarry Smith 275c07bf074SBarry Smith PetscFunctionBegin; 276a06653b4SBarry Smith ierr = PCReset_MG(pc);CHKERRQ(ierr); 277a06653b4SBarry Smith if (mglevels) { 278a06653b4SBarry Smith n = mglevels[0]->levels; 279a06653b4SBarry Smith for (i=0; i<n; i++) { 280a06653b4SBarry Smith if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 2816bf464f9SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr); 282a06653b4SBarry Smith } 2836bf464f9SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr); 284a06653b4SBarry Smith ierr = PetscFree(mglevels[i]);CHKERRQ(ierr); 285a06653b4SBarry Smith } 286a06653b4SBarry Smith ierr = PetscFree(mg->levels);CHKERRQ(ierr); 287a06653b4SBarry Smith } 288c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 289f3fbd535SBarry Smith PetscFunctionReturn(0); 290f3fbd535SBarry Smith } 291f3fbd535SBarry Smith 292f3fbd535SBarry Smith 293f3fbd535SBarry Smith 29409573ac7SBarry Smith extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**); 29509573ac7SBarry Smith extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**); 29609573ac7SBarry Smith extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**); 297f3fbd535SBarry Smith 298f3fbd535SBarry Smith /* 299f3fbd535SBarry Smith PCApply_MG - Runs either an additive, multiplicative, Kaskadic 300f3fbd535SBarry Smith or full cycle of multigrid. 301f3fbd535SBarry Smith 302f3fbd535SBarry Smith Note: 303f3fbd535SBarry Smith A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle(). 304f3fbd535SBarry Smith */ 305f3fbd535SBarry Smith #undef __FUNCT__ 306f3fbd535SBarry Smith #define __FUNCT__ "PCApply_MG" 307f3fbd535SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x) 308f3fbd535SBarry Smith { 309f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 310f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 311f3fbd535SBarry Smith PetscErrorCode ierr; 312f3fbd535SBarry Smith PetscInt levels = mglevels[0]->levels,i; 313f3fbd535SBarry Smith 314f3fbd535SBarry Smith PetscFunctionBegin; 315a9db87a2SMatthew G Knepley if (mg->stageApply) {ierr = PetscLogStagePush(mg->stageApply);CHKERRQ(ierr);} 316e1d8e5deSBarry Smith /* When the DM is supplying the matrix then it will not exist until here */ 317298cc208SBarry Smith for (i=0; i<levels; i++) { 318e1d8e5deSBarry Smith if (!mglevels[i]->A) { 319e1d8e5deSBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 320298cc208SBarry Smith ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr); 321e1d8e5deSBarry Smith } 322e1d8e5deSBarry Smith } 323e1d8e5deSBarry Smith 324f3fbd535SBarry Smith mglevels[levels-1]->b = b; 325f3fbd535SBarry Smith mglevels[levels-1]->x = x; 32631567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 327f3fbd535SBarry Smith ierr = VecSet(x,0.0);CHKERRQ(ierr); 32831567311SBarry Smith for (i=0; i<mg->cyclesperpcapply; i++) { 329f3fbd535SBarry Smith ierr = PCMGMCycle_Private(pc,mglevels+levels-1,PETSC_NULL);CHKERRQ(ierr); 330f3fbd535SBarry Smith } 331*2fa5cd67SKarl Rupp } else if (mg->am == PC_MG_ADDITIVE) { 33231567311SBarry Smith ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr); 333*2fa5cd67SKarl Rupp } else if (mg->am == PC_MG_KASKADE) { 33431567311SBarry Smith ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr); 335*2fa5cd67SKarl Rupp } else { 336f3fbd535SBarry Smith ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr); 337f3fbd535SBarry Smith } 338a9db87a2SMatthew G Knepley if (mg->stageApply) {ierr = PetscLogStagePop();CHKERRQ(ierr);} 339f3fbd535SBarry Smith PetscFunctionReturn(0); 340f3fbd535SBarry Smith } 341f3fbd535SBarry Smith 342f3fbd535SBarry Smith 343f3fbd535SBarry Smith #undef __FUNCT__ 344f3fbd535SBarry Smith #define __FUNCT__ "PCSetFromOptions_MG" 345f3fbd535SBarry Smith PetscErrorCode PCSetFromOptions_MG(PC pc) 346f3fbd535SBarry Smith { 347f3fbd535SBarry Smith PetscErrorCode ierr; 348f3fbd535SBarry Smith PetscInt m,levels = 1,cycles; 349c91913d3SJed Brown PetscBool flg,set; 350f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 351f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 352f3fbd535SBarry Smith PCMGType mgtype; 353f3fbd535SBarry Smith PCMGCycleType mgctype; 354f3fbd535SBarry Smith 355f3fbd535SBarry Smith PetscFunctionBegin; 356f3fbd535SBarry Smith ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr); 3573aeaf226SBarry Smith if (!mg->levels) { 358f3fbd535SBarry Smith ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr); 359eb3f98d2SBarry Smith if (!flg && pc->dm) { 360eb3f98d2SBarry Smith ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr); 361eb3f98d2SBarry Smith levels++; 3623aeaf226SBarry Smith mg->usedmfornumberoflevels = PETSC_TRUE; 363eb3f98d2SBarry Smith } 364f3fbd535SBarry Smith ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr); 365f3fbd535SBarry Smith } 3663aeaf226SBarry Smith mglevels = mg->levels; 3673aeaf226SBarry Smith 368f3fbd535SBarry Smith mgctype = (PCMGCycleType) mglevels[0]->cycles; 369f3fbd535SBarry Smith ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr); 370f3fbd535SBarry Smith if (flg) { 371f3fbd535SBarry Smith ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr); 372*2fa5cd67SKarl Rupp } 373f3fbd535SBarry Smith flg = PETSC_FALSE; 374c91913d3SJed Brown ierr = PetscOptionsBool("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",flg,&flg,&set);CHKERRQ(ierr); 375c91913d3SJed Brown if (set) { 376c91913d3SJed Brown ierr = PCMGSetGalerkin(pc,flg);CHKERRQ(ierr); 377f3fbd535SBarry Smith } 378f3fbd535SBarry Smith ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr); 379f3fbd535SBarry Smith if (flg) { 380f3fbd535SBarry Smith ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr); 381f3fbd535SBarry Smith } 382f3fbd535SBarry Smith ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr); 383f3fbd535SBarry Smith if (flg) { 384f3fbd535SBarry Smith ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr); 385f3fbd535SBarry Smith } 38631567311SBarry Smith mgtype = mg->am; 387f3fbd535SBarry Smith ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr); 388f3fbd535SBarry Smith if (flg) { 389f3fbd535SBarry Smith ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr); 390f3fbd535SBarry Smith } 39131567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 39231567311SBarry Smith ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGSetLevels",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr); 393f3fbd535SBarry Smith if (flg) { 394f3fbd535SBarry Smith ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr); 395f3fbd535SBarry Smith } 396f3fbd535SBarry Smith } 397f3fbd535SBarry Smith flg = PETSC_FALSE; 398acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 399f3fbd535SBarry Smith if (flg) { 400f3fbd535SBarry Smith PetscInt i; 401f3fbd535SBarry Smith char eventname[128]; 40263e6d426SJed Brown if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 403f3fbd535SBarry Smith levels = mglevels[0]->levels; 404f3fbd535SBarry Smith for (i=0; i<levels; i++) { 405f3fbd535SBarry Smith sprintf(eventname,"MGSetup Level %d",(int)i); 40663e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr); 407f3fbd535SBarry Smith sprintf(eventname,"MGSmooth Level %d",(int)i); 40863e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr); 409f3fbd535SBarry Smith if (i) { 410f3fbd535SBarry Smith sprintf(eventname,"MGResid Level %d",(int)i); 41163e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr); 412f3fbd535SBarry Smith sprintf(eventname,"MGInterp Level %d",(int)i); 41363e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr); 414f3fbd535SBarry Smith } 415f3fbd535SBarry Smith } 416eec35531SMatthew G Knepley 417ec60ca73SSatish Balay #if defined(PETSC_USE_LOG) 418eec35531SMatthew G Knepley { 419eec35531SMatthew G Knepley const char *sname = "MG Apply"; 420eec35531SMatthew G Knepley PetscStageLog stageLog; 421eec35531SMatthew G Knepley PetscInt st; 422eec35531SMatthew G Knepley 423eec35531SMatthew G Knepley PetscFunctionBegin; 424eec35531SMatthew G Knepley ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr); 425eec35531SMatthew G Knepley for (st = 0; st < stageLog->numStages; ++st) { 426eec35531SMatthew G Knepley PetscBool same; 427eec35531SMatthew G Knepley 428eec35531SMatthew G Knepley ierr = PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);CHKERRQ(ierr); 429*2fa5cd67SKarl Rupp if (same) mg->stageApply = st; 430eec35531SMatthew G Knepley } 431eec35531SMatthew G Knepley if (!mg->stageApply) { 432eec35531SMatthew G Knepley ierr = PetscLogStageRegister(sname, &mg->stageApply);CHKERRQ(ierr); 433eec35531SMatthew G Knepley } 434eec35531SMatthew G Knepley } 435ec60ca73SSatish Balay #endif 436f3fbd535SBarry Smith } 437f3fbd535SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 438f3fbd535SBarry Smith PetscFunctionReturn(0); 439f3fbd535SBarry Smith } 440f3fbd535SBarry Smith 4416a6fc655SJed Brown const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0}; 4426a6fc655SJed Brown const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0}; 443f3fbd535SBarry Smith 444f3fbd535SBarry Smith #undef __FUNCT__ 445f3fbd535SBarry Smith #define __FUNCT__ "PCView_MG" 446f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer) 447f3fbd535SBarry Smith { 448f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 449f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 450f3fbd535SBarry Smith PetscErrorCode ierr; 451f3fbd535SBarry Smith PetscInt levels = mglevels[0]->levels,i; 452219b1045SBarry Smith PetscBool iascii,isbinary,isdraw; 453f3fbd535SBarry Smith 454f3fbd535SBarry Smith PetscFunctionBegin; 455251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 4565b0b0462SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 457219b1045SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 458f3fbd535SBarry Smith if (iascii) { 45931567311SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," MG: type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,(mglevels[0]->cycles == PC_MG_CYCLE_V) ? "v" : "w");CHKERRQ(ierr); 46031567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 46131567311SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr); 462f3fbd535SBarry Smith } 463218a07d4SBarry Smith if (mg->galerkin) { 464f3fbd535SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 4654f66f45eSBarry Smith } else { 4664f66f45eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 467f3fbd535SBarry Smith } 468f3fbd535SBarry Smith for (i=0; i<levels; i++) { 469f3fbd535SBarry Smith if (!i) { 47089cce641SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);CHKERRQ(ierr); 471f3fbd535SBarry Smith } else { 47289cce641SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 473f3fbd535SBarry Smith } 474f3fbd535SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 475f3fbd535SBarry Smith ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 476f3fbd535SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 477f3fbd535SBarry Smith if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) { 478f3fbd535SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 479f3fbd535SBarry Smith } else if (i) { 48089cce641SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 481f3fbd535SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 482f3fbd535SBarry Smith ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 483f3fbd535SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 484f3fbd535SBarry Smith } 485f3fbd535SBarry Smith } 4865b0b0462SBarry Smith } else if (isbinary) { 4875b0b0462SBarry Smith for (i=levels-1; i>=0; i--) { 4885b0b0462SBarry Smith ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 4895b0b0462SBarry Smith if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) { 4905b0b0462SBarry Smith ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 4915b0b0462SBarry Smith } 4925b0b0462SBarry Smith } 493219b1045SBarry Smith } else if (isdraw) { 494219b1045SBarry Smith PetscDraw draw; 495b4375e8dSPeter Brune PetscReal x,w,y,bottom,th; 496219b1045SBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 497219b1045SBarry Smith ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr); 498219b1045SBarry Smith ierr = PetscDrawStringGetSize(draw,PETSC_NULL,&th);CHKERRQ(ierr); 499219b1045SBarry Smith bottom = y - th; 500219b1045SBarry Smith for (i=levels-1; i>=0; i--) { 501b4375e8dSPeter Brune if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) { 502219b1045SBarry Smith ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr); 503219b1045SBarry Smith ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 504219b1045SBarry Smith ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 505b4375e8dSPeter Brune } else { 506b4375e8dSPeter Brune w = 0.5*PetscMin(1.0-x,x); 507b4375e8dSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr); 508b4375e8dSPeter Brune ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 509b4375e8dSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 510b4375e8dSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr); 511b4375e8dSPeter Brune ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 512b4375e8dSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 513b4375e8dSPeter Brune } 5141cd381d6SBarry Smith ierr = PetscDrawGetBoundingBox(draw,PETSC_NULL,&bottom,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 5151cd381d6SBarry Smith bottom -= th; 516219b1045SBarry Smith } 5175b0b0462SBarry Smith } 518f3fbd535SBarry Smith PetscFunctionReturn(0); 519f3fbd535SBarry Smith } 520f3fbd535SBarry Smith 521b45d2f2cSJed Brown #include <petsc-private/dmimpl.h> 522b45d2f2cSJed Brown #include <petsc-private/kspimpl.h> 523cab2e9ccSBarry Smith 524f3fbd535SBarry Smith /* 525f3fbd535SBarry Smith Calls setup for the KSP on each level 526f3fbd535SBarry Smith */ 527f3fbd535SBarry Smith #undef __FUNCT__ 528f3fbd535SBarry Smith #define __FUNCT__ "PCSetUp_MG" 529f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc) 530f3fbd535SBarry Smith { 531f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 532f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 533f3fbd535SBarry Smith PetscErrorCode ierr; 534f3fbd535SBarry Smith PetscInt i,n = mglevels[0]->levels; 53598e04984SBarry Smith PC cpc; 536649052a6SBarry Smith PetscBool preonly,lu,redundant,cholesky,svd,dump = PETSC_FALSE,opsset; 537f3fbd535SBarry Smith Mat dA,dB; 538f3fbd535SBarry Smith MatStructure uflag; 539f3fbd535SBarry Smith Vec tvec; 540218a07d4SBarry Smith DM *dms; 541649052a6SBarry Smith PetscViewer viewer = 0; 542f3fbd535SBarry Smith 543f3fbd535SBarry Smith PetscFunctionBegin; 54401bc414fSDmitry Karpeev /* FIX: Move this to PCSetFromOptions_MG? */ 5453aeaf226SBarry Smith if (mg->usedmfornumberoflevels) { 5463aeaf226SBarry Smith PetscInt levels; 5473aeaf226SBarry Smith ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr); 5483aeaf226SBarry Smith levels++; 5493aeaf226SBarry Smith if (levels > n) { /* the problem is now being solved on a finer grid */ 5503aeaf226SBarry Smith ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr); 5513aeaf226SBarry Smith n = levels; 5523aeaf226SBarry Smith ierr = PCSetFromOptions(pc);CHKERRQ(ierr); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */ 5533aeaf226SBarry Smith mglevels = mg->levels; 5543aeaf226SBarry Smith } 5553aeaf226SBarry Smith } 55698e04984SBarry Smith ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr); 5573aeaf226SBarry Smith 558f3fbd535SBarry Smith 559f3fbd535SBarry Smith /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */ 560f3fbd535SBarry Smith /* so use those from global PC */ 561f3fbd535SBarry Smith /* Is this what we always want? What if user wants to keep old one? */ 562f3fbd535SBarry Smith ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,PETSC_NULL,&opsset);CHKERRQ(ierr); 56398e04984SBarry Smith if (opsset) { 56498e04984SBarry Smith Mat mmat; 56598e04984SBarry Smith ierr = KSPGetOperators(mglevels[n-1]->smoothd,PETSC_NULL,&mmat,PETSC_NULL);CHKERRQ(ierr); 56698e04984SBarry Smith if (mmat == pc->pmat) opsset = PETSC_FALSE; 56798e04984SBarry Smith } 56898e04984SBarry Smith if (!opsset) { 569f3fbd535SBarry 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); 570f3fbd535SBarry Smith ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr); 571f3fbd535SBarry Smith } 572f3fbd535SBarry Smith 57301bc414fSDmitry Karpeev /* Skipping this for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs? */ 574ce4cda84SJed Brown if (pc->dm && mg->galerkin != 2 && !pc->setupcalled) { 5752d2b81a6SBarry Smith /* construct the interpolation from the DMs */ 5762e499ae9SBarry Smith Mat p; 57773250ac0SBarry Smith Vec rscale; 578218a07d4SBarry Smith ierr = PetscMalloc(n*sizeof(DM),&dms);CHKERRQ(ierr); 579218a07d4SBarry Smith dms[n-1] = pc->dm; 580218a07d4SBarry Smith for (i=n-2; i>-1; i--) { 581942e3340SBarry Smith DMKSP kdm; 5822ee06e3bSJed Brown ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr); 5833c0c59f3SBarry Smith ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr); 5843c0c59f3SBarry Smith if (mg->galerkin) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);} 585942e3340SBarry Smith ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr); 586d1a3e677SJed Brown /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take 587d1a3e677SJed Brown * a bitwise OR of computing the matrix, RHS, and initial iterate. */ 5885358d0d4SBarry Smith kdm->ops->computerhs = PETSC_NULL; 589fe86f630SJed Brown kdm->rhsctx = PETSC_NULL; 59024c3aa18SBarry Smith if (!mglevels[i+1]->interpolate) { 591e727c939SJed Brown ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr); 5922d2b81a6SBarry Smith ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr); 59300ac8be1SBarry Smith if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);} 59473250ac0SBarry Smith ierr = VecDestroy(&rscale);CHKERRQ(ierr); 5956bf464f9SBarry Smith ierr = MatDestroy(&p);CHKERRQ(ierr); 596218a07d4SBarry Smith } 59724c3aa18SBarry Smith } 5982d2b81a6SBarry Smith 599218a07d4SBarry Smith for (i=n-2; i>-1; i--) { 6006bf464f9SBarry Smith ierr = DMDestroy(&dms[i]);CHKERRQ(ierr); 601218a07d4SBarry Smith } 6022d2b81a6SBarry Smith ierr = PetscFree(dms);CHKERRQ(ierr); 603ce4cda84SJed Brown } 604cab2e9ccSBarry Smith 605ce4cda84SJed Brown if (pc->dm && !pc->setupcalled) { 606ce4cda84SJed Brown /* finest smoother also gets DM but it is not active, independent of whether galerkin==2 */ 607cab2e9ccSBarry Smith ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr); 608cab2e9ccSBarry Smith ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr); 609218a07d4SBarry Smith } 610218a07d4SBarry Smith 611c91913d3SJed Brown if (mg->galerkin == 1) { 612f3fbd535SBarry Smith Mat B; 613f3fbd535SBarry Smith /* currently only handle case where mat and pmat are the same on coarser levels */ 614f3fbd535SBarry Smith ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr); 615f3fbd535SBarry Smith if (!pc->setupcalled) { 616f3fbd535SBarry Smith for (i=n-2; i>-1; i--) { 617f3fbd535SBarry Smith ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 618f3fbd535SBarry Smith ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 619f3fbd535SBarry Smith if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);} 620f3fbd535SBarry Smith dB = B; 621f3fbd535SBarry Smith } 622cd9507b2SBarry Smith if (n > 1) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);} 623f3fbd535SBarry Smith } else { 624f3fbd535SBarry Smith for (i=n-2; i>-1; i--) { 625f3fbd535SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr); 626f3fbd535SBarry Smith ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 627f3fbd535SBarry Smith ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 628f3fbd535SBarry Smith dB = B; 629f3fbd535SBarry Smith } 630f3fbd535SBarry Smith } 631ce4cda84SJed Brown } else if (!mg->galerkin && pc->dm && pc->dm->x) { 632cab2e9ccSBarry Smith /* need to restrict Jacobian location to coarser meshes for evaluation */ 633cab2e9ccSBarry Smith for (i=n-2; i>-1; i--) { 634c88c5224SJed Brown Mat R; 635c88c5224SJed Brown Vec rscale; 636cab2e9ccSBarry Smith if (!mglevels[i]->smoothd->dm->x) { 637cab2e9ccSBarry Smith Vec *vecs; 638cab2e9ccSBarry Smith ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vecs,0,PETSC_NULL);CHKERRQ(ierr); 639*2fa5cd67SKarl Rupp 640cab2e9ccSBarry Smith mglevels[i]->smoothd->dm->x = vecs[0]; 641*2fa5cd67SKarl Rupp 642cab2e9ccSBarry Smith ierr = PetscFree(vecs);CHKERRQ(ierr); 643cab2e9ccSBarry Smith } 644c88c5224SJed Brown ierr = PCMGGetRestriction(pc,i+1,&R);CHKERRQ(ierr); 645c88c5224SJed Brown ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr); 646c88c5224SJed Brown ierr = MatRestrict(R,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);CHKERRQ(ierr); 647c88c5224SJed Brown ierr = VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,rscale);CHKERRQ(ierr); 648cab2e9ccSBarry Smith } 649f3fbd535SBarry Smith } 650ccceb50cSJed Brown if (!mg->galerkin && pc->dm) { 651caa4e7f2SJed Brown for (i=n-2; i>=0; i--) { 652ccceb50cSJed Brown DM dmfine,dmcoarse; 653caa4e7f2SJed Brown Mat Restrict,Inject; 654caa4e7f2SJed Brown Vec rscale; 655ccceb50cSJed Brown ierr = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr); 656ccceb50cSJed Brown ierr = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr); 657caa4e7f2SJed Brown ierr = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr); 658caa4e7f2SJed Brown ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr); 659caa4e7f2SJed Brown Inject = PETSC_NULL; /* Callback should create it if it needs Injection */ 660caa4e7f2SJed Brown ierr = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr); 661caa4e7f2SJed Brown } 662caa4e7f2SJed Brown } 663f3fbd535SBarry Smith 664f3fbd535SBarry Smith if (!pc->setupcalled) { 665f3fbd535SBarry Smith for (i=0; i<n; i++) { 666f3fbd535SBarry Smith ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr); 667f3fbd535SBarry Smith } 668f3fbd535SBarry Smith for (i=1; i<n; i++) { 669f3fbd535SBarry Smith if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) { 670f3fbd535SBarry Smith ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr); 671f3fbd535SBarry Smith } 672f3fbd535SBarry Smith } 673f3fbd535SBarry Smith for (i=1; i<n; i++) { 674c88c5224SJed Brown ierr = PCMGGetInterpolation(pc,i,&mglevels[i]->interpolate);CHKERRQ(ierr); 675c88c5224SJed Brown ierr = PCMGGetRestriction(pc,i,&mglevels[i]->restrct);CHKERRQ(ierr); 676f3fbd535SBarry Smith } 677f3fbd535SBarry Smith for (i=0; i<n-1; i++) { 678f3fbd535SBarry Smith if (!mglevels[i]->b) { 679f3fbd535SBarry Smith Vec *vec; 680f3fbd535SBarry Smith ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr); 681f3fbd535SBarry Smith ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr); 6826bf464f9SBarry Smith ierr = VecDestroy(vec);CHKERRQ(ierr); 683f3fbd535SBarry Smith ierr = PetscFree(vec);CHKERRQ(ierr); 684f3fbd535SBarry Smith } 685f3fbd535SBarry Smith if (!mglevels[i]->r && i) { 686f3fbd535SBarry Smith ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 687f3fbd535SBarry Smith ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr); 6886bf464f9SBarry Smith ierr = VecDestroy(&tvec);CHKERRQ(ierr); 689f3fbd535SBarry Smith } 690f3fbd535SBarry Smith if (!mglevels[i]->x) { 691f3fbd535SBarry Smith ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 692f3fbd535SBarry Smith ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr); 6936bf464f9SBarry Smith ierr = VecDestroy(&tvec);CHKERRQ(ierr); 694f3fbd535SBarry Smith } 695f3fbd535SBarry Smith } 696f3fbd535SBarry Smith if (n != 1 && !mglevels[n-1]->r) { 697f3fbd535SBarry Smith /* PCMGSetR() on the finest level if user did not supply it */ 698f3fbd535SBarry Smith Vec *vec; 699f3fbd535SBarry Smith ierr = KSPGetVecs(mglevels[n-1]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr); 700f3fbd535SBarry Smith ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr); 7016bf464f9SBarry Smith ierr = VecDestroy(vec);CHKERRQ(ierr); 702f3fbd535SBarry Smith ierr = PetscFree(vec);CHKERRQ(ierr); 703f3fbd535SBarry Smith } 704f3fbd535SBarry Smith } 705f3fbd535SBarry Smith 70698e04984SBarry Smith if (pc->dm) { 70798e04984SBarry Smith /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */ 70898e04984SBarry Smith for (i=0; i<n-1; i++) { 70998e04984SBarry Smith if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX; 71098e04984SBarry Smith } 71198e04984SBarry Smith } 712f3fbd535SBarry Smith 713f3fbd535SBarry Smith for (i=1; i<n; i++) { 714f3fbd535SBarry Smith if (mglevels[i]->smoothu == mglevels[i]->smoothd) { 715f3fbd535SBarry Smith /* if doing only down then initial guess is zero */ 716f3fbd535SBarry Smith ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 717f3fbd535SBarry Smith } 71863e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 719f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr); 72063e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 721d42688cbSBarry Smith if (!mglevels[i]->residual) { 722d42688cbSBarry Smith Mat mat; 723d42688cbSBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr); 724d42688cbSBarry Smith ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr); 725d42688cbSBarry Smith } 726f3fbd535SBarry Smith } 727f3fbd535SBarry Smith for (i=1; i<n; i++) { 728f3fbd535SBarry Smith if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) { 729f3fbd535SBarry Smith Mat downmat,downpmat; 730f3fbd535SBarry Smith MatStructure matflag; 731f3fbd535SBarry Smith 732f3fbd535SBarry Smith /* check if operators have been set for up, if not use down operators to set them */ 733f3fbd535SBarry Smith ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,PETSC_NULL);CHKERRQ(ierr); 734f3fbd535SBarry Smith if (!opsset) { 735f3fbd535SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr); 736f3fbd535SBarry Smith ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr); 737f3fbd535SBarry Smith } 738f3fbd535SBarry Smith 739f3fbd535SBarry Smith ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr); 74063e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 741f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr); 74263e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 743f3fbd535SBarry Smith } 744f3fbd535SBarry Smith } 745f3fbd535SBarry Smith 746f3fbd535SBarry Smith /* 747f3fbd535SBarry Smith If coarse solver is not direct method then DO NOT USE preonly 748f3fbd535SBarry Smith */ 749251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mglevels[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr); 750f3fbd535SBarry Smith if (preonly) { 751251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr); 752251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr); 753251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr); 754251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)cpc,PCSVD,&svd);CHKERRQ(ierr); 755fd1303e9SJungho Lee if (!lu && !redundant && !cholesky && !svd) { 756f3fbd535SBarry Smith ierr = KSPSetType(mglevels[0]->smoothd,KSPGMRES);CHKERRQ(ierr); 757f3fbd535SBarry Smith } 758f3fbd535SBarry Smith } 759f3fbd535SBarry Smith 760f3fbd535SBarry Smith if (!pc->setupcalled) { 761f3fbd535SBarry Smith ierr = KSPSetFromOptions(mglevels[0]->smoothd);CHKERRQ(ierr); 762f3fbd535SBarry Smith } 763f3fbd535SBarry Smith 76463e6d426SJed Brown if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 765f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr); 76663e6d426SJed Brown if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 767f3fbd535SBarry Smith 768f3fbd535SBarry Smith /* 769f3fbd535SBarry Smith Dump the interpolation/restriction matrices plus the 770e3c5b3baSBarry Smith Jacobian/stiffness on each level. This allows MATLAB users to 771f3fbd535SBarry Smith easily check if the Galerkin condition A_c = R A_f R^T is satisfied. 772f3fbd535SBarry Smith 773f3fbd535SBarry Smith Only support one or the other at the same time. 774f3fbd535SBarry Smith */ 775f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER) 776acfcf0e5SJed Brown ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,PETSC_NULL);CHKERRQ(ierr); 777*2fa5cd67SKarl Rupp if (dump) viewer = PETSC_VIEWER_SOCKET_(((PetscObject)pc)->comm); 778f3fbd535SBarry Smith dump = PETSC_FALSE; 779f3fbd535SBarry Smith #endif 780acfcf0e5SJed Brown ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,PETSC_NULL);CHKERRQ(ierr); 781*2fa5cd67SKarl Rupp if (dump) viewer = PETSC_VIEWER_BINARY_(((PetscObject)pc)->comm); 782f3fbd535SBarry Smith 783f3fbd535SBarry Smith if (viewer) { 784f3fbd535SBarry Smith for (i=1; i<n; i++) { 785f3fbd535SBarry Smith ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr); 786f3fbd535SBarry Smith } 787f3fbd535SBarry Smith for (i=0; i<n; i++) { 788f3fbd535SBarry Smith ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr); 789f3fbd535SBarry Smith ierr = MatView(pc->mat,viewer);CHKERRQ(ierr); 790f3fbd535SBarry Smith } 791f3fbd535SBarry Smith } 792f3fbd535SBarry Smith PetscFunctionReturn(0); 793f3fbd535SBarry Smith } 794f3fbd535SBarry Smith 795f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/ 796f3fbd535SBarry Smith 797f3fbd535SBarry Smith #undef __FUNCT__ 7989dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetLevels" 7994b9ad928SBarry Smith /*@ 80097177400SBarry Smith PCMGGetLevels - Gets the number of levels to use with MG. 8014b9ad928SBarry Smith 8024b9ad928SBarry Smith Not Collective 8034b9ad928SBarry Smith 8044b9ad928SBarry Smith Input Parameter: 8054b9ad928SBarry Smith . pc - the preconditioner context 8064b9ad928SBarry Smith 8074b9ad928SBarry Smith Output parameter: 8084b9ad928SBarry Smith . levels - the number of levels 8094b9ad928SBarry Smith 8104b9ad928SBarry Smith Level: advanced 8114b9ad928SBarry Smith 8124b9ad928SBarry Smith .keywords: MG, get, levels, multigrid 8134b9ad928SBarry Smith 81497177400SBarry Smith .seealso: PCMGSetLevels() 8154b9ad928SBarry Smith @*/ 8167087cfbeSBarry Smith PetscErrorCode PCMGGetLevels(PC pc,PetscInt *levels) 8174b9ad928SBarry Smith { 818f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 8194b9ad928SBarry Smith 8204b9ad928SBarry Smith PetscFunctionBegin; 8210700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 8224482741eSBarry Smith PetscValidIntPointer(levels,2); 823f3fbd535SBarry Smith *levels = mg->nlevels; 8244b9ad928SBarry Smith PetscFunctionReturn(0); 8254b9ad928SBarry Smith } 8264b9ad928SBarry Smith 8274b9ad928SBarry Smith #undef __FUNCT__ 8289dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetType" 8294b9ad928SBarry Smith /*@ 83097177400SBarry Smith PCMGSetType - Determines the form of multigrid to use: 8314b9ad928SBarry Smith multiplicative, additive, full, or the Kaskade algorithm. 8324b9ad928SBarry Smith 833ad4df100SBarry Smith Logically Collective on PC 8344b9ad928SBarry Smith 8354b9ad928SBarry Smith Input Parameters: 8364b9ad928SBarry Smith + pc - the preconditioner context 8379dcbbd2bSBarry Smith - form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE, 8389dcbbd2bSBarry Smith PC_MG_FULL, PC_MG_KASKADE 8394b9ad928SBarry Smith 8404b9ad928SBarry Smith Options Database Key: 8414b9ad928SBarry Smith . -pc_mg_type <form> - Sets <form>, one of multiplicative, 8424b9ad928SBarry Smith additive, full, kaskade 8434b9ad928SBarry Smith 8444b9ad928SBarry Smith Level: advanced 8454b9ad928SBarry Smith 8464b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid 8474b9ad928SBarry Smith 84897177400SBarry Smith .seealso: PCMGSetLevels() 8494b9ad928SBarry Smith @*/ 8507087cfbeSBarry Smith PetscErrorCode PCMGSetType(PC pc,PCMGType form) 8514b9ad928SBarry Smith { 852f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 8534b9ad928SBarry Smith 8544b9ad928SBarry Smith PetscFunctionBegin; 8550700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 856c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc,form,2); 85731567311SBarry Smith mg->am = form; 8589dcbbd2bSBarry Smith if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG; 8594b9ad928SBarry Smith else pc->ops->applyrichardson = 0; 8604b9ad928SBarry Smith PetscFunctionReturn(0); 8614b9ad928SBarry Smith } 8624b9ad928SBarry Smith 8634b9ad928SBarry Smith #undef __FUNCT__ 8640d353602SBarry Smith #define __FUNCT__ "PCMGSetCycleType" 8654b9ad928SBarry Smith /*@ 8660d353602SBarry Smith PCMGSetCycleType - Sets the type cycles to use. Use PCMGSetCycleTypeOnLevel() for more 8674b9ad928SBarry Smith complicated cycling. 8684b9ad928SBarry Smith 869ad4df100SBarry Smith Logically Collective on PC 8704b9ad928SBarry Smith 8714b9ad928SBarry Smith Input Parameters: 872c2be2410SBarry Smith + pc - the multigrid context 8730d353602SBarry Smith - PC_MG_CYCLE_V or PC_MG_CYCLE_W 8744b9ad928SBarry Smith 8754b9ad928SBarry Smith Options Database Key: 8760d353602SBarry Smith $ -pc_mg_cycle_type v or w 8774b9ad928SBarry Smith 8784b9ad928SBarry Smith Level: advanced 8794b9ad928SBarry Smith 8804b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 8814b9ad928SBarry Smith 8820d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel() 8834b9ad928SBarry Smith @*/ 8847087cfbeSBarry Smith PetscErrorCode PCMGSetCycleType(PC pc,PCMGCycleType n) 8854b9ad928SBarry Smith { 886f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 887f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 88879416396SBarry Smith PetscInt i,levels; 8894b9ad928SBarry Smith 8904b9ad928SBarry Smith PetscFunctionBegin; 8910700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 892e7e72b3dSBarry Smith if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 893c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 894f3fbd535SBarry Smith levels = mglevels[0]->levels; 8954b9ad928SBarry Smith 896*2fa5cd67SKarl Rupp for (i=0; i<levels; i++) mglevels[i]->cycles = n; 8974b9ad928SBarry Smith PetscFunctionReturn(0); 8984b9ad928SBarry Smith } 8994b9ad928SBarry Smith 9004b9ad928SBarry Smith #undef __FUNCT__ 9018cc2d5dfSBarry Smith #define __FUNCT__ "PCMGMultiplicativeSetCycles" 9028cc2d5dfSBarry Smith /*@ 9038cc2d5dfSBarry Smith PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step 9048cc2d5dfSBarry Smith of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used 9058cc2d5dfSBarry Smith 906ad4df100SBarry Smith Logically Collective on PC 9078cc2d5dfSBarry Smith 9088cc2d5dfSBarry Smith Input Parameters: 9098cc2d5dfSBarry Smith + pc - the multigrid context 9108cc2d5dfSBarry Smith - n - number of cycles (default is 1) 9118cc2d5dfSBarry Smith 9128cc2d5dfSBarry Smith Options Database Key: 9138cc2d5dfSBarry Smith $ -pc_mg_multiplicative_cycles n 9148cc2d5dfSBarry Smith 9158cc2d5dfSBarry Smith Level: advanced 9168cc2d5dfSBarry Smith 9178cc2d5dfSBarry Smith Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType() 9188cc2d5dfSBarry Smith 9198cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 9208cc2d5dfSBarry Smith 9218cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType() 9228cc2d5dfSBarry Smith @*/ 9237087cfbeSBarry Smith PetscErrorCode PCMGMultiplicativeSetCycles(PC pc,PetscInt n) 9248cc2d5dfSBarry Smith { 925f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 926f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 9278cc2d5dfSBarry Smith PetscInt i,levels; 9288cc2d5dfSBarry Smith 9298cc2d5dfSBarry Smith PetscFunctionBegin; 9300700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 931e7e72b3dSBarry Smith if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 932c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 933f3fbd535SBarry Smith levels = mglevels[0]->levels; 9348cc2d5dfSBarry Smith 935*2fa5cd67SKarl Rupp for (i=0; i<levels; i++) mg->cyclesperpcapply = n; 9368cc2d5dfSBarry Smith PetscFunctionReturn(0); 9378cc2d5dfSBarry Smith } 9388cc2d5dfSBarry Smith 9398cc2d5dfSBarry Smith #undef __FUNCT__ 9409dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetGalerkin" 941c2be2410SBarry Smith /*@ 94297177400SBarry Smith PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the 943c2be2410SBarry Smith finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t 944c2be2410SBarry Smith 945ad4df100SBarry Smith Logically Collective on PC 946c2be2410SBarry Smith 947c2be2410SBarry Smith Input Parameters: 948c91913d3SJed Brown + pc - the multigrid context 949c91913d3SJed Brown - use - PETSC_TRUE to use the Galerkin process to compute coarse-level operators 950c2be2410SBarry Smith 951c2be2410SBarry Smith Options Database Key: 952c2be2410SBarry Smith $ -pc_mg_galerkin 953c2be2410SBarry Smith 954c2be2410SBarry Smith Level: intermediate 955c2be2410SBarry Smith 956c2be2410SBarry Smith .keywords: MG, set, Galerkin 957c2be2410SBarry Smith 9583fc8bf9cSBarry Smith .seealso: PCMGGetGalerkin() 9593fc8bf9cSBarry Smith 960c2be2410SBarry Smith @*/ 961c91913d3SJed Brown PetscErrorCode PCMGSetGalerkin(PC pc,PetscBool use) 962c2be2410SBarry Smith { 963f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 964c2be2410SBarry Smith 965c2be2410SBarry Smith PetscFunctionBegin; 9660700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 967789726d7SBarry Smith mg->galerkin = use ? 1 : 0; 968c2be2410SBarry Smith PetscFunctionReturn(0); 969c2be2410SBarry Smith } 970c2be2410SBarry Smith 971c2be2410SBarry Smith #undef __FUNCT__ 9723fc8bf9cSBarry Smith #define __FUNCT__ "PCMGGetGalerkin" 9733fc8bf9cSBarry Smith /*@ 9743fc8bf9cSBarry Smith PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. 9753fc8bf9cSBarry Smith A_i-1 = r_i * A_i * r_i^t 9763fc8bf9cSBarry Smith 9773fc8bf9cSBarry Smith Not Collective 9783fc8bf9cSBarry Smith 9793fc8bf9cSBarry Smith Input Parameter: 9803fc8bf9cSBarry Smith . pc - the multigrid context 9813fc8bf9cSBarry Smith 9823fc8bf9cSBarry Smith Output Parameter: 9833fc8bf9cSBarry Smith . gelerkin - PETSC_TRUE or PETSC_FALSE 9843fc8bf9cSBarry Smith 9853fc8bf9cSBarry Smith Options Database Key: 9863fc8bf9cSBarry Smith $ -pc_mg_galerkin 9873fc8bf9cSBarry Smith 9883fc8bf9cSBarry Smith Level: intermediate 9893fc8bf9cSBarry Smith 9903fc8bf9cSBarry Smith .keywords: MG, set, Galerkin 9913fc8bf9cSBarry Smith 9923fc8bf9cSBarry Smith .seealso: PCMGSetGalerkin() 9933fc8bf9cSBarry Smith 9943fc8bf9cSBarry Smith @*/ 9957087cfbeSBarry Smith PetscErrorCode PCMGGetGalerkin(PC pc,PetscBool *galerkin) 9963fc8bf9cSBarry Smith { 997f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 9983fc8bf9cSBarry Smith 9993fc8bf9cSBarry Smith PetscFunctionBegin; 10000700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 100113fdb3acSJose Roman *galerkin = (PetscBool)mg->galerkin; 10023fc8bf9cSBarry Smith PetscFunctionReturn(0); 10033fc8bf9cSBarry Smith } 10043fc8bf9cSBarry Smith 10053fc8bf9cSBarry Smith #undef __FUNCT__ 10069dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothDown" 10074b9ad928SBarry Smith /*@ 100897177400SBarry Smith PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to 100997177400SBarry Smith use on all levels. Use PCMGGetSmootherDown() to set different 10104b9ad928SBarry Smith pre-smoothing steps on different levels. 10114b9ad928SBarry Smith 1012ad4df100SBarry Smith Logically Collective on PC 10134b9ad928SBarry Smith 10144b9ad928SBarry Smith Input Parameters: 10154b9ad928SBarry Smith + mg - the multigrid context 10164b9ad928SBarry Smith - n - the number of smoothing steps 10174b9ad928SBarry Smith 10184b9ad928SBarry Smith Options Database Key: 10194b9ad928SBarry Smith . -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps 10204b9ad928SBarry Smith 10214b9ad928SBarry Smith Level: advanced 10224b9ad928SBarry Smith 10234b9ad928SBarry Smith .keywords: MG, smooth, down, pre-smoothing, steps, multigrid 10244b9ad928SBarry Smith 102597177400SBarry Smith .seealso: PCMGSetNumberSmoothUp() 10264b9ad928SBarry Smith @*/ 10277087cfbeSBarry Smith PetscErrorCode PCMGSetNumberSmoothDown(PC pc,PetscInt n) 10284b9ad928SBarry Smith { 1029f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1030f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 10316849ba73SBarry Smith PetscErrorCode ierr; 103279416396SBarry Smith PetscInt i,levels; 10334b9ad928SBarry Smith 10344b9ad928SBarry Smith PetscFunctionBegin; 10350700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1036e7e72b3dSBarry Smith if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 1037c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 1038f3fbd535SBarry Smith levels = mglevels[0]->levels; 10394b9ad928SBarry Smith 1040b05257ddSBarry Smith for (i=1; i<levels; i++) { 10414b9ad928SBarry Smith /* make sure smoother up and down are different */ 104297177400SBarry Smith ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 1043f3fbd535SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 1044*2fa5cd67SKarl Rupp 104531567311SBarry Smith mg->default_smoothd = n; 10464b9ad928SBarry Smith } 10474b9ad928SBarry Smith PetscFunctionReturn(0); 10484b9ad928SBarry Smith } 10494b9ad928SBarry Smith 10504b9ad928SBarry Smith #undef __FUNCT__ 10519dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothUp" 10524b9ad928SBarry Smith /*@ 105397177400SBarry Smith PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use 105497177400SBarry Smith on all levels. Use PCMGGetSmootherUp() to set different numbers of 10554b9ad928SBarry Smith post-smoothing steps on different levels. 10564b9ad928SBarry Smith 1057ad4df100SBarry Smith Logically Collective on PC 10584b9ad928SBarry Smith 10594b9ad928SBarry Smith Input Parameters: 10604b9ad928SBarry Smith + mg - the multigrid context 10614b9ad928SBarry Smith - n - the number of smoothing steps 10624b9ad928SBarry Smith 10634b9ad928SBarry Smith Options Database Key: 10644b9ad928SBarry Smith . -pc_mg_smoothup <n> - Sets number of post-smoothing steps 10654b9ad928SBarry Smith 10664b9ad928SBarry Smith Level: advanced 10674b9ad928SBarry Smith 10684b9ad928SBarry Smith Note: this does not set a value on the coarsest grid, since we assume that 1069a8c7a070SBarry Smith there is no separate smooth up on the coarsest grid. 10704b9ad928SBarry Smith 10714b9ad928SBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid 10724b9ad928SBarry Smith 107397177400SBarry Smith .seealso: PCMGSetNumberSmoothDown() 10744b9ad928SBarry Smith @*/ 10757087cfbeSBarry Smith PetscErrorCode PCMGSetNumberSmoothUp(PC pc,PetscInt n) 10764b9ad928SBarry Smith { 1077f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1078f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 10796849ba73SBarry Smith PetscErrorCode ierr; 108079416396SBarry Smith PetscInt i,levels; 10814b9ad928SBarry Smith 10824b9ad928SBarry Smith PetscFunctionBegin; 10830700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1084e7e72b3dSBarry Smith if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 1085c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 1086f3fbd535SBarry Smith levels = mglevels[0]->levels; 10874b9ad928SBarry Smith 10884b9ad928SBarry Smith for (i=1; i<levels; i++) { 10894b9ad928SBarry Smith /* make sure smoother up and down are different */ 109097177400SBarry Smith ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 1091f3fbd535SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 1092*2fa5cd67SKarl Rupp 109331567311SBarry Smith mg->default_smoothu = n; 10944b9ad928SBarry Smith } 10954b9ad928SBarry Smith PetscFunctionReturn(0); 10964b9ad928SBarry Smith } 10974b9ad928SBarry Smith 10984b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/ 10994b9ad928SBarry Smith 11003b09bd56SBarry Smith /*MC 1101ccb205f8SBarry Smith PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional 11023b09bd56SBarry Smith information about the coarser grid matrices and restriction/interpolation operators. 11033b09bd56SBarry Smith 11043b09bd56SBarry Smith Options Database Keys: 11053b09bd56SBarry Smith + -pc_mg_levels <nlevels> - number of levels including finest 11060d353602SBarry Smith . -pc_mg_cycles v or w 110779416396SBarry Smith . -pc_mg_smoothup <n> - number of smoothing steps after interpolation 11083b09bd56SBarry Smith . -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator 11098c1c2452SJed Brown . -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default 11103b09bd56SBarry Smith . -pc_mg_log - log information about time spent on each level of the solver 11113b09bd56SBarry Smith . -pc_mg_monitor - print information on the multigrid convergence 11128cf6037eSBarry Smith . -pc_mg_galerkin - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R' 11138cf6037eSBarry Smith . -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1) 11148cf6037eSBarry Smith . -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 1115e3c5b3baSBarry Smith to the Socket viewer for reading from MATLAB. 11168cf6037eSBarry Smith - -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices 11178cf6037eSBarry Smith to the binary output file called binaryoutput 11183b09bd56SBarry Smith 111924c3aa18SBarry Smith Notes: By default this uses GMRES on the fine grid smoother so this should be used with KSPFGMRES or the smoother changed to not use GMRES 11203b09bd56SBarry Smith 11218cf6037eSBarry Smith When run with a single level the smoother options are used on that level NOT the coarse grid solver options 11228cf6037eSBarry Smith 11233b09bd56SBarry Smith Level: intermediate 11243b09bd56SBarry Smith 11258f87f92bSBarry Smith Concepts: multigrid/multilevel 11263b09bd56SBarry Smith 11278cf6037eSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE 11280d353602SBarry Smith PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(), 112997177400SBarry Smith PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 113097177400SBarry Smith PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 11310d353602SBarry Smith PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 11323b09bd56SBarry Smith M*/ 11333b09bd56SBarry Smith 11344b9ad928SBarry Smith EXTERN_C_BEGIN 11354b9ad928SBarry Smith #undef __FUNCT__ 11364b9ad928SBarry Smith #define __FUNCT__ "PCCreate_MG" 11377087cfbeSBarry Smith PetscErrorCode PCCreate_MG(PC pc) 11384b9ad928SBarry Smith { 1139f3fbd535SBarry Smith PC_MG *mg; 1140f3fbd535SBarry Smith PetscErrorCode ierr; 1141f3fbd535SBarry Smith 11424b9ad928SBarry Smith PetscFunctionBegin; 1143f3fbd535SBarry Smith ierr = PetscNewLog(pc,PC_MG,&mg);CHKERRQ(ierr); 1144f3fbd535SBarry Smith pc->data = (void*)mg; 1145f3fbd535SBarry Smith mg->nlevels = -1; 1146f3fbd535SBarry Smith 11474b9ad928SBarry Smith pc->ops->apply = PCApply_MG; 11484b9ad928SBarry Smith pc->ops->setup = PCSetUp_MG; 1149a06653b4SBarry Smith pc->ops->reset = PCReset_MG; 11504b9ad928SBarry Smith pc->ops->destroy = PCDestroy_MG; 11514b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_MG; 11524b9ad928SBarry Smith pc->ops->view = PCView_MG; 11534b9ad928SBarry Smith PetscFunctionReturn(0); 11544b9ad928SBarry Smith } 11554b9ad928SBarry Smith EXTERN_C_END 1156