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); 13073250ac0SBarry Smith ierr = VecDestroy(&mglevels[i+1]->rscale);CHKERRQ(ierr); 1313aeaf226SBarry Smith } 1323aeaf226SBarry Smith 1333aeaf226SBarry Smith for (i=0; i<n; i++) { 1343aeaf226SBarry Smith ierr = MatDestroy(&mglevels[i]->A);CHKERRQ(ierr); 1353aeaf226SBarry Smith if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 1363aeaf226SBarry Smith ierr = KSPReset(mglevels[i]->smoothd);CHKERRQ(ierr); 1373aeaf226SBarry Smith } 1383aeaf226SBarry Smith ierr = KSPReset(mglevels[i]->smoothu);CHKERRQ(ierr); 1393aeaf226SBarry Smith } 1403aeaf226SBarry Smith } 1413aeaf226SBarry Smith PetscFunctionReturn(0); 1423aeaf226SBarry Smith } 1433aeaf226SBarry Smith 1443aeaf226SBarry Smith #undef __FUNCT__ 1459dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetLevels" 1464b9ad928SBarry Smith /*@C 14797177400SBarry Smith PCMGSetLevels - Sets the number of levels to use with MG. 1484b9ad928SBarry Smith Must be called before any other MG routine. 1494b9ad928SBarry Smith 150ad4df100SBarry Smith Logically Collective on PC 1514b9ad928SBarry Smith 1524b9ad928SBarry Smith Input Parameters: 1534b9ad928SBarry Smith + pc - the preconditioner context 1544b9ad928SBarry Smith . levels - the number of levels 1554b9ad928SBarry Smith - comms - optional communicators for each level; this is to allow solving the coarser problems 1564b9ad928SBarry Smith on smaller sets of processors. Use PETSC_NULL_OBJECT for default in Fortran 1574b9ad928SBarry Smith 1584b9ad928SBarry Smith Level: intermediate 1594b9ad928SBarry Smith 1604b9ad928SBarry Smith Notes: 1614b9ad928SBarry Smith If the number of levels is one then the multigrid uses the -mg_levels prefix 1624b9ad928SBarry Smith for setting the level options rather than the -mg_coarse prefix. 1634b9ad928SBarry Smith 1644b9ad928SBarry Smith .keywords: MG, set, levels, multigrid 1654b9ad928SBarry Smith 16697177400SBarry Smith .seealso: PCMGSetType(), PCMGGetLevels() 1674b9ad928SBarry Smith @*/ 1687087cfbeSBarry Smith PetscErrorCode PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms) 1694b9ad928SBarry Smith { 170dfbe8321SBarry Smith PetscErrorCode ierr; 171f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 172f3fbd535SBarry Smith MPI_Comm comm = ((PetscObject)pc)->comm; 1733aeaf226SBarry Smith PC_MG_Levels **mglevels = mg->levels; 174f3fbd535SBarry Smith PetscInt i; 175f3fbd535SBarry Smith PetscMPIInt size; 176f3fbd535SBarry Smith const char *prefix; 177f3fbd535SBarry Smith PC ipc; 1783aeaf226SBarry Smith PetscInt n; 1794b9ad928SBarry Smith 1804b9ad928SBarry Smith PetscFunctionBegin; 1810700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 182548767e0SJed Brown PetscValidLogicalCollectiveInt(pc,levels,2); 183548767e0SJed Brown if (mg->nlevels == levels) PetscFunctionReturn(0); 1843aeaf226SBarry Smith if (mglevels) { 1853aeaf226SBarry Smith /* changing the number of levels so free up the previous stuff */ 1863aeaf226SBarry Smith ierr = PCReset_MG(pc);CHKERRQ(ierr); 1873aeaf226SBarry Smith n = mglevels[0]->levels; 1883aeaf226SBarry Smith for (i=0; i<n; i++) { 1893aeaf226SBarry Smith if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 1903aeaf226SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr); 1913aeaf226SBarry Smith } 1923aeaf226SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr); 1933aeaf226SBarry Smith ierr = PetscFree(mglevels[i]);CHKERRQ(ierr); 1943aeaf226SBarry Smith } 1953aeaf226SBarry Smith ierr = PetscFree(mg->levels);CHKERRQ(ierr); 1963aeaf226SBarry Smith } 197f3fbd535SBarry Smith 198f3fbd535SBarry Smith mg->nlevels = levels; 199f3fbd535SBarry Smith 200f3fbd535SBarry Smith ierr = PetscMalloc(levels*sizeof(PC_MG*),&mglevels);CHKERRQ(ierr); 201f3fbd535SBarry Smith ierr = PetscLogObjectMemory(pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr); 202f3fbd535SBarry Smith 203f3fbd535SBarry Smith ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 204f3fbd535SBarry Smith 205a9db87a2SMatthew G Knepley mg->stageApply = 0; 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; 307a9db87a2SMatthew G Knepley if (mg->stageApply) {ierr = PetscLogStagePush(mg->stageApply);CHKERRQ(ierr);} 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 } 333a9db87a2SMatthew G Knepley if (mg->stageApply) {ierr = PetscLogStagePop();CHKERRQ(ierr);} 334f3fbd535SBarry Smith PetscFunctionReturn(0); 335f3fbd535SBarry Smith } 336f3fbd535SBarry Smith 337f3fbd535SBarry Smith 338f3fbd535SBarry Smith #undef __FUNCT__ 339f3fbd535SBarry Smith #define __FUNCT__ "PCSetFromOptions_MG" 340f3fbd535SBarry Smith PetscErrorCode PCSetFromOptions_MG(PC pc) 341f3fbd535SBarry Smith { 342f3fbd535SBarry Smith PetscErrorCode ierr; 343f3fbd535SBarry Smith PetscInt m,levels = 1,cycles; 344c91913d3SJed Brown PetscBool flg,set; 345f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 346f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 347f3fbd535SBarry Smith PCMGType mgtype; 348f3fbd535SBarry Smith PCMGCycleType mgctype; 349f3fbd535SBarry Smith 350f3fbd535SBarry Smith PetscFunctionBegin; 351f3fbd535SBarry Smith ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr); 3523aeaf226SBarry Smith if (!mg->levels) { 353f3fbd535SBarry Smith ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr); 354eb3f98d2SBarry Smith if (!flg && pc->dm) { 355eb3f98d2SBarry Smith ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr); 356eb3f98d2SBarry Smith levels++; 3573aeaf226SBarry Smith mg->usedmfornumberoflevels = PETSC_TRUE; 358eb3f98d2SBarry Smith } 359f3fbd535SBarry Smith ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr); 360f3fbd535SBarry Smith } 3613aeaf226SBarry Smith mglevels = mg->levels; 3623aeaf226SBarry Smith 363f3fbd535SBarry Smith mgctype = (PCMGCycleType) mglevels[0]->cycles; 364f3fbd535SBarry Smith ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr); 365f3fbd535SBarry Smith if (flg) { 366f3fbd535SBarry Smith ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr); 367f3fbd535SBarry Smith }; 368f3fbd535SBarry Smith flg = PETSC_FALSE; 369c91913d3SJed Brown ierr = PetscOptionsBool("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",flg,&flg,&set);CHKERRQ(ierr); 370c91913d3SJed Brown if (set) { 371c91913d3SJed Brown ierr = PCMGSetGalerkin(pc,flg);CHKERRQ(ierr); 372f3fbd535SBarry Smith } 373f3fbd535SBarry Smith ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr); 374f3fbd535SBarry Smith if (flg) { 375f3fbd535SBarry Smith ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr); 376f3fbd535SBarry Smith } 377f3fbd535SBarry Smith ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr); 378f3fbd535SBarry Smith if (flg) { 379f3fbd535SBarry Smith ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr); 380f3fbd535SBarry Smith } 38131567311SBarry Smith mgtype = mg->am; 382f3fbd535SBarry Smith ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr); 383f3fbd535SBarry Smith if (flg) { 384f3fbd535SBarry Smith ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr); 385f3fbd535SBarry Smith } 38631567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 38731567311SBarry Smith ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGSetLevels",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr); 388f3fbd535SBarry Smith if (flg) { 389f3fbd535SBarry Smith ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr); 390f3fbd535SBarry Smith } 391f3fbd535SBarry Smith } 392f3fbd535SBarry Smith flg = PETSC_FALSE; 393acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 394f3fbd535SBarry Smith if (flg) { 395f3fbd535SBarry Smith PetscInt i; 396f3fbd535SBarry Smith char eventname[128]; 39763e6d426SJed Brown if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 398f3fbd535SBarry Smith levels = mglevels[0]->levels; 399f3fbd535SBarry Smith for (i=0; i<levels; i++) { 400f3fbd535SBarry Smith sprintf(eventname,"MGSetup Level %d",(int)i); 40163e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr); 402f3fbd535SBarry Smith sprintf(eventname,"MGSmooth Level %d",(int)i); 40363e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr); 404f3fbd535SBarry Smith if (i) { 405f3fbd535SBarry Smith sprintf(eventname,"MGResid Level %d",(int)i); 40663e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr); 407f3fbd535SBarry Smith sprintf(eventname,"MGInterp Level %d",(int)i); 40863e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr); 409f3fbd535SBarry Smith } 410f3fbd535SBarry Smith } 411*eec35531SMatthew G Knepley 412*eec35531SMatthew G Knepley { 413*eec35531SMatthew G Knepley const char *sname = "MG Apply"; 414*eec35531SMatthew G Knepley PetscStageLog stageLog; 415*eec35531SMatthew G Knepley PetscInt st; 416*eec35531SMatthew G Knepley 417*eec35531SMatthew G Knepley PetscFunctionBegin; 418*eec35531SMatthew G Knepley ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr); 419*eec35531SMatthew G Knepley for(st = 0; st < stageLog->numStages; ++st) { 420*eec35531SMatthew G Knepley PetscBool same; 421*eec35531SMatthew G Knepley 422*eec35531SMatthew G Knepley ierr = PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);CHKERRQ(ierr); 423*eec35531SMatthew G Knepley if (same) {mg->stageApply = st;} 424*eec35531SMatthew G Knepley } 425*eec35531SMatthew G Knepley if (!mg->stageApply) { 426*eec35531SMatthew G Knepley ierr = PetscLogStageRegister(sname, &mg->stageApply);CHKERRQ(ierr); 427*eec35531SMatthew G Knepley } 428*eec35531SMatthew G Knepley } 429f3fbd535SBarry Smith } 430f3fbd535SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 431f3fbd535SBarry Smith PetscFunctionReturn(0); 432f3fbd535SBarry Smith } 433f3fbd535SBarry Smith 434f3fbd535SBarry Smith const char *PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0}; 435f3fbd535SBarry Smith const char *PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0}; 436f3fbd535SBarry Smith 437f3fbd535SBarry Smith #undef __FUNCT__ 438f3fbd535SBarry Smith #define __FUNCT__ "PCView_MG" 439f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer) 440f3fbd535SBarry Smith { 441f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 442f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 443f3fbd535SBarry Smith PetscErrorCode ierr; 444f3fbd535SBarry Smith PetscInt levels = mglevels[0]->levels,i; 445ace3abfcSBarry Smith PetscBool iascii; 446f3fbd535SBarry Smith 447f3fbd535SBarry Smith PetscFunctionBegin; 448251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 449f3fbd535SBarry Smith if (iascii) { 45031567311SBarry 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); 45131567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 45231567311SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr); 453f3fbd535SBarry Smith } 454218a07d4SBarry Smith if (mg->galerkin) { 455f3fbd535SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 4564f66f45eSBarry Smith } else { 4574f66f45eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 458f3fbd535SBarry Smith } 459f3fbd535SBarry Smith for (i=0; i<levels; i++) { 460f3fbd535SBarry Smith if (!i) { 46189cce641SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);CHKERRQ(ierr); 462f3fbd535SBarry Smith } else { 46389cce641SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 464f3fbd535SBarry Smith } 465f3fbd535SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 466f3fbd535SBarry Smith ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 467f3fbd535SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 468f3fbd535SBarry Smith if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) { 469f3fbd535SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 470f3fbd535SBarry Smith } else if (i){ 47189cce641SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 472f3fbd535SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 473f3fbd535SBarry Smith ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 474f3fbd535SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 475f3fbd535SBarry Smith } 476f3fbd535SBarry Smith } 477649052a6SBarry Smith } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCMG",((PetscObject)viewer)->type_name); 478f3fbd535SBarry Smith PetscFunctionReturn(0); 479f3fbd535SBarry Smith } 480f3fbd535SBarry Smith 481b45d2f2cSJed Brown #include <petsc-private/dmimpl.h> 482b45d2f2cSJed Brown #include <petsc-private/kspimpl.h> 483cab2e9ccSBarry Smith 484f3fbd535SBarry Smith /* 485f3fbd535SBarry Smith Calls setup for the KSP on each level 486f3fbd535SBarry Smith */ 487f3fbd535SBarry Smith #undef __FUNCT__ 488f3fbd535SBarry Smith #define __FUNCT__ "PCSetUp_MG" 489f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc) 490f3fbd535SBarry Smith { 491f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 492f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 493f3fbd535SBarry Smith PetscErrorCode ierr; 494f3fbd535SBarry Smith PetscInt i,n = mglevels[0]->levels; 49598e04984SBarry Smith PC cpc; 496649052a6SBarry Smith PetscBool preonly,lu,redundant,cholesky,svd,dump = PETSC_FALSE,opsset; 497f3fbd535SBarry Smith Mat dA,dB; 498f3fbd535SBarry Smith MatStructure uflag; 499f3fbd535SBarry Smith Vec tvec; 500218a07d4SBarry Smith DM *dms; 501649052a6SBarry Smith PetscViewer viewer = 0; 502f3fbd535SBarry Smith 503f3fbd535SBarry Smith PetscFunctionBegin; 50401bc414fSDmitry Karpeev /* FIX: Move this to PCSetFromOptions_MG? */ 5053aeaf226SBarry Smith if (mg->usedmfornumberoflevels) { 5063aeaf226SBarry Smith PetscInt levels; 5073aeaf226SBarry Smith ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr); 5083aeaf226SBarry Smith levels++; 5093aeaf226SBarry Smith if (levels > n) { /* the problem is now being solved on a finer grid */ 5103aeaf226SBarry Smith ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr); 5113aeaf226SBarry Smith n = levels; 5123aeaf226SBarry Smith ierr = PCSetFromOptions(pc);CHKERRQ(ierr); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */ 5133aeaf226SBarry Smith mglevels = mg->levels; 5143aeaf226SBarry Smith } 5153aeaf226SBarry Smith } 51698e04984SBarry Smith ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr); 5173aeaf226SBarry Smith 518f3fbd535SBarry Smith 519f3fbd535SBarry Smith /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */ 520f3fbd535SBarry Smith /* so use those from global PC */ 521f3fbd535SBarry Smith /* Is this what we always want? What if user wants to keep old one? */ 522f3fbd535SBarry Smith ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,PETSC_NULL,&opsset);CHKERRQ(ierr); 52398e04984SBarry Smith if (opsset) { 52498e04984SBarry Smith Mat mmat; 52598e04984SBarry Smith ierr = KSPGetOperators(mglevels[n-1]->smoothd,PETSC_NULL,&mmat,PETSC_NULL);CHKERRQ(ierr); 52698e04984SBarry Smith if (mmat == pc->pmat) opsset = PETSC_FALSE; 52798e04984SBarry Smith } 52898e04984SBarry Smith if (!opsset) { 529f3fbd535SBarry 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); 530f3fbd535SBarry Smith ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr); 531f3fbd535SBarry Smith } 532f3fbd535SBarry Smith 53301bc414fSDmitry 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? */ 534ce4cda84SJed Brown if (pc->dm && mg->galerkin != 2 && !pc->setupcalled) { 5352d2b81a6SBarry Smith /* construct the interpolation from the DMs */ 5362e499ae9SBarry Smith Mat p; 53773250ac0SBarry Smith Vec rscale; 538218a07d4SBarry Smith ierr = PetscMalloc(n*sizeof(DM),&dms);CHKERRQ(ierr); 539218a07d4SBarry Smith dms[n-1] = pc->dm; 540218a07d4SBarry Smith for (i=n-2; i>-1; i--) { 541fe86f630SJed Brown KSPDM kdm; 5422ee06e3bSJed Brown ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr); 5433c0c59f3SBarry Smith ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr); 5443c0c59f3SBarry Smith if (mg->galerkin) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);} 545fe86f630SJed Brown ierr = DMKSPGetContextWrite(dms[i],&kdm);CHKERRQ(ierr); 546fe86f630SJed Brown /* Ugly hack so that the next KSPSetUp() will use the RHS that we set */ 547fe86f630SJed Brown kdm->computerhs = PETSC_NULL; 548fe86f630SJed Brown kdm->rhsctx = PETSC_NULL; 54911629dbeSBarry Smith ierr = DMSetFunction(dms[i],0); 55011629dbeSBarry Smith ierr = DMSetInitialGuess(dms[i],0); 55124c3aa18SBarry Smith if (!mglevels[i+1]->interpolate) { 552e727c939SJed Brown ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr); 5532d2b81a6SBarry Smith ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr); 55400ac8be1SBarry Smith if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);} 55573250ac0SBarry Smith ierr = VecDestroy(&rscale);CHKERRQ(ierr); 5566bf464f9SBarry Smith ierr = MatDestroy(&p);CHKERRQ(ierr); 557218a07d4SBarry Smith } 55824c3aa18SBarry Smith } 5592d2b81a6SBarry Smith 560218a07d4SBarry Smith for (i=n-2; i>-1; i--) { 5616bf464f9SBarry Smith ierr = DMDestroy(&dms[i]);CHKERRQ(ierr); 562218a07d4SBarry Smith } 5632d2b81a6SBarry Smith ierr = PetscFree(dms);CHKERRQ(ierr); 564ce4cda84SJed Brown } 565cab2e9ccSBarry Smith 566ce4cda84SJed Brown if (pc->dm && !pc->setupcalled) { 567ce4cda84SJed Brown /* finest smoother also gets DM but it is not active, independent of whether galerkin==2 */ 568cab2e9ccSBarry Smith ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr); 569cab2e9ccSBarry Smith ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr); 570218a07d4SBarry Smith } 571218a07d4SBarry Smith 572c91913d3SJed Brown if (mg->galerkin == 1) { 573f3fbd535SBarry Smith Mat B; 574f3fbd535SBarry Smith /* currently only handle case where mat and pmat are the same on coarser levels */ 575f3fbd535SBarry Smith ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr); 576f3fbd535SBarry Smith if (!pc->setupcalled) { 577f3fbd535SBarry Smith for (i=n-2; i>-1; i--) { 578f3fbd535SBarry Smith ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 579f3fbd535SBarry Smith ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 580f3fbd535SBarry Smith if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);} 581f3fbd535SBarry Smith dB = B; 582f3fbd535SBarry Smith } 583cd9507b2SBarry Smith if (n > 1) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);} 584f3fbd535SBarry Smith } else { 585f3fbd535SBarry Smith for (i=n-2; i>-1; i--) { 586f3fbd535SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr); 587f3fbd535SBarry Smith ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 588f3fbd535SBarry Smith ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 589f3fbd535SBarry Smith dB = B; 590f3fbd535SBarry Smith } 591f3fbd535SBarry Smith } 592ce4cda84SJed Brown } else if (!mg->galerkin && pc->dm && pc->dm->x) { 593cab2e9ccSBarry Smith /* need to restrict Jacobian location to coarser meshes for evaluation */ 594cab2e9ccSBarry Smith for (i=n-2;i>-1; i--) { 595c88c5224SJed Brown Mat R; 596c88c5224SJed Brown Vec rscale; 597cab2e9ccSBarry Smith if (!mglevels[i]->smoothd->dm->x) { 598cab2e9ccSBarry Smith Vec *vecs; 599cab2e9ccSBarry Smith ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vecs,0,PETSC_NULL);CHKERRQ(ierr); 600cab2e9ccSBarry Smith mglevels[i]->smoothd->dm->x = vecs[0]; 601cab2e9ccSBarry Smith ierr = PetscFree(vecs);CHKERRQ(ierr); 602cab2e9ccSBarry Smith } 603c88c5224SJed Brown ierr = PCMGGetRestriction(pc,i+1,&R);CHKERRQ(ierr); 604c88c5224SJed Brown ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr); 605c88c5224SJed Brown ierr = MatRestrict(R,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);CHKERRQ(ierr); 606c88c5224SJed Brown ierr = VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,rscale);CHKERRQ(ierr); 607cab2e9ccSBarry Smith } 608f3fbd535SBarry Smith } 609ccceb50cSJed Brown if (!mg->galerkin && pc->dm) { 610caa4e7f2SJed Brown for (i=n-2;i>=0; i--) { 611ccceb50cSJed Brown DM dmfine,dmcoarse; 612caa4e7f2SJed Brown Mat Restrict,Inject; 613caa4e7f2SJed Brown Vec rscale; 614ccceb50cSJed Brown ierr = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr); 615ccceb50cSJed Brown ierr = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr); 616caa4e7f2SJed Brown ierr = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr); 617caa4e7f2SJed Brown ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr); 618caa4e7f2SJed Brown Inject = PETSC_NULL; /* Callback should create it if it needs Injection */ 619caa4e7f2SJed Brown ierr = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr); 620caa4e7f2SJed Brown } 621caa4e7f2SJed Brown } 622f3fbd535SBarry Smith 623f3fbd535SBarry Smith if (!pc->setupcalled) { 624f3fbd535SBarry Smith for (i=0; i<n; i++) { 625f3fbd535SBarry Smith ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr); 626f3fbd535SBarry Smith } 627f3fbd535SBarry Smith for (i=1; i<n; i++) { 628f3fbd535SBarry Smith if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) { 629f3fbd535SBarry Smith ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr); 630f3fbd535SBarry Smith } 631f3fbd535SBarry Smith } 632f3fbd535SBarry Smith for (i=1; i<n; i++) { 633c88c5224SJed Brown ierr = PCMGGetInterpolation(pc,i,&mglevels[i]->interpolate);CHKERRQ(ierr); 634c88c5224SJed Brown ierr = PCMGGetRestriction(pc,i,&mglevels[i]->restrct);CHKERRQ(ierr); 635f3fbd535SBarry Smith } 636f3fbd535SBarry Smith for (i=0; i<n-1; i++) { 637f3fbd535SBarry Smith if (!mglevels[i]->b) { 638f3fbd535SBarry Smith Vec *vec; 639f3fbd535SBarry Smith ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr); 640f3fbd535SBarry Smith ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr); 6416bf464f9SBarry Smith ierr = VecDestroy(vec);CHKERRQ(ierr); 642f3fbd535SBarry Smith ierr = PetscFree(vec);CHKERRQ(ierr); 643f3fbd535SBarry Smith } 644f3fbd535SBarry Smith if (!mglevels[i]->r && i) { 645f3fbd535SBarry Smith ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 646f3fbd535SBarry Smith ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr); 6476bf464f9SBarry Smith ierr = VecDestroy(&tvec);CHKERRQ(ierr); 648f3fbd535SBarry Smith } 649f3fbd535SBarry Smith if (!mglevels[i]->x) { 650f3fbd535SBarry Smith ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 651f3fbd535SBarry Smith ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr); 6526bf464f9SBarry Smith ierr = VecDestroy(&tvec);CHKERRQ(ierr); 653f3fbd535SBarry Smith } 654f3fbd535SBarry Smith } 655f3fbd535SBarry Smith if (n != 1 && !mglevels[n-1]->r) { 656f3fbd535SBarry Smith /* PCMGSetR() on the finest level if user did not supply it */ 657f3fbd535SBarry Smith Vec *vec; 658f3fbd535SBarry Smith ierr = KSPGetVecs(mglevels[n-1]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr); 659f3fbd535SBarry Smith ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr); 6606bf464f9SBarry Smith ierr = VecDestroy(vec);CHKERRQ(ierr); 661f3fbd535SBarry Smith ierr = PetscFree(vec);CHKERRQ(ierr); 662f3fbd535SBarry Smith } 663f3fbd535SBarry Smith } 664f3fbd535SBarry Smith 66598e04984SBarry Smith if (pc->dm) { 66698e04984SBarry Smith /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */ 66798e04984SBarry Smith for (i=0; i<n-1; i++){ 66898e04984SBarry Smith if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX; 66998e04984SBarry Smith } 67098e04984SBarry Smith } 671f3fbd535SBarry Smith 672f3fbd535SBarry Smith for (i=1; i<n; i++) { 673f3fbd535SBarry Smith if (mglevels[i]->smoothu == mglevels[i]->smoothd) { 674f3fbd535SBarry Smith /* if doing only down then initial guess is zero */ 675f3fbd535SBarry Smith ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 676f3fbd535SBarry Smith } 67763e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 678f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr); 67963e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 680d42688cbSBarry Smith if (!mglevels[i]->residual) { 681d42688cbSBarry Smith Mat mat; 682d42688cbSBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr); 683d42688cbSBarry Smith ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr); 684d42688cbSBarry Smith } 685f3fbd535SBarry Smith } 686f3fbd535SBarry Smith for (i=1; i<n; i++) { 687f3fbd535SBarry Smith if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) { 688f3fbd535SBarry Smith Mat downmat,downpmat; 689f3fbd535SBarry Smith MatStructure matflag; 690f3fbd535SBarry Smith 691f3fbd535SBarry Smith /* check if operators have been set for up, if not use down operators to set them */ 692f3fbd535SBarry Smith ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,PETSC_NULL);CHKERRQ(ierr); 693f3fbd535SBarry Smith if (!opsset) { 694f3fbd535SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr); 695f3fbd535SBarry Smith ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr); 696f3fbd535SBarry Smith } 697f3fbd535SBarry Smith 698f3fbd535SBarry Smith ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr); 69963e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 700f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr); 70163e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 702f3fbd535SBarry Smith } 703f3fbd535SBarry Smith } 704f3fbd535SBarry Smith 705f3fbd535SBarry Smith /* 706f3fbd535SBarry Smith If coarse solver is not direct method then DO NOT USE preonly 707f3fbd535SBarry Smith */ 708251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mglevels[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr); 709f3fbd535SBarry Smith if (preonly) { 710251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr); 711251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr); 712251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr); 713251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)cpc,PCSVD,&svd);CHKERRQ(ierr); 714fd1303e9SJungho Lee if (!lu && !redundant && !cholesky && !svd) { 715f3fbd535SBarry Smith ierr = KSPSetType(mglevels[0]->smoothd,KSPGMRES);CHKERRQ(ierr); 716f3fbd535SBarry Smith } 717f3fbd535SBarry Smith } 718f3fbd535SBarry Smith 719f3fbd535SBarry Smith if (!pc->setupcalled) { 720f3fbd535SBarry Smith ierr = KSPSetFromOptions(mglevels[0]->smoothd);CHKERRQ(ierr); 721f3fbd535SBarry Smith } 722f3fbd535SBarry Smith 72363e6d426SJed Brown if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 724f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr); 72563e6d426SJed Brown if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 726f3fbd535SBarry Smith 727f3fbd535SBarry Smith /* 728f3fbd535SBarry Smith Dump the interpolation/restriction matrices plus the 729e3c5b3baSBarry Smith Jacobian/stiffness on each level. This allows MATLAB users to 730f3fbd535SBarry Smith easily check if the Galerkin condition A_c = R A_f R^T is satisfied. 731f3fbd535SBarry Smith 732f3fbd535SBarry Smith Only support one or the other at the same time. 733f3fbd535SBarry Smith */ 734f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER) 735acfcf0e5SJed Brown ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,PETSC_NULL);CHKERRQ(ierr); 736f3fbd535SBarry Smith if (dump) { 737f3fbd535SBarry Smith viewer = PETSC_VIEWER_SOCKET_(((PetscObject)pc)->comm); 738f3fbd535SBarry Smith } 739f3fbd535SBarry Smith dump = PETSC_FALSE; 740f3fbd535SBarry Smith #endif 741acfcf0e5SJed Brown ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,PETSC_NULL);CHKERRQ(ierr); 742f3fbd535SBarry Smith if (dump) { 743f3fbd535SBarry Smith viewer = PETSC_VIEWER_BINARY_(((PetscObject)pc)->comm); 744f3fbd535SBarry Smith } 745f3fbd535SBarry Smith 746f3fbd535SBarry Smith if (viewer) { 747f3fbd535SBarry Smith for (i=1; i<n; i++) { 748f3fbd535SBarry Smith ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr); 749f3fbd535SBarry Smith } 750f3fbd535SBarry Smith for (i=0; i<n; i++) { 751f3fbd535SBarry Smith ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr); 752f3fbd535SBarry Smith ierr = MatView(pc->mat,viewer);CHKERRQ(ierr); 753f3fbd535SBarry Smith } 754f3fbd535SBarry Smith } 755f3fbd535SBarry Smith PetscFunctionReturn(0); 756f3fbd535SBarry Smith } 757f3fbd535SBarry Smith 758f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/ 759f3fbd535SBarry Smith 760f3fbd535SBarry Smith #undef __FUNCT__ 7619dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetLevels" 7624b9ad928SBarry Smith /*@ 76397177400SBarry Smith PCMGGetLevels - Gets the number of levels to use with MG. 7644b9ad928SBarry Smith 7654b9ad928SBarry Smith Not Collective 7664b9ad928SBarry Smith 7674b9ad928SBarry Smith Input Parameter: 7684b9ad928SBarry Smith . pc - the preconditioner context 7694b9ad928SBarry Smith 7704b9ad928SBarry Smith Output parameter: 7714b9ad928SBarry Smith . levels - the number of levels 7724b9ad928SBarry Smith 7734b9ad928SBarry Smith Level: advanced 7744b9ad928SBarry Smith 7754b9ad928SBarry Smith .keywords: MG, get, levels, multigrid 7764b9ad928SBarry Smith 77797177400SBarry Smith .seealso: PCMGSetLevels() 7784b9ad928SBarry Smith @*/ 7797087cfbeSBarry Smith PetscErrorCode PCMGGetLevels(PC pc,PetscInt *levels) 7804b9ad928SBarry Smith { 781f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 7824b9ad928SBarry Smith 7834b9ad928SBarry Smith PetscFunctionBegin; 7840700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 7854482741eSBarry Smith PetscValidIntPointer(levels,2); 786f3fbd535SBarry Smith *levels = mg->nlevels; 7874b9ad928SBarry Smith PetscFunctionReturn(0); 7884b9ad928SBarry Smith } 7894b9ad928SBarry Smith 7904b9ad928SBarry Smith #undef __FUNCT__ 7919dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetType" 7924b9ad928SBarry Smith /*@ 79397177400SBarry Smith PCMGSetType - Determines the form of multigrid to use: 7944b9ad928SBarry Smith multiplicative, additive, full, or the Kaskade algorithm. 7954b9ad928SBarry Smith 796ad4df100SBarry Smith Logically Collective on PC 7974b9ad928SBarry Smith 7984b9ad928SBarry Smith Input Parameters: 7994b9ad928SBarry Smith + pc - the preconditioner context 8009dcbbd2bSBarry Smith - form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE, 8019dcbbd2bSBarry Smith PC_MG_FULL, PC_MG_KASKADE 8024b9ad928SBarry Smith 8034b9ad928SBarry Smith Options Database Key: 8044b9ad928SBarry Smith . -pc_mg_type <form> - Sets <form>, one of multiplicative, 8054b9ad928SBarry Smith additive, full, kaskade 8064b9ad928SBarry Smith 8074b9ad928SBarry Smith Level: advanced 8084b9ad928SBarry Smith 8094b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid 8104b9ad928SBarry Smith 81197177400SBarry Smith .seealso: PCMGSetLevels() 8124b9ad928SBarry Smith @*/ 8137087cfbeSBarry Smith PetscErrorCode PCMGSetType(PC pc,PCMGType form) 8144b9ad928SBarry Smith { 815f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 8164b9ad928SBarry Smith 8174b9ad928SBarry Smith PetscFunctionBegin; 8180700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 819c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc,form,2); 82031567311SBarry Smith mg->am = form; 8219dcbbd2bSBarry Smith if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG; 8224b9ad928SBarry Smith else pc->ops->applyrichardson = 0; 8234b9ad928SBarry Smith PetscFunctionReturn(0); 8244b9ad928SBarry Smith } 8254b9ad928SBarry Smith 8264b9ad928SBarry Smith #undef __FUNCT__ 8270d353602SBarry Smith #define __FUNCT__ "PCMGSetCycleType" 8284b9ad928SBarry Smith /*@ 8290d353602SBarry Smith PCMGSetCycleType - Sets the type cycles to use. Use PCMGSetCycleTypeOnLevel() for more 8304b9ad928SBarry Smith complicated cycling. 8314b9ad928SBarry Smith 832ad4df100SBarry Smith Logically Collective on PC 8334b9ad928SBarry Smith 8344b9ad928SBarry Smith Input Parameters: 835c2be2410SBarry Smith + pc - the multigrid context 8360d353602SBarry Smith - PC_MG_CYCLE_V or PC_MG_CYCLE_W 8374b9ad928SBarry Smith 8384b9ad928SBarry Smith Options Database Key: 8390d353602SBarry Smith $ -pc_mg_cycle_type v or w 8404b9ad928SBarry Smith 8414b9ad928SBarry Smith Level: advanced 8424b9ad928SBarry Smith 8434b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 8444b9ad928SBarry Smith 8450d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel() 8464b9ad928SBarry Smith @*/ 8477087cfbeSBarry Smith PetscErrorCode PCMGSetCycleType(PC pc,PCMGCycleType n) 8484b9ad928SBarry Smith { 849f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 850f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 85179416396SBarry Smith PetscInt i,levels; 8524b9ad928SBarry Smith 8534b9ad928SBarry Smith PetscFunctionBegin; 8540700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 855e7e72b3dSBarry Smith if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 856c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 857f3fbd535SBarry Smith levels = mglevels[0]->levels; 8584b9ad928SBarry Smith 8594b9ad928SBarry Smith for (i=0; i<levels; i++) { 860f3fbd535SBarry Smith mglevels[i]->cycles = n; 8614b9ad928SBarry Smith } 8624b9ad928SBarry Smith PetscFunctionReturn(0); 8634b9ad928SBarry Smith } 8644b9ad928SBarry Smith 8654b9ad928SBarry Smith #undef __FUNCT__ 8668cc2d5dfSBarry Smith #define __FUNCT__ "PCMGMultiplicativeSetCycles" 8678cc2d5dfSBarry Smith /*@ 8688cc2d5dfSBarry Smith PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step 8698cc2d5dfSBarry Smith of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used 8708cc2d5dfSBarry Smith 871ad4df100SBarry Smith Logically Collective on PC 8728cc2d5dfSBarry Smith 8738cc2d5dfSBarry Smith Input Parameters: 8748cc2d5dfSBarry Smith + pc - the multigrid context 8758cc2d5dfSBarry Smith - n - number of cycles (default is 1) 8768cc2d5dfSBarry Smith 8778cc2d5dfSBarry Smith Options Database Key: 8788cc2d5dfSBarry Smith $ -pc_mg_multiplicative_cycles n 8798cc2d5dfSBarry Smith 8808cc2d5dfSBarry Smith Level: advanced 8818cc2d5dfSBarry Smith 8828cc2d5dfSBarry Smith Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType() 8838cc2d5dfSBarry Smith 8848cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 8858cc2d5dfSBarry Smith 8868cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType() 8878cc2d5dfSBarry Smith @*/ 8887087cfbeSBarry Smith PetscErrorCode PCMGMultiplicativeSetCycles(PC pc,PetscInt n) 8898cc2d5dfSBarry Smith { 890f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 891f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 8928cc2d5dfSBarry Smith PetscInt i,levels; 8938cc2d5dfSBarry Smith 8948cc2d5dfSBarry Smith PetscFunctionBegin; 8950700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 896e7e72b3dSBarry Smith if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 897c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 898f3fbd535SBarry Smith levels = mglevels[0]->levels; 8998cc2d5dfSBarry Smith 9008cc2d5dfSBarry Smith for (i=0; i<levels; i++) { 90131567311SBarry Smith mg->cyclesperpcapply = n; 9028cc2d5dfSBarry Smith } 9038cc2d5dfSBarry Smith PetscFunctionReturn(0); 9048cc2d5dfSBarry Smith } 9058cc2d5dfSBarry Smith 9068cc2d5dfSBarry Smith #undef __FUNCT__ 9079dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetGalerkin" 908c2be2410SBarry Smith /*@ 90997177400SBarry Smith PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the 910c2be2410SBarry Smith finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t 911c2be2410SBarry Smith 912ad4df100SBarry Smith Logically Collective on PC 913c2be2410SBarry Smith 914c2be2410SBarry Smith Input Parameters: 915c91913d3SJed Brown + pc - the multigrid context 916c91913d3SJed Brown - use - PETSC_TRUE to use the Galerkin process to compute coarse-level operators 917c2be2410SBarry Smith 918c2be2410SBarry Smith Options Database Key: 919c2be2410SBarry Smith $ -pc_mg_galerkin 920c2be2410SBarry Smith 921c2be2410SBarry Smith Level: intermediate 922c2be2410SBarry Smith 923c2be2410SBarry Smith .keywords: MG, set, Galerkin 924c2be2410SBarry Smith 9253fc8bf9cSBarry Smith .seealso: PCMGGetGalerkin() 9263fc8bf9cSBarry Smith 927c2be2410SBarry Smith @*/ 928c91913d3SJed Brown PetscErrorCode PCMGSetGalerkin(PC pc,PetscBool use) 929c2be2410SBarry Smith { 930f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 931c2be2410SBarry Smith 932c2be2410SBarry Smith PetscFunctionBegin; 9330700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 934c91913d3SJed Brown mg->galerkin = (PetscInt)use; 935c2be2410SBarry Smith PetscFunctionReturn(0); 936c2be2410SBarry Smith } 937c2be2410SBarry Smith 938c2be2410SBarry Smith #undef __FUNCT__ 9393fc8bf9cSBarry Smith #define __FUNCT__ "PCMGGetGalerkin" 9403fc8bf9cSBarry Smith /*@ 9413fc8bf9cSBarry Smith PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. 9423fc8bf9cSBarry Smith A_i-1 = r_i * A_i * r_i^t 9433fc8bf9cSBarry Smith 9443fc8bf9cSBarry Smith Not Collective 9453fc8bf9cSBarry Smith 9463fc8bf9cSBarry Smith Input Parameter: 9473fc8bf9cSBarry Smith . pc - the multigrid context 9483fc8bf9cSBarry Smith 9493fc8bf9cSBarry Smith Output Parameter: 9503fc8bf9cSBarry Smith . gelerkin - PETSC_TRUE or PETSC_FALSE 9513fc8bf9cSBarry Smith 9523fc8bf9cSBarry Smith Options Database Key: 9533fc8bf9cSBarry Smith $ -pc_mg_galerkin 9543fc8bf9cSBarry Smith 9553fc8bf9cSBarry Smith Level: intermediate 9563fc8bf9cSBarry Smith 9573fc8bf9cSBarry Smith .keywords: MG, set, Galerkin 9583fc8bf9cSBarry Smith 9593fc8bf9cSBarry Smith .seealso: PCMGSetGalerkin() 9603fc8bf9cSBarry Smith 9613fc8bf9cSBarry Smith @*/ 9627087cfbeSBarry Smith PetscErrorCode PCMGGetGalerkin(PC pc,PetscBool *galerkin) 9633fc8bf9cSBarry Smith { 964f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 9653fc8bf9cSBarry Smith 9663fc8bf9cSBarry Smith PetscFunctionBegin; 9670700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 96813fdb3acSJose Roman *galerkin = (PetscBool)mg->galerkin; 9693fc8bf9cSBarry Smith PetscFunctionReturn(0); 9703fc8bf9cSBarry Smith } 9713fc8bf9cSBarry Smith 9723fc8bf9cSBarry Smith #undef __FUNCT__ 9739dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothDown" 9744b9ad928SBarry Smith /*@ 97597177400SBarry Smith PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to 97697177400SBarry Smith use on all levels. Use PCMGGetSmootherDown() to set different 9774b9ad928SBarry Smith pre-smoothing steps on different levels. 9784b9ad928SBarry Smith 979ad4df100SBarry Smith Logically Collective on PC 9804b9ad928SBarry Smith 9814b9ad928SBarry Smith Input Parameters: 9824b9ad928SBarry Smith + mg - the multigrid context 9834b9ad928SBarry Smith - n - the number of smoothing steps 9844b9ad928SBarry Smith 9854b9ad928SBarry Smith Options Database Key: 9864b9ad928SBarry Smith . -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps 9874b9ad928SBarry Smith 9884b9ad928SBarry Smith Level: advanced 9894b9ad928SBarry Smith 9904b9ad928SBarry Smith .keywords: MG, smooth, down, pre-smoothing, steps, multigrid 9914b9ad928SBarry Smith 99297177400SBarry Smith .seealso: PCMGSetNumberSmoothUp() 9934b9ad928SBarry Smith @*/ 9947087cfbeSBarry Smith PetscErrorCode PCMGSetNumberSmoothDown(PC pc,PetscInt n) 9954b9ad928SBarry Smith { 996f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 997f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 9986849ba73SBarry Smith PetscErrorCode ierr; 99979416396SBarry Smith PetscInt i,levels; 10004b9ad928SBarry Smith 10014b9ad928SBarry Smith PetscFunctionBegin; 10020700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1003e7e72b3dSBarry Smith if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 1004c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 1005f3fbd535SBarry Smith levels = mglevels[0]->levels; 10064b9ad928SBarry Smith 1007b05257ddSBarry Smith for (i=1; i<levels; i++) { 10084b9ad928SBarry Smith /* make sure smoother up and down are different */ 100997177400SBarry Smith ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 1010f3fbd535SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 101131567311SBarry Smith mg->default_smoothd = n; 10124b9ad928SBarry Smith } 10134b9ad928SBarry Smith PetscFunctionReturn(0); 10144b9ad928SBarry Smith } 10154b9ad928SBarry Smith 10164b9ad928SBarry Smith #undef __FUNCT__ 10179dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothUp" 10184b9ad928SBarry Smith /*@ 101997177400SBarry Smith PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use 102097177400SBarry Smith on all levels. Use PCMGGetSmootherUp() to set different numbers of 10214b9ad928SBarry Smith post-smoothing steps on different levels. 10224b9ad928SBarry Smith 1023ad4df100SBarry Smith Logically Collective on PC 10244b9ad928SBarry Smith 10254b9ad928SBarry Smith Input Parameters: 10264b9ad928SBarry Smith + mg - the multigrid context 10274b9ad928SBarry Smith - n - the number of smoothing steps 10284b9ad928SBarry Smith 10294b9ad928SBarry Smith Options Database Key: 10304b9ad928SBarry Smith . -pc_mg_smoothup <n> - Sets number of post-smoothing steps 10314b9ad928SBarry Smith 10324b9ad928SBarry Smith Level: advanced 10334b9ad928SBarry Smith 10344b9ad928SBarry Smith Note: this does not set a value on the coarsest grid, since we assume that 1035a8c7a070SBarry Smith there is no separate smooth up on the coarsest grid. 10364b9ad928SBarry Smith 10374b9ad928SBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid 10384b9ad928SBarry Smith 103997177400SBarry Smith .seealso: PCMGSetNumberSmoothDown() 10404b9ad928SBarry Smith @*/ 10417087cfbeSBarry Smith PetscErrorCode PCMGSetNumberSmoothUp(PC pc,PetscInt n) 10424b9ad928SBarry Smith { 1043f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1044f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 10456849ba73SBarry Smith PetscErrorCode ierr; 104679416396SBarry Smith PetscInt i,levels; 10474b9ad928SBarry Smith 10484b9ad928SBarry Smith PetscFunctionBegin; 10490700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1050e7e72b3dSBarry Smith if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 1051c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 1052f3fbd535SBarry Smith levels = mglevels[0]->levels; 10534b9ad928SBarry Smith 10544b9ad928SBarry Smith for (i=1; i<levels; i++) { 10554b9ad928SBarry Smith /* make sure smoother up and down are different */ 105697177400SBarry Smith ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 1057f3fbd535SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 105831567311SBarry Smith mg->default_smoothu = n; 10594b9ad928SBarry Smith } 10604b9ad928SBarry Smith PetscFunctionReturn(0); 10614b9ad928SBarry Smith } 10624b9ad928SBarry Smith 10634b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/ 10644b9ad928SBarry Smith 10653b09bd56SBarry Smith /*MC 1066ccb205f8SBarry Smith PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional 10673b09bd56SBarry Smith information about the coarser grid matrices and restriction/interpolation operators. 10683b09bd56SBarry Smith 10693b09bd56SBarry Smith Options Database Keys: 10703b09bd56SBarry Smith + -pc_mg_levels <nlevels> - number of levels including finest 10710d353602SBarry Smith . -pc_mg_cycles v or w 107279416396SBarry Smith . -pc_mg_smoothup <n> - number of smoothing steps after interpolation 10733b09bd56SBarry Smith . -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator 10743b09bd56SBarry Smith . -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default 10753b09bd56SBarry Smith . -pc_mg_log - log information about time spent on each level of the solver 10763b09bd56SBarry Smith . -pc_mg_monitor - print information on the multigrid convergence 10778cf6037eSBarry Smith . -pc_mg_galerkin - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R' 10788cf6037eSBarry Smith . -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1) 10798cf6037eSBarry Smith . -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 1080e3c5b3baSBarry Smith to the Socket viewer for reading from MATLAB. 10818cf6037eSBarry Smith - -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices 10828cf6037eSBarry Smith to the binary output file called binaryoutput 10833b09bd56SBarry Smith 108424c3aa18SBarry 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 10853b09bd56SBarry Smith 10868cf6037eSBarry Smith When run with a single level the smoother options are used on that level NOT the coarse grid solver options 10878cf6037eSBarry Smith 10883b09bd56SBarry Smith Level: intermediate 10893b09bd56SBarry Smith 10908f87f92bSBarry Smith Concepts: multigrid/multilevel 10913b09bd56SBarry Smith 10928cf6037eSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE 10930d353602SBarry Smith PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(), 109497177400SBarry Smith PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 109597177400SBarry Smith PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 10960d353602SBarry Smith PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 10973b09bd56SBarry Smith M*/ 10983b09bd56SBarry Smith 10994b9ad928SBarry Smith EXTERN_C_BEGIN 11004b9ad928SBarry Smith #undef __FUNCT__ 11014b9ad928SBarry Smith #define __FUNCT__ "PCCreate_MG" 11027087cfbeSBarry Smith PetscErrorCode PCCreate_MG(PC pc) 11034b9ad928SBarry Smith { 1104f3fbd535SBarry Smith PC_MG *mg; 1105f3fbd535SBarry Smith PetscErrorCode ierr; 1106f3fbd535SBarry Smith 11074b9ad928SBarry Smith PetscFunctionBegin; 1108f3fbd535SBarry Smith ierr = PetscNewLog(pc,PC_MG,&mg);CHKERRQ(ierr); 1109f3fbd535SBarry Smith pc->data = (void*)mg; 1110f3fbd535SBarry Smith mg->nlevels = -1; 1111f3fbd535SBarry Smith 11124b9ad928SBarry Smith pc->ops->apply = PCApply_MG; 11134b9ad928SBarry Smith pc->ops->setup = PCSetUp_MG; 1114a06653b4SBarry Smith pc->ops->reset = PCReset_MG; 11154b9ad928SBarry Smith pc->ops->destroy = PCDestroy_MG; 11164b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_MG; 11174b9ad928SBarry Smith pc->ops->view = PCView_MG; 11184b9ad928SBarry Smith PetscFunctionReturn(0); 11194b9ad928SBarry Smith } 11204b9ad928SBarry Smith EXTERN_C_END 1121