1dba47a55SKris Buschelman 24b9ad928SBarry Smith /* 34b9ad928SBarry Smith Defines the multigrid preconditioner interface. 44b9ad928SBarry Smith */ 5b4876bcbSBarry Smith #include <../src/ksp/pc/impls/mg/mgimpl.h> /*I "petscksp.h" I*/ 61e25c274SJed Brown #include <petscdm.h> 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); 862fa5cd67SKarl Rupp } else if (abstol) mg->ttol = abstol; 872fa5cd67SKarl 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 1520298fd71SBarry Smith on smaller sets of processors. Use 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; 168ce94432eSBarry Smith MPI_Comm 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); 179ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 180548767e0SJed Brown if (mg->nlevels == levels) PetscFunctionReturn(0); 1813aeaf226SBarry Smith if (mglevels) { 1823aeaf226SBarry Smith /* changing the number of levels so free up the previous stuff */ 1833aeaf226SBarry Smith ierr = PCReset_MG(pc);CHKERRQ(ierr); 1843aeaf226SBarry Smith n = mglevels[0]->levels; 1853aeaf226SBarry Smith for (i=0; i<n; i++) { 1863aeaf226SBarry Smith if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 1873aeaf226SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr); 1883aeaf226SBarry Smith } 1893aeaf226SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr); 1903aeaf226SBarry Smith ierr = PetscFree(mglevels[i]);CHKERRQ(ierr); 1913aeaf226SBarry Smith } 1923aeaf226SBarry Smith ierr = PetscFree(mg->levels);CHKERRQ(ierr); 1933aeaf226SBarry Smith } 194f3fbd535SBarry Smith 195f3fbd535SBarry Smith mg->nlevels = levels; 196f3fbd535SBarry Smith 197f3fbd535SBarry Smith ierr = PetscMalloc(levels*sizeof(PC_MG*),&mglevels);CHKERRQ(ierr); 198*3bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr); 199f3fbd535SBarry Smith 200f3fbd535SBarry Smith ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 201f3fbd535SBarry Smith 202a9db87a2SMatthew G Knepley mg->stageApply = 0; 203f3fbd535SBarry Smith for (i=0; i<levels; i++) { 204f3fbd535SBarry Smith ierr = PetscNewLog(pc,PC_MG_Levels,&mglevels[i]);CHKERRQ(ierr); 2052fa5cd67SKarl Rupp 206f3fbd535SBarry Smith mglevels[i]->level = i; 207f3fbd535SBarry Smith mglevels[i]->levels = levels; 208f3fbd535SBarry Smith mglevels[i]->cycles = PC_MG_CYCLE_V; 209336babb1SJed Brown mg->default_smoothu = 2; 210336babb1SJed Brown mg->default_smoothd = 2; 21163e6d426SJed Brown mglevels[i]->eventsmoothsetup = 0; 21263e6d426SJed Brown mglevels[i]->eventsmoothsolve = 0; 21363e6d426SJed Brown mglevels[i]->eventresidual = 0; 21463e6d426SJed Brown mglevels[i]->eventinterprestrict = 0; 215f3fbd535SBarry Smith 216f3fbd535SBarry Smith if (comms) comm = comms[i]; 217f3fbd535SBarry Smith ierr = KSPCreate(comm,&mglevels[i]->smoothd);CHKERRQ(ierr); 218336babb1SJed Brown ierr = KSPSetType(mglevels[i]->smoothd,KSPCHEBYSHEV);CHKERRQ(ierr); 2190298fd71SBarry Smith ierr = KSPSetConvergenceTest(mglevels[i]->smoothd,KSPSkipConverged,NULL,NULL);CHKERRQ(ierr); 220336babb1SJed Brown ierr = KSPSetNormType(mglevels[i]->smoothd,KSP_NORM_NONE);CHKERRQ(ierr); 221336babb1SJed Brown ierr = KSPGetPC(mglevels[i]->smoothd,&ipc);CHKERRQ(ierr); 222336babb1SJed Brown ierr = PCSetType(ipc,PCSOR);CHKERRQ(ierr); 223f3fbd535SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr); 224336babb1SJed Brown ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, i ? mg->default_smoothd : 1);CHKERRQ(ierr); 225f3fbd535SBarry Smith ierr = KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);CHKERRQ(ierr); 226f3fbd535SBarry Smith 227f3fbd535SBarry Smith /* do special stuff for coarse grid */ 228f3fbd535SBarry Smith if (!i && levels > 1) { 229f3fbd535SBarry Smith ierr = KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");CHKERRQ(ierr); 230f3fbd535SBarry Smith 2319dbfc187SHong Zhang /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */ 232f3fbd535SBarry Smith ierr = KSPSetType(mglevels[0]->smoothd,KSPPREONLY);CHKERRQ(ierr); 233f3fbd535SBarry Smith ierr = KSPGetPC(mglevels[0]->smoothd,&ipc);CHKERRQ(ierr); 234f3fbd535SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 235f3fbd535SBarry Smith if (size > 1) { 2369dbfc187SHong Zhang KSP innerksp; 2379dbfc187SHong Zhang PC innerpc; 238f3fbd535SBarry Smith ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr); 2399dbfc187SHong Zhang ierr = PCRedundantGetKSP(ipc,&innerksp);CHKERRQ(ierr); 2409dbfc187SHong Zhang ierr = KSPGetPC(innerksp,&innerpc);CHKERRQ(ierr); 2419dbfc187SHong Zhang ierr = PCFactorSetShiftType(innerpc,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr); 242f3fbd535SBarry Smith } else { 243f3fbd535SBarry Smith ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr); 2449dbfc187SHong Zhang ierr = PCFactorSetShiftType(ipc,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr); 245f3fbd535SBarry Smith } 246f3fbd535SBarry Smith } else { 247f3fbd535SBarry Smith char tprefix[128]; 248f3fbd535SBarry Smith sprintf(tprefix,"mg_levels_%d_",(int)i); 249f3fbd535SBarry Smith ierr = KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);CHKERRQ(ierr); 250f3fbd535SBarry Smith } 251*3bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[i]->smoothd);CHKERRQ(ierr); 2522fa5cd67SKarl Rupp 253f3fbd535SBarry Smith mglevels[i]->smoothu = mglevels[i]->smoothd; 25431567311SBarry Smith mg->rtol = 0.0; 25531567311SBarry Smith mg->abstol = 0.0; 25631567311SBarry Smith mg->dtol = 0.0; 25731567311SBarry Smith mg->ttol = 0.0; 25831567311SBarry Smith mg->cyclesperpcapply = 1; 259f3fbd535SBarry Smith } 26031567311SBarry Smith mg->am = PC_MG_MULTIPLICATIVE; 261f3fbd535SBarry Smith mg->levels = mglevels; 2624b9ad928SBarry Smith pc->ops->applyrichardson = PCApplyRichardson_MG; 2634b9ad928SBarry Smith PetscFunctionReturn(0); 2644b9ad928SBarry Smith } 2654b9ad928SBarry Smith 266c07bf074SBarry Smith 267c07bf074SBarry Smith #undef __FUNCT__ 268c07bf074SBarry Smith #define __FUNCT__ "PCDestroy_MG" 269c07bf074SBarry Smith PetscErrorCode PCDestroy_MG(PC pc) 270c07bf074SBarry Smith { 271c07bf074SBarry Smith PetscErrorCode ierr; 272a06653b4SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 273a06653b4SBarry Smith PC_MG_Levels **mglevels = mg->levels; 274a06653b4SBarry Smith PetscInt i,n; 275c07bf074SBarry Smith 276c07bf074SBarry Smith PetscFunctionBegin; 277a06653b4SBarry Smith ierr = PCReset_MG(pc);CHKERRQ(ierr); 278a06653b4SBarry Smith if (mglevels) { 279a06653b4SBarry Smith n = mglevels[0]->levels; 280a06653b4SBarry Smith for (i=0; i<n; i++) { 281a06653b4SBarry Smith if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 2826bf464f9SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr); 283a06653b4SBarry Smith } 2846bf464f9SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr); 285a06653b4SBarry Smith ierr = PetscFree(mglevels[i]);CHKERRQ(ierr); 286a06653b4SBarry Smith } 287a06653b4SBarry Smith ierr = PetscFree(mg->levels);CHKERRQ(ierr); 288a06653b4SBarry Smith } 289c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 290f3fbd535SBarry Smith PetscFunctionReturn(0); 291f3fbd535SBarry Smith } 292f3fbd535SBarry Smith 293f3fbd535SBarry Smith 294f3fbd535SBarry Smith 29509573ac7SBarry Smith extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**); 29609573ac7SBarry Smith extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**); 29709573ac7SBarry Smith extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**); 298f3fbd535SBarry Smith 299f3fbd535SBarry Smith /* 300f3fbd535SBarry Smith PCApply_MG - Runs either an additive, multiplicative, Kaskadic 301f3fbd535SBarry Smith or full cycle of multigrid. 302f3fbd535SBarry Smith 303f3fbd535SBarry Smith Note: 304f3fbd535SBarry Smith A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle(). 305f3fbd535SBarry Smith */ 306f3fbd535SBarry Smith #undef __FUNCT__ 307f3fbd535SBarry Smith #define __FUNCT__ "PCApply_MG" 308f3fbd535SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x) 309f3fbd535SBarry Smith { 310f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 311f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 312f3fbd535SBarry Smith PetscErrorCode ierr; 313f3fbd535SBarry Smith PetscInt levels = mglevels[0]->levels,i; 314f3fbd535SBarry Smith 315f3fbd535SBarry Smith PetscFunctionBegin; 316a9db87a2SMatthew G Knepley if (mg->stageApply) {ierr = PetscLogStagePush(mg->stageApply);CHKERRQ(ierr);} 317e1d8e5deSBarry Smith /* When the DM is supplying the matrix then it will not exist until here */ 318298cc208SBarry Smith for (i=0; i<levels; i++) { 319e1d8e5deSBarry Smith if (!mglevels[i]->A) { 3200298fd71SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL,NULL);CHKERRQ(ierr); 321298cc208SBarry Smith ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr); 322e1d8e5deSBarry Smith } 323e1d8e5deSBarry Smith } 324e1d8e5deSBarry Smith 325f3fbd535SBarry Smith mglevels[levels-1]->b = b; 326f3fbd535SBarry Smith mglevels[levels-1]->x = x; 32731567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 328f3fbd535SBarry Smith ierr = VecSet(x,0.0);CHKERRQ(ierr); 32931567311SBarry Smith for (i=0; i<mg->cyclesperpcapply; i++) { 3300298fd71SBarry Smith ierr = PCMGMCycle_Private(pc,mglevels+levels-1,NULL);CHKERRQ(ierr); 331f3fbd535SBarry Smith } 3322fa5cd67SKarl Rupp } else if (mg->am == PC_MG_ADDITIVE) { 33331567311SBarry Smith ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr); 3342fa5cd67SKarl Rupp } else if (mg->am == PC_MG_KASKADE) { 33531567311SBarry Smith ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr); 3362fa5cd67SKarl Rupp } else { 337f3fbd535SBarry Smith ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr); 338f3fbd535SBarry Smith } 339a9db87a2SMatthew G Knepley if (mg->stageApply) {ierr = PetscLogStagePop();CHKERRQ(ierr);} 340f3fbd535SBarry Smith PetscFunctionReturn(0); 341f3fbd535SBarry Smith } 342f3fbd535SBarry Smith 343f3fbd535SBarry Smith 344f3fbd535SBarry Smith #undef __FUNCT__ 345f3fbd535SBarry Smith #define __FUNCT__ "PCSetFromOptions_MG" 346f3fbd535SBarry Smith PetscErrorCode PCSetFromOptions_MG(PC pc) 347f3fbd535SBarry Smith { 348f3fbd535SBarry Smith PetscErrorCode ierr; 349f3fbd535SBarry Smith PetscInt m,levels = 1,cycles; 350c91913d3SJed Brown PetscBool flg,set; 351f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 352f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 353f3fbd535SBarry Smith PCMGType mgtype; 354f3fbd535SBarry Smith PCMGCycleType mgctype; 355f3fbd535SBarry Smith 356f3fbd535SBarry Smith PetscFunctionBegin; 357f3fbd535SBarry Smith ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr); 3583aeaf226SBarry Smith if (!mg->levels) { 359f3fbd535SBarry Smith ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr); 360eb3f98d2SBarry Smith if (!flg && pc->dm) { 361eb3f98d2SBarry Smith ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr); 362eb3f98d2SBarry Smith levels++; 3633aeaf226SBarry Smith mg->usedmfornumberoflevels = PETSC_TRUE; 364eb3f98d2SBarry Smith } 3650298fd71SBarry Smith ierr = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr); 366f3fbd535SBarry Smith } 3673aeaf226SBarry Smith mglevels = mg->levels; 3683aeaf226SBarry Smith 369f3fbd535SBarry Smith mgctype = (PCMGCycleType) mglevels[0]->cycles; 370f3fbd535SBarry Smith ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr); 371f3fbd535SBarry Smith if (flg) { 372f3fbd535SBarry Smith ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr); 3732fa5cd67SKarl Rupp } 374f3fbd535SBarry Smith flg = PETSC_FALSE; 375c91913d3SJed Brown ierr = PetscOptionsBool("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",flg,&flg,&set);CHKERRQ(ierr); 376c91913d3SJed Brown if (set) { 377c91913d3SJed Brown ierr = PCMGSetGalerkin(pc,flg);CHKERRQ(ierr); 378f3fbd535SBarry Smith } 379f3fbd535SBarry Smith ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr); 380f3fbd535SBarry Smith if (flg) { 381f3fbd535SBarry Smith ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr); 382f3fbd535SBarry Smith } 383f3fbd535SBarry Smith ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr); 384f3fbd535SBarry Smith if (flg) { 385f3fbd535SBarry Smith ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr); 386f3fbd535SBarry Smith } 38731567311SBarry Smith mgtype = mg->am; 388f3fbd535SBarry Smith ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr); 389f3fbd535SBarry Smith if (flg) { 390f3fbd535SBarry Smith ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr); 391f3fbd535SBarry Smith } 39231567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 39331567311SBarry Smith ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGSetLevels",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr); 394f3fbd535SBarry Smith if (flg) { 395f3fbd535SBarry Smith ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr); 396f3fbd535SBarry Smith } 397f3fbd535SBarry Smith } 398f3fbd535SBarry Smith flg = PETSC_FALSE; 3990298fd71SBarry Smith ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,NULL);CHKERRQ(ierr); 400f3fbd535SBarry Smith if (flg) { 401f3fbd535SBarry Smith PetscInt i; 402f3fbd535SBarry Smith char eventname[128]; 403ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 404f3fbd535SBarry Smith levels = mglevels[0]->levels; 405f3fbd535SBarry Smith for (i=0; i<levels; i++) { 406f3fbd535SBarry Smith sprintf(eventname,"MGSetup Level %d",(int)i); 40763e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr); 408f3fbd535SBarry Smith sprintf(eventname,"MGSmooth Level %d",(int)i); 40963e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr); 410f3fbd535SBarry Smith if (i) { 411f3fbd535SBarry Smith sprintf(eventname,"MGResid Level %d",(int)i); 41263e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr); 413f3fbd535SBarry Smith sprintf(eventname,"MGInterp Level %d",(int)i); 41463e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr); 415f3fbd535SBarry Smith } 416f3fbd535SBarry Smith } 417eec35531SMatthew G Knepley 418ec60ca73SSatish Balay #if defined(PETSC_USE_LOG) 419eec35531SMatthew G Knepley { 420eec35531SMatthew G Knepley const char *sname = "MG Apply"; 421eec35531SMatthew G Knepley PetscStageLog stageLog; 422eec35531SMatthew G Knepley PetscInt st; 423eec35531SMatthew G Knepley 424eec35531SMatthew G Knepley PetscFunctionBegin; 425eec35531SMatthew G Knepley ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr); 426eec35531SMatthew G Knepley for (st = 0; st < stageLog->numStages; ++st) { 427eec35531SMatthew G Knepley PetscBool same; 428eec35531SMatthew G Knepley 429eec35531SMatthew G Knepley ierr = PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);CHKERRQ(ierr); 4302fa5cd67SKarl Rupp if (same) mg->stageApply = st; 431eec35531SMatthew G Knepley } 432eec35531SMatthew G Knepley if (!mg->stageApply) { 433eec35531SMatthew G Knepley ierr = PetscLogStageRegister(sname, &mg->stageApply);CHKERRQ(ierr); 434eec35531SMatthew G Knepley } 435eec35531SMatthew G Knepley } 436ec60ca73SSatish Balay #endif 437f3fbd535SBarry Smith } 438f3fbd535SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 439f3fbd535SBarry Smith PetscFunctionReturn(0); 440f3fbd535SBarry Smith } 441f3fbd535SBarry Smith 4426a6fc655SJed Brown const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0}; 4436a6fc655SJed Brown const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0}; 444f3fbd535SBarry Smith 4459804daf3SBarry Smith #include <petscdraw.h> 446f3fbd535SBarry Smith #undef __FUNCT__ 447f3fbd535SBarry Smith #define __FUNCT__ "PCView_MG" 448f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer) 449f3fbd535SBarry Smith { 450f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 451f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 452f3fbd535SBarry Smith PetscErrorCode ierr; 453e3deeeafSJed Brown PetscInt levels = mglevels ? mglevels[0]->levels : 0,i; 454219b1045SBarry Smith PetscBool iascii,isbinary,isdraw; 455f3fbd535SBarry Smith 456f3fbd535SBarry Smith PetscFunctionBegin; 457251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 4585b0b0462SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 459219b1045SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 460f3fbd535SBarry Smith if (iascii) { 461e3deeeafSJed Brown const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown"; 462e3deeeafSJed Brown ierr = PetscViewerASCIIPrintf(viewer," MG: type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,cyclename);CHKERRQ(ierr); 46331567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 46431567311SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr); 465f3fbd535SBarry Smith } 466218a07d4SBarry Smith if (mg->galerkin) { 467f3fbd535SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 4684f66f45eSBarry Smith } else { 4694f66f45eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 470f3fbd535SBarry Smith } 471f3fbd535SBarry Smith for (i=0; i<levels; i++) { 472f3fbd535SBarry Smith if (!i) { 47389cce641SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);CHKERRQ(ierr); 474f3fbd535SBarry Smith } else { 47589cce641SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 476f3fbd535SBarry Smith } 477f3fbd535SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 478f3fbd535SBarry Smith ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 479f3fbd535SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 480f3fbd535SBarry Smith if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) { 481f3fbd535SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 482f3fbd535SBarry Smith } else if (i) { 48389cce641SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 484f3fbd535SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 485f3fbd535SBarry Smith ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 486f3fbd535SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 487f3fbd535SBarry Smith } 488f3fbd535SBarry Smith } 4895b0b0462SBarry Smith } else if (isbinary) { 4905b0b0462SBarry Smith for (i=levels-1; i>=0; i--) { 4915b0b0462SBarry Smith ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 4925b0b0462SBarry Smith if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) { 4935b0b0462SBarry Smith ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 4945b0b0462SBarry Smith } 4955b0b0462SBarry Smith } 496219b1045SBarry Smith } else if (isdraw) { 497219b1045SBarry Smith PetscDraw draw; 498b4375e8dSPeter Brune PetscReal x,w,y,bottom,th; 499219b1045SBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 500219b1045SBarry Smith ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr); 5010298fd71SBarry Smith ierr = PetscDrawStringGetSize(draw,NULL,&th);CHKERRQ(ierr); 502219b1045SBarry Smith bottom = y - th; 503219b1045SBarry Smith for (i=levels-1; i>=0; i--) { 504b4375e8dSPeter Brune if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) { 505219b1045SBarry Smith ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr); 506219b1045SBarry Smith ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 507219b1045SBarry Smith ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 508b4375e8dSPeter Brune } else { 509b4375e8dSPeter Brune w = 0.5*PetscMin(1.0-x,x); 510b4375e8dSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr); 511b4375e8dSPeter Brune ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 512b4375e8dSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 513b4375e8dSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr); 514b4375e8dSPeter Brune ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 515b4375e8dSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 516b4375e8dSPeter Brune } 5170298fd71SBarry Smith ierr = PetscDrawGetBoundingBox(draw,NULL,&bottom,NULL,NULL);CHKERRQ(ierr); 5181cd381d6SBarry Smith bottom -= th; 519219b1045SBarry Smith } 5205b0b0462SBarry Smith } 521f3fbd535SBarry Smith PetscFunctionReturn(0); 522f3fbd535SBarry Smith } 523f3fbd535SBarry Smith 524b45d2f2cSJed Brown #include <petsc-private/dmimpl.h> 525b45d2f2cSJed Brown #include <petsc-private/kspimpl.h> 526cab2e9ccSBarry Smith 527f3fbd535SBarry Smith /* 528f3fbd535SBarry Smith Calls setup for the KSP on each level 529f3fbd535SBarry Smith */ 530f3fbd535SBarry Smith #undef __FUNCT__ 531f3fbd535SBarry Smith #define __FUNCT__ "PCSetUp_MG" 532f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc) 533f3fbd535SBarry Smith { 534f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 535f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 536f3fbd535SBarry Smith PetscErrorCode ierr; 537f3fbd535SBarry Smith PetscInt i,n = mglevels[0]->levels; 53898e04984SBarry Smith PC cpc; 53971b23a65SMark F. Adams PetscBool preonly,lu,redundant,cholesky,svd,dump = PETSC_FALSE,opsset,use_amat; 540f3fbd535SBarry Smith Mat dA,dB; 541f3fbd535SBarry Smith MatStructure uflag; 542f3fbd535SBarry Smith Vec tvec; 543218a07d4SBarry Smith DM *dms; 544649052a6SBarry Smith PetscViewer viewer = 0; 545f3fbd535SBarry Smith 546f3fbd535SBarry Smith PetscFunctionBegin; 54701bc414fSDmitry Karpeev /* FIX: Move this to PCSetFromOptions_MG? */ 5483aeaf226SBarry Smith if (mg->usedmfornumberoflevels) { 5493aeaf226SBarry Smith PetscInt levels; 5503aeaf226SBarry Smith ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr); 5513aeaf226SBarry Smith levels++; 5523aeaf226SBarry Smith if (levels > n) { /* the problem is now being solved on a finer grid */ 5530298fd71SBarry Smith ierr = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr); 5543aeaf226SBarry Smith n = levels; 5553aeaf226SBarry Smith ierr = PCSetFromOptions(pc);CHKERRQ(ierr); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */ 5563aeaf226SBarry Smith mglevels = mg->levels; 5573aeaf226SBarry Smith } 5583aeaf226SBarry Smith } 55998e04984SBarry Smith ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr); 5603aeaf226SBarry Smith 561f3fbd535SBarry Smith 562f3fbd535SBarry Smith /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */ 563f3fbd535SBarry Smith /* so use those from global PC */ 564f3fbd535SBarry Smith /* Is this what we always want? What if user wants to keep old one? */ 5650298fd71SBarry Smith ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset);CHKERRQ(ierr); 56698e04984SBarry Smith if (opsset) { 56798e04984SBarry Smith Mat mmat; 5680298fd71SBarry Smith ierr = KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat,NULL);CHKERRQ(ierr); 56998e04984SBarry Smith if (mmat == pc->pmat) opsset = PETSC_FALSE; 57098e04984SBarry Smith } 5714a5f07a5SMark F. Adams 5724a5f07a5SMark F. Adams if (!opsset) { 57371b23a65SMark F. Adams ierr = PCGetUseAmat(pc,&use_amat);CHKERRQ(ierr); 5744a5f07a5SMark F. Adams if(use_amat){ 575f3fbd535SBarry 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); 576f3fbd535SBarry Smith ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr); 577f3fbd535SBarry Smith } 5784a5f07a5SMark F. Adams else { 5794a5f07a5SMark F. Adams ierr = PetscInfo(pc,"Using matrix (pmat) operators to define finest grid operator \n because PCMGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n");CHKERRQ(ierr); 5804a5f07a5SMark F. Adams ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat,pc->flag);CHKERRQ(ierr); 5814a5f07a5SMark F. Adams } 5824a5f07a5SMark F. Adams } 583f3fbd535SBarry Smith 58401bc414fSDmitry 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? */ 585ce4cda84SJed Brown if (pc->dm && mg->galerkin != 2 && !pc->setupcalled) { 5862d2b81a6SBarry Smith /* construct the interpolation from the DMs */ 5872e499ae9SBarry Smith Mat p; 58873250ac0SBarry Smith Vec rscale; 589218a07d4SBarry Smith ierr = PetscMalloc(n*sizeof(DM),&dms);CHKERRQ(ierr); 590218a07d4SBarry Smith dms[n-1] = pc->dm; 591218a07d4SBarry Smith for (i=n-2; i>-1; i--) { 592942e3340SBarry Smith DMKSP kdm; 5932ee06e3bSJed Brown ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr); 5943c0c59f3SBarry Smith ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr); 5953c0c59f3SBarry Smith if (mg->galerkin) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);} 596942e3340SBarry Smith ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr); 597d1a3e677SJed Brown /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take 598d1a3e677SJed Brown * a bitwise OR of computing the matrix, RHS, and initial iterate. */ 5990298fd71SBarry Smith kdm->ops->computerhs = NULL; 6000298fd71SBarry Smith kdm->rhsctx = NULL; 60124c3aa18SBarry Smith if (!mglevels[i+1]->interpolate) { 602e727c939SJed Brown ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr); 6032d2b81a6SBarry Smith ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr); 60400ac8be1SBarry Smith if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);} 60573250ac0SBarry Smith ierr = VecDestroy(&rscale);CHKERRQ(ierr); 6066bf464f9SBarry Smith ierr = MatDestroy(&p);CHKERRQ(ierr); 607218a07d4SBarry Smith } 60824c3aa18SBarry Smith } 6092d2b81a6SBarry Smith 610218a07d4SBarry Smith for (i=n-2; i>-1; i--) { 6116bf464f9SBarry Smith ierr = DMDestroy(&dms[i]);CHKERRQ(ierr); 612218a07d4SBarry Smith } 6132d2b81a6SBarry Smith ierr = PetscFree(dms);CHKERRQ(ierr); 614ce4cda84SJed Brown } 615cab2e9ccSBarry Smith 616ce4cda84SJed Brown if (pc->dm && !pc->setupcalled) { 617ce4cda84SJed Brown /* finest smoother also gets DM but it is not active, independent of whether galerkin==2 */ 618cab2e9ccSBarry Smith ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr); 619cab2e9ccSBarry Smith ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr); 620218a07d4SBarry Smith } 621218a07d4SBarry Smith 622c91913d3SJed Brown if (mg->galerkin == 1) { 623f3fbd535SBarry Smith Mat B; 624f3fbd535SBarry Smith /* currently only handle case where mat and pmat are the same on coarser levels */ 625f3fbd535SBarry Smith ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr); 626f3fbd535SBarry Smith if (!pc->setupcalled) { 627f3fbd535SBarry Smith for (i=n-2; i>-1; i--) { 628f3fbd535SBarry Smith ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 629f3fbd535SBarry Smith ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 630f3fbd535SBarry Smith if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);} 631f3fbd535SBarry Smith dB = B; 632f3fbd535SBarry Smith } 633cd9507b2SBarry Smith if (n > 1) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);} 634f3fbd535SBarry Smith } else { 635f3fbd535SBarry Smith for (i=n-2; i>-1; i--) { 6360298fd71SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&B,NULL);CHKERRQ(ierr); 637f3fbd535SBarry Smith ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 638f3fbd535SBarry Smith ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 639f3fbd535SBarry Smith dB = B; 640f3fbd535SBarry Smith } 641f3fbd535SBarry Smith } 642ce4cda84SJed Brown } else if (!mg->galerkin && pc->dm && pc->dm->x) { 643cab2e9ccSBarry Smith /* need to restrict Jacobian location to coarser meshes for evaluation */ 644cab2e9ccSBarry Smith for (i=n-2; i>-1; i--) { 645c88c5224SJed Brown Mat R; 646c88c5224SJed Brown Vec rscale; 647cab2e9ccSBarry Smith if (!mglevels[i]->smoothd->dm->x) { 648cab2e9ccSBarry Smith Vec *vecs; 6490298fd71SBarry Smith ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vecs,0,NULL);CHKERRQ(ierr); 6502fa5cd67SKarl Rupp 651cab2e9ccSBarry Smith mglevels[i]->smoothd->dm->x = vecs[0]; 6522fa5cd67SKarl Rupp 653cab2e9ccSBarry Smith ierr = PetscFree(vecs);CHKERRQ(ierr); 654cab2e9ccSBarry Smith } 655c88c5224SJed Brown ierr = PCMGGetRestriction(pc,i+1,&R);CHKERRQ(ierr); 656c88c5224SJed Brown ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr); 657c88c5224SJed Brown ierr = MatRestrict(R,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);CHKERRQ(ierr); 658c88c5224SJed Brown ierr = VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,rscale);CHKERRQ(ierr); 659cab2e9ccSBarry Smith } 660f3fbd535SBarry Smith } 661ccceb50cSJed Brown if (!mg->galerkin && pc->dm) { 662caa4e7f2SJed Brown for (i=n-2; i>=0; i--) { 663ccceb50cSJed Brown DM dmfine,dmcoarse; 664caa4e7f2SJed Brown Mat Restrict,Inject; 665caa4e7f2SJed Brown Vec rscale; 666ccceb50cSJed Brown ierr = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr); 667ccceb50cSJed Brown ierr = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr); 668caa4e7f2SJed Brown ierr = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr); 669caa4e7f2SJed Brown ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr); 6700298fd71SBarry Smith Inject = NULL; /* Callback should create it if it needs Injection */ 671caa4e7f2SJed Brown ierr = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr); 672caa4e7f2SJed Brown } 673caa4e7f2SJed Brown } 674f3fbd535SBarry Smith 675f3fbd535SBarry Smith if (!pc->setupcalled) { 676f3fbd535SBarry Smith for (i=0; i<n; i++) { 677f3fbd535SBarry Smith ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr); 678f3fbd535SBarry Smith } 679f3fbd535SBarry Smith for (i=1; i<n; i++) { 680f3fbd535SBarry Smith if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) { 681f3fbd535SBarry Smith ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr); 682f3fbd535SBarry Smith } 683f3fbd535SBarry Smith } 684f3fbd535SBarry Smith for (i=1; i<n; i++) { 685c88c5224SJed Brown ierr = PCMGGetInterpolation(pc,i,&mglevels[i]->interpolate);CHKERRQ(ierr); 686c88c5224SJed Brown ierr = PCMGGetRestriction(pc,i,&mglevels[i]->restrct);CHKERRQ(ierr); 687f3fbd535SBarry Smith } 688f3fbd535SBarry Smith for (i=0; i<n-1; i++) { 689f3fbd535SBarry Smith if (!mglevels[i]->b) { 690f3fbd535SBarry Smith Vec *vec; 6910298fd71SBarry Smith ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr); 692f3fbd535SBarry Smith ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr); 6936bf464f9SBarry Smith ierr = VecDestroy(vec);CHKERRQ(ierr); 694f3fbd535SBarry Smith ierr = PetscFree(vec);CHKERRQ(ierr); 695f3fbd535SBarry Smith } 696f3fbd535SBarry Smith if (!mglevels[i]->r && i) { 697f3fbd535SBarry Smith ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 698f3fbd535SBarry Smith ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr); 6996bf464f9SBarry Smith ierr = VecDestroy(&tvec);CHKERRQ(ierr); 700f3fbd535SBarry Smith } 701f3fbd535SBarry Smith if (!mglevels[i]->x) { 702f3fbd535SBarry Smith ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 703f3fbd535SBarry Smith ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr); 7046bf464f9SBarry Smith ierr = VecDestroy(&tvec);CHKERRQ(ierr); 705f3fbd535SBarry Smith } 706f3fbd535SBarry Smith } 707f3fbd535SBarry Smith if (n != 1 && !mglevels[n-1]->r) { 708f3fbd535SBarry Smith /* PCMGSetR() on the finest level if user did not supply it */ 709f3fbd535SBarry Smith Vec *vec; 7100298fd71SBarry Smith ierr = KSPGetVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr); 711f3fbd535SBarry Smith ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr); 7126bf464f9SBarry Smith ierr = VecDestroy(vec);CHKERRQ(ierr); 713f3fbd535SBarry Smith ierr = PetscFree(vec);CHKERRQ(ierr); 714f3fbd535SBarry Smith } 715f3fbd535SBarry Smith } 716f3fbd535SBarry Smith 71798e04984SBarry Smith if (pc->dm) { 71898e04984SBarry Smith /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */ 71998e04984SBarry Smith for (i=0; i<n-1; i++) { 72098e04984SBarry Smith if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX; 72198e04984SBarry Smith } 72298e04984SBarry Smith } 723f3fbd535SBarry Smith 724f3fbd535SBarry Smith for (i=1; i<n; i++) { 725f3fbd535SBarry Smith if (mglevels[i]->smoothu == mglevels[i]->smoothd) { 726f3fbd535SBarry Smith /* if doing only down then initial guess is zero */ 727f3fbd535SBarry Smith ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 728f3fbd535SBarry Smith } 72963e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 730f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr); 73163e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 732d42688cbSBarry Smith if (!mglevels[i]->residual) { 733d42688cbSBarry Smith Mat mat; 7340298fd71SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&mat,NULL);CHKERRQ(ierr); 735d0e4de75SBarry Smith ierr = PCMGSetResidual(pc,i,PCMGResidual_Default,mat);CHKERRQ(ierr); 736d42688cbSBarry Smith } 737f3fbd535SBarry Smith } 738f3fbd535SBarry Smith for (i=1; i<n; i++) { 739f3fbd535SBarry Smith if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) { 740f3fbd535SBarry Smith Mat downmat,downpmat; 741f3fbd535SBarry Smith MatStructure matflag; 742f3fbd535SBarry Smith 743f3fbd535SBarry Smith /* check if operators have been set for up, if not use down operators to set them */ 7440298fd71SBarry Smith ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);CHKERRQ(ierr); 745f3fbd535SBarry Smith if (!opsset) { 746f3fbd535SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr); 747f3fbd535SBarry Smith ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr); 748f3fbd535SBarry Smith } 749f3fbd535SBarry Smith 750f3fbd535SBarry Smith ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr); 75163e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 752f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr); 75363e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 754f3fbd535SBarry Smith } 755f3fbd535SBarry Smith } 756f3fbd535SBarry Smith 757f3fbd535SBarry Smith /* 758f3fbd535SBarry Smith If coarse solver is not direct method then DO NOT USE preonly 759f3fbd535SBarry Smith */ 760251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mglevels[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr); 761f3fbd535SBarry Smith if (preonly) { 762251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr); 763251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr); 764251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr); 765251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)cpc,PCSVD,&svd);CHKERRQ(ierr); 766fd1303e9SJungho Lee if (!lu && !redundant && !cholesky && !svd) { 767f3fbd535SBarry Smith ierr = KSPSetType(mglevels[0]->smoothd,KSPGMRES);CHKERRQ(ierr); 768f3fbd535SBarry Smith } 769f3fbd535SBarry Smith } 770f3fbd535SBarry Smith 771f3fbd535SBarry Smith if (!pc->setupcalled) { 772f3fbd535SBarry Smith ierr = KSPSetFromOptions(mglevels[0]->smoothd);CHKERRQ(ierr); 773f3fbd535SBarry Smith } 774f3fbd535SBarry Smith 77563e6d426SJed Brown if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 776f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr); 77763e6d426SJed Brown if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 778f3fbd535SBarry Smith 779f3fbd535SBarry Smith /* 780f3fbd535SBarry Smith Dump the interpolation/restriction matrices plus the 781e3c5b3baSBarry Smith Jacobian/stiffness on each level. This allows MATLAB users to 782f3fbd535SBarry Smith easily check if the Galerkin condition A_c = R A_f R^T is satisfied. 783f3fbd535SBarry Smith 784f3fbd535SBarry Smith Only support one or the other at the same time. 785f3fbd535SBarry Smith */ 786f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER) 7870298fd71SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);CHKERRQ(ierr); 788ce94432eSBarry Smith if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc)); 789f3fbd535SBarry Smith dump = PETSC_FALSE; 790f3fbd535SBarry Smith #endif 7910298fd71SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);CHKERRQ(ierr); 792ce94432eSBarry Smith if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc)); 793f3fbd535SBarry Smith 794f3fbd535SBarry Smith if (viewer) { 795f3fbd535SBarry Smith for (i=1; i<n; i++) { 796f3fbd535SBarry Smith ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr); 797f3fbd535SBarry Smith } 798f3fbd535SBarry Smith for (i=0; i<n; i++) { 799f3fbd535SBarry Smith ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr); 800f3fbd535SBarry Smith ierr = MatView(pc->mat,viewer);CHKERRQ(ierr); 801f3fbd535SBarry Smith } 802f3fbd535SBarry Smith } 803f3fbd535SBarry Smith PetscFunctionReturn(0); 804f3fbd535SBarry Smith } 805f3fbd535SBarry Smith 806f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/ 807f3fbd535SBarry Smith 808f3fbd535SBarry Smith #undef __FUNCT__ 8099dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetLevels" 8104b9ad928SBarry Smith /*@ 81197177400SBarry Smith PCMGGetLevels - Gets the number of levels to use with MG. 8124b9ad928SBarry Smith 8134b9ad928SBarry Smith Not Collective 8144b9ad928SBarry Smith 8154b9ad928SBarry Smith Input Parameter: 8164b9ad928SBarry Smith . pc - the preconditioner context 8174b9ad928SBarry Smith 8184b9ad928SBarry Smith Output parameter: 8194b9ad928SBarry Smith . levels - the number of levels 8204b9ad928SBarry Smith 8214b9ad928SBarry Smith Level: advanced 8224b9ad928SBarry Smith 8234b9ad928SBarry Smith .keywords: MG, get, levels, multigrid 8244b9ad928SBarry Smith 82597177400SBarry Smith .seealso: PCMGSetLevels() 8264b9ad928SBarry Smith @*/ 8277087cfbeSBarry Smith PetscErrorCode PCMGGetLevels(PC pc,PetscInt *levels) 8284b9ad928SBarry Smith { 829f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 8304b9ad928SBarry Smith 8314b9ad928SBarry Smith PetscFunctionBegin; 8320700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 8334482741eSBarry Smith PetscValidIntPointer(levels,2); 834f3fbd535SBarry Smith *levels = mg->nlevels; 8354b9ad928SBarry Smith PetscFunctionReturn(0); 8364b9ad928SBarry Smith } 8374b9ad928SBarry Smith 8384b9ad928SBarry Smith #undef __FUNCT__ 8399dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetType" 8404b9ad928SBarry Smith /*@ 84197177400SBarry Smith PCMGSetType - Determines the form of multigrid to use: 8424b9ad928SBarry Smith multiplicative, additive, full, or the Kaskade algorithm. 8434b9ad928SBarry Smith 844ad4df100SBarry Smith Logically Collective on PC 8454b9ad928SBarry Smith 8464b9ad928SBarry Smith Input Parameters: 8474b9ad928SBarry Smith + pc - the preconditioner context 8489dcbbd2bSBarry Smith - form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE, 8499dcbbd2bSBarry Smith PC_MG_FULL, PC_MG_KASKADE 8504b9ad928SBarry Smith 8514b9ad928SBarry Smith Options Database Key: 8524b9ad928SBarry Smith . -pc_mg_type <form> - Sets <form>, one of multiplicative, 8534b9ad928SBarry Smith additive, full, kaskade 8544b9ad928SBarry Smith 8554b9ad928SBarry Smith Level: advanced 8564b9ad928SBarry Smith 8574b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid 8584b9ad928SBarry Smith 85997177400SBarry Smith .seealso: PCMGSetLevels() 8604b9ad928SBarry Smith @*/ 8617087cfbeSBarry Smith PetscErrorCode PCMGSetType(PC pc,PCMGType form) 8624b9ad928SBarry Smith { 863f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 8644b9ad928SBarry Smith 8654b9ad928SBarry Smith PetscFunctionBegin; 8660700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 867c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc,form,2); 86831567311SBarry Smith mg->am = form; 8699dcbbd2bSBarry Smith if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG; 8704b9ad928SBarry Smith else pc->ops->applyrichardson = 0; 8714b9ad928SBarry Smith PetscFunctionReturn(0); 8724b9ad928SBarry Smith } 8734b9ad928SBarry Smith 8744b9ad928SBarry Smith #undef __FUNCT__ 8750d353602SBarry Smith #define __FUNCT__ "PCMGSetCycleType" 8764b9ad928SBarry Smith /*@ 8770d353602SBarry Smith PCMGSetCycleType - Sets the type cycles to use. Use PCMGSetCycleTypeOnLevel() for more 8784b9ad928SBarry Smith complicated cycling. 8794b9ad928SBarry Smith 880ad4df100SBarry Smith Logically Collective on PC 8814b9ad928SBarry Smith 8824b9ad928SBarry Smith Input Parameters: 883c2be2410SBarry Smith + pc - the multigrid context 8840d353602SBarry Smith - PC_MG_CYCLE_V or PC_MG_CYCLE_W 8854b9ad928SBarry Smith 8864b9ad928SBarry Smith Options Database Key: 8870d353602SBarry Smith $ -pc_mg_cycle_type v or w 8884b9ad928SBarry Smith 8894b9ad928SBarry Smith Level: advanced 8904b9ad928SBarry Smith 8914b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 8924b9ad928SBarry Smith 8930d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel() 8944b9ad928SBarry Smith @*/ 8957087cfbeSBarry Smith PetscErrorCode PCMGSetCycleType(PC pc,PCMGCycleType n) 8964b9ad928SBarry Smith { 897f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 898f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 89979416396SBarry Smith PetscInt i,levels; 9004b9ad928SBarry Smith 9014b9ad928SBarry Smith PetscFunctionBegin; 9020700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 903ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 904c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 905f3fbd535SBarry Smith levels = mglevels[0]->levels; 9064b9ad928SBarry Smith 9072fa5cd67SKarl Rupp for (i=0; i<levels; i++) mglevels[i]->cycles = n; 9084b9ad928SBarry Smith PetscFunctionReturn(0); 9094b9ad928SBarry Smith } 9104b9ad928SBarry Smith 9114b9ad928SBarry Smith #undef __FUNCT__ 9128cc2d5dfSBarry Smith #define __FUNCT__ "PCMGMultiplicativeSetCycles" 9138cc2d5dfSBarry Smith /*@ 9148cc2d5dfSBarry Smith PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step 9158cc2d5dfSBarry Smith of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used 9168cc2d5dfSBarry Smith 917ad4df100SBarry Smith Logically Collective on PC 9188cc2d5dfSBarry Smith 9198cc2d5dfSBarry Smith Input Parameters: 9208cc2d5dfSBarry Smith + pc - the multigrid context 9218cc2d5dfSBarry Smith - n - number of cycles (default is 1) 9228cc2d5dfSBarry Smith 9238cc2d5dfSBarry Smith Options Database Key: 9248cc2d5dfSBarry Smith $ -pc_mg_multiplicative_cycles n 9258cc2d5dfSBarry Smith 9268cc2d5dfSBarry Smith Level: advanced 9278cc2d5dfSBarry Smith 9288cc2d5dfSBarry Smith Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType() 9298cc2d5dfSBarry Smith 9308cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 9318cc2d5dfSBarry Smith 9328cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType() 9338cc2d5dfSBarry Smith @*/ 9347087cfbeSBarry Smith PetscErrorCode PCMGMultiplicativeSetCycles(PC pc,PetscInt n) 9358cc2d5dfSBarry Smith { 936f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 937f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 9388cc2d5dfSBarry Smith PetscInt i,levels; 9398cc2d5dfSBarry Smith 9408cc2d5dfSBarry Smith PetscFunctionBegin; 9410700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 942ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 943c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 944f3fbd535SBarry Smith levels = mglevels[0]->levels; 9458cc2d5dfSBarry Smith 9462fa5cd67SKarl Rupp for (i=0; i<levels; i++) mg->cyclesperpcapply = n; 9478cc2d5dfSBarry Smith PetscFunctionReturn(0); 9488cc2d5dfSBarry Smith } 9498cc2d5dfSBarry Smith 9508cc2d5dfSBarry Smith #undef __FUNCT__ 9519dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetGalerkin" 952c2be2410SBarry Smith /*@ 95397177400SBarry Smith PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the 954c2be2410SBarry Smith finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t 955c2be2410SBarry Smith 956ad4df100SBarry Smith Logically Collective on PC 957c2be2410SBarry Smith 958c2be2410SBarry Smith Input Parameters: 959c91913d3SJed Brown + pc - the multigrid context 960c91913d3SJed Brown - use - PETSC_TRUE to use the Galerkin process to compute coarse-level operators 961c2be2410SBarry Smith 962c2be2410SBarry Smith Options Database Key: 963c2be2410SBarry Smith $ -pc_mg_galerkin 964c2be2410SBarry Smith 965c2be2410SBarry Smith Level: intermediate 966c2be2410SBarry Smith 967c2be2410SBarry Smith .keywords: MG, set, Galerkin 968c2be2410SBarry Smith 9693fc8bf9cSBarry Smith .seealso: PCMGGetGalerkin() 9703fc8bf9cSBarry Smith 971c2be2410SBarry Smith @*/ 972c91913d3SJed Brown PetscErrorCode PCMGSetGalerkin(PC pc,PetscBool use) 973c2be2410SBarry Smith { 974f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 975c2be2410SBarry Smith 976c2be2410SBarry Smith PetscFunctionBegin; 9770700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 978789726d7SBarry Smith mg->galerkin = use ? 1 : 0; 979c2be2410SBarry Smith PetscFunctionReturn(0); 980c2be2410SBarry Smith } 981c2be2410SBarry Smith 982c2be2410SBarry Smith #undef __FUNCT__ 9833fc8bf9cSBarry Smith #define __FUNCT__ "PCMGGetGalerkin" 9843fc8bf9cSBarry Smith /*@ 9853fc8bf9cSBarry Smith PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. 9863fc8bf9cSBarry Smith A_i-1 = r_i * A_i * r_i^t 9873fc8bf9cSBarry Smith 9883fc8bf9cSBarry Smith Not Collective 9893fc8bf9cSBarry Smith 9903fc8bf9cSBarry Smith Input Parameter: 9913fc8bf9cSBarry Smith . pc - the multigrid context 9923fc8bf9cSBarry Smith 9933fc8bf9cSBarry Smith Output Parameter: 9943fc8bf9cSBarry Smith . gelerkin - PETSC_TRUE or PETSC_FALSE 9953fc8bf9cSBarry Smith 9963fc8bf9cSBarry Smith Options Database Key: 9973fc8bf9cSBarry Smith $ -pc_mg_galerkin 9983fc8bf9cSBarry Smith 9993fc8bf9cSBarry Smith Level: intermediate 10003fc8bf9cSBarry Smith 10013fc8bf9cSBarry Smith .keywords: MG, set, Galerkin 10023fc8bf9cSBarry Smith 10033fc8bf9cSBarry Smith .seealso: PCMGSetGalerkin() 10043fc8bf9cSBarry Smith 10053fc8bf9cSBarry Smith @*/ 10067087cfbeSBarry Smith PetscErrorCode PCMGGetGalerkin(PC pc,PetscBool *galerkin) 10073fc8bf9cSBarry Smith { 1008f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 10093fc8bf9cSBarry Smith 10103fc8bf9cSBarry Smith PetscFunctionBegin; 10110700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 101213fdb3acSJose Roman *galerkin = (PetscBool)mg->galerkin; 10133fc8bf9cSBarry Smith PetscFunctionReturn(0); 10143fc8bf9cSBarry Smith } 10153fc8bf9cSBarry Smith 10163fc8bf9cSBarry Smith #undef __FUNCT__ 10179dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothDown" 10184b9ad928SBarry Smith /*@ 101997177400SBarry Smith PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to 102097177400SBarry Smith use on all levels. Use PCMGGetSmootherDown() to set different 10214b9ad928SBarry Smith pre-smoothing steps on different levels. 10224b9ad928SBarry Smith 1023ad4df100SBarry Smith Logically Collective on PC 10244b9ad928SBarry Smith 10254b9ad928SBarry Smith Input Parameters: 10264b9ad928SBarry Smith + mg - the multigrid context 10274b9ad928SBarry Smith - n - the number of smoothing steps 10284b9ad928SBarry Smith 10294b9ad928SBarry Smith Options Database Key: 10304b9ad928SBarry Smith . -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps 10314b9ad928SBarry Smith 10324b9ad928SBarry Smith Level: advanced 10334b9ad928SBarry Smith 10344b9ad928SBarry Smith .keywords: MG, smooth, down, pre-smoothing, steps, multigrid 10354b9ad928SBarry Smith 103697177400SBarry Smith .seealso: PCMGSetNumberSmoothUp() 10374b9ad928SBarry Smith @*/ 10387087cfbeSBarry Smith PetscErrorCode PCMGSetNumberSmoothDown(PC pc,PetscInt n) 10394b9ad928SBarry Smith { 1040f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1041f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 10426849ba73SBarry Smith PetscErrorCode ierr; 104379416396SBarry Smith PetscInt i,levels; 10444b9ad928SBarry Smith 10454b9ad928SBarry Smith PetscFunctionBegin; 10460700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1047ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 1048c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 1049f3fbd535SBarry Smith levels = mglevels[0]->levels; 10504b9ad928SBarry Smith 1051b05257ddSBarry Smith for (i=1; i<levels; i++) { 10524b9ad928SBarry Smith /* make sure smoother up and down are different */ 10530298fd71SBarry Smith ierr = PCMGGetSmootherUp(pc,i,NULL);CHKERRQ(ierr); 1054f3fbd535SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 10552fa5cd67SKarl Rupp 105631567311SBarry Smith mg->default_smoothd = n; 10574b9ad928SBarry Smith } 10584b9ad928SBarry Smith PetscFunctionReturn(0); 10594b9ad928SBarry Smith } 10604b9ad928SBarry Smith 10614b9ad928SBarry Smith #undef __FUNCT__ 10629dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothUp" 10634b9ad928SBarry Smith /*@ 106497177400SBarry Smith PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use 106597177400SBarry Smith on all levels. Use PCMGGetSmootherUp() to set different numbers of 10664b9ad928SBarry Smith post-smoothing steps on different levels. 10674b9ad928SBarry Smith 1068ad4df100SBarry Smith Logically Collective on PC 10694b9ad928SBarry Smith 10704b9ad928SBarry Smith Input Parameters: 10714b9ad928SBarry Smith + mg - the multigrid context 10724b9ad928SBarry Smith - n - the number of smoothing steps 10734b9ad928SBarry Smith 10744b9ad928SBarry Smith Options Database Key: 10754b9ad928SBarry Smith . -pc_mg_smoothup <n> - Sets number of post-smoothing steps 10764b9ad928SBarry Smith 10774b9ad928SBarry Smith Level: advanced 10784b9ad928SBarry Smith 10794b9ad928SBarry Smith Note: this does not set a value on the coarsest grid, since we assume that 1080a8c7a070SBarry Smith there is no separate smooth up on the coarsest grid. 10814b9ad928SBarry Smith 10824b9ad928SBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid 10834b9ad928SBarry Smith 108497177400SBarry Smith .seealso: PCMGSetNumberSmoothDown() 10854b9ad928SBarry Smith @*/ 10867087cfbeSBarry Smith PetscErrorCode PCMGSetNumberSmoothUp(PC pc,PetscInt n) 10874b9ad928SBarry Smith { 1088f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1089f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 10906849ba73SBarry Smith PetscErrorCode ierr; 109179416396SBarry Smith PetscInt i,levels; 10924b9ad928SBarry Smith 10934b9ad928SBarry Smith PetscFunctionBegin; 10940700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1095ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 1096c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 1097f3fbd535SBarry Smith levels = mglevels[0]->levels; 10984b9ad928SBarry Smith 10994b9ad928SBarry Smith for (i=1; i<levels; i++) { 11004b9ad928SBarry Smith /* make sure smoother up and down are different */ 11010298fd71SBarry Smith ierr = PCMGGetSmootherUp(pc,i,NULL);CHKERRQ(ierr); 1102f3fbd535SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 11032fa5cd67SKarl Rupp 110431567311SBarry Smith mg->default_smoothu = n; 11054b9ad928SBarry Smith } 11064b9ad928SBarry Smith PetscFunctionReturn(0); 11074b9ad928SBarry Smith } 11084b9ad928SBarry Smith 11094b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/ 11104b9ad928SBarry Smith 11113b09bd56SBarry Smith /*MC 1112ccb205f8SBarry Smith PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional 11133b09bd56SBarry Smith information about the coarser grid matrices and restriction/interpolation operators. 11143b09bd56SBarry Smith 11153b09bd56SBarry Smith Options Database Keys: 11163b09bd56SBarry Smith + -pc_mg_levels <nlevels> - number of levels including finest 11170d353602SBarry Smith . -pc_mg_cycles v or w 111879416396SBarry Smith . -pc_mg_smoothup <n> - number of smoothing steps after interpolation 11193b09bd56SBarry Smith . -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator 11208c1c2452SJed Brown . -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default 11213b09bd56SBarry Smith . -pc_mg_log - log information about time spent on each level of the solver 11223b09bd56SBarry Smith . -pc_mg_monitor - print information on the multigrid convergence 11238cf6037eSBarry Smith . -pc_mg_galerkin - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R' 11248cf6037eSBarry Smith . -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1) 11258cf6037eSBarry Smith . -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 1126e3c5b3baSBarry Smith to the Socket viewer for reading from MATLAB. 11278cf6037eSBarry Smith - -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices 11288cf6037eSBarry Smith to the binary output file called binaryoutput 11293b09bd56SBarry Smith 113024c3aa18SBarry 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 11313b09bd56SBarry Smith 11328cf6037eSBarry Smith When run with a single level the smoother options are used on that level NOT the coarse grid solver options 11338cf6037eSBarry Smith 11343b09bd56SBarry Smith Level: intermediate 11353b09bd56SBarry Smith 11368f87f92bSBarry Smith Concepts: multigrid/multilevel 11373b09bd56SBarry Smith 11388cf6037eSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE 11390d353602SBarry Smith PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(), 114097177400SBarry Smith PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 114197177400SBarry Smith PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 11420d353602SBarry Smith PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 11433b09bd56SBarry Smith M*/ 11443b09bd56SBarry Smith 11454b9ad928SBarry Smith #undef __FUNCT__ 11464b9ad928SBarry Smith #define __FUNCT__ "PCCreate_MG" 11478cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc) 11484b9ad928SBarry Smith { 1149f3fbd535SBarry Smith PC_MG *mg; 1150f3fbd535SBarry Smith PetscErrorCode ierr; 1151f3fbd535SBarry Smith 11524b9ad928SBarry Smith PetscFunctionBegin; 1153f3fbd535SBarry Smith ierr = PetscNewLog(pc,PC_MG,&mg);CHKERRQ(ierr); 1154f3fbd535SBarry Smith pc->data = (void*)mg; 1155f3fbd535SBarry Smith mg->nlevels = -1; 1156f3fbd535SBarry Smith 115737a44384SMark Adams pc->useAmat = PETSC_TRUE; 115837a44384SMark Adams 11594b9ad928SBarry Smith pc->ops->apply = PCApply_MG; 11604b9ad928SBarry Smith pc->ops->setup = PCSetUp_MG; 1161a06653b4SBarry Smith pc->ops->reset = PCReset_MG; 11624b9ad928SBarry Smith pc->ops->destroy = PCDestroy_MG; 11634b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_MG; 11644b9ad928SBarry Smith pc->ops->view = PCView_MG; 11654b9ad928SBarry Smith PetscFunctionReturn(0); 11664b9ad928SBarry Smith } 1167