1dba47a55SKris Buschelman 24b9ad928SBarry Smith /* 34b9ad928SBarry Smith Defines the multigrid preconditioner interface. 44b9ad928SBarry Smith */ 5af0996ceSBarry Smith #include <petsc/private/pcmgimpl.h> /*I "petscksp.h" I*/ 61e25c274SJed Brown #include <petscdm.h> 74b9ad928SBarry Smith 84b9ad928SBarry Smith #undef __FUNCT__ 99dcbbd2bSBarry Smith #define __FUNCT__ "PCMGMCycle_Private" 1031567311SBarry Smith PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG_Levels **mglevelsin,PCRichardsonConvergedReason *reason) 114b9ad928SBarry Smith { 1231567311SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1331567311SBarry Smith PC_MG_Levels *mgc,*mglevels = *mglevelsin; 146849ba73SBarry Smith PetscErrorCode ierr; 1531567311SBarry Smith PetscInt cycles = (mglevels->level == 1) ? 1 : (PetscInt) mglevels->cycles; 16*a0808db4SHong Zhang PC subpc; 17*a0808db4SHong Zhang PCFailedReason pcreason; 184b9ad928SBarry Smith 194b9ad928SBarry Smith PetscFunctionBegin; 2063e6d426SJed Brown if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 2131567311SBarry Smith ierr = KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x);CHKERRQ(ierr); /* pre-smooth */ 22*a0808db4SHong Zhang ierr = KSPGetPC(mglevels->smoothd,&subpc);CHKERRQ(ierr); 23*a0808db4SHong Zhang ierr = PCGetSetUpFailedReason(subpc,&pcreason);CHKERRQ(ierr); 24*a0808db4SHong Zhang if (pcreason) { 25*a0808db4SHong Zhang pc->failedreason = PC_SUBPC_ERROR; 26*a0808db4SHong Zhang } 2763e6d426SJed Brown if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 2831567311SBarry Smith if (mglevels->level) { /* not the coarsest grid */ 2963e6d426SJed Brown if (mglevels->eventresidual) {ierr = PetscLogEventBegin(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);} 3031567311SBarry Smith ierr = (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);CHKERRQ(ierr); 3163e6d426SJed Brown if (mglevels->eventresidual) {ierr = PetscLogEventEnd(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);} 324b9ad928SBarry Smith 334b9ad928SBarry Smith /* if on finest level and have convergence criteria set */ 3431567311SBarry Smith if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) { 354b9ad928SBarry Smith PetscReal rnorm; 3631567311SBarry Smith ierr = VecNorm(mglevels->r,NORM_2,&rnorm);CHKERRQ(ierr); 374b9ad928SBarry Smith if (rnorm <= mg->ttol) { 3870441072SBarry Smith if (rnorm < mg->abstol) { 394d0a8057SBarry Smith *reason = PCRICHARDSON_CONVERGED_ATOL; 4057622a8eSBarry Smith ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %g is less than absolute tolerance %g\n",(double)rnorm,(double)mg->abstol);CHKERRQ(ierr); 414b9ad928SBarry Smith } else { 424d0a8057SBarry Smith *reason = PCRICHARDSON_CONVERGED_RTOL; 4357622a8eSBarry Smith ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %g is less than relative tolerance times initial residual norm %g\n",(double)rnorm,(double)mg->ttol);CHKERRQ(ierr); 444b9ad928SBarry Smith } 454b9ad928SBarry Smith PetscFunctionReturn(0); 464b9ad928SBarry Smith } 474b9ad928SBarry Smith } 484b9ad928SBarry Smith 4931567311SBarry Smith mgc = *(mglevelsin - 1); 5063e6d426SJed Brown if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 5131567311SBarry Smith ierr = MatRestrict(mglevels->restrct,mglevels->r,mgc->b);CHKERRQ(ierr); 5263e6d426SJed Brown if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 53efb30889SBarry Smith ierr = VecSet(mgc->x,0.0);CHKERRQ(ierr); 544b9ad928SBarry Smith while (cycles--) { 5531567311SBarry Smith ierr = PCMGMCycle_Private(pc,mglevelsin-1,reason);CHKERRQ(ierr); 564b9ad928SBarry Smith } 5763e6d426SJed Brown if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 5831567311SBarry Smith ierr = MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x);CHKERRQ(ierr); 5963e6d426SJed Brown if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 6063e6d426SJed Brown if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 6131567311SBarry Smith ierr = KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x);CHKERRQ(ierr); /* post smooth */ 6263e6d426SJed Brown if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 634b9ad928SBarry Smith } 644b9ad928SBarry Smith PetscFunctionReturn(0); 654b9ad928SBarry Smith } 664b9ad928SBarry Smith 674b9ad928SBarry Smith #undef __FUNCT__ 684b9ad928SBarry Smith #define __FUNCT__ "PCApplyRichardson_MG" 69ace3abfcSBarry 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) 704b9ad928SBarry Smith { 71f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 72f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 73dfbe8321SBarry Smith PetscErrorCode ierr; 74f3fbd535SBarry Smith PetscInt levels = mglevels[0]->levels,i; 754b9ad928SBarry Smith 764b9ad928SBarry Smith PetscFunctionBegin; 77a762d673SBarry Smith /* When the DM is supplying the matrix then it will not exist until here */ 78a762d673SBarry Smith for (i=0; i<levels; i++) { 79a762d673SBarry Smith if (!mglevels[i]->A) { 80a762d673SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr); 81a762d673SBarry Smith ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr); 82a762d673SBarry Smith } 83a762d673SBarry Smith } 84f3fbd535SBarry Smith mglevels[levels-1]->b = b; 85f3fbd535SBarry Smith mglevels[levels-1]->x = x; 864b9ad928SBarry Smith 8731567311SBarry Smith mg->rtol = rtol; 8831567311SBarry Smith mg->abstol = abstol; 8931567311SBarry Smith mg->dtol = dtol; 904b9ad928SBarry Smith if (rtol) { 914b9ad928SBarry Smith /* compute initial residual norm for relative convergence test */ 924b9ad928SBarry Smith PetscReal rnorm; 937319c654SBarry Smith if (zeroguess) { 947319c654SBarry Smith ierr = VecNorm(b,NORM_2,&rnorm);CHKERRQ(ierr); 957319c654SBarry Smith } else { 96f3fbd535SBarry Smith ierr = (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);CHKERRQ(ierr); 974b9ad928SBarry Smith ierr = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr); 987319c654SBarry Smith } 9931567311SBarry Smith mg->ttol = PetscMax(rtol*rnorm,abstol); 1002fa5cd67SKarl Rupp } else if (abstol) mg->ttol = abstol; 1012fa5cd67SKarl Rupp else mg->ttol = 0.0; 1024b9ad928SBarry Smith 1034d0a8057SBarry Smith /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't 104336babb1SJed Brown stop prematurely due to small residual */ 1054d0a8057SBarry Smith for (i=1; i<levels; i++) { 106f3fbd535SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 107f3fbd535SBarry Smith if (mglevels[i]->smoothu != mglevels[i]->smoothd) { 108f3fbd535SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 1094b9ad928SBarry Smith } 1104d0a8057SBarry Smith } 1114d0a8057SBarry Smith 1124d0a8057SBarry Smith *reason = (PCRichardsonConvergedReason)0; 1134d0a8057SBarry Smith for (i=0; i<its; i++) { 114f3fbd535SBarry Smith ierr = PCMGMCycle_Private(pc,mglevels+levels-1,reason);CHKERRQ(ierr); 1154d0a8057SBarry Smith if (*reason) break; 1164d0a8057SBarry Smith } 1174d0a8057SBarry Smith if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS; 1184d0a8057SBarry Smith *outits = i; 1194b9ad928SBarry Smith PetscFunctionReturn(0); 1204b9ad928SBarry Smith } 1214b9ad928SBarry Smith 1224b9ad928SBarry Smith #undef __FUNCT__ 1233aeaf226SBarry Smith #define __FUNCT__ "PCReset_MG" 1243aeaf226SBarry Smith PetscErrorCode PCReset_MG(PC pc) 1253aeaf226SBarry Smith { 1263aeaf226SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1273aeaf226SBarry Smith PC_MG_Levels **mglevels = mg->levels; 1283aeaf226SBarry Smith PetscErrorCode ierr; 1293aeaf226SBarry Smith PetscInt i,n; 1303aeaf226SBarry Smith 1313aeaf226SBarry Smith PetscFunctionBegin; 1323aeaf226SBarry Smith if (mglevels) { 1333aeaf226SBarry Smith n = mglevels[0]->levels; 1343aeaf226SBarry Smith for (i=0; i<n-1; i++) { 1353aeaf226SBarry Smith ierr = VecDestroy(&mglevels[i+1]->r);CHKERRQ(ierr); 1363aeaf226SBarry Smith ierr = VecDestroy(&mglevels[i]->b);CHKERRQ(ierr); 1373aeaf226SBarry Smith ierr = VecDestroy(&mglevels[i]->x);CHKERRQ(ierr); 1383aeaf226SBarry Smith ierr = MatDestroy(&mglevels[i+1]->restrct);CHKERRQ(ierr); 1393aeaf226SBarry Smith ierr = MatDestroy(&mglevels[i+1]->interpolate);CHKERRQ(ierr); 14073250ac0SBarry Smith ierr = VecDestroy(&mglevels[i+1]->rscale);CHKERRQ(ierr); 1413aeaf226SBarry Smith } 1423aeaf226SBarry Smith 1433aeaf226SBarry Smith for (i=0; i<n; i++) { 1443aeaf226SBarry Smith ierr = MatDestroy(&mglevels[i]->A);CHKERRQ(ierr); 1453aeaf226SBarry Smith if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 1463aeaf226SBarry Smith ierr = KSPReset(mglevels[i]->smoothd);CHKERRQ(ierr); 1473aeaf226SBarry Smith } 1483aeaf226SBarry Smith ierr = KSPReset(mglevels[i]->smoothu);CHKERRQ(ierr); 1493aeaf226SBarry Smith } 1503aeaf226SBarry Smith } 1513aeaf226SBarry Smith PetscFunctionReturn(0); 1523aeaf226SBarry Smith } 1533aeaf226SBarry Smith 1543aeaf226SBarry Smith #undef __FUNCT__ 1559dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetLevels" 1564b9ad928SBarry Smith /*@C 15797177400SBarry Smith PCMGSetLevels - Sets the number of levels to use with MG. 1584b9ad928SBarry Smith Must be called before any other MG routine. 1594b9ad928SBarry Smith 160ad4df100SBarry Smith Logically Collective on PC 1614b9ad928SBarry Smith 1624b9ad928SBarry Smith Input Parameters: 1634b9ad928SBarry Smith + pc - the preconditioner context 1644b9ad928SBarry Smith . levels - the number of levels 1654b9ad928SBarry Smith - comms - optional communicators for each level; this is to allow solving the coarser problems 1660298fd71SBarry Smith on smaller sets of processors. Use NULL_OBJECT for default in Fortran 1674b9ad928SBarry Smith 1684b9ad928SBarry Smith Level: intermediate 1694b9ad928SBarry Smith 1704b9ad928SBarry Smith Notes: 1714b9ad928SBarry Smith If the number of levels is one then the multigrid uses the -mg_levels prefix 1724b9ad928SBarry Smith for setting the level options rather than the -mg_coarse prefix. 1734b9ad928SBarry Smith 1744b9ad928SBarry Smith .keywords: MG, set, levels, multigrid 1754b9ad928SBarry Smith 17697177400SBarry Smith .seealso: PCMGSetType(), PCMGGetLevels() 1774b9ad928SBarry Smith @*/ 1787087cfbeSBarry Smith PetscErrorCode PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms) 1794b9ad928SBarry Smith { 180dfbe8321SBarry Smith PetscErrorCode ierr; 181f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 182ce94432eSBarry Smith MPI_Comm comm; 1833aeaf226SBarry Smith PC_MG_Levels **mglevels = mg->levels; 18410eca3edSLisandro Dalcin PCMGType mgtype = mg->am; 18510167fecSBarry Smith PetscInt mgctype = (PetscInt) PC_MG_CYCLE_V; 186f3fbd535SBarry Smith PetscInt i; 187f3fbd535SBarry Smith PetscMPIInt size; 188f3fbd535SBarry Smith const char *prefix; 189f3fbd535SBarry Smith PC ipc; 1903aeaf226SBarry Smith PetscInt n; 1914b9ad928SBarry Smith 1924b9ad928SBarry Smith PetscFunctionBegin; 1930700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 194548767e0SJed Brown PetscValidLogicalCollectiveInt(pc,levels,2); 195ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 196548767e0SJed Brown if (mg->nlevels == levels) PetscFunctionReturn(0); 1973aeaf226SBarry Smith if (mglevels) { 19810eca3edSLisandro Dalcin mgctype = mglevels[0]->cycles; 1993aeaf226SBarry Smith /* changing the number of levels so free up the previous stuff */ 2003aeaf226SBarry Smith ierr = PCReset_MG(pc);CHKERRQ(ierr); 2013aeaf226SBarry Smith n = mglevels[0]->levels; 2023aeaf226SBarry Smith for (i=0; i<n; i++) { 2033aeaf226SBarry Smith if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 2043aeaf226SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr); 2053aeaf226SBarry Smith } 2063aeaf226SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr); 2073aeaf226SBarry Smith ierr = PetscFree(mglevels[i]);CHKERRQ(ierr); 2083aeaf226SBarry Smith } 2093aeaf226SBarry Smith ierr = PetscFree(mg->levels);CHKERRQ(ierr); 2103aeaf226SBarry Smith } 211f3fbd535SBarry Smith 212f3fbd535SBarry Smith mg->nlevels = levels; 213f3fbd535SBarry Smith 214785e854fSJed Brown ierr = PetscMalloc1(levels,&mglevels);CHKERRQ(ierr); 2153bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr); 216f3fbd535SBarry Smith 217f3fbd535SBarry Smith ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 218f3fbd535SBarry Smith 219a9db87a2SMatthew G Knepley mg->stageApply = 0; 220f3fbd535SBarry Smith for (i=0; i<levels; i++) { 221b00a9115SJed Brown ierr = PetscNewLog(pc,&mglevels[i]);CHKERRQ(ierr); 2222fa5cd67SKarl Rupp 223f3fbd535SBarry Smith mglevels[i]->level = i; 224f3fbd535SBarry Smith mglevels[i]->levels = levels; 22510eca3edSLisandro Dalcin mglevels[i]->cycles = mgctype; 226336babb1SJed Brown mg->default_smoothu = 2; 227336babb1SJed Brown mg->default_smoothd = 2; 22863e6d426SJed Brown mglevels[i]->eventsmoothsetup = 0; 22963e6d426SJed Brown mglevels[i]->eventsmoothsolve = 0; 23063e6d426SJed Brown mglevels[i]->eventresidual = 0; 23163e6d426SJed Brown mglevels[i]->eventinterprestrict = 0; 232f3fbd535SBarry Smith 233f3fbd535SBarry Smith if (comms) comm = comms[i]; 234f3fbd535SBarry Smith ierr = KSPCreate(comm,&mglevels[i]->smoothd);CHKERRQ(ierr); 235422a814eSBarry Smith ierr = KSPSetErrorIfNotConverged(mglevels[i]->smoothd,pc->erroriffailure);CHKERRQ(ierr); 236336babb1SJed Brown ierr = KSPSetType(mglevels[i]->smoothd,KSPCHEBYSHEV);CHKERRQ(ierr); 2370059c7bdSJed Brown ierr = KSPSetConvergenceTest(mglevels[i]->smoothd,KSPConvergedSkip,NULL,NULL);CHKERRQ(ierr); 238336babb1SJed Brown ierr = KSPSetNormType(mglevels[i]->smoothd,KSP_NORM_NONE);CHKERRQ(ierr); 239336babb1SJed Brown ierr = KSPGetPC(mglevels[i]->smoothd,&ipc);CHKERRQ(ierr); 240336babb1SJed Brown ierr = PCSetType(ipc,PCSOR);CHKERRQ(ierr); 241f3fbd535SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr); 242336babb1SJed Brown ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, i ? mg->default_smoothd : 1);CHKERRQ(ierr); 243f3fbd535SBarry Smith ierr = KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);CHKERRQ(ierr); 244ab83eea4SMatthew G. Knepley ierr = PetscObjectComposedDataSetInt((PetscObject) mglevels[i]->smoothd, PetscMGLevelId, mglevels[i]->level);CHKERRQ(ierr); 245f3fbd535SBarry Smith 246f3fbd535SBarry Smith /* do special stuff for coarse grid */ 247f3fbd535SBarry Smith if (!i && levels > 1) { 248f3fbd535SBarry Smith ierr = KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");CHKERRQ(ierr); 249f3fbd535SBarry Smith 2509dbfc187SHong Zhang /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */ 251f3fbd535SBarry Smith ierr = KSPSetType(mglevels[0]->smoothd,KSPPREONLY);CHKERRQ(ierr); 252f3fbd535SBarry Smith ierr = KSPGetPC(mglevels[0]->smoothd,&ipc);CHKERRQ(ierr); 253f3fbd535SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 254f3fbd535SBarry Smith if (size > 1) { 2559dbfc187SHong Zhang KSP innerksp; 2569dbfc187SHong Zhang PC innerpc; 257f3fbd535SBarry Smith ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr); 2589dbfc187SHong Zhang ierr = PCRedundantGetKSP(ipc,&innerksp);CHKERRQ(ierr); 2599dbfc187SHong Zhang ierr = KSPGetPC(innerksp,&innerpc);CHKERRQ(ierr); 2609dbfc187SHong Zhang ierr = PCFactorSetShiftType(innerpc,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr); 261f3fbd535SBarry Smith } else { 262f3fbd535SBarry Smith ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr); 2639dbfc187SHong Zhang ierr = PCFactorSetShiftType(ipc,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr); 264f3fbd535SBarry Smith } 265f3fbd535SBarry Smith } else { 266f3fbd535SBarry Smith char tprefix[128]; 267f3fbd535SBarry Smith sprintf(tprefix,"mg_levels_%d_",(int)i); 268f3fbd535SBarry Smith ierr = KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);CHKERRQ(ierr); 269f3fbd535SBarry Smith } 2703bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[i]->smoothd);CHKERRQ(ierr); 2712fa5cd67SKarl Rupp 272f3fbd535SBarry Smith mglevels[i]->smoothu = mglevels[i]->smoothd; 27331567311SBarry Smith mg->rtol = 0.0; 27431567311SBarry Smith mg->abstol = 0.0; 27531567311SBarry Smith mg->dtol = 0.0; 27631567311SBarry Smith mg->ttol = 0.0; 27731567311SBarry Smith mg->cyclesperpcapply = 1; 278f3fbd535SBarry Smith } 279f3fbd535SBarry Smith mg->levels = mglevels; 28010eca3edSLisandro Dalcin ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr); 2814b9ad928SBarry Smith PetscFunctionReturn(0); 2824b9ad928SBarry Smith } 2834b9ad928SBarry Smith 284c07bf074SBarry Smith 285c07bf074SBarry Smith #undef __FUNCT__ 286c07bf074SBarry Smith #define __FUNCT__ "PCDestroy_MG" 287c07bf074SBarry Smith PetscErrorCode PCDestroy_MG(PC pc) 288c07bf074SBarry Smith { 289c07bf074SBarry Smith PetscErrorCode ierr; 290a06653b4SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 291a06653b4SBarry Smith PC_MG_Levels **mglevels = mg->levels; 292a06653b4SBarry Smith PetscInt i,n; 293c07bf074SBarry Smith 294c07bf074SBarry Smith PetscFunctionBegin; 295a06653b4SBarry Smith ierr = PCReset_MG(pc);CHKERRQ(ierr); 296a06653b4SBarry Smith if (mglevels) { 297a06653b4SBarry Smith n = mglevels[0]->levels; 298a06653b4SBarry Smith for (i=0; i<n; i++) { 299a06653b4SBarry Smith if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 3006bf464f9SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr); 301a06653b4SBarry Smith } 3026bf464f9SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr); 303a06653b4SBarry Smith ierr = PetscFree(mglevels[i]);CHKERRQ(ierr); 304a06653b4SBarry Smith } 305a06653b4SBarry Smith ierr = PetscFree(mg->levels);CHKERRQ(ierr); 306a06653b4SBarry Smith } 307c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 308f3fbd535SBarry Smith PetscFunctionReturn(0); 309f3fbd535SBarry Smith } 310f3fbd535SBarry Smith 311f3fbd535SBarry Smith 312f3fbd535SBarry Smith 31309573ac7SBarry Smith extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**); 31409573ac7SBarry Smith extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**); 31509573ac7SBarry Smith extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**); 316f3fbd535SBarry Smith 317f3fbd535SBarry Smith /* 318f3fbd535SBarry Smith PCApply_MG - Runs either an additive, multiplicative, Kaskadic 319f3fbd535SBarry Smith or full cycle of multigrid. 320f3fbd535SBarry Smith 321f3fbd535SBarry Smith Note: 322f3fbd535SBarry Smith A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle(). 323f3fbd535SBarry Smith */ 324f3fbd535SBarry Smith #undef __FUNCT__ 325f3fbd535SBarry Smith #define __FUNCT__ "PCApply_MG" 326f3fbd535SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x) 327f3fbd535SBarry Smith { 328f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 329f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 330f3fbd535SBarry Smith PetscErrorCode ierr; 331f3fbd535SBarry Smith PetscInt levels = mglevels[0]->levels,i; 332f3fbd535SBarry Smith 333f3fbd535SBarry Smith PetscFunctionBegin; 334a9db87a2SMatthew G Knepley if (mg->stageApply) {ierr = PetscLogStagePush(mg->stageApply);CHKERRQ(ierr);} 335e1d8e5deSBarry Smith /* When the DM is supplying the matrix then it will not exist until here */ 336298cc208SBarry Smith for (i=0; i<levels; i++) { 337e1d8e5deSBarry Smith if (!mglevels[i]->A) { 33823ee1639SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr); 339298cc208SBarry Smith ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr); 340e1d8e5deSBarry Smith } 341e1d8e5deSBarry Smith } 342e1d8e5deSBarry Smith 343f3fbd535SBarry Smith mglevels[levels-1]->b = b; 344f3fbd535SBarry Smith mglevels[levels-1]->x = x; 34531567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 346f3fbd535SBarry Smith ierr = VecSet(x,0.0);CHKERRQ(ierr); 34731567311SBarry Smith for (i=0; i<mg->cyclesperpcapply; i++) { 3480298fd71SBarry Smith ierr = PCMGMCycle_Private(pc,mglevels+levels-1,NULL);CHKERRQ(ierr); 349f3fbd535SBarry Smith } 3502fa5cd67SKarl Rupp } else if (mg->am == PC_MG_ADDITIVE) { 35131567311SBarry Smith ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr); 3522fa5cd67SKarl Rupp } else if (mg->am == PC_MG_KASKADE) { 35331567311SBarry Smith ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr); 3542fa5cd67SKarl Rupp } else { 355f3fbd535SBarry Smith ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr); 356f3fbd535SBarry Smith } 357a9db87a2SMatthew G Knepley if (mg->stageApply) {ierr = PetscLogStagePop();CHKERRQ(ierr);} 358f3fbd535SBarry Smith PetscFunctionReturn(0); 359f3fbd535SBarry Smith } 360f3fbd535SBarry Smith 361f3fbd535SBarry Smith 362f3fbd535SBarry Smith #undef __FUNCT__ 363f3fbd535SBarry Smith #define __FUNCT__ "PCSetFromOptions_MG" 3644416b707SBarry Smith PetscErrorCode PCSetFromOptions_MG(PetscOptionItems *PetscOptionsObject,PC pc) 365f3fbd535SBarry Smith { 366f3fbd535SBarry Smith PetscErrorCode ierr; 367f3fbd535SBarry Smith PetscInt m,levels = 1,cycles; 368c91913d3SJed Brown PetscBool flg,set; 369f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 370f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 371f3fbd535SBarry Smith PCMGType mgtype; 372f3fbd535SBarry Smith PCMGCycleType mgctype; 373f3fbd535SBarry Smith 374f3fbd535SBarry Smith PetscFunctionBegin; 375e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Multigrid options");CHKERRQ(ierr); 3763aeaf226SBarry Smith if (!mg->levels) { 377f3fbd535SBarry Smith ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr); 378eb3f98d2SBarry Smith if (!flg && pc->dm) { 379eb3f98d2SBarry Smith ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr); 380eb3f98d2SBarry Smith levels++; 3813aeaf226SBarry Smith mg->usedmfornumberoflevels = PETSC_TRUE; 382eb3f98d2SBarry Smith } 3830298fd71SBarry Smith ierr = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr); 384f3fbd535SBarry Smith } 3853aeaf226SBarry Smith mglevels = mg->levels; 3863aeaf226SBarry Smith 387f3fbd535SBarry Smith mgctype = (PCMGCycleType) mglevels[0]->cycles; 388f3fbd535SBarry Smith ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr); 389f3fbd535SBarry Smith if (flg) { 390f3fbd535SBarry Smith ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr); 3912fa5cd67SKarl Rupp } 392f3fbd535SBarry Smith flg = PETSC_FALSE; 393c91913d3SJed Brown ierr = PetscOptionsBool("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",flg,&flg,&set);CHKERRQ(ierr); 394c91913d3SJed Brown if (set) { 395c91913d3SJed Brown ierr = PCMGSetGalerkin(pc,flg);CHKERRQ(ierr); 396f3fbd535SBarry Smith } 39794ae4db5SBarry Smith ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",mg->default_smoothu,&m,&flg);CHKERRQ(ierr); 398f3fbd535SBarry Smith if (flg) { 399f3fbd535SBarry Smith ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr); 400f3fbd535SBarry Smith } 40194ae4db5SBarry Smith ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",mg->default_smoothd,&m,&flg);CHKERRQ(ierr); 402f3fbd535SBarry Smith if (flg) { 403f3fbd535SBarry Smith ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr); 404f3fbd535SBarry Smith } 40531567311SBarry Smith mgtype = mg->am; 406f3fbd535SBarry Smith ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr); 407f3fbd535SBarry Smith if (flg) { 408f3fbd535SBarry Smith ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr); 409f3fbd535SBarry Smith } 41031567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 4115363c1e0SLisandro Dalcin ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGMultiplicativeSetCycles",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr); 412f3fbd535SBarry Smith if (flg) { 413f3fbd535SBarry Smith ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr); 414f3fbd535SBarry Smith } 415f3fbd535SBarry Smith } 416f3fbd535SBarry Smith flg = PETSC_FALSE; 4170298fd71SBarry Smith ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,NULL);CHKERRQ(ierr); 418f3fbd535SBarry Smith if (flg) { 419f3fbd535SBarry Smith PetscInt i; 420f3fbd535SBarry Smith char eventname[128]; 421ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 422f3fbd535SBarry Smith levels = mglevels[0]->levels; 423f3fbd535SBarry Smith for (i=0; i<levels; i++) { 424f3fbd535SBarry Smith sprintf(eventname,"MGSetup Level %d",(int)i); 42563e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr); 426f3fbd535SBarry Smith sprintf(eventname,"MGSmooth Level %d",(int)i); 42763e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr); 428f3fbd535SBarry Smith if (i) { 429f3fbd535SBarry Smith sprintf(eventname,"MGResid Level %d",(int)i); 43063e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr); 431f3fbd535SBarry Smith sprintf(eventname,"MGInterp Level %d",(int)i); 43263e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr); 433f3fbd535SBarry Smith } 434f3fbd535SBarry Smith } 435eec35531SMatthew G Knepley 436ec60ca73SSatish Balay #if defined(PETSC_USE_LOG) 437eec35531SMatthew G Knepley { 438eec35531SMatthew G Knepley const char *sname = "MG Apply"; 439eec35531SMatthew G Knepley PetscStageLog stageLog; 440eec35531SMatthew G Knepley PetscInt st; 441eec35531SMatthew G Knepley 442eec35531SMatthew G Knepley PetscFunctionBegin; 443eec35531SMatthew G Knepley ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr); 444eec35531SMatthew G Knepley for (st = 0; st < stageLog->numStages; ++st) { 445eec35531SMatthew G Knepley PetscBool same; 446eec35531SMatthew G Knepley 447eec35531SMatthew G Knepley ierr = PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);CHKERRQ(ierr); 4482fa5cd67SKarl Rupp if (same) mg->stageApply = st; 449eec35531SMatthew G Knepley } 450eec35531SMatthew G Knepley if (!mg->stageApply) { 451eec35531SMatthew G Knepley ierr = PetscLogStageRegister(sname, &mg->stageApply);CHKERRQ(ierr); 452eec35531SMatthew G Knepley } 453eec35531SMatthew G Knepley } 454ec60ca73SSatish Balay #endif 455f3fbd535SBarry Smith } 456f3fbd535SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 457f3fbd535SBarry Smith PetscFunctionReturn(0); 458f3fbd535SBarry Smith } 459f3fbd535SBarry Smith 4606a6fc655SJed Brown const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0}; 4616a6fc655SJed Brown const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0}; 462f3fbd535SBarry Smith 4639804daf3SBarry Smith #include <petscdraw.h> 464f3fbd535SBarry Smith #undef __FUNCT__ 465f3fbd535SBarry Smith #define __FUNCT__ "PCView_MG" 466f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer) 467f3fbd535SBarry Smith { 468f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 469f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 470f3fbd535SBarry Smith PetscErrorCode ierr; 471e3deeeafSJed Brown PetscInt levels = mglevels ? mglevels[0]->levels : 0,i; 472219b1045SBarry Smith PetscBool iascii,isbinary,isdraw; 473f3fbd535SBarry Smith 474f3fbd535SBarry Smith PetscFunctionBegin; 475251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 4765b0b0462SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 477219b1045SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 478f3fbd535SBarry Smith if (iascii) { 479e3deeeafSJed Brown const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown"; 480e3deeeafSJed Brown ierr = PetscViewerASCIIPrintf(viewer," MG: type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,cyclename);CHKERRQ(ierr); 48131567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 48231567311SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr); 483f3fbd535SBarry Smith } 484218a07d4SBarry Smith if (mg->galerkin) { 485f3fbd535SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 4864f66f45eSBarry Smith } else { 4874f66f45eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 488f3fbd535SBarry Smith } 4895adeb434SBarry Smith if (mg->view){ 4905adeb434SBarry Smith ierr = (*mg->view)(pc,viewer);CHKERRQ(ierr); 4915adeb434SBarry Smith } 492f3fbd535SBarry Smith for (i=0; i<levels; i++) { 493f3fbd535SBarry Smith if (!i) { 49489cce641SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);CHKERRQ(ierr); 495f3fbd535SBarry Smith } else { 49689cce641SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 497f3fbd535SBarry Smith } 498f3fbd535SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 499f3fbd535SBarry Smith ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 500f3fbd535SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 501f3fbd535SBarry Smith if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) { 502f3fbd535SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 503f3fbd535SBarry Smith } else if (i) { 50489cce641SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 505f3fbd535SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 506f3fbd535SBarry Smith ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 507f3fbd535SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 508f3fbd535SBarry Smith } 509f3fbd535SBarry Smith } 5105b0b0462SBarry Smith } else if (isbinary) { 5115b0b0462SBarry Smith for (i=levels-1; i>=0; i--) { 5125b0b0462SBarry Smith ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 5135b0b0462SBarry Smith if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) { 5145b0b0462SBarry Smith ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 5155b0b0462SBarry Smith } 5165b0b0462SBarry Smith } 517219b1045SBarry Smith } else if (isdraw) { 518219b1045SBarry Smith PetscDraw draw; 519b4375e8dSPeter Brune PetscReal x,w,y,bottom,th; 520219b1045SBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 521219b1045SBarry Smith ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr); 5220298fd71SBarry Smith ierr = PetscDrawStringGetSize(draw,NULL,&th);CHKERRQ(ierr); 523219b1045SBarry Smith bottom = y - th; 524219b1045SBarry Smith for (i=levels-1; i>=0; i--) { 525b4375e8dSPeter Brune if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) { 526219b1045SBarry Smith ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr); 527219b1045SBarry Smith ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 528219b1045SBarry Smith ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 529b4375e8dSPeter Brune } else { 530b4375e8dSPeter Brune w = 0.5*PetscMin(1.0-x,x); 531b4375e8dSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr); 532b4375e8dSPeter Brune ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 533b4375e8dSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 534b4375e8dSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr); 535b4375e8dSPeter Brune ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 536b4375e8dSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 537b4375e8dSPeter Brune } 5380298fd71SBarry Smith ierr = PetscDrawGetBoundingBox(draw,NULL,&bottom,NULL,NULL);CHKERRQ(ierr); 5391cd381d6SBarry Smith bottom -= th; 540219b1045SBarry Smith } 5415b0b0462SBarry Smith } 542f3fbd535SBarry Smith PetscFunctionReturn(0); 543f3fbd535SBarry Smith } 544f3fbd535SBarry Smith 545af0996ceSBarry Smith #include <petsc/private/dmimpl.h> 546af0996ceSBarry Smith #include <petsc/private/kspimpl.h> 547cab2e9ccSBarry Smith 548f3fbd535SBarry Smith /* 549f3fbd535SBarry Smith Calls setup for the KSP on each level 550f3fbd535SBarry Smith */ 551f3fbd535SBarry Smith #undef __FUNCT__ 552f3fbd535SBarry Smith #define __FUNCT__ "PCSetUp_MG" 553f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc) 554f3fbd535SBarry Smith { 555f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 556f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 557f3fbd535SBarry Smith PetscErrorCode ierr; 558f3fbd535SBarry Smith PetscInt i,n = mglevels[0]->levels; 55998e04984SBarry Smith PC cpc; 5603db492dfSBarry Smith PetscBool dump = PETSC_FALSE,opsset,use_amat,missinginterpolate = PETSC_FALSE; 561f3fbd535SBarry Smith Mat dA,dB; 562f3fbd535SBarry Smith Vec tvec; 563218a07d4SBarry Smith DM *dms; 564649052a6SBarry Smith PetscViewer viewer = 0; 565f3fbd535SBarry Smith 566f3fbd535SBarry Smith PetscFunctionBegin; 56701bc414fSDmitry Karpeev /* FIX: Move this to PCSetFromOptions_MG? */ 5683aeaf226SBarry Smith if (mg->usedmfornumberoflevels) { 5693aeaf226SBarry Smith PetscInt levels; 5703aeaf226SBarry Smith ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr); 5713aeaf226SBarry Smith levels++; 5723aeaf226SBarry Smith if (levels > n) { /* the problem is now being solved on a finer grid */ 5730298fd71SBarry Smith ierr = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr); 5743aeaf226SBarry Smith n = levels; 5753aeaf226SBarry Smith ierr = PCSetFromOptions(pc);CHKERRQ(ierr); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */ 5763aeaf226SBarry Smith mglevels = mg->levels; 5773aeaf226SBarry Smith } 5783aeaf226SBarry Smith } 57998e04984SBarry Smith ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr); 5803aeaf226SBarry Smith 581f3fbd535SBarry Smith 582f3fbd535SBarry Smith /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */ 583f3fbd535SBarry Smith /* so use those from global PC */ 584f3fbd535SBarry Smith /* Is this what we always want? What if user wants to keep old one? */ 5850298fd71SBarry Smith ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset);CHKERRQ(ierr); 58698e04984SBarry Smith if (opsset) { 58798e04984SBarry Smith Mat mmat; 58823ee1639SBarry Smith ierr = KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat);CHKERRQ(ierr); 58998e04984SBarry Smith if (mmat == pc->pmat) opsset = PETSC_FALSE; 59098e04984SBarry Smith } 5914a5f07a5SMark F. Adams 5924a5f07a5SMark F. Adams if (!opsset) { 59371b23a65SMark F. Adams ierr = PCGetUseAmat(pc,&use_amat);CHKERRQ(ierr); 5944a5f07a5SMark F. Adams if(use_amat){ 595f3fbd535SBarry 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); 59623ee1639SBarry Smith ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat);CHKERRQ(ierr); 597f3fbd535SBarry Smith } 5984a5f07a5SMark F. Adams else { 5994a5f07a5SMark F. Adams ierr = PetscInfo(pc,"Using matrix (pmat) operators to define finest grid operator \n because PCMGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n");CHKERRQ(ierr); 60023ee1639SBarry Smith ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat);CHKERRQ(ierr); 6014a5f07a5SMark F. Adams } 6024a5f07a5SMark F. Adams } 603f3fbd535SBarry Smith 6049c7ffb39SBarry Smith for (i=n-1; i>0; i--) { 6059c7ffb39SBarry Smith if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) { 6069c7ffb39SBarry Smith missinginterpolate = PETSC_TRUE; 6079c7ffb39SBarry Smith continue; 6089c7ffb39SBarry Smith } 6099c7ffb39SBarry Smith } 6109c7ffb39SBarry Smith /* 6119c7ffb39SBarry Smith Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS) 6129c7ffb39SBarry Smith Skipping for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs? 6139c7ffb39SBarry Smith */ 6149c7ffb39SBarry Smith if (missinginterpolate && pc->dm && mg->galerkin != 2 && !pc->setupcalled) { 6152d2b81a6SBarry Smith /* construct the interpolation from the DMs */ 6162e499ae9SBarry Smith Mat p; 61773250ac0SBarry Smith Vec rscale; 618785e854fSJed Brown ierr = PetscMalloc1(n,&dms);CHKERRQ(ierr); 619218a07d4SBarry Smith dms[n-1] = pc->dm; 620ef1267afSMatthew G. Knepley /* Separately create them so we do not get DMKSP interference between levels */ 621ef1267afSMatthew G. Knepley for (i=n-2; i>-1; i--) {ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr);} 622218a07d4SBarry Smith for (i=n-2; i>-1; i--) { 623942e3340SBarry Smith DMKSP kdm; 6243c0c59f3SBarry Smith ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr); 6253c0c59f3SBarry Smith if (mg->galerkin) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);} 626942e3340SBarry Smith ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr); 627d1a3e677SJed Brown /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take 628d1a3e677SJed Brown * a bitwise OR of computing the matrix, RHS, and initial iterate. */ 6290298fd71SBarry Smith kdm->ops->computerhs = NULL; 6300298fd71SBarry Smith kdm->rhsctx = NULL; 63124c3aa18SBarry Smith if (!mglevels[i+1]->interpolate) { 632e727c939SJed Brown ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr); 6332d2b81a6SBarry Smith ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr); 63400ac8be1SBarry Smith if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);} 63573250ac0SBarry Smith ierr = VecDestroy(&rscale);CHKERRQ(ierr); 6366bf464f9SBarry Smith ierr = MatDestroy(&p);CHKERRQ(ierr); 637218a07d4SBarry Smith } 63824c3aa18SBarry Smith } 6392d2b81a6SBarry Smith 640ef1267afSMatthew G. Knepley for (i=n-2; i>-1; i--) {ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);} 6412d2b81a6SBarry Smith ierr = PetscFree(dms);CHKERRQ(ierr); 642ce4cda84SJed Brown } 643cab2e9ccSBarry Smith 644ce4cda84SJed Brown if (pc->dm && !pc->setupcalled) { 645ce4cda84SJed Brown /* finest smoother also gets DM but it is not active, independent of whether galerkin==2 */ 646cab2e9ccSBarry Smith ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr); 647cab2e9ccSBarry Smith ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr); 648218a07d4SBarry Smith } 649218a07d4SBarry Smith 650c91913d3SJed Brown if (mg->galerkin == 1) { 651f3fbd535SBarry Smith Mat B; 652f3fbd535SBarry Smith /* currently only handle case where mat and pmat are the same on coarser levels */ 65323ee1639SBarry Smith ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB);CHKERRQ(ierr); 654f3fbd535SBarry Smith if (!pc->setupcalled) { 655f3fbd535SBarry Smith for (i=n-2; i>-1; i--) { 656b9367914SBarry Smith if (!mglevels[i+1]->restrct && !mglevels[i+1]->interpolate) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must provide interpolation or restriction for each MG level except level 0"); 657b9367914SBarry Smith if (!mglevels[i+1]->interpolate) { 658b9367914SBarry Smith ierr = PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct);CHKERRQ(ierr); 659b9367914SBarry Smith } 660b9367914SBarry Smith if (!mglevels[i+1]->restrct) { 661b9367914SBarry Smith ierr = PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate);CHKERRQ(ierr); 662b9367914SBarry Smith } 663b9367914SBarry Smith if (mglevels[i+1]->interpolate == mglevels[i+1]->restrct) { 664f3fbd535SBarry Smith ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 6657d537d92SJed Brown } else { 6667d537d92SJed Brown ierr = MatMatMatMult(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 6677d537d92SJed Brown } 66823ee1639SBarry Smith ierr = KSPSetOperators(mglevels[i]->smoothd,B,B);CHKERRQ(ierr); 669f3fbd535SBarry Smith if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);} 670f3fbd535SBarry Smith dB = B; 671f3fbd535SBarry Smith } 672cd9507b2SBarry Smith if (n > 1) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);} 673f3fbd535SBarry Smith } else { 674f3fbd535SBarry Smith for (i=n-2; i>-1; i--) { 675b9367914SBarry Smith if (!mglevels[i+1]->restrct && !mglevels[i+1]->interpolate) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must provide interpolation or restriction for each MG level except level 0"); 676b9367914SBarry Smith if (!mglevels[i+1]->interpolate) { 677b9367914SBarry Smith ierr = PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct);CHKERRQ(ierr); 678b9367914SBarry Smith } 679b9367914SBarry Smith if (!mglevels[i+1]->restrct) { 680b9367914SBarry Smith ierr = PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate);CHKERRQ(ierr); 681b9367914SBarry Smith } 68223ee1639SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&B);CHKERRQ(ierr); 683b9367914SBarry Smith if (mglevels[i+1]->interpolate == mglevels[i+1]->restrct) { 684f3fbd535SBarry Smith ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 6857d537d92SJed Brown } else { 6867d537d92SJed Brown ierr = MatMatMatMult(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 6877d537d92SJed Brown } 68823ee1639SBarry Smith ierr = KSPSetOperators(mglevels[i]->smoothd,B,B);CHKERRQ(ierr); 689f3fbd535SBarry Smith dB = B; 690f3fbd535SBarry Smith } 691f3fbd535SBarry Smith } 692ce4cda84SJed Brown } else if (!mg->galerkin && pc->dm && pc->dm->x) { 693cab2e9ccSBarry Smith /* need to restrict Jacobian location to coarser meshes for evaluation */ 694cab2e9ccSBarry Smith for (i=n-2; i>-1; i--) { 695c88c5224SJed Brown Mat R; 696c88c5224SJed Brown Vec rscale; 697cab2e9ccSBarry Smith if (!mglevels[i]->smoothd->dm->x) { 698cab2e9ccSBarry Smith Vec *vecs; 6992a7a6963SBarry Smith ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vecs,0,NULL);CHKERRQ(ierr); 7002fa5cd67SKarl Rupp 701cab2e9ccSBarry Smith mglevels[i]->smoothd->dm->x = vecs[0]; 7022fa5cd67SKarl Rupp 703cab2e9ccSBarry Smith ierr = PetscFree(vecs);CHKERRQ(ierr); 704cab2e9ccSBarry Smith } 705c88c5224SJed Brown ierr = PCMGGetRestriction(pc,i+1,&R);CHKERRQ(ierr); 706c88c5224SJed Brown ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr); 707c88c5224SJed Brown ierr = MatRestrict(R,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);CHKERRQ(ierr); 708c88c5224SJed Brown ierr = VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,rscale);CHKERRQ(ierr); 709cab2e9ccSBarry Smith } 710f3fbd535SBarry Smith } 711ccceb50cSJed Brown if (!mg->galerkin && pc->dm) { 712caa4e7f2SJed Brown for (i=n-2; i>=0; i--) { 713ccceb50cSJed Brown DM dmfine,dmcoarse; 714caa4e7f2SJed Brown Mat Restrict,Inject; 715caa4e7f2SJed Brown Vec rscale; 716ccceb50cSJed Brown ierr = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr); 717ccceb50cSJed Brown ierr = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr); 718caa4e7f2SJed Brown ierr = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr); 719caa4e7f2SJed Brown ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr); 7200298fd71SBarry Smith Inject = NULL; /* Callback should create it if it needs Injection */ 721caa4e7f2SJed Brown ierr = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr); 722caa4e7f2SJed Brown } 723caa4e7f2SJed Brown } 724f3fbd535SBarry Smith 725f3fbd535SBarry Smith if (!pc->setupcalled) { 726f3fbd535SBarry Smith for (i=0; i<n; i++) { 727f3fbd535SBarry Smith ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr); 728f3fbd535SBarry Smith } 729f3fbd535SBarry Smith for (i=1; i<n; i++) { 730f3fbd535SBarry Smith if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) { 731f3fbd535SBarry Smith ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr); 732f3fbd535SBarry Smith } 733f3fbd535SBarry Smith } 734f3fbd535SBarry Smith for (i=1; i<n; i++) { 735c88c5224SJed Brown ierr = PCMGGetInterpolation(pc,i,&mglevels[i]->interpolate);CHKERRQ(ierr); 736c88c5224SJed Brown ierr = PCMGGetRestriction(pc,i,&mglevels[i]->restrct);CHKERRQ(ierr); 737f3fbd535SBarry Smith } 738f3fbd535SBarry Smith for (i=0; i<n-1; i++) { 739f3fbd535SBarry Smith if (!mglevels[i]->b) { 740f3fbd535SBarry Smith Vec *vec; 7412a7a6963SBarry Smith ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr); 742f3fbd535SBarry Smith ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr); 7436bf464f9SBarry Smith ierr = VecDestroy(vec);CHKERRQ(ierr); 744f3fbd535SBarry Smith ierr = PetscFree(vec);CHKERRQ(ierr); 745f3fbd535SBarry Smith } 746f3fbd535SBarry Smith if (!mglevels[i]->r && i) { 747f3fbd535SBarry Smith ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 748f3fbd535SBarry Smith ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr); 7496bf464f9SBarry Smith ierr = VecDestroy(&tvec);CHKERRQ(ierr); 750f3fbd535SBarry Smith } 751f3fbd535SBarry Smith if (!mglevels[i]->x) { 752f3fbd535SBarry Smith ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 753f3fbd535SBarry Smith ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr); 7546bf464f9SBarry Smith ierr = VecDestroy(&tvec);CHKERRQ(ierr); 755f3fbd535SBarry Smith } 756f3fbd535SBarry Smith } 757f3fbd535SBarry Smith if (n != 1 && !mglevels[n-1]->r) { 758f3fbd535SBarry Smith /* PCMGSetR() on the finest level if user did not supply it */ 759f3fbd535SBarry Smith Vec *vec; 7602a7a6963SBarry Smith ierr = KSPCreateVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr); 761f3fbd535SBarry Smith ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr); 7626bf464f9SBarry Smith ierr = VecDestroy(vec);CHKERRQ(ierr); 763f3fbd535SBarry Smith ierr = PetscFree(vec);CHKERRQ(ierr); 764f3fbd535SBarry Smith } 765f3fbd535SBarry Smith } 766f3fbd535SBarry Smith 76798e04984SBarry Smith if (pc->dm) { 76898e04984SBarry Smith /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */ 76998e04984SBarry Smith for (i=0; i<n-1; i++) { 77098e04984SBarry Smith if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX; 77198e04984SBarry Smith } 77298e04984SBarry Smith } 773f3fbd535SBarry Smith 774f3fbd535SBarry Smith for (i=1; i<n; i++) { 7752a21e185SBarry Smith if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1){ 776f3fbd535SBarry Smith /* if doing only down then initial guess is zero */ 777f3fbd535SBarry Smith ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 778f3fbd535SBarry Smith } 77963e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 780f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr); 781899639b0SHong Zhang if (mglevels[i]->smoothd->reason == KSP_DIVERGED_PCSETUP_FAILED) { 782899639b0SHong Zhang pc->failedreason = PC_SUBPC_ERROR; 783899639b0SHong Zhang } 78463e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 785d42688cbSBarry Smith if (!mglevels[i]->residual) { 786d42688cbSBarry Smith Mat mat; 78723ee1639SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&mat);CHKERRQ(ierr); 78854b2cd4bSJed Brown ierr = PCMGSetResidual(pc,i,PCMGResidualDefault,mat);CHKERRQ(ierr); 789d42688cbSBarry Smith } 790f3fbd535SBarry Smith } 791f3fbd535SBarry Smith for (i=1; i<n; i++) { 792f3fbd535SBarry Smith if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) { 793f3fbd535SBarry Smith Mat downmat,downpmat; 794f3fbd535SBarry Smith 795f3fbd535SBarry Smith /* check if operators have been set for up, if not use down operators to set them */ 7960298fd71SBarry Smith ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);CHKERRQ(ierr); 797f3fbd535SBarry Smith if (!opsset) { 79823ee1639SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);CHKERRQ(ierr); 79923ee1639SBarry Smith ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat);CHKERRQ(ierr); 800f3fbd535SBarry Smith } 801f3fbd535SBarry Smith 802f3fbd535SBarry Smith ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr); 80363e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 804f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr); 805899639b0SHong Zhang if (mglevels[i]->smoothu->reason == KSP_DIVERGED_PCSETUP_FAILED) { 806899639b0SHong Zhang pc->failedreason = PC_SUBPC_ERROR; 807899639b0SHong Zhang } 80863e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 809f3fbd535SBarry Smith } 810f3fbd535SBarry Smith } 811f3fbd535SBarry Smith 81263e6d426SJed Brown if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 813f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr); 814899639b0SHong Zhang if (mglevels[0]->smoothd->reason == KSP_DIVERGED_PCSETUP_FAILED) { 815899639b0SHong Zhang pc->failedreason = PC_SUBPC_ERROR; 816899639b0SHong Zhang } 81763e6d426SJed Brown if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 818f3fbd535SBarry Smith 819f3fbd535SBarry Smith /* 820f3fbd535SBarry Smith Dump the interpolation/restriction matrices plus the 821e3c5b3baSBarry Smith Jacobian/stiffness on each level. This allows MATLAB users to 822f3fbd535SBarry Smith easily check if the Galerkin condition A_c = R A_f R^T is satisfied. 823f3fbd535SBarry Smith 824f3fbd535SBarry Smith Only support one or the other at the same time. 825f3fbd535SBarry Smith */ 826f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER) 827c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);CHKERRQ(ierr); 828ce94432eSBarry Smith if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc)); 829f3fbd535SBarry Smith dump = PETSC_FALSE; 830f3fbd535SBarry Smith #endif 831c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);CHKERRQ(ierr); 832ce94432eSBarry Smith if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc)); 833f3fbd535SBarry Smith 834f3fbd535SBarry Smith if (viewer) { 835f3fbd535SBarry Smith for (i=1; i<n; i++) { 836f3fbd535SBarry Smith ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr); 837f3fbd535SBarry Smith } 838f3fbd535SBarry Smith for (i=0; i<n; i++) { 839f3fbd535SBarry Smith ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr); 840f3fbd535SBarry Smith ierr = MatView(pc->mat,viewer);CHKERRQ(ierr); 841f3fbd535SBarry Smith } 842f3fbd535SBarry Smith } 843f3fbd535SBarry Smith PetscFunctionReturn(0); 844f3fbd535SBarry Smith } 845f3fbd535SBarry Smith 846f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/ 847f3fbd535SBarry Smith 848f3fbd535SBarry Smith #undef __FUNCT__ 8499dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetLevels" 8504b9ad928SBarry Smith /*@ 85197177400SBarry Smith PCMGGetLevels - Gets the number of levels to use with MG. 8524b9ad928SBarry Smith 8534b9ad928SBarry Smith Not Collective 8544b9ad928SBarry Smith 8554b9ad928SBarry Smith Input Parameter: 8564b9ad928SBarry Smith . pc - the preconditioner context 8574b9ad928SBarry Smith 8584b9ad928SBarry Smith Output parameter: 8594b9ad928SBarry Smith . levels - the number of levels 8604b9ad928SBarry Smith 8614b9ad928SBarry Smith Level: advanced 8624b9ad928SBarry Smith 8634b9ad928SBarry Smith .keywords: MG, get, levels, multigrid 8644b9ad928SBarry Smith 86597177400SBarry Smith .seealso: PCMGSetLevels() 8664b9ad928SBarry Smith @*/ 8677087cfbeSBarry Smith PetscErrorCode PCMGGetLevels(PC pc,PetscInt *levels) 8684b9ad928SBarry Smith { 869f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 8704b9ad928SBarry Smith 8714b9ad928SBarry Smith PetscFunctionBegin; 8720700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 8734482741eSBarry Smith PetscValidIntPointer(levels,2); 874f3fbd535SBarry Smith *levels = mg->nlevels; 8754b9ad928SBarry Smith PetscFunctionReturn(0); 8764b9ad928SBarry Smith } 8774b9ad928SBarry Smith 8784b9ad928SBarry Smith #undef __FUNCT__ 8799dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetType" 8804b9ad928SBarry Smith /*@ 88197177400SBarry Smith PCMGSetType - Determines the form of multigrid to use: 8824b9ad928SBarry Smith multiplicative, additive, full, or the Kaskade algorithm. 8834b9ad928SBarry Smith 884ad4df100SBarry Smith Logically Collective on PC 8854b9ad928SBarry Smith 8864b9ad928SBarry Smith Input Parameters: 8874b9ad928SBarry Smith + pc - the preconditioner context 8889dcbbd2bSBarry Smith - form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE, 8899dcbbd2bSBarry Smith PC_MG_FULL, PC_MG_KASKADE 8904b9ad928SBarry Smith 8914b9ad928SBarry Smith Options Database Key: 8924b9ad928SBarry Smith . -pc_mg_type <form> - Sets <form>, one of multiplicative, 8934b9ad928SBarry Smith additive, full, kaskade 8944b9ad928SBarry Smith 8954b9ad928SBarry Smith Level: advanced 8964b9ad928SBarry Smith 8974b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid 8984b9ad928SBarry Smith 89997177400SBarry Smith .seealso: PCMGSetLevels() 9004b9ad928SBarry Smith @*/ 9017087cfbeSBarry Smith PetscErrorCode PCMGSetType(PC pc,PCMGType form) 9024b9ad928SBarry Smith { 903f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 9044b9ad928SBarry Smith 9054b9ad928SBarry Smith PetscFunctionBegin; 9060700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 907c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc,form,2); 90831567311SBarry Smith mg->am = form; 9099dcbbd2bSBarry Smith if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG; 910c60c7ad4SBarry Smith else pc->ops->applyrichardson = NULL; 911c60c7ad4SBarry Smith PetscFunctionReturn(0); 912c60c7ad4SBarry Smith } 913c60c7ad4SBarry Smith 914f0af28e8SLisandro Dalcin #undef __FUNCT__ 915f0af28e8SLisandro Dalcin #define __FUNCT__ "PCMGGetType" 916c60c7ad4SBarry Smith /*@ 917c60c7ad4SBarry Smith PCMGGetType - Determines the form of multigrid to use: 918c60c7ad4SBarry Smith multiplicative, additive, full, or the Kaskade algorithm. 919c60c7ad4SBarry Smith 920c60c7ad4SBarry Smith Logically Collective on PC 921c60c7ad4SBarry Smith 922c60c7ad4SBarry Smith Input Parameter: 923c60c7ad4SBarry Smith . pc - the preconditioner context 924c60c7ad4SBarry Smith 925c60c7ad4SBarry Smith Output Parameter: 926c60c7ad4SBarry Smith . type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE 927c60c7ad4SBarry Smith 928c60c7ad4SBarry Smith 929c60c7ad4SBarry Smith Level: advanced 930c60c7ad4SBarry Smith 931c60c7ad4SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid 932c60c7ad4SBarry Smith 933c60c7ad4SBarry Smith .seealso: PCMGSetLevels() 934c60c7ad4SBarry Smith @*/ 935c60c7ad4SBarry Smith PetscErrorCode PCMGGetType(PC pc,PCMGType *type) 936c60c7ad4SBarry Smith { 937c60c7ad4SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 938c60c7ad4SBarry Smith 939c60c7ad4SBarry Smith PetscFunctionBegin; 940c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 941c60c7ad4SBarry Smith *type = mg->am; 9424b9ad928SBarry Smith PetscFunctionReturn(0); 9434b9ad928SBarry Smith } 9444b9ad928SBarry Smith 9454b9ad928SBarry Smith #undef __FUNCT__ 9460d353602SBarry Smith #define __FUNCT__ "PCMGSetCycleType" 9474b9ad928SBarry Smith /*@ 9480d353602SBarry Smith PCMGSetCycleType - Sets the type cycles to use. Use PCMGSetCycleTypeOnLevel() for more 9494b9ad928SBarry Smith complicated cycling. 9504b9ad928SBarry Smith 951ad4df100SBarry Smith Logically Collective on PC 9524b9ad928SBarry Smith 9534b9ad928SBarry Smith Input Parameters: 954c2be2410SBarry Smith + pc - the multigrid context 9550d353602SBarry Smith - PC_MG_CYCLE_V or PC_MG_CYCLE_W 9564b9ad928SBarry Smith 9574b9ad928SBarry Smith Options Database Key: 9584446f3b4SBarry Smith . -pc_mg_cycle_type <v,w> 9594b9ad928SBarry Smith 9604b9ad928SBarry Smith Level: advanced 9614b9ad928SBarry Smith 9624b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 9634b9ad928SBarry Smith 9640d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel() 9654b9ad928SBarry Smith @*/ 9667087cfbeSBarry Smith PetscErrorCode PCMGSetCycleType(PC pc,PCMGCycleType n) 9674b9ad928SBarry Smith { 968f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 969f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 97079416396SBarry Smith PetscInt i,levels; 9714b9ad928SBarry Smith 9724b9ad928SBarry Smith PetscFunctionBegin; 9730700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 974ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 97569fbec6eSBarry Smith PetscValidLogicalCollectiveEnum(pc,n,2); 976f3fbd535SBarry Smith levels = mglevels[0]->levels; 9774b9ad928SBarry Smith 9782fa5cd67SKarl Rupp for (i=0; i<levels; i++) mglevels[i]->cycles = n; 9794b9ad928SBarry Smith PetscFunctionReturn(0); 9804b9ad928SBarry Smith } 9814b9ad928SBarry Smith 9824b9ad928SBarry Smith #undef __FUNCT__ 9838cc2d5dfSBarry Smith #define __FUNCT__ "PCMGMultiplicativeSetCycles" 9848cc2d5dfSBarry Smith /*@ 9858cc2d5dfSBarry Smith PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step 9868cc2d5dfSBarry Smith of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used 9878cc2d5dfSBarry Smith 988ad4df100SBarry Smith Logically Collective on PC 9898cc2d5dfSBarry Smith 9908cc2d5dfSBarry Smith Input Parameters: 9918cc2d5dfSBarry Smith + pc - the multigrid context 9928cc2d5dfSBarry Smith - n - number of cycles (default is 1) 9938cc2d5dfSBarry Smith 9948cc2d5dfSBarry Smith Options Database Key: 995e1bc860dSBarry Smith . -pc_mg_multiplicative_cycles n 9968cc2d5dfSBarry Smith 9978cc2d5dfSBarry Smith Level: advanced 9988cc2d5dfSBarry Smith 9998cc2d5dfSBarry Smith Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType() 10008cc2d5dfSBarry Smith 10018cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 10028cc2d5dfSBarry Smith 10038cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType() 10048cc2d5dfSBarry Smith @*/ 10057087cfbeSBarry Smith PetscErrorCode PCMGMultiplicativeSetCycles(PC pc,PetscInt n) 10068cc2d5dfSBarry Smith { 1007f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 10088cc2d5dfSBarry Smith 10098cc2d5dfSBarry Smith PetscFunctionBegin; 10100700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1011c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 10122a21e185SBarry Smith mg->cyclesperpcapply = n; 10138cc2d5dfSBarry Smith PetscFunctionReturn(0); 10148cc2d5dfSBarry Smith } 10158cc2d5dfSBarry Smith 10168cc2d5dfSBarry Smith #undef __FUNCT__ 1017fb15c04eSBarry Smith #define __FUNCT__ "PCMGSetGalerkin_MG" 1018fb15c04eSBarry Smith PetscErrorCode PCMGSetGalerkin_MG(PC pc,PetscBool use) 1019fb15c04eSBarry Smith { 1020fb15c04eSBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1021fb15c04eSBarry Smith 1022fb15c04eSBarry Smith PetscFunctionBegin; 1023fb15c04eSBarry Smith mg->galerkin = use ? 1 : 0; 1024fb15c04eSBarry Smith PetscFunctionReturn(0); 1025fb15c04eSBarry Smith } 1026fb15c04eSBarry Smith 1027fb15c04eSBarry Smith #undef __FUNCT__ 10289dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetGalerkin" 1029c2be2410SBarry Smith /*@ 103097177400SBarry Smith PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the 103182b87d87SMatthew G. Knepley finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i 1032c2be2410SBarry Smith 1033ad4df100SBarry Smith Logically Collective on PC 1034c2be2410SBarry Smith 1035c2be2410SBarry Smith Input Parameters: 1036c91913d3SJed Brown + pc - the multigrid context 1037c91913d3SJed Brown - use - PETSC_TRUE to use the Galerkin process to compute coarse-level operators 1038c2be2410SBarry Smith 1039c2be2410SBarry Smith Options Database Key: 10401c1aac46SBarry Smith . -pc_mg_galerkin <true,false> 1041c2be2410SBarry Smith 1042c2be2410SBarry Smith Level: intermediate 1043c2be2410SBarry Smith 10441c1aac46SBarry Smith Notes: Some codes that use PCMG such as PCGAMG use Galerkin internally while constructing the hierarchy and thus do not 10451c1aac46SBarry Smith use the PCMG construction of the coarser grids. 10461c1aac46SBarry Smith 1047c2be2410SBarry Smith .keywords: MG, set, Galerkin 1048c2be2410SBarry Smith 10493fc8bf9cSBarry Smith .seealso: PCMGGetGalerkin() 10503fc8bf9cSBarry Smith 1051c2be2410SBarry Smith @*/ 1052c91913d3SJed Brown PetscErrorCode PCMGSetGalerkin(PC pc,PetscBool use) 1053c2be2410SBarry Smith { 1054fb15c04eSBarry Smith PetscErrorCode ierr; 1055c2be2410SBarry Smith 1056c2be2410SBarry Smith PetscFunctionBegin; 10570700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 10587a4c0024SBarry Smith ierr = PetscTryMethod(pc,"PCMGSetGalerkin_C",(PC,PetscBool),(pc,use));CHKERRQ(ierr); 1059c2be2410SBarry Smith PetscFunctionReturn(0); 1060c2be2410SBarry Smith } 1061c2be2410SBarry Smith 1062c2be2410SBarry Smith #undef __FUNCT__ 10633fc8bf9cSBarry Smith #define __FUNCT__ "PCMGGetGalerkin" 10643fc8bf9cSBarry Smith /*@ 10653fc8bf9cSBarry Smith PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. 106682b87d87SMatthew G. Knepley A_i-1 = r_i * A_i * p_i 10673fc8bf9cSBarry Smith 10683fc8bf9cSBarry Smith Not Collective 10693fc8bf9cSBarry Smith 10703fc8bf9cSBarry Smith Input Parameter: 10713fc8bf9cSBarry Smith . pc - the multigrid context 10723fc8bf9cSBarry Smith 10733fc8bf9cSBarry Smith Output Parameter: 1074e1bc860dSBarry Smith . galerkin - PETSC_TRUE or PETSC_FALSE 10753fc8bf9cSBarry Smith 10763fc8bf9cSBarry Smith Options Database Key: 1077e1bc860dSBarry Smith . -pc_mg_galerkin 10783fc8bf9cSBarry Smith 10793fc8bf9cSBarry Smith Level: intermediate 10803fc8bf9cSBarry Smith 10813fc8bf9cSBarry Smith .keywords: MG, set, Galerkin 10823fc8bf9cSBarry Smith 10833fc8bf9cSBarry Smith .seealso: PCMGSetGalerkin() 10843fc8bf9cSBarry Smith 10853fc8bf9cSBarry Smith @*/ 10867087cfbeSBarry Smith PetscErrorCode PCMGGetGalerkin(PC pc,PetscBool *galerkin) 10873fc8bf9cSBarry Smith { 1088f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 10893fc8bf9cSBarry Smith 10903fc8bf9cSBarry Smith PetscFunctionBegin; 10910700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 109213fdb3acSJose Roman *galerkin = (PetscBool)mg->galerkin; 10933fc8bf9cSBarry Smith PetscFunctionReturn(0); 10943fc8bf9cSBarry Smith } 10953fc8bf9cSBarry Smith 10963fc8bf9cSBarry Smith #undef __FUNCT__ 10979dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothDown" 10984b9ad928SBarry Smith /*@ 109997177400SBarry Smith PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to 110097177400SBarry Smith use on all levels. Use PCMGGetSmootherDown() to set different 11014b9ad928SBarry Smith pre-smoothing steps on different levels. 11024b9ad928SBarry Smith 1103ad4df100SBarry Smith Logically Collective on PC 11044b9ad928SBarry Smith 11054b9ad928SBarry Smith Input Parameters: 11064b9ad928SBarry Smith + mg - the multigrid context 11074b9ad928SBarry Smith - n - the number of smoothing steps 11084b9ad928SBarry Smith 11094b9ad928SBarry Smith Options Database Key: 11104b9ad928SBarry Smith . -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps 11114b9ad928SBarry Smith 11124b9ad928SBarry Smith Level: advanced 11134b9ad928SBarry Smith 11144b9ad928SBarry Smith .keywords: MG, smooth, down, pre-smoothing, steps, multigrid 11154b9ad928SBarry Smith 111697177400SBarry Smith .seealso: PCMGSetNumberSmoothUp() 11174b9ad928SBarry Smith @*/ 11187087cfbeSBarry Smith PetscErrorCode PCMGSetNumberSmoothDown(PC pc,PetscInt n) 11194b9ad928SBarry Smith { 1120f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1121f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 11226849ba73SBarry Smith PetscErrorCode ierr; 112379416396SBarry Smith PetscInt i,levels; 11244b9ad928SBarry Smith 11254b9ad928SBarry Smith PetscFunctionBegin; 11260700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1127ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 1128c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 1129f3fbd535SBarry Smith levels = mglevels[0]->levels; 11304b9ad928SBarry Smith 1131b05257ddSBarry Smith for (i=1; i<levels; i++) { 11324b9ad928SBarry Smith /* make sure smoother up and down are different */ 11330298fd71SBarry Smith ierr = PCMGGetSmootherUp(pc,i,NULL);CHKERRQ(ierr); 1134f3fbd535SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 11352fa5cd67SKarl Rupp 113631567311SBarry Smith mg->default_smoothd = n; 11374b9ad928SBarry Smith } 11384b9ad928SBarry Smith PetscFunctionReturn(0); 11394b9ad928SBarry Smith } 11404b9ad928SBarry Smith 11414b9ad928SBarry Smith #undef __FUNCT__ 11429dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothUp" 11434b9ad928SBarry Smith /*@ 114497177400SBarry Smith PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use 114597177400SBarry Smith on all levels. Use PCMGGetSmootherUp() to set different numbers of 11464b9ad928SBarry Smith post-smoothing steps on different levels. 11474b9ad928SBarry Smith 1148ad4df100SBarry Smith Logically Collective on PC 11494b9ad928SBarry Smith 11504b9ad928SBarry Smith Input Parameters: 11514b9ad928SBarry Smith + mg - the multigrid context 11524b9ad928SBarry Smith - n - the number of smoothing steps 11534b9ad928SBarry Smith 11544b9ad928SBarry Smith Options Database Key: 11554b9ad928SBarry Smith . -pc_mg_smoothup <n> - Sets number of post-smoothing steps 11564b9ad928SBarry Smith 11574b9ad928SBarry Smith Level: advanced 11584b9ad928SBarry Smith 11594b9ad928SBarry Smith Note: this does not set a value on the coarsest grid, since we assume that 1160a8c7a070SBarry Smith there is no separate smooth up on the coarsest grid. 11614b9ad928SBarry Smith 11624b9ad928SBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid 11634b9ad928SBarry Smith 116497177400SBarry Smith .seealso: PCMGSetNumberSmoothDown() 11654b9ad928SBarry Smith @*/ 11667087cfbeSBarry Smith PetscErrorCode PCMGSetNumberSmoothUp(PC pc,PetscInt n) 11674b9ad928SBarry Smith { 1168f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1169f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 11706849ba73SBarry Smith PetscErrorCode ierr; 117179416396SBarry Smith PetscInt i,levels; 11724b9ad928SBarry Smith 11734b9ad928SBarry Smith PetscFunctionBegin; 11740700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1175ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 1176c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 1177f3fbd535SBarry Smith levels = mglevels[0]->levels; 11784b9ad928SBarry Smith 11794b9ad928SBarry Smith for (i=1; i<levels; i++) { 11804b9ad928SBarry Smith /* make sure smoother up and down are different */ 11810298fd71SBarry Smith ierr = PCMGGetSmootherUp(pc,i,NULL);CHKERRQ(ierr); 1182f3fbd535SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 11832fa5cd67SKarl Rupp 118431567311SBarry Smith mg->default_smoothu = n; 11854b9ad928SBarry Smith } 11864b9ad928SBarry Smith PetscFunctionReturn(0); 11874b9ad928SBarry Smith } 11884b9ad928SBarry Smith 11894b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/ 11904b9ad928SBarry Smith 11913b09bd56SBarry Smith /*MC 1192ccb205f8SBarry Smith PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional 11933b09bd56SBarry Smith information about the coarser grid matrices and restriction/interpolation operators. 11943b09bd56SBarry Smith 11953b09bd56SBarry Smith Options Database Keys: 11963b09bd56SBarry Smith + -pc_mg_levels <nlevels> - number of levels including finest 1197b4519db0SPatrick Sanan . -pc_mg_cycles <v,w> - 119879416396SBarry Smith . -pc_mg_smoothup <n> - number of smoothing steps after interpolation 11993b09bd56SBarry Smith . -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator 12008c1c2452SJed Brown . -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default 12013b09bd56SBarry Smith . -pc_mg_log - log information about time spent on each level of the solver 12028cf6037eSBarry Smith . -pc_mg_galerkin - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R' 12038cf6037eSBarry Smith . -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1) 12048cf6037eSBarry Smith . -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 1205e3c5b3baSBarry Smith to the Socket viewer for reading from MATLAB. 12068cf6037eSBarry Smith - -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices 12078cf6037eSBarry Smith to the binary output file called binaryoutput 12083b09bd56SBarry Smith 120924c3aa18SBarry 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 12103b09bd56SBarry Smith 12118cf6037eSBarry Smith When run with a single level the smoother options are used on that level NOT the coarse grid solver options 12128cf6037eSBarry Smith 12133b09bd56SBarry Smith Level: intermediate 12143b09bd56SBarry Smith 12158f87f92bSBarry Smith Concepts: multigrid/multilevel 12163b09bd56SBarry Smith 12178cf6037eSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE 12180d353602SBarry Smith PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(), 121997177400SBarry Smith PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 122097177400SBarry Smith PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 12210d353602SBarry Smith PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 12223b09bd56SBarry Smith M*/ 12233b09bd56SBarry Smith 12244b9ad928SBarry Smith #undef __FUNCT__ 12254b9ad928SBarry Smith #define __FUNCT__ "PCCreate_MG" 12268cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc) 12274b9ad928SBarry Smith { 1228f3fbd535SBarry Smith PC_MG *mg; 1229f3fbd535SBarry Smith PetscErrorCode ierr; 1230f3fbd535SBarry Smith 12314b9ad928SBarry Smith PetscFunctionBegin; 1232b00a9115SJed Brown ierr = PetscNewLog(pc,&mg);CHKERRQ(ierr); 1233f3fbd535SBarry Smith pc->data = (void*)mg; 1234f3fbd535SBarry Smith mg->nlevels = -1; 123510eca3edSLisandro Dalcin mg->am = PC_MG_MULTIPLICATIVE; 1236f3fbd535SBarry Smith 123737a44384SMark Adams pc->useAmat = PETSC_TRUE; 123837a44384SMark Adams 12394b9ad928SBarry Smith pc->ops->apply = PCApply_MG; 12404b9ad928SBarry Smith pc->ops->setup = PCSetUp_MG; 1241a06653b4SBarry Smith pc->ops->reset = PCReset_MG; 12424b9ad928SBarry Smith pc->ops->destroy = PCDestroy_MG; 12434b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_MG; 12444b9ad928SBarry Smith pc->ops->view = PCView_MG; 1245fb15c04eSBarry Smith 1246fb15c04eSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",PCMGSetGalerkin_MG);CHKERRQ(ierr); 12474b9ad928SBarry Smith PetscFunctionReturn(0); 12484b9ad928SBarry Smith } 1249