1dba47a55SKris Buschelman 24b9ad928SBarry Smith /* 34b9ad928SBarry Smith Defines the multigrid preconditioner interface. 44b9ad928SBarry Smith */ 5c6db04a5SJed Brown #include <../src/ksp/pc/impls/mg/mgimpl.h> /*I "petscpcmg.h" I*/ 64b9ad928SBarry Smith 74b9ad928SBarry Smith 84b9ad928SBarry Smith #undef __FUNCT__ 99dcbbd2bSBarry Smith #define __FUNCT__ "PCMGMCycle_Private" 1031567311SBarry Smith PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG_Levels **mglevelsin,PCRichardsonConvergedReason *reason) 114b9ad928SBarry Smith { 1231567311SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1331567311SBarry Smith PC_MG_Levels *mgc,*mglevels = *mglevelsin; 146849ba73SBarry Smith PetscErrorCode ierr; 1531567311SBarry Smith PetscInt cycles = (mglevels->level == 1) ? 1 : (PetscInt) mglevels->cycles; 164b9ad928SBarry Smith 174b9ad928SBarry Smith PetscFunctionBegin; 184b9ad928SBarry Smith 1963e6d426SJed Brown if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 2031567311SBarry Smith ierr = KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x);CHKERRQ(ierr); /* pre-smooth */ 2163e6d426SJed Brown if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 2231567311SBarry Smith if (mglevels->level) { /* not the coarsest grid */ 2363e6d426SJed Brown if (mglevels->eventresidual) {ierr = PetscLogEventBegin(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);} 2431567311SBarry Smith ierr = (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);CHKERRQ(ierr); 2563e6d426SJed Brown if (mglevels->eventresidual) {ierr = PetscLogEventEnd(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);} 264b9ad928SBarry Smith 274b9ad928SBarry Smith /* if on finest level and have convergence criteria set */ 2831567311SBarry Smith if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) { 294b9ad928SBarry Smith PetscReal rnorm; 3031567311SBarry Smith ierr = VecNorm(mglevels->r,NORM_2,&rnorm);CHKERRQ(ierr); 314b9ad928SBarry Smith if (rnorm <= mg->ttol) { 3270441072SBarry Smith if (rnorm < mg->abstol) { 334d0a8057SBarry Smith *reason = PCRICHARDSON_CONVERGED_ATOL; 341e2582c4SBarry Smith ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %G is less than absolute tolerance %G\n",rnorm,mg->abstol);CHKERRQ(ierr); 354b9ad928SBarry Smith } else { 364d0a8057SBarry Smith *reason = PCRICHARDSON_CONVERGED_RTOL; 371e2582c4SBarry 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); 384b9ad928SBarry Smith } 394b9ad928SBarry Smith PetscFunctionReturn(0); 404b9ad928SBarry Smith } 414b9ad928SBarry Smith } 424b9ad928SBarry Smith 4331567311SBarry Smith mgc = *(mglevelsin - 1); 4463e6d426SJed Brown if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 4531567311SBarry Smith ierr = MatRestrict(mglevels->restrct,mglevels->r,mgc->b);CHKERRQ(ierr); 4663e6d426SJed Brown if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 47efb30889SBarry Smith ierr = VecSet(mgc->x,0.0);CHKERRQ(ierr); 484b9ad928SBarry Smith while (cycles--) { 4931567311SBarry Smith ierr = PCMGMCycle_Private(pc,mglevelsin-1,reason);CHKERRQ(ierr); 504b9ad928SBarry Smith } 5163e6d426SJed Brown if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 5231567311SBarry Smith ierr = MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x);CHKERRQ(ierr); 5363e6d426SJed Brown if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 5463e6d426SJed Brown if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 5531567311SBarry Smith ierr = KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x);CHKERRQ(ierr); /* post smooth */ 5663e6d426SJed Brown if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 574b9ad928SBarry Smith } 584b9ad928SBarry Smith PetscFunctionReturn(0); 594b9ad928SBarry Smith } 604b9ad928SBarry Smith 614b9ad928SBarry Smith #undef __FUNCT__ 624b9ad928SBarry Smith #define __FUNCT__ "PCApplyRichardson_MG" 63ace3abfcSBarry 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) 644b9ad928SBarry Smith { 65f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 66f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 67dfbe8321SBarry Smith PetscErrorCode ierr; 68f3fbd535SBarry Smith PetscInt levels = mglevels[0]->levels,i; 694b9ad928SBarry Smith 704b9ad928SBarry Smith PetscFunctionBegin; 71f3fbd535SBarry Smith mglevels[levels-1]->b = b; 72f3fbd535SBarry Smith mglevels[levels-1]->x = x; 734b9ad928SBarry Smith 7431567311SBarry Smith mg->rtol = rtol; 7531567311SBarry Smith mg->abstol = abstol; 7631567311SBarry Smith mg->dtol = dtol; 774b9ad928SBarry Smith if (rtol) { 784b9ad928SBarry Smith /* compute initial residual norm for relative convergence test */ 794b9ad928SBarry Smith PetscReal rnorm; 807319c654SBarry Smith if (zeroguess) { 817319c654SBarry Smith ierr = VecNorm(b,NORM_2,&rnorm);CHKERRQ(ierr); 827319c654SBarry Smith } else { 83f3fbd535SBarry Smith ierr = (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);CHKERRQ(ierr); 844b9ad928SBarry Smith ierr = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr); 857319c654SBarry Smith } 8631567311SBarry Smith mg->ttol = PetscMax(rtol*rnorm,abstol); 8770441072SBarry Smith } else if (abstol) { 8831567311SBarry Smith mg->ttol = abstol; 894b9ad928SBarry Smith } else { 9031567311SBarry Smith mg->ttol = 0.0; 914b9ad928SBarry Smith } 924b9ad928SBarry Smith 934d0a8057SBarry Smith /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't 944d0a8057SBarry Smith stop prematurely do to small residual */ 954d0a8057SBarry Smith for (i=1; i<levels; i++) { 96f3fbd535SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 97f3fbd535SBarry Smith if (mglevels[i]->smoothu != mglevels[i]->smoothd) { 98f3fbd535SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 994b9ad928SBarry Smith } 1004d0a8057SBarry Smith } 1014d0a8057SBarry Smith 1024d0a8057SBarry Smith *reason = (PCRichardsonConvergedReason)0; 1034d0a8057SBarry Smith for (i=0; i<its; i++) { 104f3fbd535SBarry Smith ierr = PCMGMCycle_Private(pc,mglevels+levels-1,reason);CHKERRQ(ierr); 1054d0a8057SBarry Smith if (*reason) break; 1064d0a8057SBarry Smith } 1074d0a8057SBarry Smith if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS; 1084d0a8057SBarry Smith *outits = i; 1094b9ad928SBarry Smith PetscFunctionReturn(0); 1104b9ad928SBarry Smith } 1114b9ad928SBarry Smith 1124b9ad928SBarry Smith #undef __FUNCT__ 1133aeaf226SBarry Smith #define __FUNCT__ "PCReset_MG" 1143aeaf226SBarry Smith PetscErrorCode PCReset_MG(PC pc) 1153aeaf226SBarry Smith { 1163aeaf226SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1173aeaf226SBarry Smith PC_MG_Levels **mglevels = mg->levels; 1183aeaf226SBarry Smith PetscErrorCode ierr; 1193aeaf226SBarry Smith PetscInt i,n; 1203aeaf226SBarry Smith 1213aeaf226SBarry Smith PetscFunctionBegin; 1223aeaf226SBarry Smith if (mglevels) { 1233aeaf226SBarry Smith n = mglevels[0]->levels; 1243aeaf226SBarry Smith for (i=0; i<n-1; i++) { 1253aeaf226SBarry Smith ierr = VecDestroy(&mglevels[i+1]->r);CHKERRQ(ierr); 1263aeaf226SBarry Smith ierr = VecDestroy(&mglevels[i]->b);CHKERRQ(ierr); 1273aeaf226SBarry Smith ierr = VecDestroy(&mglevels[i]->x);CHKERRQ(ierr); 1283aeaf226SBarry Smith ierr = MatDestroy(&mglevels[i+1]->restrct);CHKERRQ(ierr); 1293aeaf226SBarry Smith ierr = MatDestroy(&mglevels[i+1]->interpolate);CHKERRQ(ierr); 1303aeaf226SBarry Smith } 1313aeaf226SBarry Smith 1323aeaf226SBarry Smith for (i=0; i<n; i++) { 1333aeaf226SBarry Smith ierr = MatDestroy(&mglevels[i]->A);CHKERRQ(ierr); 1343aeaf226SBarry Smith if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 1353aeaf226SBarry Smith ierr = KSPReset(mglevels[i]->smoothd);CHKERRQ(ierr); 1363aeaf226SBarry Smith } 1373aeaf226SBarry Smith ierr = KSPReset(mglevels[i]->smoothu);CHKERRQ(ierr); 1383aeaf226SBarry Smith } 1393aeaf226SBarry Smith } 1403aeaf226SBarry Smith PetscFunctionReturn(0); 1413aeaf226SBarry Smith } 1423aeaf226SBarry Smith 1433aeaf226SBarry Smith #undef __FUNCT__ 1449dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetLevels" 1454b9ad928SBarry Smith /*@C 14697177400SBarry Smith PCMGSetLevels - Sets the number of levels to use with MG. 1474b9ad928SBarry Smith Must be called before any other MG routine. 1484b9ad928SBarry Smith 149ad4df100SBarry Smith Logically Collective on PC 1504b9ad928SBarry Smith 1514b9ad928SBarry Smith Input Parameters: 1524b9ad928SBarry Smith + pc - the preconditioner context 1534b9ad928SBarry Smith . levels - the number of levels 1544b9ad928SBarry Smith - comms - optional communicators for each level; this is to allow solving the coarser problems 1554b9ad928SBarry Smith on smaller sets of processors. Use PETSC_NULL_OBJECT for default in Fortran 1564b9ad928SBarry Smith 1574b9ad928SBarry Smith Level: intermediate 1584b9ad928SBarry Smith 1594b9ad928SBarry Smith Notes: 1604b9ad928SBarry Smith If the number of levels is one then the multigrid uses the -mg_levels prefix 1614b9ad928SBarry Smith for setting the level options rather than the -mg_coarse prefix. 1624b9ad928SBarry Smith 1634b9ad928SBarry Smith .keywords: MG, set, levels, multigrid 1644b9ad928SBarry Smith 16597177400SBarry Smith .seealso: PCMGSetType(), PCMGGetLevels() 1664b9ad928SBarry Smith @*/ 1677087cfbeSBarry Smith PetscErrorCode PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms) 1684b9ad928SBarry Smith { 169dfbe8321SBarry Smith PetscErrorCode ierr; 170f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 171f3fbd535SBarry Smith MPI_Comm comm = ((PetscObject)pc)->comm; 1723aeaf226SBarry Smith PC_MG_Levels **mglevels = mg->levels; 173f3fbd535SBarry Smith PetscInt i; 174f3fbd535SBarry Smith PetscMPIInt size; 175f3fbd535SBarry Smith const char *prefix; 176f3fbd535SBarry Smith PC ipc; 1773aeaf226SBarry Smith PetscInt n; 1784b9ad928SBarry Smith 1794b9ad928SBarry Smith PetscFunctionBegin; 1800700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 181548767e0SJed Brown PetscValidLogicalCollectiveInt(pc,levels,2); 182548767e0SJed Brown if (mg->nlevels == levels) PetscFunctionReturn(0); 1833aeaf226SBarry Smith if (mglevels) { 1843aeaf226SBarry Smith /* changing the number of levels so free up the previous stuff */ 1853aeaf226SBarry Smith ierr = PCReset_MG(pc);CHKERRQ(ierr); 1863aeaf226SBarry Smith n = mglevels[0]->levels; 1873aeaf226SBarry Smith for (i=0; i<n; i++) { 1883aeaf226SBarry Smith if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 1893aeaf226SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr); 1903aeaf226SBarry Smith } 1913aeaf226SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr); 1923aeaf226SBarry Smith ierr = PetscFree(mglevels[i]);CHKERRQ(ierr); 1933aeaf226SBarry Smith } 1943aeaf226SBarry Smith ierr = PetscFree(mg->levels);CHKERRQ(ierr); 1953aeaf226SBarry Smith } 196f3fbd535SBarry Smith 197f3fbd535SBarry Smith mg->nlevels = levels; 198218a07d4SBarry Smith mg->galerkin = PETSC_FALSE; 199218a07d4SBarry Smith mg->galerkinused = PETSC_FALSE; 200f3fbd535SBarry Smith 201f3fbd535SBarry Smith ierr = PetscMalloc(levels*sizeof(PC_MG*),&mglevels);CHKERRQ(ierr); 202f3fbd535SBarry Smith ierr = PetscLogObjectMemory(pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr); 203f3fbd535SBarry Smith 204f3fbd535SBarry Smith ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 205f3fbd535SBarry Smith 206f3fbd535SBarry Smith for (i=0; i<levels; i++) { 207f3fbd535SBarry Smith ierr = PetscNewLog(pc,PC_MG_Levels,&mglevels[i]);CHKERRQ(ierr); 208f3fbd535SBarry Smith mglevels[i]->level = i; 209f3fbd535SBarry Smith mglevels[i]->levels = levels; 210f3fbd535SBarry Smith mglevels[i]->cycles = PC_MG_CYCLE_V; 21131567311SBarry Smith mg->default_smoothu = 1; 21231567311SBarry Smith mg->default_smoothd = 1; 21363e6d426SJed Brown mglevels[i]->eventsmoothsetup = 0; 21463e6d426SJed Brown mglevels[i]->eventsmoothsolve = 0; 21563e6d426SJed Brown mglevels[i]->eventresidual = 0; 21663e6d426SJed Brown mglevels[i]->eventinterprestrict = 0; 217f3fbd535SBarry Smith 218f3fbd535SBarry Smith if (comms) comm = comms[i]; 219f3fbd535SBarry Smith ierr = KSPCreate(comm,&mglevels[i]->smoothd);CHKERRQ(ierr); 220f3fbd535SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr); 22131567311SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg->default_smoothd);CHKERRQ(ierr); 222f3fbd535SBarry Smith ierr = KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);CHKERRQ(ierr); 223f3fbd535SBarry Smith 224f3fbd535SBarry Smith /* do special stuff for coarse grid */ 225f3fbd535SBarry Smith if (!i && levels > 1) { 226f3fbd535SBarry Smith ierr = KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");CHKERRQ(ierr); 227f3fbd535SBarry Smith 228f3fbd535SBarry Smith /* coarse solve is (redundant) LU by default */ 229f3fbd535SBarry Smith ierr = KSPSetType(mglevels[0]->smoothd,KSPPREONLY);CHKERRQ(ierr); 230f3fbd535SBarry Smith ierr = KSPGetPC(mglevels[0]->smoothd,&ipc);CHKERRQ(ierr); 231f3fbd535SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 232f3fbd535SBarry Smith if (size > 1) { 233f3fbd535SBarry Smith ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr); 234f3fbd535SBarry Smith } else { 235f3fbd535SBarry Smith ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr); 236f3fbd535SBarry Smith } 237f3fbd535SBarry Smith 238f3fbd535SBarry Smith } else { 239f3fbd535SBarry Smith char tprefix[128]; 240f3fbd535SBarry Smith sprintf(tprefix,"mg_levels_%d_",(int)i); 241f3fbd535SBarry Smith ierr = KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);CHKERRQ(ierr); 242f3fbd535SBarry Smith } 243f3fbd535SBarry Smith ierr = PetscLogObjectParent(pc,mglevels[i]->smoothd);CHKERRQ(ierr); 244f3fbd535SBarry Smith mglevels[i]->smoothu = mglevels[i]->smoothd; 24531567311SBarry Smith mg->rtol = 0.0; 24631567311SBarry Smith mg->abstol = 0.0; 24731567311SBarry Smith mg->dtol = 0.0; 24831567311SBarry Smith mg->ttol = 0.0; 24931567311SBarry Smith mg->cyclesperpcapply = 1; 250f3fbd535SBarry Smith } 25131567311SBarry Smith mg->am = PC_MG_MULTIPLICATIVE; 252f3fbd535SBarry Smith mg->levels = mglevels; 2534b9ad928SBarry Smith pc->ops->applyrichardson = PCApplyRichardson_MG; 2544b9ad928SBarry Smith PetscFunctionReturn(0); 2554b9ad928SBarry Smith } 2564b9ad928SBarry Smith 257c07bf074SBarry Smith 258c07bf074SBarry Smith #undef __FUNCT__ 259c07bf074SBarry Smith #define __FUNCT__ "PCDestroy_MG" 260c07bf074SBarry Smith PetscErrorCode PCDestroy_MG(PC pc) 261c07bf074SBarry Smith { 262c07bf074SBarry Smith PetscErrorCode ierr; 263a06653b4SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 264a06653b4SBarry Smith PC_MG_Levels **mglevels = mg->levels; 265a06653b4SBarry Smith PetscInt i,n; 266c07bf074SBarry Smith 267c07bf074SBarry Smith PetscFunctionBegin; 268a06653b4SBarry Smith ierr = PCReset_MG(pc);CHKERRQ(ierr); 269a06653b4SBarry Smith if (mglevels) { 270a06653b4SBarry Smith n = mglevels[0]->levels; 271a06653b4SBarry Smith for (i=0; i<n; i++) { 272a06653b4SBarry Smith if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 2736bf464f9SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr); 274a06653b4SBarry Smith } 2756bf464f9SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr); 276a06653b4SBarry Smith ierr = PetscFree(mglevels[i]);CHKERRQ(ierr); 277a06653b4SBarry Smith } 278a06653b4SBarry Smith ierr = PetscFree(mg->levels);CHKERRQ(ierr); 279a06653b4SBarry Smith } 280c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 281f3fbd535SBarry Smith PetscFunctionReturn(0); 282f3fbd535SBarry Smith } 283f3fbd535SBarry Smith 284f3fbd535SBarry Smith 285f3fbd535SBarry Smith 28609573ac7SBarry Smith extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**); 28709573ac7SBarry Smith extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**); 28809573ac7SBarry Smith extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**); 289f3fbd535SBarry Smith 290f3fbd535SBarry Smith /* 291f3fbd535SBarry Smith PCApply_MG - Runs either an additive, multiplicative, Kaskadic 292f3fbd535SBarry Smith or full cycle of multigrid. 293f3fbd535SBarry Smith 294f3fbd535SBarry Smith Note: 295f3fbd535SBarry Smith A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle(). 296f3fbd535SBarry Smith */ 297f3fbd535SBarry Smith #undef __FUNCT__ 298f3fbd535SBarry Smith #define __FUNCT__ "PCApply_MG" 299f3fbd535SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x) 300f3fbd535SBarry Smith { 301f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 302f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 303f3fbd535SBarry Smith PetscErrorCode ierr; 304f3fbd535SBarry Smith PetscInt levels = mglevels[0]->levels,i; 305f3fbd535SBarry Smith 306f3fbd535SBarry Smith PetscFunctionBegin; 307e1d8e5deSBarry Smith 308e1d8e5deSBarry Smith /* When the DM is supplying the matrix then it will not exist until here */ 309298cc208SBarry Smith for (i=0; i<levels; i++) { 310e1d8e5deSBarry Smith if (!mglevels[i]->A) { 311e1d8e5deSBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 312298cc208SBarry Smith ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr); 313e1d8e5deSBarry Smith } 314e1d8e5deSBarry Smith } 315e1d8e5deSBarry Smith 316f3fbd535SBarry Smith mglevels[levels-1]->b = b; 317f3fbd535SBarry Smith mglevels[levels-1]->x = x; 31831567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 319f3fbd535SBarry Smith ierr = VecSet(x,0.0);CHKERRQ(ierr); 32031567311SBarry Smith for (i=0; i<mg->cyclesperpcapply; i++) { 321f3fbd535SBarry Smith ierr = PCMGMCycle_Private(pc,mglevels+levels-1,PETSC_NULL);CHKERRQ(ierr); 322f3fbd535SBarry Smith } 323f3fbd535SBarry Smith } 32431567311SBarry Smith else if (mg->am == PC_MG_ADDITIVE) { 32531567311SBarry Smith ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr); 326f3fbd535SBarry Smith } 32731567311SBarry Smith else if (mg->am == PC_MG_KASKADE) { 32831567311SBarry Smith ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr); 329f3fbd535SBarry Smith } 330f3fbd535SBarry Smith else { 331f3fbd535SBarry Smith ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr); 332f3fbd535SBarry Smith } 333f3fbd535SBarry Smith PetscFunctionReturn(0); 334f3fbd535SBarry Smith } 335f3fbd535SBarry Smith 336f3fbd535SBarry Smith 337f3fbd535SBarry Smith #undef __FUNCT__ 338f3fbd535SBarry Smith #define __FUNCT__ "PCSetFromOptions_MG" 339f3fbd535SBarry Smith PetscErrorCode PCSetFromOptions_MG(PC pc) 340f3fbd535SBarry Smith { 341f3fbd535SBarry Smith PetscErrorCode ierr; 342f3fbd535SBarry Smith PetscInt m,levels = 1,cycles; 343ace3abfcSBarry Smith PetscBool flg; 344f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 345f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 346f3fbd535SBarry Smith PCMGType mgtype; 347f3fbd535SBarry Smith PCMGCycleType mgctype; 348f3fbd535SBarry Smith 349f3fbd535SBarry Smith PetscFunctionBegin; 350f3fbd535SBarry Smith ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr); 3513aeaf226SBarry Smith if (!mg->levels) { 352f3fbd535SBarry Smith ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr); 353eb3f98d2SBarry Smith if (!flg && pc->dm) { 354eb3f98d2SBarry Smith ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr); 355eb3f98d2SBarry Smith levels++; 3563aeaf226SBarry Smith mg->usedmfornumberoflevels = PETSC_TRUE; 357eb3f98d2SBarry Smith } 358f3fbd535SBarry Smith ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr); 359f3fbd535SBarry Smith } 3603aeaf226SBarry Smith mglevels = mg->levels; 3613aeaf226SBarry Smith 362f3fbd535SBarry Smith mgctype = (PCMGCycleType) mglevels[0]->cycles; 363f3fbd535SBarry Smith ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr); 364f3fbd535SBarry Smith if (flg) { 365f3fbd535SBarry Smith ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr); 366f3fbd535SBarry Smith }; 367f3fbd535SBarry Smith flg = PETSC_FALSE; 368acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 369f3fbd535SBarry Smith if (flg) { 370f3fbd535SBarry Smith ierr = PCMGSetGalerkin(pc);CHKERRQ(ierr); 371f3fbd535SBarry Smith } 372f3fbd535SBarry Smith ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr); 373f3fbd535SBarry Smith if (flg) { 374f3fbd535SBarry Smith ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr); 375f3fbd535SBarry Smith } 376f3fbd535SBarry Smith ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr); 377f3fbd535SBarry Smith if (flg) { 378f3fbd535SBarry Smith ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr); 379f3fbd535SBarry Smith } 38031567311SBarry Smith mgtype = mg->am; 381f3fbd535SBarry Smith ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr); 382f3fbd535SBarry Smith if (flg) { 383f3fbd535SBarry Smith ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr); 384f3fbd535SBarry Smith } 38531567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 38631567311SBarry Smith ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGSetLevels",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr); 387f3fbd535SBarry Smith if (flg) { 388f3fbd535SBarry Smith ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr); 389f3fbd535SBarry Smith } 390f3fbd535SBarry Smith } 391f3fbd535SBarry Smith flg = PETSC_FALSE; 392acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 393f3fbd535SBarry Smith if (flg) { 394f3fbd535SBarry Smith PetscInt i; 395f3fbd535SBarry Smith char eventname[128]; 39663e6d426SJed Brown if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 397f3fbd535SBarry Smith levels = mglevels[0]->levels; 398f3fbd535SBarry Smith for (i=0; i<levels; i++) { 399f3fbd535SBarry Smith sprintf(eventname,"MGSetup Level %d",(int)i); 40063e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr); 401f3fbd535SBarry Smith sprintf(eventname,"MGSmooth Level %d",(int)i); 40263e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr); 403f3fbd535SBarry Smith if (i) { 404f3fbd535SBarry Smith sprintf(eventname,"MGResid Level %d",(int)i); 40563e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr); 406f3fbd535SBarry Smith sprintf(eventname,"MGInterp Level %d",(int)i); 40763e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr); 408f3fbd535SBarry Smith } 409f3fbd535SBarry Smith } 410f3fbd535SBarry Smith } 411f3fbd535SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 412f3fbd535SBarry Smith PetscFunctionReturn(0); 413f3fbd535SBarry Smith } 414f3fbd535SBarry Smith 415f3fbd535SBarry Smith const char *PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0}; 416f3fbd535SBarry Smith const char *PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0}; 417f3fbd535SBarry Smith 418f3fbd535SBarry Smith #undef __FUNCT__ 419f3fbd535SBarry Smith #define __FUNCT__ "PCView_MG" 420f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer) 421f3fbd535SBarry Smith { 422f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 423f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 424f3fbd535SBarry Smith PetscErrorCode ierr; 425f3fbd535SBarry Smith PetscInt levels = mglevels[0]->levels,i; 426ace3abfcSBarry Smith PetscBool iascii; 427f3fbd535SBarry Smith 428f3fbd535SBarry Smith PetscFunctionBegin; 4292692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 430f3fbd535SBarry Smith if (iascii) { 43131567311SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," MG: type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,(mglevels[0]->cycles == PC_MG_CYCLE_V) ? "v" : "w");CHKERRQ(ierr); 43231567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 43331567311SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr); 434f3fbd535SBarry Smith } 435218a07d4SBarry Smith if (mg->galerkin) { 436f3fbd535SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 4374f66f45eSBarry Smith } else { 4384f66f45eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 439f3fbd535SBarry Smith } 440f3fbd535SBarry Smith for (i=0; i<levels; i++) { 441f3fbd535SBarry Smith if (!i) { 4422d19a89fSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %D smooths=%D --------------------\n",i,mg->default_smoothd);CHKERRQ(ierr); 443f3fbd535SBarry Smith } else { 44431567311SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D smooths=%D --------------------\n",i,mg->default_smoothd);CHKERRQ(ierr); 445f3fbd535SBarry Smith } 446f3fbd535SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 447f3fbd535SBarry Smith ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 448f3fbd535SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 449f3fbd535SBarry Smith if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) { 450f3fbd535SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 451f3fbd535SBarry Smith } else if (i){ 45231567311SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D smooths=%D --------------------\n",i,mg->default_smoothu);CHKERRQ(ierr); 453f3fbd535SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 454f3fbd535SBarry Smith ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 455f3fbd535SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 456f3fbd535SBarry Smith } 457f3fbd535SBarry Smith } 458f3fbd535SBarry Smith } else { 45965e19b50SBarry Smith SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCMG",((PetscObject)viewer)->type_name); 460f3fbd535SBarry Smith } 461f3fbd535SBarry Smith PetscFunctionReturn(0); 462f3fbd535SBarry Smith } 463f3fbd535SBarry Smith 464*cab2e9ccSBarry Smith #include <private/dmimpl.h> 465*cab2e9ccSBarry Smith #include <private/kspimpl.h> 466*cab2e9ccSBarry Smith 467f3fbd535SBarry Smith /* 468f3fbd535SBarry Smith Calls setup for the KSP on each level 469f3fbd535SBarry Smith */ 470f3fbd535SBarry Smith #undef __FUNCT__ 471f3fbd535SBarry Smith #define __FUNCT__ "PCSetUp_MG" 472f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc) 473f3fbd535SBarry Smith { 474f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 475f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 476f3fbd535SBarry Smith PetscErrorCode ierr; 477f3fbd535SBarry Smith PetscInt i,n = mglevels[0]->levels; 478f3fbd535SBarry Smith PC cpc,mpc; 479fd1303e9SJungho Lee PetscBool preonly,lu,redundant,cholesky,svd,monitor = PETSC_FALSE,dump = PETSC_FALSE,opsset; 480f3fbd535SBarry Smith PetscViewerASCIIMonitor ascii; 481f3fbd535SBarry Smith PetscViewer viewer = PETSC_NULL; 482f3fbd535SBarry Smith MPI_Comm comm; 483f3fbd535SBarry Smith Mat dA,dB; 484f3fbd535SBarry Smith MatStructure uflag; 485f3fbd535SBarry Smith Vec tvec; 486218a07d4SBarry Smith DM *dms; 487f3fbd535SBarry Smith 488f3fbd535SBarry Smith PetscFunctionBegin; 4893aeaf226SBarry Smith if (mg->usedmfornumberoflevels) { 4903aeaf226SBarry Smith PetscInt levels; 4913aeaf226SBarry Smith ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr); 4923aeaf226SBarry Smith levels++; 4933aeaf226SBarry Smith if (levels > n) { /* the problem is now being solved on a finer grid */ 4943aeaf226SBarry Smith ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr); 4953aeaf226SBarry Smith n = levels; 4963aeaf226SBarry Smith ierr = PCSetFromOptions(pc);CHKERRQ(ierr); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */ 4973aeaf226SBarry Smith mglevels = mg->levels; 4983aeaf226SBarry Smith } 4993aeaf226SBarry Smith } 5003aeaf226SBarry Smith 501f3fbd535SBarry Smith 502f3fbd535SBarry Smith /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */ 503f3fbd535SBarry Smith /* so use those from global PC */ 504f3fbd535SBarry Smith /* Is this what we always want? What if user wants to keep old one? */ 505f3fbd535SBarry Smith ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,PETSC_NULL,&opsset);CHKERRQ(ierr); 506f3fbd535SBarry Smith ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr); 507f3fbd535SBarry Smith ierr = KSPGetPC(mglevels[n-1]->smoothd,&mpc);CHKERRQ(ierr); 5081338a6b9SJed Brown if (!opsset || ((cpc->setupcalled == 1) && (mpc->setupcalled == 2)) || ((mpc == cpc) && (mpc->setupcalled == 2))) { 509f3fbd535SBarry 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); 510f3fbd535SBarry Smith ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr); 511f3fbd535SBarry Smith } 512f3fbd535SBarry Smith 5132d2b81a6SBarry Smith if (pc->dm && !pc->setupcalled) { 5142d2b81a6SBarry Smith /* construct the interpolation from the DMs */ 5152e499ae9SBarry Smith Mat p; 516218a07d4SBarry Smith ierr = PetscMalloc(n*sizeof(DM),&dms);CHKERRQ(ierr); 517218a07d4SBarry Smith dms[n-1] = pc->dm; 518218a07d4SBarry Smith for (i=n-2; i>-1; i--) { 519218a07d4SBarry Smith ierr = DMCoarsen(dms[i+1],PETSC_NULL,&dms[i]);CHKERRQ(ierr); 5203c0c59f3SBarry Smith ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr); 5213c0c59f3SBarry Smith if (mg->galerkin) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);} 52211629dbeSBarry Smith ierr = DMSetFunction(dms[i],0); 52311629dbeSBarry Smith ierr = DMSetInitialGuess(dms[i],0); 52424c3aa18SBarry Smith if (!mglevels[i+1]->interpolate) { 5252d2b81a6SBarry Smith ierr = DMGetInterpolation(dms[i],dms[i+1],&p,PETSC_NULL);CHKERRQ(ierr); 5262d2b81a6SBarry Smith ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr); 5276bf464f9SBarry Smith ierr = MatDestroy(&p);CHKERRQ(ierr); 528218a07d4SBarry Smith } 52924c3aa18SBarry Smith } 5302d2b81a6SBarry Smith 531218a07d4SBarry Smith for (i=n-2; i>-1; i--) { 5326bf464f9SBarry Smith ierr = DMDestroy(&dms[i]);CHKERRQ(ierr); 533218a07d4SBarry Smith } 5342d2b81a6SBarry Smith ierr = PetscFree(dms);CHKERRQ(ierr); 535*cab2e9ccSBarry Smith 536*cab2e9ccSBarry Smith /* finest smoother also gets DM but it is not active */ 537*cab2e9ccSBarry Smith ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr); 538*cab2e9ccSBarry Smith ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr); 539218a07d4SBarry Smith } 540218a07d4SBarry Smith 541218a07d4SBarry Smith if (mg->galerkin) { 542f3fbd535SBarry Smith Mat B; 543218a07d4SBarry Smith mg->galerkinused = PETSC_TRUE; 544f3fbd535SBarry Smith /* currently only handle case where mat and pmat are the same on coarser levels */ 545f3fbd535SBarry Smith ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr); 546f3fbd535SBarry Smith if (!pc->setupcalled) { 547f3fbd535SBarry Smith for (i=n-2; i>-1; i--) { 548f3fbd535SBarry Smith ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 549f3fbd535SBarry Smith ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 550f3fbd535SBarry Smith if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);} 551f3fbd535SBarry Smith dB = B; 552f3fbd535SBarry Smith } 553cd9507b2SBarry Smith if (n > 1) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);} 554f3fbd535SBarry Smith } else { 555f3fbd535SBarry Smith for (i=n-2; i>-1; i--) { 556f3fbd535SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr); 557f3fbd535SBarry Smith ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 558f3fbd535SBarry Smith ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 559f3fbd535SBarry Smith dB = B; 560f3fbd535SBarry Smith } 561f3fbd535SBarry Smith } 562*cab2e9ccSBarry Smith } else if (pc->dm && pc->dm->x) { 563*cab2e9ccSBarry Smith /* need to restrict Jacobian location to coarser meshes for evaluation */ 564*cab2e9ccSBarry Smith for (i=n-2;i>-1; i--) { 565*cab2e9ccSBarry Smith if (!mglevels[i]->smoothd->dm->x) { 566*cab2e9ccSBarry Smith Vec *vecs; 567*cab2e9ccSBarry Smith ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vecs,0,PETSC_NULL);CHKERRQ(ierr); 568*cab2e9ccSBarry Smith mglevels[i]->smoothd->dm->x = vecs[0]; 569*cab2e9ccSBarry Smith ierr = PetscFree(vecs);CHKERRQ(ierr); 570*cab2e9ccSBarry Smith } 571*cab2e9ccSBarry Smith ierr = MatRestrict(mglevels[i+1]->interpolate,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);CHKERRQ(ierr); 572*cab2e9ccSBarry Smith } 573f3fbd535SBarry Smith } 574f3fbd535SBarry Smith 575f3fbd535SBarry Smith if (!pc->setupcalled) { 576acfcf0e5SJed Brown ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_monitor",&monitor,PETSC_NULL);CHKERRQ(ierr); 577f3fbd535SBarry Smith 578f3fbd535SBarry Smith for (i=0; i<n; i++) { 579f3fbd535SBarry Smith if (monitor) { 580f3fbd535SBarry Smith ierr = PetscObjectGetComm((PetscObject)mglevels[i]->smoothd,&comm);CHKERRQ(ierr); 581f3fbd535SBarry Smith ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr); 582c2efdce3SBarry Smith ierr = KSPMonitorSet(mglevels[i]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void**))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 583f3fbd535SBarry Smith } 584f3fbd535SBarry Smith ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr); 585f3fbd535SBarry Smith } 586f3fbd535SBarry Smith for (i=1; i<n; i++) { 587f3fbd535SBarry Smith if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) { 588f3fbd535SBarry Smith if (monitor) { 589f3fbd535SBarry Smith ierr = PetscObjectGetComm((PetscObject)mglevels[i]->smoothu,&comm);CHKERRQ(ierr); 590f3fbd535SBarry Smith ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr); 591c2efdce3SBarry Smith ierr = KSPMonitorSet(mglevels[i]->smoothu,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void**))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 592f3fbd535SBarry Smith } 593f3fbd535SBarry Smith ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr); 594f3fbd535SBarry Smith } 595f3fbd535SBarry Smith } 596f3fbd535SBarry Smith for (i=1; i<n; i++) { 597f3fbd535SBarry Smith if (mglevels[i]->restrct && !mglevels[i]->interpolate) { 598f3fbd535SBarry Smith ierr = PCMGSetInterpolation(pc,i,mglevels[i]->restrct);CHKERRQ(ierr); 599f3fbd535SBarry Smith } 600f3fbd535SBarry Smith if (!mglevels[i]->restrct && mglevels[i]->interpolate) { 601f3fbd535SBarry Smith ierr = PCMGSetRestriction(pc,i,mglevels[i]->interpolate);CHKERRQ(ierr); 602f3fbd535SBarry Smith } 603f3fbd535SBarry Smith #if defined(PETSC_USE_DEBUG) 604f3fbd535SBarry Smith if (!mglevels[i]->restrct || !mglevels[i]->interpolate) { 60565e19b50SBarry Smith SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to set restriction or interpolation on level %d",(int)i); 606f3fbd535SBarry Smith } 607f3fbd535SBarry Smith #endif 608f3fbd535SBarry Smith } 609f3fbd535SBarry Smith for (i=0; i<n-1; i++) { 610f3fbd535SBarry Smith if (!mglevels[i]->b) { 611f3fbd535SBarry Smith Vec *vec; 612f3fbd535SBarry Smith ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr); 613f3fbd535SBarry Smith ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr); 6146bf464f9SBarry Smith ierr = VecDestroy(vec);CHKERRQ(ierr); 615f3fbd535SBarry Smith ierr = PetscFree(vec);CHKERRQ(ierr); 616f3fbd535SBarry Smith } 617f3fbd535SBarry Smith if (!mglevels[i]->r && i) { 618f3fbd535SBarry Smith ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 619f3fbd535SBarry Smith ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr); 6206bf464f9SBarry Smith ierr = VecDestroy(&tvec);CHKERRQ(ierr); 621f3fbd535SBarry Smith } 622f3fbd535SBarry Smith if (!mglevels[i]->x) { 623f3fbd535SBarry Smith ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 624f3fbd535SBarry Smith ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr); 6256bf464f9SBarry Smith ierr = VecDestroy(&tvec);CHKERRQ(ierr); 626f3fbd535SBarry Smith } 627f3fbd535SBarry Smith } 628f3fbd535SBarry Smith if (n != 1 && !mglevels[n-1]->r) { 629f3fbd535SBarry Smith /* PCMGSetR() on the finest level if user did not supply it */ 630f3fbd535SBarry Smith Vec *vec; 631f3fbd535SBarry Smith ierr = KSPGetVecs(mglevels[n-1]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr); 632f3fbd535SBarry Smith ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr); 6336bf464f9SBarry Smith ierr = VecDestroy(vec);CHKERRQ(ierr); 634f3fbd535SBarry Smith ierr = PetscFree(vec);CHKERRQ(ierr); 635f3fbd535SBarry Smith } 636f3fbd535SBarry Smith } 637f3fbd535SBarry Smith 638f3fbd535SBarry Smith 639f3fbd535SBarry Smith for (i=1; i<n; i++) { 640f3fbd535SBarry Smith if (mglevels[i]->smoothu == mglevels[i]->smoothd) { 641f3fbd535SBarry Smith /* if doing only down then initial guess is zero */ 642f3fbd535SBarry Smith ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 643f3fbd535SBarry Smith } 64463e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 645f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr); 64663e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 647d42688cbSBarry Smith if (!mglevels[i]->residual) { 648d42688cbSBarry Smith Mat mat; 649d42688cbSBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr); 650d42688cbSBarry Smith ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr); 651d42688cbSBarry Smith } 652f3fbd535SBarry Smith } 653f3fbd535SBarry Smith for (i=1; i<n; i++) { 654f3fbd535SBarry Smith if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) { 655f3fbd535SBarry Smith Mat downmat,downpmat; 656f3fbd535SBarry Smith MatStructure matflag; 657ace3abfcSBarry Smith PetscBool opsset; 658f3fbd535SBarry Smith 659f3fbd535SBarry Smith /* check if operators have been set for up, if not use down operators to set them */ 660f3fbd535SBarry Smith ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,PETSC_NULL);CHKERRQ(ierr); 661f3fbd535SBarry Smith if (!opsset) { 662f3fbd535SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr); 663f3fbd535SBarry Smith ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr); 664f3fbd535SBarry Smith } 665f3fbd535SBarry Smith 666f3fbd535SBarry Smith ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr); 66763e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 668f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr); 66963e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 670f3fbd535SBarry Smith } 671f3fbd535SBarry Smith } 672f3fbd535SBarry Smith 673f3fbd535SBarry Smith /* 674f3fbd535SBarry Smith If coarse solver is not direct method then DO NOT USE preonly 675f3fbd535SBarry Smith */ 676f3fbd535SBarry Smith ierr = PetscTypeCompare((PetscObject)mglevels[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr); 677f3fbd535SBarry Smith if (preonly) { 678f3fbd535SBarry Smith ierr = PetscTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr); 679f3fbd535SBarry Smith ierr = PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr); 680f3fbd535SBarry Smith ierr = PetscTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr); 681fd1303e9SJungho Lee ierr = PetscTypeCompare((PetscObject)cpc,PCSVD,&svd);CHKERRQ(ierr); 682fd1303e9SJungho Lee if (!lu && !redundant && !cholesky && !svd) { 683f3fbd535SBarry Smith ierr = KSPSetType(mglevels[0]->smoothd,KSPGMRES);CHKERRQ(ierr); 684f3fbd535SBarry Smith } 685f3fbd535SBarry Smith } 686f3fbd535SBarry Smith 687f3fbd535SBarry Smith if (!pc->setupcalled) { 688f3fbd535SBarry Smith if (monitor) { 689f3fbd535SBarry Smith ierr = PetscObjectGetComm((PetscObject)mglevels[0]->smoothd,&comm);CHKERRQ(ierr); 690f3fbd535SBarry Smith ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n,&ascii);CHKERRQ(ierr); 691c2efdce3SBarry Smith ierr = KSPMonitorSet(mglevels[0]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void**))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 692f3fbd535SBarry Smith } 693f3fbd535SBarry Smith ierr = KSPSetFromOptions(mglevels[0]->smoothd);CHKERRQ(ierr); 694f3fbd535SBarry Smith } 695f3fbd535SBarry Smith 69663e6d426SJed Brown if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 697f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr); 69863e6d426SJed Brown if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 699f3fbd535SBarry Smith 700f3fbd535SBarry Smith /* 701f3fbd535SBarry Smith Dump the interpolation/restriction matrices plus the 702e3c5b3baSBarry Smith Jacobian/stiffness on each level. This allows MATLAB users to 703f3fbd535SBarry Smith easily check if the Galerkin condition A_c = R A_f R^T is satisfied. 704f3fbd535SBarry Smith 705f3fbd535SBarry Smith Only support one or the other at the same time. 706f3fbd535SBarry Smith */ 707f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER) 708acfcf0e5SJed Brown ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,PETSC_NULL);CHKERRQ(ierr); 709f3fbd535SBarry Smith if (dump) { 710f3fbd535SBarry Smith viewer = PETSC_VIEWER_SOCKET_(((PetscObject)pc)->comm); 711f3fbd535SBarry Smith } 712f3fbd535SBarry Smith dump = PETSC_FALSE; 713f3fbd535SBarry Smith #endif 714acfcf0e5SJed Brown ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,PETSC_NULL);CHKERRQ(ierr); 715f3fbd535SBarry Smith if (dump) { 716f3fbd535SBarry Smith viewer = PETSC_VIEWER_BINARY_(((PetscObject)pc)->comm); 717f3fbd535SBarry Smith } 718f3fbd535SBarry Smith 719f3fbd535SBarry Smith if (viewer) { 720f3fbd535SBarry Smith for (i=1; i<n; i++) { 721f3fbd535SBarry Smith ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr); 722f3fbd535SBarry Smith } 723f3fbd535SBarry Smith for (i=0; i<n; i++) { 724f3fbd535SBarry Smith ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr); 725f3fbd535SBarry Smith ierr = MatView(pc->mat,viewer);CHKERRQ(ierr); 726f3fbd535SBarry Smith } 727f3fbd535SBarry Smith } 728f3fbd535SBarry Smith PetscFunctionReturn(0); 729f3fbd535SBarry Smith } 730f3fbd535SBarry Smith 731f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/ 732f3fbd535SBarry Smith 733f3fbd535SBarry Smith #undef __FUNCT__ 7349dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetLevels" 7354b9ad928SBarry Smith /*@ 73697177400SBarry Smith PCMGGetLevels - Gets the number of levels to use with MG. 7374b9ad928SBarry Smith 7384b9ad928SBarry Smith Not Collective 7394b9ad928SBarry Smith 7404b9ad928SBarry Smith Input Parameter: 7414b9ad928SBarry Smith . pc - the preconditioner context 7424b9ad928SBarry Smith 7434b9ad928SBarry Smith Output parameter: 7444b9ad928SBarry Smith . levels - the number of levels 7454b9ad928SBarry Smith 7464b9ad928SBarry Smith Level: advanced 7474b9ad928SBarry Smith 7484b9ad928SBarry Smith .keywords: MG, get, levels, multigrid 7494b9ad928SBarry Smith 75097177400SBarry Smith .seealso: PCMGSetLevels() 7514b9ad928SBarry Smith @*/ 7527087cfbeSBarry Smith PetscErrorCode PCMGGetLevels(PC pc,PetscInt *levels) 7534b9ad928SBarry Smith { 754f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 7554b9ad928SBarry Smith 7564b9ad928SBarry Smith PetscFunctionBegin; 7570700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 7584482741eSBarry Smith PetscValidIntPointer(levels,2); 759f3fbd535SBarry Smith *levels = mg->nlevels; 7604b9ad928SBarry Smith PetscFunctionReturn(0); 7614b9ad928SBarry Smith } 7624b9ad928SBarry Smith 7634b9ad928SBarry Smith #undef __FUNCT__ 7649dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetType" 7654b9ad928SBarry Smith /*@ 76697177400SBarry Smith PCMGSetType - Determines the form of multigrid to use: 7674b9ad928SBarry Smith multiplicative, additive, full, or the Kaskade algorithm. 7684b9ad928SBarry Smith 769ad4df100SBarry Smith Logically Collective on PC 7704b9ad928SBarry Smith 7714b9ad928SBarry Smith Input Parameters: 7724b9ad928SBarry Smith + pc - the preconditioner context 7739dcbbd2bSBarry Smith - form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE, 7749dcbbd2bSBarry Smith PC_MG_FULL, PC_MG_KASKADE 7754b9ad928SBarry Smith 7764b9ad928SBarry Smith Options Database Key: 7774b9ad928SBarry Smith . -pc_mg_type <form> - Sets <form>, one of multiplicative, 7784b9ad928SBarry Smith additive, full, kaskade 7794b9ad928SBarry Smith 7804b9ad928SBarry Smith Level: advanced 7814b9ad928SBarry Smith 7824b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid 7834b9ad928SBarry Smith 78497177400SBarry Smith .seealso: PCMGSetLevels() 7854b9ad928SBarry Smith @*/ 7867087cfbeSBarry Smith PetscErrorCode PCMGSetType(PC pc,PCMGType form) 7874b9ad928SBarry Smith { 788f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 7894b9ad928SBarry Smith 7904b9ad928SBarry Smith PetscFunctionBegin; 7910700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 792c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc,form,2); 79331567311SBarry Smith mg->am = form; 7949dcbbd2bSBarry Smith if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG; 7954b9ad928SBarry Smith else pc->ops->applyrichardson = 0; 7964b9ad928SBarry Smith PetscFunctionReturn(0); 7974b9ad928SBarry Smith } 7984b9ad928SBarry Smith 7994b9ad928SBarry Smith #undef __FUNCT__ 8000d353602SBarry Smith #define __FUNCT__ "PCMGSetCycleType" 8014b9ad928SBarry Smith /*@ 8020d353602SBarry Smith PCMGSetCycleType - Sets the type cycles to use. Use PCMGSetCycleTypeOnLevel() for more 8034b9ad928SBarry Smith complicated cycling. 8044b9ad928SBarry Smith 805ad4df100SBarry Smith Logically Collective on PC 8064b9ad928SBarry Smith 8074b9ad928SBarry Smith Input Parameters: 808c2be2410SBarry Smith + pc - the multigrid context 8090d353602SBarry Smith - PC_MG_CYCLE_V or PC_MG_CYCLE_W 8104b9ad928SBarry Smith 8114b9ad928SBarry Smith Options Database Key: 8120d353602SBarry Smith $ -pc_mg_cycle_type v or w 8134b9ad928SBarry Smith 8144b9ad928SBarry Smith Level: advanced 8154b9ad928SBarry Smith 8164b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 8174b9ad928SBarry Smith 8180d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel() 8194b9ad928SBarry Smith @*/ 8207087cfbeSBarry Smith PetscErrorCode PCMGSetCycleType(PC pc,PCMGCycleType n) 8214b9ad928SBarry Smith { 822f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 823f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 82479416396SBarry Smith PetscInt i,levels; 8254b9ad928SBarry Smith 8264b9ad928SBarry Smith PetscFunctionBegin; 8270700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 828e7e72b3dSBarry Smith if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 829c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 830f3fbd535SBarry Smith levels = mglevels[0]->levels; 8314b9ad928SBarry Smith 8324b9ad928SBarry Smith for (i=0; i<levels; i++) { 833f3fbd535SBarry Smith mglevels[i]->cycles = n; 8344b9ad928SBarry Smith } 8354b9ad928SBarry Smith PetscFunctionReturn(0); 8364b9ad928SBarry Smith } 8374b9ad928SBarry Smith 8384b9ad928SBarry Smith #undef __FUNCT__ 8398cc2d5dfSBarry Smith #define __FUNCT__ "PCMGMultiplicativeSetCycles" 8408cc2d5dfSBarry Smith /*@ 8418cc2d5dfSBarry Smith PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step 8428cc2d5dfSBarry Smith of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used 8438cc2d5dfSBarry Smith 844ad4df100SBarry Smith Logically Collective on PC 8458cc2d5dfSBarry Smith 8468cc2d5dfSBarry Smith Input Parameters: 8478cc2d5dfSBarry Smith + pc - the multigrid context 8488cc2d5dfSBarry Smith - n - number of cycles (default is 1) 8498cc2d5dfSBarry Smith 8508cc2d5dfSBarry Smith Options Database Key: 8518cc2d5dfSBarry Smith $ -pc_mg_multiplicative_cycles n 8528cc2d5dfSBarry Smith 8538cc2d5dfSBarry Smith Level: advanced 8548cc2d5dfSBarry Smith 8558cc2d5dfSBarry Smith Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType() 8568cc2d5dfSBarry Smith 8578cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 8588cc2d5dfSBarry Smith 8598cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType() 8608cc2d5dfSBarry Smith @*/ 8617087cfbeSBarry Smith PetscErrorCode PCMGMultiplicativeSetCycles(PC pc,PetscInt n) 8628cc2d5dfSBarry Smith { 863f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 864f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 8658cc2d5dfSBarry Smith PetscInt i,levels; 8668cc2d5dfSBarry Smith 8678cc2d5dfSBarry Smith PetscFunctionBegin; 8680700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 869e7e72b3dSBarry Smith if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 870c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 871f3fbd535SBarry Smith levels = mglevels[0]->levels; 8728cc2d5dfSBarry Smith 8738cc2d5dfSBarry Smith for (i=0; i<levels; i++) { 87431567311SBarry Smith mg->cyclesperpcapply = n; 8758cc2d5dfSBarry Smith } 8768cc2d5dfSBarry Smith PetscFunctionReturn(0); 8778cc2d5dfSBarry Smith } 8788cc2d5dfSBarry Smith 8798cc2d5dfSBarry Smith #undef __FUNCT__ 8809dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetGalerkin" 881c2be2410SBarry Smith /*@ 88297177400SBarry Smith PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the 883c2be2410SBarry Smith finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t 884c2be2410SBarry Smith 885ad4df100SBarry Smith Logically Collective on PC 886c2be2410SBarry Smith 887c2be2410SBarry Smith Input Parameters: 8883fc8bf9cSBarry Smith . pc - the multigrid context 889c2be2410SBarry Smith 890c2be2410SBarry Smith Options Database Key: 891c2be2410SBarry Smith $ -pc_mg_galerkin 892c2be2410SBarry Smith 893c2be2410SBarry Smith Level: intermediate 894c2be2410SBarry Smith 895c2be2410SBarry Smith .keywords: MG, set, Galerkin 896c2be2410SBarry Smith 8973fc8bf9cSBarry Smith .seealso: PCMGGetGalerkin() 8983fc8bf9cSBarry Smith 899c2be2410SBarry Smith @*/ 9007087cfbeSBarry Smith PetscErrorCode PCMGSetGalerkin(PC pc) 901c2be2410SBarry Smith { 902f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 903c2be2410SBarry Smith 904c2be2410SBarry Smith PetscFunctionBegin; 9050700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 906218a07d4SBarry Smith mg->galerkin = PETSC_TRUE; 907c2be2410SBarry Smith PetscFunctionReturn(0); 908c2be2410SBarry Smith } 909c2be2410SBarry Smith 910c2be2410SBarry Smith #undef __FUNCT__ 9113fc8bf9cSBarry Smith #define __FUNCT__ "PCMGGetGalerkin" 9123fc8bf9cSBarry Smith /*@ 9133fc8bf9cSBarry Smith PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. 9143fc8bf9cSBarry Smith A_i-1 = r_i * A_i * r_i^t 9153fc8bf9cSBarry Smith 9163fc8bf9cSBarry Smith Not Collective 9173fc8bf9cSBarry Smith 9183fc8bf9cSBarry Smith Input Parameter: 9193fc8bf9cSBarry Smith . pc - the multigrid context 9203fc8bf9cSBarry Smith 9213fc8bf9cSBarry Smith Output Parameter: 9223fc8bf9cSBarry Smith . gelerkin - PETSC_TRUE or PETSC_FALSE 9233fc8bf9cSBarry Smith 9243fc8bf9cSBarry Smith Options Database Key: 9253fc8bf9cSBarry Smith $ -pc_mg_galerkin 9263fc8bf9cSBarry Smith 9273fc8bf9cSBarry Smith Level: intermediate 9283fc8bf9cSBarry Smith 9293fc8bf9cSBarry Smith .keywords: MG, set, Galerkin 9303fc8bf9cSBarry Smith 9313fc8bf9cSBarry Smith .seealso: PCMGSetGalerkin() 9323fc8bf9cSBarry Smith 9333fc8bf9cSBarry Smith @*/ 9347087cfbeSBarry Smith PetscErrorCode PCMGGetGalerkin(PC pc,PetscBool *galerkin) 9353fc8bf9cSBarry Smith { 936f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 9373fc8bf9cSBarry Smith 9383fc8bf9cSBarry Smith PetscFunctionBegin; 9390700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 940218a07d4SBarry Smith *galerkin = mg->galerkin; 9413fc8bf9cSBarry Smith PetscFunctionReturn(0); 9423fc8bf9cSBarry Smith } 9433fc8bf9cSBarry Smith 9443fc8bf9cSBarry Smith #undef __FUNCT__ 9459dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothDown" 9464b9ad928SBarry Smith /*@ 94797177400SBarry Smith PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to 94897177400SBarry Smith use on all levels. Use PCMGGetSmootherDown() to set different 9494b9ad928SBarry Smith pre-smoothing steps on different levels. 9504b9ad928SBarry Smith 951ad4df100SBarry Smith Logically Collective on PC 9524b9ad928SBarry Smith 9534b9ad928SBarry Smith Input Parameters: 9544b9ad928SBarry Smith + mg - the multigrid context 9554b9ad928SBarry Smith - n - the number of smoothing steps 9564b9ad928SBarry Smith 9574b9ad928SBarry Smith Options Database Key: 9584b9ad928SBarry Smith . -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps 9594b9ad928SBarry Smith 9604b9ad928SBarry Smith Level: advanced 9614b9ad928SBarry Smith 9624b9ad928SBarry Smith .keywords: MG, smooth, down, pre-smoothing, steps, multigrid 9634b9ad928SBarry Smith 96497177400SBarry Smith .seealso: PCMGSetNumberSmoothUp() 9654b9ad928SBarry Smith @*/ 9667087cfbeSBarry Smith PetscErrorCode PCMGSetNumberSmoothDown(PC pc,PetscInt n) 9674b9ad928SBarry Smith { 968f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 969f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 9706849ba73SBarry Smith PetscErrorCode ierr; 97179416396SBarry Smith PetscInt i,levels; 9724b9ad928SBarry Smith 9734b9ad928SBarry Smith PetscFunctionBegin; 9740700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 975e7e72b3dSBarry Smith if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 976c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 977f3fbd535SBarry Smith levels = mglevels[0]->levels; 9784b9ad928SBarry Smith 979b05257ddSBarry Smith for (i=1; i<levels; i++) { 9804b9ad928SBarry Smith /* make sure smoother up and down are different */ 98197177400SBarry Smith ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 982f3fbd535SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 98331567311SBarry Smith mg->default_smoothd = n; 9844b9ad928SBarry Smith } 9854b9ad928SBarry Smith PetscFunctionReturn(0); 9864b9ad928SBarry Smith } 9874b9ad928SBarry Smith 9884b9ad928SBarry Smith #undef __FUNCT__ 9899dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothUp" 9904b9ad928SBarry Smith /*@ 99197177400SBarry Smith PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use 99297177400SBarry Smith on all levels. Use PCMGGetSmootherUp() to set different numbers of 9934b9ad928SBarry Smith post-smoothing steps on different levels. 9944b9ad928SBarry Smith 995ad4df100SBarry Smith Logically Collective on PC 9964b9ad928SBarry Smith 9974b9ad928SBarry Smith Input Parameters: 9984b9ad928SBarry Smith + mg - the multigrid context 9994b9ad928SBarry Smith - n - the number of smoothing steps 10004b9ad928SBarry Smith 10014b9ad928SBarry Smith Options Database Key: 10024b9ad928SBarry Smith . -pc_mg_smoothup <n> - Sets number of post-smoothing steps 10034b9ad928SBarry Smith 10044b9ad928SBarry Smith Level: advanced 10054b9ad928SBarry Smith 10064b9ad928SBarry Smith Note: this does not set a value on the coarsest grid, since we assume that 1007a8c7a070SBarry Smith there is no separate smooth up on the coarsest grid. 10084b9ad928SBarry Smith 10094b9ad928SBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid 10104b9ad928SBarry Smith 101197177400SBarry Smith .seealso: PCMGSetNumberSmoothDown() 10124b9ad928SBarry Smith @*/ 10137087cfbeSBarry Smith PetscErrorCode PCMGSetNumberSmoothUp(PC pc,PetscInt n) 10144b9ad928SBarry Smith { 1015f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1016f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 10176849ba73SBarry Smith PetscErrorCode ierr; 101879416396SBarry Smith PetscInt i,levels; 10194b9ad928SBarry Smith 10204b9ad928SBarry Smith PetscFunctionBegin; 10210700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1022e7e72b3dSBarry Smith if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 1023c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 1024f3fbd535SBarry Smith levels = mglevels[0]->levels; 10254b9ad928SBarry Smith 10264b9ad928SBarry Smith for (i=1; i<levels; i++) { 10274b9ad928SBarry Smith /* make sure smoother up and down are different */ 102897177400SBarry Smith ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 1029f3fbd535SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 103031567311SBarry Smith mg->default_smoothu = n; 10314b9ad928SBarry Smith } 10324b9ad928SBarry Smith PetscFunctionReturn(0); 10334b9ad928SBarry Smith } 10344b9ad928SBarry Smith 10354b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/ 10364b9ad928SBarry Smith 10373b09bd56SBarry Smith /*MC 1038ccb205f8SBarry Smith PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional 10393b09bd56SBarry Smith information about the coarser grid matrices and restriction/interpolation operators. 10403b09bd56SBarry Smith 10413b09bd56SBarry Smith Options Database Keys: 10423b09bd56SBarry Smith + -pc_mg_levels <nlevels> - number of levels including finest 10430d353602SBarry Smith . -pc_mg_cycles v or w 104479416396SBarry Smith . -pc_mg_smoothup <n> - number of smoothing steps after interpolation 10453b09bd56SBarry Smith . -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator 10463b09bd56SBarry Smith . -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default 10473b09bd56SBarry Smith . -pc_mg_log - log information about time spent on each level of the solver 10483b09bd56SBarry Smith . -pc_mg_monitor - print information on the multigrid convergence 104968eff7e6SBarry Smith . -pc_mg_galerkin - use Galerkin process to compute coarser operators 10503b09bd56SBarry Smith - -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 1051e3c5b3baSBarry Smith to the Socket viewer for reading from MATLAB. 10523b09bd56SBarry Smith 105324c3aa18SBarry 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 10543b09bd56SBarry Smith 10553b09bd56SBarry Smith Level: intermediate 10563b09bd56SBarry Smith 10578f87f92bSBarry Smith Concepts: multigrid/multilevel 10583b09bd56SBarry Smith 105924c3aa18SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC 10600d353602SBarry Smith PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(), 106197177400SBarry Smith PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 106297177400SBarry Smith PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 10630d353602SBarry Smith PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 10643b09bd56SBarry Smith M*/ 10653b09bd56SBarry Smith 10664b9ad928SBarry Smith EXTERN_C_BEGIN 10674b9ad928SBarry Smith #undef __FUNCT__ 10684b9ad928SBarry Smith #define __FUNCT__ "PCCreate_MG" 10697087cfbeSBarry Smith PetscErrorCode PCCreate_MG(PC pc) 10704b9ad928SBarry Smith { 1071f3fbd535SBarry Smith PC_MG *mg; 1072f3fbd535SBarry Smith PetscErrorCode ierr; 1073f3fbd535SBarry Smith 10744b9ad928SBarry Smith PetscFunctionBegin; 1075f3fbd535SBarry Smith ierr = PetscNewLog(pc,PC_MG,&mg);CHKERRQ(ierr); 1076f3fbd535SBarry Smith pc->data = (void*)mg; 1077f3fbd535SBarry Smith mg->nlevels = -1; 1078f3fbd535SBarry Smith 10794b9ad928SBarry Smith pc->ops->apply = PCApply_MG; 10804b9ad928SBarry Smith pc->ops->setup = PCSetUp_MG; 1081a06653b4SBarry Smith pc->ops->reset = PCReset_MG; 10824b9ad928SBarry Smith pc->ops->destroy = PCDestroy_MG; 10834b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_MG; 10844b9ad928SBarry Smith pc->ops->view = PCView_MG; 10854b9ad928SBarry Smith PetscFunctionReturn(0); 10864b9ad928SBarry Smith } 10874b9ad928SBarry Smith EXTERN_C_END 1088