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; 164b9ad928SBarry Smith 174b9ad928SBarry Smith PetscFunctionBegin; 1863e6d426SJed Brown if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 1931567311SBarry Smith ierr = KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x);CHKERRQ(ierr); /* pre-smooth */ 2063e6d426SJed Brown if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 2131567311SBarry Smith if (mglevels->level) { /* not the coarsest grid */ 2263e6d426SJed Brown if (mglevels->eventresidual) {ierr = PetscLogEventBegin(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);} 2331567311SBarry Smith ierr = (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);CHKERRQ(ierr); 2463e6d426SJed Brown if (mglevels->eventresidual) {ierr = PetscLogEventEnd(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);} 254b9ad928SBarry Smith 264b9ad928SBarry Smith /* if on finest level and have convergence criteria set */ 2731567311SBarry Smith if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) { 284b9ad928SBarry Smith PetscReal rnorm; 2931567311SBarry Smith ierr = VecNorm(mglevels->r,NORM_2,&rnorm);CHKERRQ(ierr); 304b9ad928SBarry Smith if (rnorm <= mg->ttol) { 3170441072SBarry Smith if (rnorm < mg->abstol) { 324d0a8057SBarry Smith *reason = PCRICHARDSON_CONVERGED_ATOL; 3357622a8eSBarry Smith ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %g is less than absolute tolerance %g\n",(double)rnorm,(double)mg->abstol);CHKERRQ(ierr); 344b9ad928SBarry Smith } else { 354d0a8057SBarry Smith *reason = PCRICHARDSON_CONVERGED_RTOL; 3657622a8eSBarry Smith ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %g is less than relative tolerance times initial residual norm %g\n",(double)rnorm,(double)mg->ttol);CHKERRQ(ierr); 374b9ad928SBarry Smith } 384b9ad928SBarry Smith PetscFunctionReturn(0); 394b9ad928SBarry Smith } 404b9ad928SBarry Smith } 414b9ad928SBarry Smith 4231567311SBarry Smith mgc = *(mglevelsin - 1); 4363e6d426SJed Brown if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 4431567311SBarry Smith ierr = MatRestrict(mglevels->restrct,mglevels->r,mgc->b);CHKERRQ(ierr); 4563e6d426SJed Brown if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 46efb30889SBarry Smith ierr = VecSet(mgc->x,0.0);CHKERRQ(ierr); 474b9ad928SBarry Smith while (cycles--) { 4831567311SBarry Smith ierr = PCMGMCycle_Private(pc,mglevelsin-1,reason);CHKERRQ(ierr); 494b9ad928SBarry Smith } 5063e6d426SJed Brown if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 5131567311SBarry Smith ierr = MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x);CHKERRQ(ierr); 5263e6d426SJed Brown if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 5363e6d426SJed Brown if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 5431567311SBarry Smith ierr = KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x);CHKERRQ(ierr); /* post smooth */ 5563e6d426SJed Brown if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 564b9ad928SBarry Smith } 574b9ad928SBarry Smith PetscFunctionReturn(0); 584b9ad928SBarry Smith } 594b9ad928SBarry Smith 604b9ad928SBarry Smith #undef __FUNCT__ 614b9ad928SBarry Smith #define __FUNCT__ "PCApplyRichardson_MG" 62ace3abfcSBarry Smith static PetscErrorCode PCApplyRichardson_MG(PC pc,Vec b,Vec x,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscBool zeroguess,PetscInt *outits,PCRichardsonConvergedReason *reason) 634b9ad928SBarry Smith { 64f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 65f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 66dfbe8321SBarry Smith PetscErrorCode ierr; 67f3fbd535SBarry Smith PetscInt levels = mglevels[0]->levels,i; 684b9ad928SBarry Smith 694b9ad928SBarry Smith PetscFunctionBegin; 70a762d673SBarry Smith /* When the DM is supplying the matrix then it will not exist until here */ 71a762d673SBarry Smith for (i=0; i<levels; i++) { 72a762d673SBarry Smith if (!mglevels[i]->A) { 73a762d673SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr); 74a762d673SBarry Smith ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr); 75a762d673SBarry Smith } 76a762d673SBarry Smith } 77f3fbd535SBarry Smith mglevels[levels-1]->b = b; 78f3fbd535SBarry Smith mglevels[levels-1]->x = x; 794b9ad928SBarry Smith 8031567311SBarry Smith mg->rtol = rtol; 8131567311SBarry Smith mg->abstol = abstol; 8231567311SBarry Smith mg->dtol = dtol; 834b9ad928SBarry Smith if (rtol) { 844b9ad928SBarry Smith /* compute initial residual norm for relative convergence test */ 854b9ad928SBarry Smith PetscReal rnorm; 867319c654SBarry Smith if (zeroguess) { 877319c654SBarry Smith ierr = VecNorm(b,NORM_2,&rnorm);CHKERRQ(ierr); 887319c654SBarry Smith } else { 89f3fbd535SBarry Smith ierr = (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);CHKERRQ(ierr); 904b9ad928SBarry Smith ierr = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr); 917319c654SBarry Smith } 9231567311SBarry Smith mg->ttol = PetscMax(rtol*rnorm,abstol); 932fa5cd67SKarl Rupp } else if (abstol) mg->ttol = abstol; 942fa5cd67SKarl Rupp else mg->ttol = 0.0; 954b9ad928SBarry Smith 964d0a8057SBarry Smith /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't 97336babb1SJed Brown stop prematurely due to small residual */ 984d0a8057SBarry Smith for (i=1; i<levels; i++) { 99f3fbd535SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 100f3fbd535SBarry Smith if (mglevels[i]->smoothu != mglevels[i]->smoothd) { 101f3fbd535SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 1024b9ad928SBarry Smith } 1034d0a8057SBarry Smith } 1044d0a8057SBarry Smith 1054d0a8057SBarry Smith *reason = (PCRichardsonConvergedReason)0; 1064d0a8057SBarry Smith for (i=0; i<its; i++) { 107f3fbd535SBarry Smith ierr = PCMGMCycle_Private(pc,mglevels+levels-1,reason);CHKERRQ(ierr); 1084d0a8057SBarry Smith if (*reason) break; 1094d0a8057SBarry Smith } 1104d0a8057SBarry Smith if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS; 1114d0a8057SBarry Smith *outits = i; 1124b9ad928SBarry Smith PetscFunctionReturn(0); 1134b9ad928SBarry Smith } 1144b9ad928SBarry Smith 1154b9ad928SBarry Smith #undef __FUNCT__ 1163aeaf226SBarry Smith #define __FUNCT__ "PCReset_MG" 1173aeaf226SBarry Smith PetscErrorCode PCReset_MG(PC pc) 1183aeaf226SBarry Smith { 1193aeaf226SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1203aeaf226SBarry Smith PC_MG_Levels **mglevels = mg->levels; 1213aeaf226SBarry Smith PetscErrorCode ierr; 1223aeaf226SBarry Smith PetscInt i,n; 1233aeaf226SBarry Smith 1243aeaf226SBarry Smith PetscFunctionBegin; 1253aeaf226SBarry Smith if (mglevels) { 1263aeaf226SBarry Smith n = mglevels[0]->levels; 1273aeaf226SBarry Smith for (i=0; i<n-1; i++) { 1283aeaf226SBarry Smith ierr = VecDestroy(&mglevels[i+1]->r);CHKERRQ(ierr); 1293aeaf226SBarry Smith ierr = VecDestroy(&mglevels[i]->b);CHKERRQ(ierr); 1303aeaf226SBarry Smith ierr = VecDestroy(&mglevels[i]->x);CHKERRQ(ierr); 1313aeaf226SBarry Smith ierr = MatDestroy(&mglevels[i+1]->restrct);CHKERRQ(ierr); 1323aeaf226SBarry Smith ierr = MatDestroy(&mglevels[i+1]->interpolate);CHKERRQ(ierr); 13373250ac0SBarry Smith ierr = VecDestroy(&mglevels[i+1]->rscale);CHKERRQ(ierr); 1343aeaf226SBarry Smith } 1353aeaf226SBarry Smith 1363aeaf226SBarry Smith for (i=0; i<n; i++) { 1373aeaf226SBarry Smith ierr = MatDestroy(&mglevels[i]->A);CHKERRQ(ierr); 1383aeaf226SBarry Smith if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 1393aeaf226SBarry Smith ierr = KSPReset(mglevels[i]->smoothd);CHKERRQ(ierr); 1403aeaf226SBarry Smith } 1413aeaf226SBarry Smith ierr = KSPReset(mglevels[i]->smoothu);CHKERRQ(ierr); 1423aeaf226SBarry Smith } 1433aeaf226SBarry Smith } 1443aeaf226SBarry Smith PetscFunctionReturn(0); 1453aeaf226SBarry Smith } 1463aeaf226SBarry Smith 1473aeaf226SBarry Smith #undef __FUNCT__ 1489dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetLevels" 1494b9ad928SBarry Smith /*@C 15097177400SBarry Smith PCMGSetLevels - Sets the number of levels to use with MG. 1514b9ad928SBarry Smith Must be called before any other MG routine. 1524b9ad928SBarry Smith 153ad4df100SBarry Smith Logically Collective on PC 1544b9ad928SBarry Smith 1554b9ad928SBarry Smith Input Parameters: 1564b9ad928SBarry Smith + pc - the preconditioner context 1574b9ad928SBarry Smith . levels - the number of levels 1584b9ad928SBarry Smith - comms - optional communicators for each level; this is to allow solving the coarser problems 1590298fd71SBarry Smith on smaller sets of processors. Use NULL_OBJECT for default in Fortran 1604b9ad928SBarry Smith 1614b9ad928SBarry Smith Level: intermediate 1624b9ad928SBarry Smith 1634b9ad928SBarry Smith Notes: 1644b9ad928SBarry Smith If the number of levels is one then the multigrid uses the -mg_levels prefix 1654b9ad928SBarry Smith for setting the level options rather than the -mg_coarse prefix. 1664b9ad928SBarry Smith 1674b9ad928SBarry Smith .keywords: MG, set, levels, multigrid 1684b9ad928SBarry Smith 16997177400SBarry Smith .seealso: PCMGSetType(), PCMGGetLevels() 1704b9ad928SBarry Smith @*/ 1717087cfbeSBarry Smith PetscErrorCode PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms) 1724b9ad928SBarry Smith { 173dfbe8321SBarry Smith PetscErrorCode ierr; 174f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 175ce94432eSBarry Smith MPI_Comm comm; 1763aeaf226SBarry Smith PC_MG_Levels **mglevels = mg->levels; 177f3fbd535SBarry Smith PetscInt i; 178f3fbd535SBarry Smith PetscMPIInt size; 179f3fbd535SBarry Smith const char *prefix; 180f3fbd535SBarry Smith PC ipc; 1813aeaf226SBarry Smith PetscInt n; 1824b9ad928SBarry Smith 1834b9ad928SBarry Smith PetscFunctionBegin; 1840700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 185548767e0SJed Brown PetscValidLogicalCollectiveInt(pc,levels,2); 186ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 187548767e0SJed Brown if (mg->nlevels == levels) PetscFunctionReturn(0); 1883aeaf226SBarry Smith if (mglevels) { 1893aeaf226SBarry Smith /* changing the number of levels so free up the previous stuff */ 1903aeaf226SBarry Smith ierr = PCReset_MG(pc);CHKERRQ(ierr); 1913aeaf226SBarry Smith n = mglevels[0]->levels; 1923aeaf226SBarry Smith for (i=0; i<n; i++) { 1933aeaf226SBarry Smith if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 1943aeaf226SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr); 1953aeaf226SBarry Smith } 1963aeaf226SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr); 1973aeaf226SBarry Smith ierr = PetscFree(mglevels[i]);CHKERRQ(ierr); 1983aeaf226SBarry Smith } 1993aeaf226SBarry Smith ierr = PetscFree(mg->levels);CHKERRQ(ierr); 2003aeaf226SBarry Smith } 201f3fbd535SBarry Smith 202f3fbd535SBarry Smith mg->nlevels = levels; 203f3fbd535SBarry Smith 204785e854fSJed Brown ierr = PetscMalloc1(levels,&mglevels);CHKERRQ(ierr); 2053bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr); 206f3fbd535SBarry Smith 207f3fbd535SBarry Smith ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 208f3fbd535SBarry Smith 209a9db87a2SMatthew G Knepley mg->stageApply = 0; 210f3fbd535SBarry Smith for (i=0; i<levels; i++) { 211b00a9115SJed Brown ierr = PetscNewLog(pc,&mglevels[i]);CHKERRQ(ierr); 2122fa5cd67SKarl Rupp 213f3fbd535SBarry Smith mglevels[i]->level = i; 214f3fbd535SBarry Smith mglevels[i]->levels = levels; 215f3fbd535SBarry Smith mglevels[i]->cycles = PC_MG_CYCLE_V; 216336babb1SJed Brown mg->default_smoothu = 2; 217336babb1SJed Brown mg->default_smoothd = 2; 21863e6d426SJed Brown mglevels[i]->eventsmoothsetup = 0; 21963e6d426SJed Brown mglevels[i]->eventsmoothsolve = 0; 22063e6d426SJed Brown mglevels[i]->eventresidual = 0; 22163e6d426SJed Brown mglevels[i]->eventinterprestrict = 0; 222f3fbd535SBarry Smith 223f3fbd535SBarry Smith if (comms) comm = comms[i]; 224f3fbd535SBarry Smith ierr = KSPCreate(comm,&mglevels[i]->smoothd);CHKERRQ(ierr); 225422a814eSBarry Smith ierr = KSPSetErrorIfNotConverged(mglevels[i]->smoothd,pc->erroriffailure);CHKERRQ(ierr); 226336babb1SJed Brown ierr = KSPSetType(mglevels[i]->smoothd,KSPCHEBYSHEV);CHKERRQ(ierr); 2270059c7bdSJed Brown ierr = KSPSetConvergenceTest(mglevels[i]->smoothd,KSPConvergedSkip,NULL,NULL);CHKERRQ(ierr); 228336babb1SJed Brown ierr = KSPSetNormType(mglevels[i]->smoothd,KSP_NORM_NONE);CHKERRQ(ierr); 229336babb1SJed Brown ierr = KSPGetPC(mglevels[i]->smoothd,&ipc);CHKERRQ(ierr); 230336babb1SJed Brown ierr = PCSetType(ipc,PCSOR);CHKERRQ(ierr); 231f3fbd535SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr); 232336babb1SJed Brown ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, i ? mg->default_smoothd : 1);CHKERRQ(ierr); 233f3fbd535SBarry Smith ierr = KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);CHKERRQ(ierr); 234f3fbd535SBarry Smith 235f3fbd535SBarry Smith /* do special stuff for coarse grid */ 236f3fbd535SBarry Smith if (!i && levels > 1) { 237f3fbd535SBarry Smith ierr = KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");CHKERRQ(ierr); 238f3fbd535SBarry Smith 2399dbfc187SHong Zhang /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */ 240f3fbd535SBarry Smith ierr = KSPSetType(mglevels[0]->smoothd,KSPPREONLY);CHKERRQ(ierr); 241f3fbd535SBarry Smith ierr = KSPGetPC(mglevels[0]->smoothd,&ipc);CHKERRQ(ierr); 242f3fbd535SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 243f3fbd535SBarry Smith if (size > 1) { 2449dbfc187SHong Zhang KSP innerksp; 2459dbfc187SHong Zhang PC innerpc; 246f3fbd535SBarry Smith ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr); 2479dbfc187SHong Zhang ierr = PCRedundantGetKSP(ipc,&innerksp);CHKERRQ(ierr); 2489dbfc187SHong Zhang ierr = KSPGetPC(innerksp,&innerpc);CHKERRQ(ierr); 2499dbfc187SHong Zhang ierr = PCFactorSetShiftType(innerpc,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr); 250f3fbd535SBarry Smith } else { 251f3fbd535SBarry Smith ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr); 2529dbfc187SHong Zhang ierr = PCFactorSetShiftType(ipc,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr); 253f3fbd535SBarry Smith } 254f3fbd535SBarry Smith } else { 255f3fbd535SBarry Smith char tprefix[128]; 256f3fbd535SBarry Smith sprintf(tprefix,"mg_levels_%d_",(int)i); 257f3fbd535SBarry Smith ierr = KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);CHKERRQ(ierr); 258f3fbd535SBarry Smith } 2593bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[i]->smoothd);CHKERRQ(ierr); 2602fa5cd67SKarl Rupp 261f3fbd535SBarry Smith mglevels[i]->smoothu = mglevels[i]->smoothd; 26231567311SBarry Smith mg->rtol = 0.0; 26331567311SBarry Smith mg->abstol = 0.0; 26431567311SBarry Smith mg->dtol = 0.0; 26531567311SBarry Smith mg->ttol = 0.0; 26631567311SBarry Smith mg->cyclesperpcapply = 1; 267f3fbd535SBarry Smith } 26831567311SBarry Smith mg->am = PC_MG_MULTIPLICATIVE; 269f3fbd535SBarry Smith mg->levels = mglevels; 2704b9ad928SBarry Smith pc->ops->applyrichardson = PCApplyRichardson_MG; 2714b9ad928SBarry Smith PetscFunctionReturn(0); 2724b9ad928SBarry Smith } 2734b9ad928SBarry Smith 274c07bf074SBarry Smith 275c07bf074SBarry Smith #undef __FUNCT__ 276c07bf074SBarry Smith #define __FUNCT__ "PCDestroy_MG" 277c07bf074SBarry Smith PetscErrorCode PCDestroy_MG(PC pc) 278c07bf074SBarry Smith { 279c07bf074SBarry Smith PetscErrorCode ierr; 280a06653b4SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 281a06653b4SBarry Smith PC_MG_Levels **mglevels = mg->levels; 282a06653b4SBarry Smith PetscInt i,n; 283c07bf074SBarry Smith 284c07bf074SBarry Smith PetscFunctionBegin; 285a06653b4SBarry Smith ierr = PCReset_MG(pc);CHKERRQ(ierr); 286a06653b4SBarry Smith if (mglevels) { 287a06653b4SBarry Smith n = mglevels[0]->levels; 288a06653b4SBarry Smith for (i=0; i<n; i++) { 289a06653b4SBarry Smith if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 2906bf464f9SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr); 291a06653b4SBarry Smith } 2926bf464f9SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr); 293a06653b4SBarry Smith ierr = PetscFree(mglevels[i]);CHKERRQ(ierr); 294a06653b4SBarry Smith } 295a06653b4SBarry Smith ierr = PetscFree(mg->levels);CHKERRQ(ierr); 296a06653b4SBarry Smith } 297c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 298f3fbd535SBarry Smith PetscFunctionReturn(0); 299f3fbd535SBarry Smith } 300f3fbd535SBarry Smith 301f3fbd535SBarry Smith 302f3fbd535SBarry Smith 30309573ac7SBarry Smith extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**); 30409573ac7SBarry Smith extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**); 30509573ac7SBarry Smith extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**); 306f3fbd535SBarry Smith 307f3fbd535SBarry Smith /* 308f3fbd535SBarry Smith PCApply_MG - Runs either an additive, multiplicative, Kaskadic 309f3fbd535SBarry Smith or full cycle of multigrid. 310f3fbd535SBarry Smith 311f3fbd535SBarry Smith Note: 312f3fbd535SBarry Smith A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle(). 313f3fbd535SBarry Smith */ 314f3fbd535SBarry Smith #undef __FUNCT__ 315f3fbd535SBarry Smith #define __FUNCT__ "PCApply_MG" 316f3fbd535SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x) 317f3fbd535SBarry Smith { 318f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 319f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 320f3fbd535SBarry Smith PetscErrorCode ierr; 321f3fbd535SBarry Smith PetscInt levels = mglevels[0]->levels,i; 322f3fbd535SBarry Smith 323f3fbd535SBarry Smith PetscFunctionBegin; 324a9db87a2SMatthew G Knepley if (mg->stageApply) {ierr = PetscLogStagePush(mg->stageApply);CHKERRQ(ierr);} 325e1d8e5deSBarry Smith /* When the DM is supplying the matrix then it will not exist until here */ 326298cc208SBarry Smith for (i=0; i<levels; i++) { 327e1d8e5deSBarry Smith if (!mglevels[i]->A) { 32823ee1639SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr); 329298cc208SBarry Smith ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr); 330e1d8e5deSBarry Smith } 331e1d8e5deSBarry Smith } 332e1d8e5deSBarry Smith 333f3fbd535SBarry Smith mglevels[levels-1]->b = b; 334f3fbd535SBarry Smith mglevels[levels-1]->x = x; 33531567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 336f3fbd535SBarry Smith ierr = VecSet(x,0.0);CHKERRQ(ierr); 33731567311SBarry Smith for (i=0; i<mg->cyclesperpcapply; i++) { 3380298fd71SBarry Smith ierr = PCMGMCycle_Private(pc,mglevels+levels-1,NULL);CHKERRQ(ierr); 339f3fbd535SBarry Smith } 3402fa5cd67SKarl Rupp } else if (mg->am == PC_MG_ADDITIVE) { 34131567311SBarry Smith ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr); 3422fa5cd67SKarl Rupp } else if (mg->am == PC_MG_KASKADE) { 34331567311SBarry Smith ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr); 3442fa5cd67SKarl Rupp } else { 345f3fbd535SBarry Smith ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr); 346f3fbd535SBarry Smith } 347a9db87a2SMatthew G Knepley if (mg->stageApply) {ierr = PetscLogStagePop();CHKERRQ(ierr);} 348f3fbd535SBarry Smith PetscFunctionReturn(0); 349f3fbd535SBarry Smith } 350f3fbd535SBarry Smith 351f3fbd535SBarry Smith 352f3fbd535SBarry Smith #undef __FUNCT__ 353f3fbd535SBarry Smith #define __FUNCT__ "PCSetFromOptions_MG" 3548c34d3f5SBarry Smith PetscErrorCode PCSetFromOptions_MG(PetscOptions *PetscOptionsObject,PC pc) 355f3fbd535SBarry Smith { 356f3fbd535SBarry Smith PetscErrorCode ierr; 357f3fbd535SBarry Smith PetscInt m,levels = 1,cycles; 358c91913d3SJed Brown PetscBool flg,set; 359f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 360f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 361f3fbd535SBarry Smith PCMGType mgtype; 362f3fbd535SBarry Smith PCMGCycleType mgctype; 363f3fbd535SBarry Smith 364f3fbd535SBarry Smith PetscFunctionBegin; 365e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Multigrid options");CHKERRQ(ierr); 3663aeaf226SBarry Smith if (!mg->levels) { 367f3fbd535SBarry Smith ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr); 368eb3f98d2SBarry Smith if (!flg && pc->dm) { 369eb3f98d2SBarry Smith ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr); 370eb3f98d2SBarry Smith levels++; 3713aeaf226SBarry Smith mg->usedmfornumberoflevels = PETSC_TRUE; 372eb3f98d2SBarry Smith } 3730298fd71SBarry Smith ierr = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr); 374f3fbd535SBarry Smith } 3753aeaf226SBarry Smith mglevels = mg->levels; 3763aeaf226SBarry Smith 377f3fbd535SBarry Smith mgctype = (PCMGCycleType) mglevels[0]->cycles; 378f3fbd535SBarry Smith ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr); 379f3fbd535SBarry Smith if (flg) { 380f3fbd535SBarry Smith ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr); 3812fa5cd67SKarl Rupp } 382f3fbd535SBarry Smith flg = PETSC_FALSE; 383c91913d3SJed Brown ierr = PetscOptionsBool("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",flg,&flg,&set);CHKERRQ(ierr); 384c91913d3SJed Brown if (set) { 385c91913d3SJed Brown ierr = PCMGSetGalerkin(pc,flg);CHKERRQ(ierr); 386f3fbd535SBarry Smith } 38794ae4db5SBarry Smith ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",mg->default_smoothu,&m,&flg);CHKERRQ(ierr); 388f3fbd535SBarry Smith if (flg) { 389f3fbd535SBarry Smith ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr); 390f3fbd535SBarry Smith } 39194ae4db5SBarry Smith ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",mg->default_smoothd,&m,&flg);CHKERRQ(ierr); 392f3fbd535SBarry Smith if (flg) { 393f3fbd535SBarry Smith ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr); 394f3fbd535SBarry Smith } 39531567311SBarry Smith mgtype = mg->am; 396f3fbd535SBarry Smith ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr); 397f3fbd535SBarry Smith if (flg) { 398f3fbd535SBarry Smith ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr); 399f3fbd535SBarry Smith } 40031567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 40131567311SBarry Smith ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGSetLevels",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr); 402f3fbd535SBarry Smith if (flg) { 403f3fbd535SBarry Smith ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr); 404f3fbd535SBarry Smith } 405f3fbd535SBarry Smith } 406f3fbd535SBarry Smith flg = PETSC_FALSE; 4070298fd71SBarry Smith ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,NULL);CHKERRQ(ierr); 408f3fbd535SBarry Smith if (flg) { 409f3fbd535SBarry Smith PetscInt i; 410f3fbd535SBarry Smith char eventname[128]; 411ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 412f3fbd535SBarry Smith levels = mglevels[0]->levels; 413f3fbd535SBarry Smith for (i=0; i<levels; i++) { 414f3fbd535SBarry Smith sprintf(eventname,"MGSetup Level %d",(int)i); 41563e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr); 416f3fbd535SBarry Smith sprintf(eventname,"MGSmooth Level %d",(int)i); 41763e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr); 418f3fbd535SBarry Smith if (i) { 419f3fbd535SBarry Smith sprintf(eventname,"MGResid Level %d",(int)i); 42063e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr); 421f3fbd535SBarry Smith sprintf(eventname,"MGInterp Level %d",(int)i); 42263e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr); 423f3fbd535SBarry Smith } 424f3fbd535SBarry Smith } 425eec35531SMatthew G Knepley 426ec60ca73SSatish Balay #if defined(PETSC_USE_LOG) 427eec35531SMatthew G Knepley { 428eec35531SMatthew G Knepley const char *sname = "MG Apply"; 429eec35531SMatthew G Knepley PetscStageLog stageLog; 430eec35531SMatthew G Knepley PetscInt st; 431eec35531SMatthew G Knepley 432eec35531SMatthew G Knepley PetscFunctionBegin; 433eec35531SMatthew G Knepley ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr); 434eec35531SMatthew G Knepley for (st = 0; st < stageLog->numStages; ++st) { 435eec35531SMatthew G Knepley PetscBool same; 436eec35531SMatthew G Knepley 437eec35531SMatthew G Knepley ierr = PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);CHKERRQ(ierr); 4382fa5cd67SKarl Rupp if (same) mg->stageApply = st; 439eec35531SMatthew G Knepley } 440eec35531SMatthew G Knepley if (!mg->stageApply) { 441eec35531SMatthew G Knepley ierr = PetscLogStageRegister(sname, &mg->stageApply);CHKERRQ(ierr); 442eec35531SMatthew G Knepley } 443eec35531SMatthew G Knepley } 444ec60ca73SSatish Balay #endif 445f3fbd535SBarry Smith } 446f3fbd535SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 447f3fbd535SBarry Smith PetscFunctionReturn(0); 448f3fbd535SBarry Smith } 449f3fbd535SBarry Smith 4506a6fc655SJed Brown const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0}; 4516a6fc655SJed Brown const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0}; 452f3fbd535SBarry Smith 4539804daf3SBarry Smith #include <petscdraw.h> 454f3fbd535SBarry Smith #undef __FUNCT__ 455f3fbd535SBarry Smith #define __FUNCT__ "PCView_MG" 456f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer) 457f3fbd535SBarry Smith { 458f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 459f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 460f3fbd535SBarry Smith PetscErrorCode ierr; 461e3deeeafSJed Brown PetscInt levels = mglevels ? mglevels[0]->levels : 0,i; 462219b1045SBarry Smith PetscBool iascii,isbinary,isdraw; 463f3fbd535SBarry Smith 464f3fbd535SBarry Smith PetscFunctionBegin; 465251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 4665b0b0462SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 467219b1045SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 468f3fbd535SBarry Smith if (iascii) { 469e3deeeafSJed Brown const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown"; 470e3deeeafSJed Brown ierr = PetscViewerASCIIPrintf(viewer," MG: type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,cyclename);CHKERRQ(ierr); 47131567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 47231567311SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr); 473f3fbd535SBarry Smith } 474218a07d4SBarry Smith if (mg->galerkin) { 475f3fbd535SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 4764f66f45eSBarry Smith } else { 4774f66f45eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 478f3fbd535SBarry Smith } 4795adeb434SBarry Smith if (mg->view){ 4805adeb434SBarry Smith ierr = (*mg->view)(pc,viewer);CHKERRQ(ierr); 4815adeb434SBarry Smith } 482f3fbd535SBarry Smith for (i=0; i<levels; i++) { 483f3fbd535SBarry Smith if (!i) { 48489cce641SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);CHKERRQ(ierr); 485f3fbd535SBarry Smith } else { 48689cce641SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 487f3fbd535SBarry Smith } 488f3fbd535SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 489f3fbd535SBarry Smith ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 490f3fbd535SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 491f3fbd535SBarry Smith if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) { 492f3fbd535SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 493f3fbd535SBarry Smith } else if (i) { 49489cce641SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 495f3fbd535SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 496f3fbd535SBarry Smith ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 497f3fbd535SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 498f3fbd535SBarry Smith } 499f3fbd535SBarry Smith } 5005b0b0462SBarry Smith } else if (isbinary) { 5015b0b0462SBarry Smith for (i=levels-1; i>=0; i--) { 5025b0b0462SBarry Smith ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 5035b0b0462SBarry Smith if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) { 5045b0b0462SBarry Smith ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 5055b0b0462SBarry Smith } 5065b0b0462SBarry Smith } 507219b1045SBarry Smith } else if (isdraw) { 508219b1045SBarry Smith PetscDraw draw; 509b4375e8dSPeter Brune PetscReal x,w,y,bottom,th; 510219b1045SBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 511219b1045SBarry Smith ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr); 5120298fd71SBarry Smith ierr = PetscDrawStringGetSize(draw,NULL,&th);CHKERRQ(ierr); 513219b1045SBarry Smith bottom = y - th; 514219b1045SBarry Smith for (i=levels-1; i>=0; i--) { 515b4375e8dSPeter Brune if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) { 516219b1045SBarry Smith ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr); 517219b1045SBarry Smith ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 518219b1045SBarry Smith ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 519b4375e8dSPeter Brune } else { 520b4375e8dSPeter Brune w = 0.5*PetscMin(1.0-x,x); 521b4375e8dSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr); 522b4375e8dSPeter Brune ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 523b4375e8dSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 524b4375e8dSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr); 525b4375e8dSPeter Brune ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 526b4375e8dSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 527b4375e8dSPeter Brune } 5280298fd71SBarry Smith ierr = PetscDrawGetBoundingBox(draw,NULL,&bottom,NULL,NULL);CHKERRQ(ierr); 5291cd381d6SBarry Smith bottom -= th; 530219b1045SBarry Smith } 5315b0b0462SBarry Smith } 532f3fbd535SBarry Smith PetscFunctionReturn(0); 533f3fbd535SBarry Smith } 534f3fbd535SBarry Smith 535af0996ceSBarry Smith #include <petsc/private/dmimpl.h> 536af0996ceSBarry Smith #include <petsc/private/kspimpl.h> 537cab2e9ccSBarry Smith 538f3fbd535SBarry Smith /* 539f3fbd535SBarry Smith Calls setup for the KSP on each level 540f3fbd535SBarry Smith */ 541f3fbd535SBarry Smith #undef __FUNCT__ 542f3fbd535SBarry Smith #define __FUNCT__ "PCSetUp_MG" 543f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc) 544f3fbd535SBarry Smith { 545f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 546f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 547f3fbd535SBarry Smith PetscErrorCode ierr; 548f3fbd535SBarry Smith PetscInt i,n = mglevels[0]->levels; 54998e04984SBarry Smith PC cpc; 5509c7ffb39SBarry Smith PetscBool preonly,lu,redundant,cholesky,svd,dump = PETSC_FALSE,opsset,use_amat,missinginterpolate = PETSC_FALSE; 551f3fbd535SBarry Smith Mat dA,dB; 552f3fbd535SBarry Smith Vec tvec; 553218a07d4SBarry Smith DM *dms; 554649052a6SBarry Smith PetscViewer viewer = 0; 555f3fbd535SBarry Smith 556f3fbd535SBarry Smith PetscFunctionBegin; 55701bc414fSDmitry Karpeev /* FIX: Move this to PCSetFromOptions_MG? */ 5583aeaf226SBarry Smith if (mg->usedmfornumberoflevels) { 5593aeaf226SBarry Smith PetscInt levels; 5603aeaf226SBarry Smith ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr); 5613aeaf226SBarry Smith levels++; 5623aeaf226SBarry Smith if (levels > n) { /* the problem is now being solved on a finer grid */ 5630298fd71SBarry Smith ierr = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr); 5643aeaf226SBarry Smith n = levels; 5653aeaf226SBarry Smith ierr = PCSetFromOptions(pc);CHKERRQ(ierr); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */ 5663aeaf226SBarry Smith mglevels = mg->levels; 5673aeaf226SBarry Smith } 5683aeaf226SBarry Smith } 56998e04984SBarry Smith ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr); 5703aeaf226SBarry Smith 571f3fbd535SBarry Smith 572f3fbd535SBarry Smith /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */ 573f3fbd535SBarry Smith /* so use those from global PC */ 574f3fbd535SBarry Smith /* Is this what we always want? What if user wants to keep old one? */ 5750298fd71SBarry Smith ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset);CHKERRQ(ierr); 57698e04984SBarry Smith if (opsset) { 57798e04984SBarry Smith Mat mmat; 57823ee1639SBarry Smith ierr = KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat);CHKERRQ(ierr); 57998e04984SBarry Smith if (mmat == pc->pmat) opsset = PETSC_FALSE; 58098e04984SBarry Smith } 5814a5f07a5SMark F. Adams 5824a5f07a5SMark F. Adams if (!opsset) { 58371b23a65SMark F. Adams ierr = PCGetUseAmat(pc,&use_amat);CHKERRQ(ierr); 5844a5f07a5SMark F. Adams if(use_amat){ 585f3fbd535SBarry 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); 58623ee1639SBarry Smith ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat);CHKERRQ(ierr); 587f3fbd535SBarry Smith } 5884a5f07a5SMark F. Adams else { 5894a5f07a5SMark 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); 59023ee1639SBarry Smith ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat);CHKERRQ(ierr); 5914a5f07a5SMark F. Adams } 5924a5f07a5SMark F. Adams } 593f3fbd535SBarry Smith 5949c7ffb39SBarry Smith for (i=n-1; i>0; i--) { 5959c7ffb39SBarry Smith if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) { 5969c7ffb39SBarry Smith missinginterpolate = PETSC_TRUE; 5979c7ffb39SBarry Smith continue; 5989c7ffb39SBarry Smith } 5999c7ffb39SBarry Smith } 6009c7ffb39SBarry Smith /* 6019c7ffb39SBarry Smith Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS) 6029c7ffb39SBarry Smith Skipping for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs? 6039c7ffb39SBarry Smith */ 6049c7ffb39SBarry Smith if (missinginterpolate && pc->dm && mg->galerkin != 2 && !pc->setupcalled) { 6052d2b81a6SBarry Smith /* construct the interpolation from the DMs */ 6062e499ae9SBarry Smith Mat p; 60773250ac0SBarry Smith Vec rscale; 608785e854fSJed Brown ierr = PetscMalloc1(n,&dms);CHKERRQ(ierr); 609218a07d4SBarry Smith dms[n-1] = pc->dm; 610ef1267afSMatthew G. Knepley /* Separately create them so we do not get DMKSP interference between levels */ 611ef1267afSMatthew G. Knepley for (i=n-2; i>-1; i--) {ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr);} 612218a07d4SBarry Smith for (i=n-2; i>-1; i--) { 613942e3340SBarry Smith DMKSP kdm; 6143c0c59f3SBarry Smith ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr); 6153c0c59f3SBarry Smith if (mg->galerkin) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);} 616942e3340SBarry Smith ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr); 617d1a3e677SJed Brown /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take 618d1a3e677SJed Brown * a bitwise OR of computing the matrix, RHS, and initial iterate. */ 6190298fd71SBarry Smith kdm->ops->computerhs = NULL; 6200298fd71SBarry Smith kdm->rhsctx = NULL; 62124c3aa18SBarry Smith if (!mglevels[i+1]->interpolate) { 622e727c939SJed Brown ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr); 6232d2b81a6SBarry Smith ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr); 62400ac8be1SBarry Smith if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);} 62573250ac0SBarry Smith ierr = VecDestroy(&rscale);CHKERRQ(ierr); 6266bf464f9SBarry Smith ierr = MatDestroy(&p);CHKERRQ(ierr); 627218a07d4SBarry Smith } 62824c3aa18SBarry Smith } 6292d2b81a6SBarry Smith 630ef1267afSMatthew G. Knepley for (i=n-2; i>-1; i--) {ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);} 6312d2b81a6SBarry Smith ierr = PetscFree(dms);CHKERRQ(ierr); 632ce4cda84SJed Brown } 633cab2e9ccSBarry Smith 634ce4cda84SJed Brown if (pc->dm && !pc->setupcalled) { 635ce4cda84SJed Brown /* finest smoother also gets DM but it is not active, independent of whether galerkin==2 */ 636cab2e9ccSBarry Smith ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr); 637cab2e9ccSBarry Smith ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr); 638218a07d4SBarry Smith } 639218a07d4SBarry Smith 640c91913d3SJed Brown if (mg->galerkin == 1) { 641f3fbd535SBarry Smith Mat B; 642f3fbd535SBarry Smith /* currently only handle case where mat and pmat are the same on coarser levels */ 64323ee1639SBarry Smith ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB);CHKERRQ(ierr); 644f3fbd535SBarry Smith if (!pc->setupcalled) { 645f3fbd535SBarry Smith for (i=n-2; i>-1; i--) { 646b9367914SBarry 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"); 647b9367914SBarry Smith if (!mglevels[i+1]->interpolate) { 648b9367914SBarry Smith ierr = PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct);CHKERRQ(ierr); 649b9367914SBarry Smith } 650b9367914SBarry Smith if (!mglevels[i+1]->restrct) { 651b9367914SBarry Smith ierr = PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate);CHKERRQ(ierr); 652b9367914SBarry Smith } 653b9367914SBarry Smith if (mglevels[i+1]->interpolate == mglevels[i+1]->restrct) { 654f3fbd535SBarry Smith ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 6557d537d92SJed Brown } else { 6567d537d92SJed Brown ierr = MatMatMatMult(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 6577d537d92SJed Brown } 65823ee1639SBarry Smith ierr = KSPSetOperators(mglevels[i]->smoothd,B,B);CHKERRQ(ierr); 659f3fbd535SBarry Smith if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);} 660f3fbd535SBarry Smith dB = B; 661f3fbd535SBarry Smith } 662cd9507b2SBarry Smith if (n > 1) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);} 663f3fbd535SBarry Smith } else { 664f3fbd535SBarry Smith for (i=n-2; i>-1; i--) { 665b9367914SBarry 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"); 666b9367914SBarry Smith if (!mglevels[i+1]->interpolate) { 667b9367914SBarry Smith ierr = PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct);CHKERRQ(ierr); 668b9367914SBarry Smith } 669b9367914SBarry Smith if (!mglevels[i+1]->restrct) { 670b9367914SBarry Smith ierr = PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate);CHKERRQ(ierr); 671b9367914SBarry Smith } 67223ee1639SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&B);CHKERRQ(ierr); 673b9367914SBarry Smith if (mglevels[i+1]->interpolate == mglevels[i+1]->restrct) { 674f3fbd535SBarry Smith ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 6757d537d92SJed Brown } else { 6767d537d92SJed Brown ierr = MatMatMatMult(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 6777d537d92SJed Brown } 67823ee1639SBarry Smith ierr = KSPSetOperators(mglevels[i]->smoothd,B,B);CHKERRQ(ierr); 679f3fbd535SBarry Smith dB = B; 680f3fbd535SBarry Smith } 681f3fbd535SBarry Smith } 682ce4cda84SJed Brown } else if (!mg->galerkin && pc->dm && pc->dm->x) { 683cab2e9ccSBarry Smith /* need to restrict Jacobian location to coarser meshes for evaluation */ 684cab2e9ccSBarry Smith for (i=n-2; i>-1; i--) { 685c88c5224SJed Brown Mat R; 686c88c5224SJed Brown Vec rscale; 687cab2e9ccSBarry Smith if (!mglevels[i]->smoothd->dm->x) { 688cab2e9ccSBarry Smith Vec *vecs; 6892a7a6963SBarry Smith ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vecs,0,NULL);CHKERRQ(ierr); 6902fa5cd67SKarl Rupp 691cab2e9ccSBarry Smith mglevels[i]->smoothd->dm->x = vecs[0]; 6922fa5cd67SKarl Rupp 693cab2e9ccSBarry Smith ierr = PetscFree(vecs);CHKERRQ(ierr); 694cab2e9ccSBarry Smith } 695c88c5224SJed Brown ierr = PCMGGetRestriction(pc,i+1,&R);CHKERRQ(ierr); 696c88c5224SJed Brown ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr); 697c88c5224SJed Brown ierr = MatRestrict(R,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);CHKERRQ(ierr); 698c88c5224SJed Brown ierr = VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,rscale);CHKERRQ(ierr); 699cab2e9ccSBarry Smith } 700f3fbd535SBarry Smith } 701ccceb50cSJed Brown if (!mg->galerkin && pc->dm) { 702caa4e7f2SJed Brown for (i=n-2; i>=0; i--) { 703ccceb50cSJed Brown DM dmfine,dmcoarse; 704caa4e7f2SJed Brown Mat Restrict,Inject; 705caa4e7f2SJed Brown Vec rscale; 706ccceb50cSJed Brown ierr = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr); 707ccceb50cSJed Brown ierr = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr); 708caa4e7f2SJed Brown ierr = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr); 709caa4e7f2SJed Brown ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr); 7100298fd71SBarry Smith Inject = NULL; /* Callback should create it if it needs Injection */ 711caa4e7f2SJed Brown ierr = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr); 712caa4e7f2SJed Brown } 713caa4e7f2SJed Brown } 714f3fbd535SBarry Smith 715f3fbd535SBarry Smith if (!pc->setupcalled) { 716f3fbd535SBarry Smith for (i=0; i<n; i++) { 717f3fbd535SBarry Smith ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr); 718f3fbd535SBarry Smith } 719f3fbd535SBarry Smith for (i=1; i<n; i++) { 720f3fbd535SBarry Smith if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) { 721f3fbd535SBarry Smith ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr); 722f3fbd535SBarry Smith } 723f3fbd535SBarry Smith } 724f3fbd535SBarry Smith for (i=1; i<n; i++) { 725c88c5224SJed Brown ierr = PCMGGetInterpolation(pc,i,&mglevels[i]->interpolate);CHKERRQ(ierr); 726c88c5224SJed Brown ierr = PCMGGetRestriction(pc,i,&mglevels[i]->restrct);CHKERRQ(ierr); 727f3fbd535SBarry Smith } 728f3fbd535SBarry Smith for (i=0; i<n-1; i++) { 729f3fbd535SBarry Smith if (!mglevels[i]->b) { 730f3fbd535SBarry Smith Vec *vec; 7312a7a6963SBarry Smith ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr); 732f3fbd535SBarry Smith ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr); 7336bf464f9SBarry Smith ierr = VecDestroy(vec);CHKERRQ(ierr); 734f3fbd535SBarry Smith ierr = PetscFree(vec);CHKERRQ(ierr); 735f3fbd535SBarry Smith } 736f3fbd535SBarry Smith if (!mglevels[i]->r && i) { 737f3fbd535SBarry Smith ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 738f3fbd535SBarry Smith ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr); 7396bf464f9SBarry Smith ierr = VecDestroy(&tvec);CHKERRQ(ierr); 740f3fbd535SBarry Smith } 741f3fbd535SBarry Smith if (!mglevels[i]->x) { 742f3fbd535SBarry Smith ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 743f3fbd535SBarry Smith ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr); 7446bf464f9SBarry Smith ierr = VecDestroy(&tvec);CHKERRQ(ierr); 745f3fbd535SBarry Smith } 746f3fbd535SBarry Smith } 747f3fbd535SBarry Smith if (n != 1 && !mglevels[n-1]->r) { 748f3fbd535SBarry Smith /* PCMGSetR() on the finest level if user did not supply it */ 749f3fbd535SBarry Smith Vec *vec; 7502a7a6963SBarry Smith ierr = KSPCreateVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr); 751f3fbd535SBarry Smith ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr); 7526bf464f9SBarry Smith ierr = VecDestroy(vec);CHKERRQ(ierr); 753f3fbd535SBarry Smith ierr = PetscFree(vec);CHKERRQ(ierr); 754f3fbd535SBarry Smith } 755f3fbd535SBarry Smith } 756f3fbd535SBarry Smith 75798e04984SBarry Smith if (pc->dm) { 75898e04984SBarry Smith /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */ 75998e04984SBarry Smith for (i=0; i<n-1; i++) { 76098e04984SBarry Smith if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX; 76198e04984SBarry Smith } 76298e04984SBarry Smith } 763f3fbd535SBarry Smith 764f3fbd535SBarry Smith for (i=1; i<n; i++) { 7652a21e185SBarry Smith if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1){ 766f3fbd535SBarry Smith /* if doing only down then initial guess is zero */ 767f3fbd535SBarry Smith ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 768f3fbd535SBarry Smith } 76963e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 770f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr); 77163e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 772d42688cbSBarry Smith if (!mglevels[i]->residual) { 773d42688cbSBarry Smith Mat mat; 77423ee1639SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&mat);CHKERRQ(ierr); 77554b2cd4bSJed Brown ierr = PCMGSetResidual(pc,i,PCMGResidualDefault,mat);CHKERRQ(ierr); 776d42688cbSBarry Smith } 777f3fbd535SBarry Smith } 778f3fbd535SBarry Smith for (i=1; i<n; i++) { 779f3fbd535SBarry Smith if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) { 780f3fbd535SBarry Smith Mat downmat,downpmat; 781f3fbd535SBarry Smith 782f3fbd535SBarry Smith /* check if operators have been set for up, if not use down operators to set them */ 7830298fd71SBarry Smith ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);CHKERRQ(ierr); 784f3fbd535SBarry Smith if (!opsset) { 78523ee1639SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);CHKERRQ(ierr); 78623ee1639SBarry Smith ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat);CHKERRQ(ierr); 787f3fbd535SBarry Smith } 788f3fbd535SBarry Smith 789f3fbd535SBarry Smith ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr); 79063e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 791f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr); 79263e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 793f3fbd535SBarry Smith } 794f3fbd535SBarry Smith } 795f3fbd535SBarry Smith 796f3fbd535SBarry Smith /* 797f3fbd535SBarry Smith If coarse solver is not direct method then DO NOT USE preonly 798f3fbd535SBarry Smith */ 799251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mglevels[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr); 800f3fbd535SBarry Smith if (preonly) { 801251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr); 802251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr); 803251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr); 804251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)cpc,PCSVD,&svd);CHKERRQ(ierr); 805fd1303e9SJungho Lee if (!lu && !redundant && !cholesky && !svd) { 806f3fbd535SBarry Smith ierr = KSPSetType(mglevels[0]->smoothd,KSPGMRES);CHKERRQ(ierr); 807f3fbd535SBarry Smith } 808f3fbd535SBarry Smith } 809f3fbd535SBarry Smith 81063e6d426SJed Brown if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 811f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr); 81263e6d426SJed Brown if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 813f3fbd535SBarry Smith 814f3fbd535SBarry Smith /* 815f3fbd535SBarry Smith Dump the interpolation/restriction matrices plus the 816e3c5b3baSBarry Smith Jacobian/stiffness on each level. This allows MATLAB users to 817f3fbd535SBarry Smith easily check if the Galerkin condition A_c = R A_f R^T is satisfied. 818f3fbd535SBarry Smith 819f3fbd535SBarry Smith Only support one or the other at the same time. 820f3fbd535SBarry Smith */ 821f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER) 8220298fd71SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);CHKERRQ(ierr); 823ce94432eSBarry Smith if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc)); 824f3fbd535SBarry Smith dump = PETSC_FALSE; 825f3fbd535SBarry Smith #endif 8260298fd71SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);CHKERRQ(ierr); 827ce94432eSBarry Smith if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc)); 828f3fbd535SBarry Smith 829f3fbd535SBarry Smith if (viewer) { 830f3fbd535SBarry Smith for (i=1; i<n; i++) { 831f3fbd535SBarry Smith ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr); 832f3fbd535SBarry Smith } 833f3fbd535SBarry Smith for (i=0; i<n; i++) { 834f3fbd535SBarry Smith ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr); 835f3fbd535SBarry Smith ierr = MatView(pc->mat,viewer);CHKERRQ(ierr); 836f3fbd535SBarry Smith } 837f3fbd535SBarry Smith } 838f3fbd535SBarry Smith PetscFunctionReturn(0); 839f3fbd535SBarry Smith } 840f3fbd535SBarry Smith 841f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/ 842f3fbd535SBarry Smith 843f3fbd535SBarry Smith #undef __FUNCT__ 8449dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetLevels" 8454b9ad928SBarry Smith /*@ 84697177400SBarry Smith PCMGGetLevels - Gets the number of levels to use with MG. 8474b9ad928SBarry Smith 8484b9ad928SBarry Smith Not Collective 8494b9ad928SBarry Smith 8504b9ad928SBarry Smith Input Parameter: 8514b9ad928SBarry Smith . pc - the preconditioner context 8524b9ad928SBarry Smith 8534b9ad928SBarry Smith Output parameter: 8544b9ad928SBarry Smith . levels - the number of levels 8554b9ad928SBarry Smith 8564b9ad928SBarry Smith Level: advanced 8574b9ad928SBarry Smith 8584b9ad928SBarry Smith .keywords: MG, get, levels, multigrid 8594b9ad928SBarry Smith 86097177400SBarry Smith .seealso: PCMGSetLevels() 8614b9ad928SBarry Smith @*/ 8627087cfbeSBarry Smith PetscErrorCode PCMGGetLevels(PC pc,PetscInt *levels) 8634b9ad928SBarry Smith { 864f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 8654b9ad928SBarry Smith 8664b9ad928SBarry Smith PetscFunctionBegin; 8670700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 8684482741eSBarry Smith PetscValidIntPointer(levels,2); 869f3fbd535SBarry Smith *levels = mg->nlevels; 8704b9ad928SBarry Smith PetscFunctionReturn(0); 8714b9ad928SBarry Smith } 8724b9ad928SBarry Smith 8734b9ad928SBarry Smith #undef __FUNCT__ 8749dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetType" 8754b9ad928SBarry Smith /*@ 87697177400SBarry Smith PCMGSetType - Determines the form of multigrid to use: 8774b9ad928SBarry Smith multiplicative, additive, full, or the Kaskade algorithm. 8784b9ad928SBarry Smith 879ad4df100SBarry Smith Logically Collective on PC 8804b9ad928SBarry Smith 8814b9ad928SBarry Smith Input Parameters: 8824b9ad928SBarry Smith + pc - the preconditioner context 8839dcbbd2bSBarry Smith - form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE, 8849dcbbd2bSBarry Smith PC_MG_FULL, PC_MG_KASKADE 8854b9ad928SBarry Smith 8864b9ad928SBarry Smith Options Database Key: 8874b9ad928SBarry Smith . -pc_mg_type <form> - Sets <form>, one of multiplicative, 8884b9ad928SBarry Smith additive, full, kaskade 8894b9ad928SBarry Smith 8904b9ad928SBarry Smith Level: advanced 8914b9ad928SBarry Smith 8924b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid 8934b9ad928SBarry Smith 89497177400SBarry Smith .seealso: PCMGSetLevels() 8954b9ad928SBarry Smith @*/ 8967087cfbeSBarry Smith PetscErrorCode PCMGSetType(PC pc,PCMGType form) 8974b9ad928SBarry Smith { 898f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 8994b9ad928SBarry Smith 9004b9ad928SBarry Smith PetscFunctionBegin; 9010700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 902c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc,form,2); 90331567311SBarry Smith mg->am = form; 9049dcbbd2bSBarry Smith if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG; 905c60c7ad4SBarry Smith else pc->ops->applyrichardson = NULL; 906c60c7ad4SBarry Smith PetscFunctionReturn(0); 907c60c7ad4SBarry Smith } 908c60c7ad4SBarry Smith 909c60c7ad4SBarry Smith /*@ 910c60c7ad4SBarry Smith PCMGGetType - Determines the form of multigrid to use: 911c60c7ad4SBarry Smith multiplicative, additive, full, or the Kaskade algorithm. 912c60c7ad4SBarry Smith 913c60c7ad4SBarry Smith Logically Collective on PC 914c60c7ad4SBarry Smith 915c60c7ad4SBarry Smith Input Parameter: 916c60c7ad4SBarry Smith . pc - the preconditioner context 917c60c7ad4SBarry Smith 918c60c7ad4SBarry Smith Output Parameter: 919c60c7ad4SBarry Smith . type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE 920c60c7ad4SBarry Smith 921c60c7ad4SBarry Smith 922c60c7ad4SBarry Smith Level: advanced 923c60c7ad4SBarry Smith 924c60c7ad4SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid 925c60c7ad4SBarry Smith 926c60c7ad4SBarry Smith .seealso: PCMGSetLevels() 927c60c7ad4SBarry Smith @*/ 928c60c7ad4SBarry Smith PetscErrorCode PCMGGetType(PC pc,PCMGType *type) 929c60c7ad4SBarry Smith { 930c60c7ad4SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 931c60c7ad4SBarry Smith 932c60c7ad4SBarry Smith PetscFunctionBegin; 933c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 934c60c7ad4SBarry Smith *type = mg->am; 9354b9ad928SBarry Smith PetscFunctionReturn(0); 9364b9ad928SBarry Smith } 9374b9ad928SBarry Smith 9384b9ad928SBarry Smith #undef __FUNCT__ 9390d353602SBarry Smith #define __FUNCT__ "PCMGSetCycleType" 9404b9ad928SBarry Smith /*@ 9410d353602SBarry Smith PCMGSetCycleType - Sets the type cycles to use. Use PCMGSetCycleTypeOnLevel() for more 9424b9ad928SBarry Smith complicated cycling. 9434b9ad928SBarry Smith 944ad4df100SBarry Smith Logically Collective on PC 9454b9ad928SBarry Smith 9464b9ad928SBarry Smith Input Parameters: 947c2be2410SBarry Smith + pc - the multigrid context 9480d353602SBarry Smith - PC_MG_CYCLE_V or PC_MG_CYCLE_W 9494b9ad928SBarry Smith 9504b9ad928SBarry Smith Options Database Key: 9514446f3b4SBarry Smith . -pc_mg_cycle_type <v,w> 9524b9ad928SBarry Smith 9534b9ad928SBarry Smith Level: advanced 9544b9ad928SBarry Smith 9554b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 9564b9ad928SBarry Smith 9570d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel() 9584b9ad928SBarry Smith @*/ 9597087cfbeSBarry Smith PetscErrorCode PCMGSetCycleType(PC pc,PCMGCycleType n) 9604b9ad928SBarry Smith { 961f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 962f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 96379416396SBarry Smith PetscInt i,levels; 9644b9ad928SBarry Smith 9654b9ad928SBarry Smith PetscFunctionBegin; 9660700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 967ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 968*69fbec6eSBarry Smith PetscValidLogicalCollectiveEnum(pc,n,2); 969f3fbd535SBarry Smith levels = mglevels[0]->levels; 9704b9ad928SBarry Smith 9712fa5cd67SKarl Rupp for (i=0; i<levels; i++) mglevels[i]->cycles = n; 9724b9ad928SBarry Smith PetscFunctionReturn(0); 9734b9ad928SBarry Smith } 9744b9ad928SBarry Smith 9754b9ad928SBarry Smith #undef __FUNCT__ 9768cc2d5dfSBarry Smith #define __FUNCT__ "PCMGMultiplicativeSetCycles" 9778cc2d5dfSBarry Smith /*@ 9788cc2d5dfSBarry Smith PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step 9798cc2d5dfSBarry Smith of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used 9808cc2d5dfSBarry Smith 981ad4df100SBarry Smith Logically Collective on PC 9828cc2d5dfSBarry Smith 9838cc2d5dfSBarry Smith Input Parameters: 9848cc2d5dfSBarry Smith + pc - the multigrid context 9858cc2d5dfSBarry Smith - n - number of cycles (default is 1) 9868cc2d5dfSBarry Smith 9878cc2d5dfSBarry Smith Options Database Key: 988e1bc860dSBarry Smith . -pc_mg_multiplicative_cycles n 9898cc2d5dfSBarry Smith 9908cc2d5dfSBarry Smith Level: advanced 9918cc2d5dfSBarry Smith 9928cc2d5dfSBarry Smith Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType() 9938cc2d5dfSBarry Smith 9948cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 9958cc2d5dfSBarry Smith 9968cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType() 9978cc2d5dfSBarry Smith @*/ 9987087cfbeSBarry Smith PetscErrorCode PCMGMultiplicativeSetCycles(PC pc,PetscInt n) 9998cc2d5dfSBarry Smith { 1000f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 10018cc2d5dfSBarry Smith 10028cc2d5dfSBarry Smith PetscFunctionBegin; 10030700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1004c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 10052a21e185SBarry Smith mg->cyclesperpcapply = n; 10068cc2d5dfSBarry Smith PetscFunctionReturn(0); 10078cc2d5dfSBarry Smith } 10088cc2d5dfSBarry Smith 10098cc2d5dfSBarry Smith #undef __FUNCT__ 1010fb15c04eSBarry Smith #define __FUNCT__ "PCMGSetGalerkin_MG" 1011fb15c04eSBarry Smith PetscErrorCode PCMGSetGalerkin_MG(PC pc,PetscBool use) 1012fb15c04eSBarry Smith { 1013fb15c04eSBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1014fb15c04eSBarry Smith 1015fb15c04eSBarry Smith PetscFunctionBegin; 1016fb15c04eSBarry Smith mg->galerkin = use ? 1 : 0; 1017fb15c04eSBarry Smith PetscFunctionReturn(0); 1018fb15c04eSBarry Smith } 1019fb15c04eSBarry Smith 1020fb15c04eSBarry Smith #undef __FUNCT__ 10219dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetGalerkin" 1022c2be2410SBarry Smith /*@ 102397177400SBarry Smith PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the 102482b87d87SMatthew G. Knepley finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i 1025c2be2410SBarry Smith 1026ad4df100SBarry Smith Logically Collective on PC 1027c2be2410SBarry Smith 1028c2be2410SBarry Smith Input Parameters: 1029c91913d3SJed Brown + pc - the multigrid context 1030c91913d3SJed Brown - use - PETSC_TRUE to use the Galerkin process to compute coarse-level operators 1031c2be2410SBarry Smith 1032c2be2410SBarry Smith Options Database Key: 10331c1aac46SBarry Smith . -pc_mg_galerkin <true,false> 1034c2be2410SBarry Smith 1035c2be2410SBarry Smith Level: intermediate 1036c2be2410SBarry Smith 10371c1aac46SBarry Smith Notes: Some codes that use PCMG such as PCGAMG use Galerkin internally while constructing the hierarchy and thus do not 10381c1aac46SBarry Smith use the PCMG construction of the coarser grids. 10391c1aac46SBarry Smith 1040c2be2410SBarry Smith .keywords: MG, set, Galerkin 1041c2be2410SBarry Smith 10423fc8bf9cSBarry Smith .seealso: PCMGGetGalerkin() 10433fc8bf9cSBarry Smith 1044c2be2410SBarry Smith @*/ 1045c91913d3SJed Brown PetscErrorCode PCMGSetGalerkin(PC pc,PetscBool use) 1046c2be2410SBarry Smith { 1047fb15c04eSBarry Smith PetscErrorCode ierr; 1048c2be2410SBarry Smith 1049c2be2410SBarry Smith PetscFunctionBegin; 10500700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 10517a4c0024SBarry Smith ierr = PetscTryMethod(pc,"PCMGSetGalerkin_C",(PC,PetscBool),(pc,use));CHKERRQ(ierr); 1052c2be2410SBarry Smith PetscFunctionReturn(0); 1053c2be2410SBarry Smith } 1054c2be2410SBarry Smith 1055c2be2410SBarry Smith #undef __FUNCT__ 10563fc8bf9cSBarry Smith #define __FUNCT__ "PCMGGetGalerkin" 10573fc8bf9cSBarry Smith /*@ 10583fc8bf9cSBarry Smith PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. 105982b87d87SMatthew G. Knepley A_i-1 = r_i * A_i * p_i 10603fc8bf9cSBarry Smith 10613fc8bf9cSBarry Smith Not Collective 10623fc8bf9cSBarry Smith 10633fc8bf9cSBarry Smith Input Parameter: 10643fc8bf9cSBarry Smith . pc - the multigrid context 10653fc8bf9cSBarry Smith 10663fc8bf9cSBarry Smith Output Parameter: 1067e1bc860dSBarry Smith . galerkin - PETSC_TRUE or PETSC_FALSE 10683fc8bf9cSBarry Smith 10693fc8bf9cSBarry Smith Options Database Key: 1070e1bc860dSBarry Smith . -pc_mg_galerkin 10713fc8bf9cSBarry Smith 10723fc8bf9cSBarry Smith Level: intermediate 10733fc8bf9cSBarry Smith 10743fc8bf9cSBarry Smith .keywords: MG, set, Galerkin 10753fc8bf9cSBarry Smith 10763fc8bf9cSBarry Smith .seealso: PCMGSetGalerkin() 10773fc8bf9cSBarry Smith 10783fc8bf9cSBarry Smith @*/ 10797087cfbeSBarry Smith PetscErrorCode PCMGGetGalerkin(PC pc,PetscBool *galerkin) 10803fc8bf9cSBarry Smith { 1081f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 10823fc8bf9cSBarry Smith 10833fc8bf9cSBarry Smith PetscFunctionBegin; 10840700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 108513fdb3acSJose Roman *galerkin = (PetscBool)mg->galerkin; 10863fc8bf9cSBarry Smith PetscFunctionReturn(0); 10873fc8bf9cSBarry Smith } 10883fc8bf9cSBarry Smith 10893fc8bf9cSBarry Smith #undef __FUNCT__ 10909dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothDown" 10914b9ad928SBarry Smith /*@ 109297177400SBarry Smith PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to 109397177400SBarry Smith use on all levels. Use PCMGGetSmootherDown() to set different 10944b9ad928SBarry Smith pre-smoothing steps on different levels. 10954b9ad928SBarry Smith 1096ad4df100SBarry Smith Logically Collective on PC 10974b9ad928SBarry Smith 10984b9ad928SBarry Smith Input Parameters: 10994b9ad928SBarry Smith + mg - the multigrid context 11004b9ad928SBarry Smith - n - the number of smoothing steps 11014b9ad928SBarry Smith 11024b9ad928SBarry Smith Options Database Key: 11034b9ad928SBarry Smith . -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps 11044b9ad928SBarry Smith 11054b9ad928SBarry Smith Level: advanced 11064b9ad928SBarry Smith 11074b9ad928SBarry Smith .keywords: MG, smooth, down, pre-smoothing, steps, multigrid 11084b9ad928SBarry Smith 110997177400SBarry Smith .seealso: PCMGSetNumberSmoothUp() 11104b9ad928SBarry Smith @*/ 11117087cfbeSBarry Smith PetscErrorCode PCMGSetNumberSmoothDown(PC pc,PetscInt n) 11124b9ad928SBarry Smith { 1113f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1114f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 11156849ba73SBarry Smith PetscErrorCode ierr; 111679416396SBarry Smith PetscInt i,levels; 11174b9ad928SBarry Smith 11184b9ad928SBarry Smith PetscFunctionBegin; 11190700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1120ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 1121c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 1122f3fbd535SBarry Smith levels = mglevels[0]->levels; 11234b9ad928SBarry Smith 1124b05257ddSBarry Smith for (i=1; i<levels; i++) { 11254b9ad928SBarry Smith /* make sure smoother up and down are different */ 11260298fd71SBarry Smith ierr = PCMGGetSmootherUp(pc,i,NULL);CHKERRQ(ierr); 1127f3fbd535SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 11282fa5cd67SKarl Rupp 112931567311SBarry Smith mg->default_smoothd = n; 11304b9ad928SBarry Smith } 11314b9ad928SBarry Smith PetscFunctionReturn(0); 11324b9ad928SBarry Smith } 11334b9ad928SBarry Smith 11344b9ad928SBarry Smith #undef __FUNCT__ 11359dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothUp" 11364b9ad928SBarry Smith /*@ 113797177400SBarry Smith PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use 113897177400SBarry Smith on all levels. Use PCMGGetSmootherUp() to set different numbers of 11394b9ad928SBarry Smith post-smoothing steps on different levels. 11404b9ad928SBarry Smith 1141ad4df100SBarry Smith Logically Collective on PC 11424b9ad928SBarry Smith 11434b9ad928SBarry Smith Input Parameters: 11444b9ad928SBarry Smith + mg - the multigrid context 11454b9ad928SBarry Smith - n - the number of smoothing steps 11464b9ad928SBarry Smith 11474b9ad928SBarry Smith Options Database Key: 11484b9ad928SBarry Smith . -pc_mg_smoothup <n> - Sets number of post-smoothing steps 11494b9ad928SBarry Smith 11504b9ad928SBarry Smith Level: advanced 11514b9ad928SBarry Smith 11524b9ad928SBarry Smith Note: this does not set a value on the coarsest grid, since we assume that 1153a8c7a070SBarry Smith there is no separate smooth up on the coarsest grid. 11544b9ad928SBarry Smith 11554b9ad928SBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid 11564b9ad928SBarry Smith 115797177400SBarry Smith .seealso: PCMGSetNumberSmoothDown() 11584b9ad928SBarry Smith @*/ 11597087cfbeSBarry Smith PetscErrorCode PCMGSetNumberSmoothUp(PC pc,PetscInt n) 11604b9ad928SBarry Smith { 1161f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1162f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 11636849ba73SBarry Smith PetscErrorCode ierr; 116479416396SBarry Smith PetscInt i,levels; 11654b9ad928SBarry Smith 11664b9ad928SBarry Smith PetscFunctionBegin; 11670700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1168ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 1169c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 1170f3fbd535SBarry Smith levels = mglevels[0]->levels; 11714b9ad928SBarry Smith 11724b9ad928SBarry Smith for (i=1; i<levels; i++) { 11734b9ad928SBarry Smith /* make sure smoother up and down are different */ 11740298fd71SBarry Smith ierr = PCMGGetSmootherUp(pc,i,NULL);CHKERRQ(ierr); 1175f3fbd535SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 11762fa5cd67SKarl Rupp 117731567311SBarry Smith mg->default_smoothu = n; 11784b9ad928SBarry Smith } 11794b9ad928SBarry Smith PetscFunctionReturn(0); 11804b9ad928SBarry Smith } 11814b9ad928SBarry Smith 11824b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/ 11834b9ad928SBarry Smith 11843b09bd56SBarry Smith /*MC 1185ccb205f8SBarry Smith PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional 11863b09bd56SBarry Smith information about the coarser grid matrices and restriction/interpolation operators. 11873b09bd56SBarry Smith 11883b09bd56SBarry Smith Options Database Keys: 11893b09bd56SBarry Smith + -pc_mg_levels <nlevels> - number of levels including finest 1190b4519db0SPatrick Sanan . -pc_mg_cycles <v,w> - 119179416396SBarry Smith . -pc_mg_smoothup <n> - number of smoothing steps after interpolation 11923b09bd56SBarry Smith . -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator 11938c1c2452SJed Brown . -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default 11943b09bd56SBarry Smith . -pc_mg_log - log information about time spent on each level of the solver 11958cf6037eSBarry Smith . -pc_mg_galerkin - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R' 11968cf6037eSBarry Smith . -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1) 11978cf6037eSBarry Smith . -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 1198e3c5b3baSBarry Smith to the Socket viewer for reading from MATLAB. 11998cf6037eSBarry Smith - -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices 12008cf6037eSBarry Smith to the binary output file called binaryoutput 12013b09bd56SBarry Smith 120224c3aa18SBarry 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 12033b09bd56SBarry Smith 12048cf6037eSBarry Smith When run with a single level the smoother options are used on that level NOT the coarse grid solver options 12058cf6037eSBarry Smith 12063b09bd56SBarry Smith Level: intermediate 12073b09bd56SBarry Smith 12088f87f92bSBarry Smith Concepts: multigrid/multilevel 12093b09bd56SBarry Smith 12108cf6037eSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE 12110d353602SBarry Smith PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(), 121297177400SBarry Smith PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 121397177400SBarry Smith PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 12140d353602SBarry Smith PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 12153b09bd56SBarry Smith M*/ 12163b09bd56SBarry Smith 12174b9ad928SBarry Smith #undef __FUNCT__ 12184b9ad928SBarry Smith #define __FUNCT__ "PCCreate_MG" 12198cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc) 12204b9ad928SBarry Smith { 1221f3fbd535SBarry Smith PC_MG *mg; 1222f3fbd535SBarry Smith PetscErrorCode ierr; 1223f3fbd535SBarry Smith 12244b9ad928SBarry Smith PetscFunctionBegin; 1225b00a9115SJed Brown ierr = PetscNewLog(pc,&mg);CHKERRQ(ierr); 1226f3fbd535SBarry Smith pc->data = (void*)mg; 1227f3fbd535SBarry Smith mg->nlevels = -1; 1228f3fbd535SBarry Smith 122937a44384SMark Adams pc->useAmat = PETSC_TRUE; 123037a44384SMark Adams 12314b9ad928SBarry Smith pc->ops->apply = PCApply_MG; 12324b9ad928SBarry Smith pc->ops->setup = PCSetUp_MG; 1233a06653b4SBarry Smith pc->ops->reset = PCReset_MG; 12344b9ad928SBarry Smith pc->ops->destroy = PCDestroy_MG; 12354b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_MG; 12364b9ad928SBarry Smith pc->ops->view = PCView_MG; 1237fb15c04eSBarry Smith 1238fb15c04eSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",PCMGSetGalerkin_MG);CHKERRQ(ierr); 12394b9ad928SBarry Smith PetscFunctionReturn(0); 12404b9ad928SBarry Smith } 1241