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; 3357622a8eSBarry Smith ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %g is less than absolute tolerance %g\n",(double)rnorm,(double)mg->abstol);CHKERRQ(ierr); 344b9ad928SBarry Smith } else { 354d0a8057SBarry Smith *reason = PCRICHARDSON_CONVERGED_RTOL; 3657622a8eSBarry Smith ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %g is less than relative tolerance times initial residual norm %g\n",(double)rnorm,(double)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 197785e854fSJed Brown ierr = PetscMalloc1(levels,&mglevels);CHKERRQ(ierr); 1983bb1ff40SBarry 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++) { 204b00a9115SJed Brown ierr = PetscNewLog(pc,&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); 2190059c7bdSJed Brown ierr = KSPSetConvergenceTest(mglevels[i]->smoothd,KSPConvergedSkip,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 } 2513bb1ff40SBarry 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) { 320*23ee1639SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,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 Vec tvec; 542218a07d4SBarry Smith DM *dms; 543649052a6SBarry Smith PetscViewer viewer = 0; 544f3fbd535SBarry Smith 545f3fbd535SBarry Smith PetscFunctionBegin; 54601bc414fSDmitry Karpeev /* FIX: Move this to PCSetFromOptions_MG? */ 5473aeaf226SBarry Smith if (mg->usedmfornumberoflevels) { 5483aeaf226SBarry Smith PetscInt levels; 5493aeaf226SBarry Smith ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr); 5503aeaf226SBarry Smith levels++; 5513aeaf226SBarry Smith if (levels > n) { /* the problem is now being solved on a finer grid */ 5520298fd71SBarry Smith ierr = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr); 5533aeaf226SBarry Smith n = levels; 5543aeaf226SBarry Smith ierr = PCSetFromOptions(pc);CHKERRQ(ierr); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */ 5553aeaf226SBarry Smith mglevels = mg->levels; 5563aeaf226SBarry Smith } 5573aeaf226SBarry Smith } 55898e04984SBarry Smith ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr); 5593aeaf226SBarry Smith 560f3fbd535SBarry Smith 561f3fbd535SBarry Smith /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */ 562f3fbd535SBarry Smith /* so use those from global PC */ 563f3fbd535SBarry Smith /* Is this what we always want? What if user wants to keep old one? */ 5640298fd71SBarry Smith ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset);CHKERRQ(ierr); 56598e04984SBarry Smith if (opsset) { 56698e04984SBarry Smith Mat mmat; 567*23ee1639SBarry Smith ierr = KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat);CHKERRQ(ierr); 56898e04984SBarry Smith if (mmat == pc->pmat) opsset = PETSC_FALSE; 56998e04984SBarry Smith } 5704a5f07a5SMark F. Adams 5714a5f07a5SMark F. Adams if (!opsset) { 57271b23a65SMark F. Adams ierr = PCGetUseAmat(pc,&use_amat);CHKERRQ(ierr); 5734a5f07a5SMark F. Adams if(use_amat){ 574f3fbd535SBarry 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); 575*23ee1639SBarry Smith ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat);CHKERRQ(ierr); 576f3fbd535SBarry Smith } 5774a5f07a5SMark F. Adams else { 5784a5f07a5SMark 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); 579*23ee1639SBarry Smith ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat);CHKERRQ(ierr); 5804a5f07a5SMark F. Adams } 5814a5f07a5SMark F. Adams } 582f3fbd535SBarry Smith 58301bc414fSDmitry 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? */ 584ce4cda84SJed Brown if (pc->dm && mg->galerkin != 2 && !pc->setupcalled) { 5852d2b81a6SBarry Smith /* construct the interpolation from the DMs */ 5862e499ae9SBarry Smith Mat p; 58773250ac0SBarry Smith Vec rscale; 588785e854fSJed Brown ierr = PetscMalloc1(n,&dms);CHKERRQ(ierr); 589218a07d4SBarry Smith dms[n-1] = pc->dm; 590218a07d4SBarry Smith for (i=n-2; i>-1; i--) { 591942e3340SBarry Smith DMKSP kdm; 5922ee06e3bSJed Brown ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr); 5933c0c59f3SBarry Smith ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr); 5943c0c59f3SBarry Smith if (mg->galerkin) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);} 595942e3340SBarry Smith ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr); 596d1a3e677SJed Brown /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take 597d1a3e677SJed Brown * a bitwise OR of computing the matrix, RHS, and initial iterate. */ 5980298fd71SBarry Smith kdm->ops->computerhs = NULL; 5990298fd71SBarry Smith kdm->rhsctx = NULL; 60024c3aa18SBarry Smith if (!mglevels[i+1]->interpolate) { 601e727c939SJed Brown ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr); 6022d2b81a6SBarry Smith ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr); 60300ac8be1SBarry Smith if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);} 60473250ac0SBarry Smith ierr = VecDestroy(&rscale);CHKERRQ(ierr); 6056bf464f9SBarry Smith ierr = MatDestroy(&p);CHKERRQ(ierr); 606218a07d4SBarry Smith } 60724c3aa18SBarry Smith } 6082d2b81a6SBarry Smith 609218a07d4SBarry Smith for (i=n-2; i>-1; i--) { 6106bf464f9SBarry Smith ierr = DMDestroy(&dms[i]);CHKERRQ(ierr); 611218a07d4SBarry Smith } 6122d2b81a6SBarry Smith ierr = PetscFree(dms);CHKERRQ(ierr); 613ce4cda84SJed Brown } 614cab2e9ccSBarry Smith 615ce4cda84SJed Brown if (pc->dm && !pc->setupcalled) { 616ce4cda84SJed Brown /* finest smoother also gets DM but it is not active, independent of whether galerkin==2 */ 617cab2e9ccSBarry Smith ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr); 618cab2e9ccSBarry Smith ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr); 619218a07d4SBarry Smith } 620218a07d4SBarry Smith 621c91913d3SJed Brown if (mg->galerkin == 1) { 622f3fbd535SBarry Smith Mat B; 623f3fbd535SBarry Smith /* currently only handle case where mat and pmat are the same on coarser levels */ 624*23ee1639SBarry Smith ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB);CHKERRQ(ierr); 625f3fbd535SBarry Smith if (!pc->setupcalled) { 626f3fbd535SBarry Smith for (i=n-2; i>-1; i--) { 6279f3d9e6bSJed Brown if (mglevels[i+1]->interpolate == mglevels[i+1]->restrct || !mglevels[i+1]->restrct) { 628f3fbd535SBarry Smith ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 6297d537d92SJed Brown } else { 6307d537d92SJed Brown ierr = MatMatMatMult(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 6317d537d92SJed Brown } 632*23ee1639SBarry Smith ierr = KSPSetOperators(mglevels[i]->smoothd,B,B);CHKERRQ(ierr); 633f3fbd535SBarry Smith if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);} 634f3fbd535SBarry Smith dB = B; 635f3fbd535SBarry Smith } 636cd9507b2SBarry Smith if (n > 1) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);} 637f3fbd535SBarry Smith } else { 638f3fbd535SBarry Smith for (i=n-2; i>-1; i--) { 639*23ee1639SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&B);CHKERRQ(ierr); 6409f3d9e6bSJed Brown if (mglevels[i+1]->interpolate == mglevels[i+1]->restrct || !mglevels[i+1]->restrct) { 641f3fbd535SBarry Smith ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 6427d537d92SJed Brown } else { 6437d537d92SJed Brown ierr = MatMatMatMult(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 6447d537d92SJed Brown } 645*23ee1639SBarry Smith ierr = KSPSetOperators(mglevels[i]->smoothd,B,B);CHKERRQ(ierr); 646f3fbd535SBarry Smith dB = B; 647f3fbd535SBarry Smith } 648f3fbd535SBarry Smith } 649ce4cda84SJed Brown } else if (!mg->galerkin && pc->dm && pc->dm->x) { 650cab2e9ccSBarry Smith /* need to restrict Jacobian location to coarser meshes for evaluation */ 651cab2e9ccSBarry Smith for (i=n-2; i>-1; i--) { 652c88c5224SJed Brown Mat R; 653c88c5224SJed Brown Vec rscale; 654cab2e9ccSBarry Smith if (!mglevels[i]->smoothd->dm->x) { 655cab2e9ccSBarry Smith Vec *vecs; 6560298fd71SBarry Smith ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vecs,0,NULL);CHKERRQ(ierr); 6572fa5cd67SKarl Rupp 658cab2e9ccSBarry Smith mglevels[i]->smoothd->dm->x = vecs[0]; 6592fa5cd67SKarl Rupp 660cab2e9ccSBarry Smith ierr = PetscFree(vecs);CHKERRQ(ierr); 661cab2e9ccSBarry Smith } 662c88c5224SJed Brown ierr = PCMGGetRestriction(pc,i+1,&R);CHKERRQ(ierr); 663c88c5224SJed Brown ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr); 664c88c5224SJed Brown ierr = MatRestrict(R,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);CHKERRQ(ierr); 665c88c5224SJed Brown ierr = VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,rscale);CHKERRQ(ierr); 666cab2e9ccSBarry Smith } 667f3fbd535SBarry Smith } 668ccceb50cSJed Brown if (!mg->galerkin && pc->dm) { 669caa4e7f2SJed Brown for (i=n-2; i>=0; i--) { 670ccceb50cSJed Brown DM dmfine,dmcoarse; 671caa4e7f2SJed Brown Mat Restrict,Inject; 672caa4e7f2SJed Brown Vec rscale; 673ccceb50cSJed Brown ierr = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr); 674ccceb50cSJed Brown ierr = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr); 675caa4e7f2SJed Brown ierr = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr); 676caa4e7f2SJed Brown ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr); 6770298fd71SBarry Smith Inject = NULL; /* Callback should create it if it needs Injection */ 678caa4e7f2SJed Brown ierr = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr); 679caa4e7f2SJed Brown } 680caa4e7f2SJed Brown } 681f3fbd535SBarry Smith 682f3fbd535SBarry Smith if (!pc->setupcalled) { 683f3fbd535SBarry Smith for (i=0; i<n; i++) { 684f3fbd535SBarry Smith ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr); 685f3fbd535SBarry Smith } 686f3fbd535SBarry Smith for (i=1; i<n; i++) { 687f3fbd535SBarry Smith if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) { 688f3fbd535SBarry Smith ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr); 689f3fbd535SBarry Smith } 690f3fbd535SBarry Smith } 691f3fbd535SBarry Smith for (i=1; i<n; i++) { 692c88c5224SJed Brown ierr = PCMGGetInterpolation(pc,i,&mglevels[i]->interpolate);CHKERRQ(ierr); 693c88c5224SJed Brown ierr = PCMGGetRestriction(pc,i,&mglevels[i]->restrct);CHKERRQ(ierr); 694f3fbd535SBarry Smith } 695f3fbd535SBarry Smith for (i=0; i<n-1; i++) { 696f3fbd535SBarry Smith if (!mglevels[i]->b) { 697f3fbd535SBarry Smith Vec *vec; 6980298fd71SBarry Smith ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr); 699f3fbd535SBarry Smith ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr); 7006bf464f9SBarry Smith ierr = VecDestroy(vec);CHKERRQ(ierr); 701f3fbd535SBarry Smith ierr = PetscFree(vec);CHKERRQ(ierr); 702f3fbd535SBarry Smith } 703f3fbd535SBarry Smith if (!mglevels[i]->r && i) { 704f3fbd535SBarry Smith ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 705f3fbd535SBarry Smith ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr); 7066bf464f9SBarry Smith ierr = VecDestroy(&tvec);CHKERRQ(ierr); 707f3fbd535SBarry Smith } 708f3fbd535SBarry Smith if (!mglevels[i]->x) { 709f3fbd535SBarry Smith ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 710f3fbd535SBarry Smith ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr); 7116bf464f9SBarry Smith ierr = VecDestroy(&tvec);CHKERRQ(ierr); 712f3fbd535SBarry Smith } 713f3fbd535SBarry Smith } 714f3fbd535SBarry Smith if (n != 1 && !mglevels[n-1]->r) { 715f3fbd535SBarry Smith /* PCMGSetR() on the finest level if user did not supply it */ 716f3fbd535SBarry Smith Vec *vec; 7170298fd71SBarry Smith ierr = KSPGetVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr); 718f3fbd535SBarry Smith ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr); 7196bf464f9SBarry Smith ierr = VecDestroy(vec);CHKERRQ(ierr); 720f3fbd535SBarry Smith ierr = PetscFree(vec);CHKERRQ(ierr); 721f3fbd535SBarry Smith } 722f3fbd535SBarry Smith } 723f3fbd535SBarry Smith 72498e04984SBarry Smith if (pc->dm) { 72598e04984SBarry Smith /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */ 72698e04984SBarry Smith for (i=0; i<n-1; i++) { 72798e04984SBarry Smith if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX; 72898e04984SBarry Smith } 72998e04984SBarry Smith } 730f3fbd535SBarry Smith 731f3fbd535SBarry Smith for (i=1; i<n; i++) { 732f3fbd535SBarry Smith if (mglevels[i]->smoothu == mglevels[i]->smoothd) { 733f3fbd535SBarry Smith /* if doing only down then initial guess is zero */ 734f3fbd535SBarry Smith ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 735f3fbd535SBarry Smith } 73663e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 737f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr); 73863e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 739d42688cbSBarry Smith if (!mglevels[i]->residual) { 740d42688cbSBarry Smith Mat mat; 741*23ee1639SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&mat);CHKERRQ(ierr); 74254b2cd4bSJed Brown ierr = PCMGSetResidual(pc,i,PCMGResidualDefault,mat);CHKERRQ(ierr); 743d42688cbSBarry Smith } 744f3fbd535SBarry Smith } 745f3fbd535SBarry Smith for (i=1; i<n; i++) { 746f3fbd535SBarry Smith if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) { 747f3fbd535SBarry Smith Mat downmat,downpmat; 748f3fbd535SBarry Smith 749f3fbd535SBarry Smith /* check if operators have been set for up, if not use down operators to set them */ 7500298fd71SBarry Smith ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);CHKERRQ(ierr); 751f3fbd535SBarry Smith if (!opsset) { 752*23ee1639SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);CHKERRQ(ierr); 753*23ee1639SBarry Smith ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat);CHKERRQ(ierr); 754f3fbd535SBarry Smith } 755f3fbd535SBarry Smith 756f3fbd535SBarry Smith ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr); 75763e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 758f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr); 75963e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 760f3fbd535SBarry Smith } 761f3fbd535SBarry Smith } 762f3fbd535SBarry Smith 763f3fbd535SBarry Smith /* 764f3fbd535SBarry Smith If coarse solver is not direct method then DO NOT USE preonly 765f3fbd535SBarry Smith */ 766251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mglevels[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr); 767f3fbd535SBarry Smith if (preonly) { 768251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr); 769251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr); 770251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr); 771251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)cpc,PCSVD,&svd);CHKERRQ(ierr); 772fd1303e9SJungho Lee if (!lu && !redundant && !cholesky && !svd) { 773f3fbd535SBarry Smith ierr = KSPSetType(mglevels[0]->smoothd,KSPGMRES);CHKERRQ(ierr); 774f3fbd535SBarry Smith } 775f3fbd535SBarry Smith } 776f3fbd535SBarry Smith 777f3fbd535SBarry Smith if (!pc->setupcalled) { 778f3fbd535SBarry Smith ierr = KSPSetFromOptions(mglevels[0]->smoothd);CHKERRQ(ierr); 779f3fbd535SBarry Smith } 780f3fbd535SBarry Smith 78163e6d426SJed Brown if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 782f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr); 78363e6d426SJed Brown if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 784f3fbd535SBarry Smith 785f3fbd535SBarry Smith /* 786f3fbd535SBarry Smith Dump the interpolation/restriction matrices plus the 787e3c5b3baSBarry Smith Jacobian/stiffness on each level. This allows MATLAB users to 788f3fbd535SBarry Smith easily check if the Galerkin condition A_c = R A_f R^T is satisfied. 789f3fbd535SBarry Smith 790f3fbd535SBarry Smith Only support one or the other at the same time. 791f3fbd535SBarry Smith */ 792f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER) 7930298fd71SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);CHKERRQ(ierr); 794ce94432eSBarry Smith if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc)); 795f3fbd535SBarry Smith dump = PETSC_FALSE; 796f3fbd535SBarry Smith #endif 7970298fd71SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);CHKERRQ(ierr); 798ce94432eSBarry Smith if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc)); 799f3fbd535SBarry Smith 800f3fbd535SBarry Smith if (viewer) { 801f3fbd535SBarry Smith for (i=1; i<n; i++) { 802f3fbd535SBarry Smith ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr); 803f3fbd535SBarry Smith } 804f3fbd535SBarry Smith for (i=0; i<n; i++) { 805f3fbd535SBarry Smith ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr); 806f3fbd535SBarry Smith ierr = MatView(pc->mat,viewer);CHKERRQ(ierr); 807f3fbd535SBarry Smith } 808f3fbd535SBarry Smith } 809f3fbd535SBarry Smith PetscFunctionReturn(0); 810f3fbd535SBarry Smith } 811f3fbd535SBarry Smith 812f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/ 813f3fbd535SBarry Smith 814f3fbd535SBarry Smith #undef __FUNCT__ 8159dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetLevels" 8164b9ad928SBarry Smith /*@ 81797177400SBarry Smith PCMGGetLevels - Gets the number of levels to use with MG. 8184b9ad928SBarry Smith 8194b9ad928SBarry Smith Not Collective 8204b9ad928SBarry Smith 8214b9ad928SBarry Smith Input Parameter: 8224b9ad928SBarry Smith . pc - the preconditioner context 8234b9ad928SBarry Smith 8244b9ad928SBarry Smith Output parameter: 8254b9ad928SBarry Smith . levels - the number of levels 8264b9ad928SBarry Smith 8274b9ad928SBarry Smith Level: advanced 8284b9ad928SBarry Smith 8294b9ad928SBarry Smith .keywords: MG, get, levels, multigrid 8304b9ad928SBarry Smith 83197177400SBarry Smith .seealso: PCMGSetLevels() 8324b9ad928SBarry Smith @*/ 8337087cfbeSBarry Smith PetscErrorCode PCMGGetLevels(PC pc,PetscInt *levels) 8344b9ad928SBarry Smith { 835f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 8364b9ad928SBarry Smith 8374b9ad928SBarry Smith PetscFunctionBegin; 8380700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 8394482741eSBarry Smith PetscValidIntPointer(levels,2); 840f3fbd535SBarry Smith *levels = mg->nlevels; 8414b9ad928SBarry Smith PetscFunctionReturn(0); 8424b9ad928SBarry Smith } 8434b9ad928SBarry Smith 8444b9ad928SBarry Smith #undef __FUNCT__ 8459dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetType" 8464b9ad928SBarry Smith /*@ 84797177400SBarry Smith PCMGSetType - Determines the form of multigrid to use: 8484b9ad928SBarry Smith multiplicative, additive, full, or the Kaskade algorithm. 8494b9ad928SBarry Smith 850ad4df100SBarry Smith Logically Collective on PC 8514b9ad928SBarry Smith 8524b9ad928SBarry Smith Input Parameters: 8534b9ad928SBarry Smith + pc - the preconditioner context 8549dcbbd2bSBarry Smith - form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE, 8559dcbbd2bSBarry Smith PC_MG_FULL, PC_MG_KASKADE 8564b9ad928SBarry Smith 8574b9ad928SBarry Smith Options Database Key: 8584b9ad928SBarry Smith . -pc_mg_type <form> - Sets <form>, one of multiplicative, 8594b9ad928SBarry Smith additive, full, kaskade 8604b9ad928SBarry Smith 8614b9ad928SBarry Smith Level: advanced 8624b9ad928SBarry Smith 8634b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid 8644b9ad928SBarry Smith 86597177400SBarry Smith .seealso: PCMGSetLevels() 8664b9ad928SBarry Smith @*/ 8677087cfbeSBarry Smith PetscErrorCode PCMGSetType(PC pc,PCMGType form) 8684b9ad928SBarry Smith { 869f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 8704b9ad928SBarry Smith 8714b9ad928SBarry Smith PetscFunctionBegin; 8720700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 873c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc,form,2); 87431567311SBarry Smith mg->am = form; 8759dcbbd2bSBarry Smith if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG; 8764b9ad928SBarry Smith else pc->ops->applyrichardson = 0; 8774b9ad928SBarry Smith PetscFunctionReturn(0); 8784b9ad928SBarry Smith } 8794b9ad928SBarry Smith 8804b9ad928SBarry Smith #undef __FUNCT__ 8810d353602SBarry Smith #define __FUNCT__ "PCMGSetCycleType" 8824b9ad928SBarry Smith /*@ 8830d353602SBarry Smith PCMGSetCycleType - Sets the type cycles to use. Use PCMGSetCycleTypeOnLevel() for more 8844b9ad928SBarry Smith complicated cycling. 8854b9ad928SBarry Smith 886ad4df100SBarry Smith Logically Collective on PC 8874b9ad928SBarry Smith 8884b9ad928SBarry Smith Input Parameters: 889c2be2410SBarry Smith + pc - the multigrid context 8900d353602SBarry Smith - PC_MG_CYCLE_V or PC_MG_CYCLE_W 8914b9ad928SBarry Smith 8924b9ad928SBarry Smith Options Database Key: 8930d353602SBarry Smith $ -pc_mg_cycle_type v or w 8944b9ad928SBarry Smith 8954b9ad928SBarry Smith Level: advanced 8964b9ad928SBarry Smith 8974b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 8984b9ad928SBarry Smith 8990d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel() 9004b9ad928SBarry Smith @*/ 9017087cfbeSBarry Smith PetscErrorCode PCMGSetCycleType(PC pc,PCMGCycleType n) 9024b9ad928SBarry Smith { 903f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 904f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 90579416396SBarry Smith PetscInt i,levels; 9064b9ad928SBarry Smith 9074b9ad928SBarry Smith PetscFunctionBegin; 9080700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 909ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 910c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 911f3fbd535SBarry Smith levels = mglevels[0]->levels; 9124b9ad928SBarry Smith 9132fa5cd67SKarl Rupp for (i=0; i<levels; i++) mglevels[i]->cycles = n; 9144b9ad928SBarry Smith PetscFunctionReturn(0); 9154b9ad928SBarry Smith } 9164b9ad928SBarry Smith 9174b9ad928SBarry Smith #undef __FUNCT__ 9188cc2d5dfSBarry Smith #define __FUNCT__ "PCMGMultiplicativeSetCycles" 9198cc2d5dfSBarry Smith /*@ 9208cc2d5dfSBarry Smith PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step 9218cc2d5dfSBarry Smith of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used 9228cc2d5dfSBarry Smith 923ad4df100SBarry Smith Logically Collective on PC 9248cc2d5dfSBarry Smith 9258cc2d5dfSBarry Smith Input Parameters: 9268cc2d5dfSBarry Smith + pc - the multigrid context 9278cc2d5dfSBarry Smith - n - number of cycles (default is 1) 9288cc2d5dfSBarry Smith 9298cc2d5dfSBarry Smith Options Database Key: 9308cc2d5dfSBarry Smith $ -pc_mg_multiplicative_cycles n 9318cc2d5dfSBarry Smith 9328cc2d5dfSBarry Smith Level: advanced 9338cc2d5dfSBarry Smith 9348cc2d5dfSBarry Smith Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType() 9358cc2d5dfSBarry Smith 9368cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 9378cc2d5dfSBarry Smith 9388cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType() 9398cc2d5dfSBarry Smith @*/ 9407087cfbeSBarry Smith PetscErrorCode PCMGMultiplicativeSetCycles(PC pc,PetscInt n) 9418cc2d5dfSBarry Smith { 942f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 943f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 9448cc2d5dfSBarry Smith PetscInt i,levels; 9458cc2d5dfSBarry Smith 9468cc2d5dfSBarry Smith PetscFunctionBegin; 9470700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 948ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 949c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 950f3fbd535SBarry Smith levels = mglevels[0]->levels; 9518cc2d5dfSBarry Smith 9522fa5cd67SKarl Rupp for (i=0; i<levels; i++) mg->cyclesperpcapply = n; 9538cc2d5dfSBarry Smith PetscFunctionReturn(0); 9548cc2d5dfSBarry Smith } 9558cc2d5dfSBarry Smith 9568cc2d5dfSBarry Smith #undef __FUNCT__ 9579dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetGalerkin" 958c2be2410SBarry Smith /*@ 95997177400SBarry Smith PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the 96082b87d87SMatthew G. Knepley finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i 961c2be2410SBarry Smith 962ad4df100SBarry Smith Logically Collective on PC 963c2be2410SBarry Smith 964c2be2410SBarry Smith Input Parameters: 965c91913d3SJed Brown + pc - the multigrid context 966c91913d3SJed Brown - use - PETSC_TRUE to use the Galerkin process to compute coarse-level operators 967c2be2410SBarry Smith 968c2be2410SBarry Smith Options Database Key: 969c2be2410SBarry Smith $ -pc_mg_galerkin 970c2be2410SBarry Smith 971c2be2410SBarry Smith Level: intermediate 972c2be2410SBarry Smith 973c2be2410SBarry Smith .keywords: MG, set, Galerkin 974c2be2410SBarry Smith 9753fc8bf9cSBarry Smith .seealso: PCMGGetGalerkin() 9763fc8bf9cSBarry Smith 977c2be2410SBarry Smith @*/ 978c91913d3SJed Brown PetscErrorCode PCMGSetGalerkin(PC pc,PetscBool use) 979c2be2410SBarry Smith { 980f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 981c2be2410SBarry Smith 982c2be2410SBarry Smith PetscFunctionBegin; 9830700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 984789726d7SBarry Smith mg->galerkin = use ? 1 : 0; 985c2be2410SBarry Smith PetscFunctionReturn(0); 986c2be2410SBarry Smith } 987c2be2410SBarry Smith 988c2be2410SBarry Smith #undef __FUNCT__ 9893fc8bf9cSBarry Smith #define __FUNCT__ "PCMGGetGalerkin" 9903fc8bf9cSBarry Smith /*@ 9913fc8bf9cSBarry Smith PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. 99282b87d87SMatthew G. Knepley A_i-1 = r_i * A_i * p_i 9933fc8bf9cSBarry Smith 9943fc8bf9cSBarry Smith Not Collective 9953fc8bf9cSBarry Smith 9963fc8bf9cSBarry Smith Input Parameter: 9973fc8bf9cSBarry Smith . pc - the multigrid context 9983fc8bf9cSBarry Smith 9993fc8bf9cSBarry Smith Output Parameter: 10003fc8bf9cSBarry Smith . gelerkin - PETSC_TRUE or PETSC_FALSE 10013fc8bf9cSBarry Smith 10023fc8bf9cSBarry Smith Options Database Key: 10033fc8bf9cSBarry Smith $ -pc_mg_galerkin 10043fc8bf9cSBarry Smith 10053fc8bf9cSBarry Smith Level: intermediate 10063fc8bf9cSBarry Smith 10073fc8bf9cSBarry Smith .keywords: MG, set, Galerkin 10083fc8bf9cSBarry Smith 10093fc8bf9cSBarry Smith .seealso: PCMGSetGalerkin() 10103fc8bf9cSBarry Smith 10113fc8bf9cSBarry Smith @*/ 10127087cfbeSBarry Smith PetscErrorCode PCMGGetGalerkin(PC pc,PetscBool *galerkin) 10133fc8bf9cSBarry Smith { 1014f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 10153fc8bf9cSBarry Smith 10163fc8bf9cSBarry Smith PetscFunctionBegin; 10170700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 101813fdb3acSJose Roman *galerkin = (PetscBool)mg->galerkin; 10193fc8bf9cSBarry Smith PetscFunctionReturn(0); 10203fc8bf9cSBarry Smith } 10213fc8bf9cSBarry Smith 10223fc8bf9cSBarry Smith #undef __FUNCT__ 10239dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothDown" 10244b9ad928SBarry Smith /*@ 102597177400SBarry Smith PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to 102697177400SBarry Smith use on all levels. Use PCMGGetSmootherDown() to set different 10274b9ad928SBarry Smith pre-smoothing steps on different levels. 10284b9ad928SBarry Smith 1029ad4df100SBarry Smith Logically Collective on PC 10304b9ad928SBarry Smith 10314b9ad928SBarry Smith Input Parameters: 10324b9ad928SBarry Smith + mg - the multigrid context 10334b9ad928SBarry Smith - n - the number of smoothing steps 10344b9ad928SBarry Smith 10354b9ad928SBarry Smith Options Database Key: 10364b9ad928SBarry Smith . -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps 10374b9ad928SBarry Smith 10384b9ad928SBarry Smith Level: advanced 10394b9ad928SBarry Smith 10404b9ad928SBarry Smith .keywords: MG, smooth, down, pre-smoothing, steps, multigrid 10414b9ad928SBarry Smith 104297177400SBarry Smith .seealso: PCMGSetNumberSmoothUp() 10434b9ad928SBarry Smith @*/ 10447087cfbeSBarry Smith PetscErrorCode PCMGSetNumberSmoothDown(PC pc,PetscInt n) 10454b9ad928SBarry Smith { 1046f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1047f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 10486849ba73SBarry Smith PetscErrorCode ierr; 104979416396SBarry Smith PetscInt i,levels; 10504b9ad928SBarry Smith 10514b9ad928SBarry Smith PetscFunctionBegin; 10520700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1053ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 1054c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 1055f3fbd535SBarry Smith levels = mglevels[0]->levels; 10564b9ad928SBarry Smith 1057b05257ddSBarry Smith for (i=1; i<levels; i++) { 10584b9ad928SBarry Smith /* make sure smoother up and down are different */ 10590298fd71SBarry Smith ierr = PCMGGetSmootherUp(pc,i,NULL);CHKERRQ(ierr); 1060f3fbd535SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 10612fa5cd67SKarl Rupp 106231567311SBarry Smith mg->default_smoothd = n; 10634b9ad928SBarry Smith } 10644b9ad928SBarry Smith PetscFunctionReturn(0); 10654b9ad928SBarry Smith } 10664b9ad928SBarry Smith 10674b9ad928SBarry Smith #undef __FUNCT__ 10689dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothUp" 10694b9ad928SBarry Smith /*@ 107097177400SBarry Smith PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use 107197177400SBarry Smith on all levels. Use PCMGGetSmootherUp() to set different numbers of 10724b9ad928SBarry Smith post-smoothing steps on different levels. 10734b9ad928SBarry Smith 1074ad4df100SBarry Smith Logically Collective on PC 10754b9ad928SBarry Smith 10764b9ad928SBarry Smith Input Parameters: 10774b9ad928SBarry Smith + mg - the multigrid context 10784b9ad928SBarry Smith - n - the number of smoothing steps 10794b9ad928SBarry Smith 10804b9ad928SBarry Smith Options Database Key: 10814b9ad928SBarry Smith . -pc_mg_smoothup <n> - Sets number of post-smoothing steps 10824b9ad928SBarry Smith 10834b9ad928SBarry Smith Level: advanced 10844b9ad928SBarry Smith 10854b9ad928SBarry Smith Note: this does not set a value on the coarsest grid, since we assume that 1086a8c7a070SBarry Smith there is no separate smooth up on the coarsest grid. 10874b9ad928SBarry Smith 10884b9ad928SBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid 10894b9ad928SBarry Smith 109097177400SBarry Smith .seealso: PCMGSetNumberSmoothDown() 10914b9ad928SBarry Smith @*/ 10927087cfbeSBarry Smith PetscErrorCode PCMGSetNumberSmoothUp(PC pc,PetscInt n) 10934b9ad928SBarry Smith { 1094f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1095f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 10966849ba73SBarry Smith PetscErrorCode ierr; 109779416396SBarry Smith PetscInt i,levels; 10984b9ad928SBarry Smith 10994b9ad928SBarry Smith PetscFunctionBegin; 11000700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1101ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 1102c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 1103f3fbd535SBarry Smith levels = mglevels[0]->levels; 11044b9ad928SBarry Smith 11054b9ad928SBarry Smith for (i=1; i<levels; i++) { 11064b9ad928SBarry Smith /* make sure smoother up and down are different */ 11070298fd71SBarry Smith ierr = PCMGGetSmootherUp(pc,i,NULL);CHKERRQ(ierr); 1108f3fbd535SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 11092fa5cd67SKarl Rupp 111031567311SBarry Smith mg->default_smoothu = n; 11114b9ad928SBarry Smith } 11124b9ad928SBarry Smith PetscFunctionReturn(0); 11134b9ad928SBarry Smith } 11144b9ad928SBarry Smith 11154b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/ 11164b9ad928SBarry Smith 11173b09bd56SBarry Smith /*MC 1118ccb205f8SBarry Smith PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional 11193b09bd56SBarry Smith information about the coarser grid matrices and restriction/interpolation operators. 11203b09bd56SBarry Smith 11213b09bd56SBarry Smith Options Database Keys: 11223b09bd56SBarry Smith + -pc_mg_levels <nlevels> - number of levels including finest 11230d353602SBarry Smith . -pc_mg_cycles v or w 112479416396SBarry Smith . -pc_mg_smoothup <n> - number of smoothing steps after interpolation 11253b09bd56SBarry Smith . -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator 11268c1c2452SJed Brown . -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default 11273b09bd56SBarry Smith . -pc_mg_log - log information about time spent on each level of the solver 11283b09bd56SBarry Smith . -pc_mg_monitor - print information on the multigrid convergence 11298cf6037eSBarry Smith . -pc_mg_galerkin - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R' 11308cf6037eSBarry Smith . -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1) 11318cf6037eSBarry Smith . -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 1132e3c5b3baSBarry Smith to the Socket viewer for reading from MATLAB. 11338cf6037eSBarry Smith - -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices 11348cf6037eSBarry Smith to the binary output file called binaryoutput 11353b09bd56SBarry Smith 113624c3aa18SBarry 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 11373b09bd56SBarry Smith 11388cf6037eSBarry Smith When run with a single level the smoother options are used on that level NOT the coarse grid solver options 11398cf6037eSBarry Smith 11403b09bd56SBarry Smith Level: intermediate 11413b09bd56SBarry Smith 11428f87f92bSBarry Smith Concepts: multigrid/multilevel 11433b09bd56SBarry Smith 11448cf6037eSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE 11450d353602SBarry Smith PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(), 114697177400SBarry Smith PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 114797177400SBarry Smith PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 11480d353602SBarry Smith PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 11493b09bd56SBarry Smith M*/ 11503b09bd56SBarry Smith 11514b9ad928SBarry Smith #undef __FUNCT__ 11524b9ad928SBarry Smith #define __FUNCT__ "PCCreate_MG" 11538cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc) 11544b9ad928SBarry Smith { 1155f3fbd535SBarry Smith PC_MG *mg; 1156f3fbd535SBarry Smith PetscErrorCode ierr; 1157f3fbd535SBarry Smith 11584b9ad928SBarry Smith PetscFunctionBegin; 1159b00a9115SJed Brown ierr = PetscNewLog(pc,&mg);CHKERRQ(ierr); 1160f3fbd535SBarry Smith pc->data = (void*)mg; 1161f3fbd535SBarry Smith mg->nlevels = -1; 1162f3fbd535SBarry Smith 116337a44384SMark Adams pc->useAmat = PETSC_TRUE; 116437a44384SMark Adams 11654b9ad928SBarry Smith pc->ops->apply = PCApply_MG; 11664b9ad928SBarry Smith pc->ops->setup = PCSetUp_MG; 1167a06653b4SBarry Smith pc->ops->reset = PCReset_MG; 11684b9ad928SBarry Smith pc->ops->destroy = PCDestroy_MG; 11694b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_MG; 11704b9ad928SBarry Smith pc->ops->view = PCView_MG; 11714b9ad928SBarry Smith PetscFunctionReturn(0); 11724b9ad928SBarry Smith } 1173