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> 7391689abSStefano Zampini PETSC_INTERN PetscErrorCode PCPreSolveChangeRHS(PC,PetscBool*); 84b9ad928SBarry Smith 931567311SBarry Smith PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG_Levels **mglevelsin,PCRichardsonConvergedReason *reason) 104b9ad928SBarry Smith { 1131567311SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1231567311SBarry Smith PC_MG_Levels *mgc,*mglevels = *mglevelsin; 136849ba73SBarry Smith PetscErrorCode ierr; 1431567311SBarry Smith PetscInt cycles = (mglevels->level == 1) ? 1 : (PetscInt) mglevels->cycles; 15a0808db4SHong Zhang PC subpc; 16a0808db4SHong Zhang PCFailedReason pcreason; 174b9ad928SBarry Smith 184b9ad928SBarry Smith PetscFunctionBegin; 1963e6d426SJed Brown if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 2031567311SBarry Smith ierr = KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x);CHKERRQ(ierr); /* pre-smooth */ 21a0808db4SHong Zhang ierr = KSPGetPC(mglevels->smoothd,&subpc);CHKERRQ(ierr); 22a0808db4SHong Zhang ierr = PCGetSetUpFailedReason(subpc,&pcreason);CHKERRQ(ierr); 23a0808db4SHong Zhang if (pcreason) { 24a0808db4SHong Zhang pc->failedreason = PC_SUBPC_ERROR; 25a0808db4SHong Zhang } 2663e6d426SJed Brown if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 2731567311SBarry Smith if (mglevels->level) { /* not the coarsest grid */ 2863e6d426SJed Brown if (mglevels->eventresidual) {ierr = PetscLogEventBegin(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);} 2931567311SBarry Smith ierr = (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);CHKERRQ(ierr); 3063e6d426SJed Brown if (mglevels->eventresidual) {ierr = PetscLogEventEnd(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);} 314b9ad928SBarry Smith 324b9ad928SBarry Smith /* if on finest level and have convergence criteria set */ 3331567311SBarry Smith if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) { 344b9ad928SBarry Smith PetscReal rnorm; 3531567311SBarry Smith ierr = VecNorm(mglevels->r,NORM_2,&rnorm);CHKERRQ(ierr); 364b9ad928SBarry Smith if (rnorm <= mg->ttol) { 3770441072SBarry Smith if (rnorm < mg->abstol) { 384d0a8057SBarry Smith *reason = PCRICHARDSON_CONVERGED_ATOL; 3957622a8eSBarry 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); 404b9ad928SBarry Smith } else { 414d0a8057SBarry Smith *reason = PCRICHARDSON_CONVERGED_RTOL; 4257622a8eSBarry 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); 434b9ad928SBarry Smith } 444b9ad928SBarry Smith PetscFunctionReturn(0); 454b9ad928SBarry Smith } 464b9ad928SBarry Smith } 474b9ad928SBarry Smith 4831567311SBarry Smith mgc = *(mglevelsin - 1); 4963e6d426SJed Brown if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 5031567311SBarry Smith ierr = MatRestrict(mglevels->restrct,mglevels->r,mgc->b);CHKERRQ(ierr); 5163e6d426SJed Brown if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 52efb30889SBarry Smith ierr = VecSet(mgc->x,0.0);CHKERRQ(ierr); 534b9ad928SBarry Smith while (cycles--) { 5431567311SBarry Smith ierr = PCMGMCycle_Private(pc,mglevelsin-1,reason);CHKERRQ(ierr); 554b9ad928SBarry Smith } 5663e6d426SJed Brown if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 5731567311SBarry Smith ierr = MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x);CHKERRQ(ierr); 5863e6d426SJed Brown if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 5963e6d426SJed Brown if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 6031567311SBarry Smith ierr = KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x);CHKERRQ(ierr); /* post smooth */ 6163e6d426SJed Brown if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 624b9ad928SBarry Smith } 634b9ad928SBarry Smith PetscFunctionReturn(0); 644b9ad928SBarry Smith } 654b9ad928SBarry Smith 66ace3abfcSBarry 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) 674b9ad928SBarry Smith { 68f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 69f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 70dfbe8321SBarry Smith PetscErrorCode ierr; 71391689abSStefano Zampini PC tpc; 72391689abSStefano Zampini PetscBool changeu,changed; 73f3fbd535SBarry Smith PetscInt levels = mglevels[0]->levels,i; 744b9ad928SBarry Smith 754b9ad928SBarry Smith PetscFunctionBegin; 76a762d673SBarry Smith /* When the DM is supplying the matrix then it will not exist until here */ 77a762d673SBarry Smith for (i=0; i<levels; i++) { 78a762d673SBarry Smith if (!mglevels[i]->A) { 79a762d673SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr); 80a762d673SBarry Smith ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr); 81a762d673SBarry Smith } 82a762d673SBarry Smith } 83391689abSStefano Zampini 84391689abSStefano Zampini ierr = KSPGetPC(mglevels[levels-1]->smoothd,&tpc);CHKERRQ(ierr); 85391689abSStefano Zampini ierr = PCPreSolveChangeRHS(tpc,&changed);CHKERRQ(ierr); 86391689abSStefano Zampini ierr = KSPGetPC(mglevels[levels-1]->smoothu,&tpc);CHKERRQ(ierr); 87391689abSStefano Zampini ierr = PCPreSolveChangeRHS(tpc,&changeu);CHKERRQ(ierr); 88391689abSStefano Zampini if (!changed && !changeu) { 89391689abSStefano Zampini ierr = VecDestroy(&mglevels[levels-1]->b);CHKERRQ(ierr); 90f3fbd535SBarry Smith mglevels[levels-1]->b = b; 91391689abSStefano Zampini } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */ 92391689abSStefano Zampini if (!mglevels[levels-1]->b) { 93391689abSStefano Zampini Vec *vec; 94391689abSStefano Zampini 95391689abSStefano Zampini ierr = KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr); 96391689abSStefano Zampini mglevels[levels-1]->b = *vec; 97*7ae21d43SStefano Zampini ierr = PetscFree(vec);CHKERRQ(ierr); 98391689abSStefano Zampini } 99391689abSStefano Zampini ierr = VecCopy(b,mglevels[levels-1]->b);CHKERRQ(ierr); 100391689abSStefano Zampini } 101f3fbd535SBarry Smith mglevels[levels-1]->x = x; 1024b9ad928SBarry Smith 10331567311SBarry Smith mg->rtol = rtol; 10431567311SBarry Smith mg->abstol = abstol; 10531567311SBarry Smith mg->dtol = dtol; 1064b9ad928SBarry Smith if (rtol) { 1074b9ad928SBarry Smith /* compute initial residual norm for relative convergence test */ 1084b9ad928SBarry Smith PetscReal rnorm; 1097319c654SBarry Smith if (zeroguess) { 1107319c654SBarry Smith ierr = VecNorm(b,NORM_2,&rnorm);CHKERRQ(ierr); 1117319c654SBarry Smith } else { 112f3fbd535SBarry Smith ierr = (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);CHKERRQ(ierr); 1134b9ad928SBarry Smith ierr = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr); 1147319c654SBarry Smith } 11531567311SBarry Smith mg->ttol = PetscMax(rtol*rnorm,abstol); 1162fa5cd67SKarl Rupp } else if (abstol) mg->ttol = abstol; 1172fa5cd67SKarl Rupp else mg->ttol = 0.0; 1184b9ad928SBarry Smith 1194d0a8057SBarry Smith /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't 120336babb1SJed Brown stop prematurely due to small residual */ 1214d0a8057SBarry Smith for (i=1; i<levels; i++) { 122f3fbd535SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 123f3fbd535SBarry Smith if (mglevels[i]->smoothu != mglevels[i]->smoothd) { 12423067569SBarry Smith /* For Richardson the initial guess is nonzero since it is solving in each cycle the original system not just applying as a preconditioner */ 12523067569SBarry Smith ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 126f3fbd535SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 1274b9ad928SBarry Smith } 1284d0a8057SBarry Smith } 1294d0a8057SBarry Smith 1304d0a8057SBarry Smith *reason = (PCRichardsonConvergedReason)0; 1314d0a8057SBarry Smith for (i=0; i<its; i++) { 132f3fbd535SBarry Smith ierr = PCMGMCycle_Private(pc,mglevels+levels-1,reason);CHKERRQ(ierr); 1334d0a8057SBarry Smith if (*reason) break; 1344d0a8057SBarry Smith } 1354d0a8057SBarry Smith if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS; 1364d0a8057SBarry Smith *outits = i; 137391689abSStefano Zampini if (!changed && !changeu) mglevels[levels-1]->b = NULL; 1384b9ad928SBarry Smith PetscFunctionReturn(0); 1394b9ad928SBarry Smith } 1404b9ad928SBarry Smith 1413aeaf226SBarry Smith PetscErrorCode PCReset_MG(PC pc) 1423aeaf226SBarry Smith { 1433aeaf226SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1443aeaf226SBarry Smith PC_MG_Levels **mglevels = mg->levels; 1453aeaf226SBarry Smith PetscErrorCode ierr; 1463aeaf226SBarry Smith PetscInt i,n; 1473aeaf226SBarry Smith 1483aeaf226SBarry Smith PetscFunctionBegin; 1493aeaf226SBarry Smith if (mglevels) { 1503aeaf226SBarry Smith n = mglevels[0]->levels; 1513aeaf226SBarry Smith for (i=0; i<n-1; i++) { 1523aeaf226SBarry Smith ierr = VecDestroy(&mglevels[i+1]->r);CHKERRQ(ierr); 1533aeaf226SBarry Smith ierr = VecDestroy(&mglevels[i]->b);CHKERRQ(ierr); 1543aeaf226SBarry Smith ierr = VecDestroy(&mglevels[i]->x);CHKERRQ(ierr); 1553aeaf226SBarry Smith ierr = MatDestroy(&mglevels[i+1]->restrct);CHKERRQ(ierr); 1563aeaf226SBarry Smith ierr = MatDestroy(&mglevels[i+1]->interpolate);CHKERRQ(ierr); 1570de8493bSLawrence Mitchell ierr = MatDestroy(&mglevels[i+1]->inject);CHKERRQ(ierr); 15873250ac0SBarry Smith ierr = VecDestroy(&mglevels[i+1]->rscale);CHKERRQ(ierr); 1593aeaf226SBarry Smith } 160391689abSStefano Zampini /* this is not null only if the smoother on the finest level 161391689abSStefano Zampini changes the rhs during PreSolve */ 162391689abSStefano Zampini ierr = VecDestroy(&mglevels[n-1]->b);CHKERRQ(ierr); 1633aeaf226SBarry Smith 1643aeaf226SBarry Smith for (i=0; i<n; i++) { 1653aeaf226SBarry Smith ierr = MatDestroy(&mglevels[i]->A);CHKERRQ(ierr); 1663aeaf226SBarry Smith if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 1673aeaf226SBarry Smith ierr = KSPReset(mglevels[i]->smoothd);CHKERRQ(ierr); 1683aeaf226SBarry Smith } 1693aeaf226SBarry Smith ierr = KSPReset(mglevels[i]->smoothu);CHKERRQ(ierr); 1703aeaf226SBarry Smith } 1713aeaf226SBarry Smith } 1723aeaf226SBarry Smith PetscFunctionReturn(0); 1733aeaf226SBarry Smith } 1743aeaf226SBarry Smith 1754b9ad928SBarry Smith /*@C 17697177400SBarry Smith PCMGSetLevels - Sets the number of levels to use with MG. 1774b9ad928SBarry Smith Must be called before any other MG routine. 1784b9ad928SBarry Smith 179ad4df100SBarry Smith Logically Collective on PC 1804b9ad928SBarry Smith 1814b9ad928SBarry Smith Input Parameters: 1824b9ad928SBarry Smith + pc - the preconditioner context 1834b9ad928SBarry Smith . levels - the number of levels 1844b9ad928SBarry Smith - comms - optional communicators for each level; this is to allow solving the coarser problems 1852bf68e3eSBarry Smith on smaller sets of processors. 1864b9ad928SBarry Smith 1874b9ad928SBarry Smith Level: intermediate 1884b9ad928SBarry Smith 1894b9ad928SBarry Smith Notes: 1904b9ad928SBarry Smith If the number of levels is one then the multigrid uses the -mg_levels prefix 1914b9ad928SBarry Smith for setting the level options rather than the -mg_coarse prefix. 1924b9ad928SBarry Smith 1934b9ad928SBarry Smith .keywords: MG, set, levels, multigrid 1944b9ad928SBarry Smith 19597177400SBarry Smith .seealso: PCMGSetType(), PCMGGetLevels() 1964b9ad928SBarry Smith @*/ 1977087cfbeSBarry Smith PetscErrorCode PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms) 1984b9ad928SBarry Smith { 199dfbe8321SBarry Smith PetscErrorCode ierr; 200f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 201ce94432eSBarry Smith MPI_Comm comm; 2023aeaf226SBarry Smith PC_MG_Levels **mglevels = mg->levels; 20310eca3edSLisandro Dalcin PCMGType mgtype = mg->am; 20410167fecSBarry Smith PetscInt mgctype = (PetscInt) PC_MG_CYCLE_V; 205f3fbd535SBarry Smith PetscInt i; 206f3fbd535SBarry Smith PetscMPIInt size; 207f3fbd535SBarry Smith const char *prefix; 208f3fbd535SBarry Smith PC ipc; 2093aeaf226SBarry Smith PetscInt n; 2104b9ad928SBarry Smith 2114b9ad928SBarry Smith PetscFunctionBegin; 2120700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 213548767e0SJed Brown PetscValidLogicalCollectiveInt(pc,levels,2); 214548767e0SJed Brown if (mg->nlevels == levels) PetscFunctionReturn(0); 2151a97d7b6SStefano Zampini ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 2163aeaf226SBarry Smith if (mglevels) { 21710eca3edSLisandro Dalcin mgctype = mglevels[0]->cycles; 2183aeaf226SBarry Smith /* changing the number of levels so free up the previous stuff */ 2193aeaf226SBarry Smith ierr = PCReset_MG(pc);CHKERRQ(ierr); 2203aeaf226SBarry Smith n = mglevels[0]->levels; 2213aeaf226SBarry Smith for (i=0; i<n; i++) { 2223aeaf226SBarry Smith if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 2233aeaf226SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr); 2243aeaf226SBarry Smith } 2253aeaf226SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr); 2263aeaf226SBarry Smith ierr = PetscFree(mglevels[i]);CHKERRQ(ierr); 2273aeaf226SBarry Smith } 2283aeaf226SBarry Smith ierr = PetscFree(mg->levels);CHKERRQ(ierr); 2293aeaf226SBarry Smith } 230f3fbd535SBarry Smith 231f3fbd535SBarry Smith mg->nlevels = levels; 232f3fbd535SBarry Smith 233785e854fSJed Brown ierr = PetscMalloc1(levels,&mglevels);CHKERRQ(ierr); 2343bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr); 235f3fbd535SBarry Smith 236f3fbd535SBarry Smith ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 237f3fbd535SBarry Smith 238a9db87a2SMatthew G Knepley mg->stageApply = 0; 239f3fbd535SBarry Smith for (i=0; i<levels; i++) { 240b00a9115SJed Brown ierr = PetscNewLog(pc,&mglevels[i]);CHKERRQ(ierr); 2412fa5cd67SKarl Rupp 242f3fbd535SBarry Smith mglevels[i]->level = i; 243f3fbd535SBarry Smith mglevels[i]->levels = levels; 24410eca3edSLisandro Dalcin mglevels[i]->cycles = mgctype; 245336babb1SJed Brown mg->default_smoothu = 2; 246336babb1SJed Brown mg->default_smoothd = 2; 24763e6d426SJed Brown mglevels[i]->eventsmoothsetup = 0; 24863e6d426SJed Brown mglevels[i]->eventsmoothsolve = 0; 24963e6d426SJed Brown mglevels[i]->eventresidual = 0; 25063e6d426SJed Brown mglevels[i]->eventinterprestrict = 0; 251f3fbd535SBarry Smith 252f3fbd535SBarry Smith if (comms) comm = comms[i]; 253f3fbd535SBarry Smith ierr = KSPCreate(comm,&mglevels[i]->smoothd);CHKERRQ(ierr); 254422a814eSBarry Smith ierr = KSPSetErrorIfNotConverged(mglevels[i]->smoothd,pc->erroriffailure);CHKERRQ(ierr); 2550ee9a668SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr); 2560ee9a668SBarry Smith ierr = KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);CHKERRQ(ierr); 2570ee9a668SBarry Smith ierr = PetscObjectComposedDataSetInt((PetscObject) mglevels[i]->smoothd, PetscMGLevelId, mglevels[i]->level);CHKERRQ(ierr); 2580ee9a668SBarry Smith if (i || levels == 1) { 2590ee9a668SBarry Smith char tprefix[128]; 2600ee9a668SBarry Smith 261336babb1SJed Brown ierr = KSPSetType(mglevels[i]->smoothd,KSPCHEBYSHEV);CHKERRQ(ierr); 2620059c7bdSJed Brown ierr = KSPSetConvergenceTest(mglevels[i]->smoothd,KSPConvergedSkip,NULL,NULL);CHKERRQ(ierr); 263336babb1SJed Brown ierr = KSPSetNormType(mglevels[i]->smoothd,KSP_NORM_NONE);CHKERRQ(ierr); 264336babb1SJed Brown ierr = KSPGetPC(mglevels[i]->smoothd,&ipc);CHKERRQ(ierr); 265336babb1SJed Brown ierr = PCSetType(ipc,PCSOR);CHKERRQ(ierr); 2660ee9a668SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg->default_smoothd);CHKERRQ(ierr); 267f3fbd535SBarry Smith 2680ee9a668SBarry Smith sprintf(tprefix,"mg_levels_%d_",(int)i); 2690ee9a668SBarry Smith ierr = KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);CHKERRQ(ierr); 2700ee9a668SBarry Smith } else { 271f3fbd535SBarry Smith ierr = KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");CHKERRQ(ierr); 272f3fbd535SBarry Smith 2739dbfc187SHong Zhang /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */ 274f3fbd535SBarry Smith ierr = KSPSetType(mglevels[0]->smoothd,KSPPREONLY);CHKERRQ(ierr); 275f3fbd535SBarry Smith ierr = KSPGetPC(mglevels[0]->smoothd,&ipc);CHKERRQ(ierr); 276f3fbd535SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 277f3fbd535SBarry Smith if (size > 1) { 278f3fbd535SBarry Smith ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr); 279f3fbd535SBarry Smith } else { 280f3fbd535SBarry Smith ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr); 281f3fbd535SBarry Smith } 282753b7fb9SBarry Smith ierr = PCFactorSetShiftType(ipc,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr); 283f3fbd535SBarry Smith } 2843bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[i]->smoothd);CHKERRQ(ierr); 2852fa5cd67SKarl Rupp 286f3fbd535SBarry Smith mglevels[i]->smoothu = mglevels[i]->smoothd; 28731567311SBarry Smith mg->rtol = 0.0; 28831567311SBarry Smith mg->abstol = 0.0; 28931567311SBarry Smith mg->dtol = 0.0; 29031567311SBarry Smith mg->ttol = 0.0; 29131567311SBarry Smith mg->cyclesperpcapply = 1; 292f3fbd535SBarry Smith } 293f3fbd535SBarry Smith mg->levels = mglevels; 29410eca3edSLisandro Dalcin ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr); 2954b9ad928SBarry Smith PetscFunctionReturn(0); 2964b9ad928SBarry Smith } 2974b9ad928SBarry Smith 298c07bf074SBarry Smith 299c07bf074SBarry Smith PetscErrorCode PCDestroy_MG(PC pc) 300c07bf074SBarry Smith { 301c07bf074SBarry Smith PetscErrorCode ierr; 302a06653b4SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 303a06653b4SBarry Smith PC_MG_Levels **mglevels = mg->levels; 304a06653b4SBarry Smith PetscInt i,n; 305c07bf074SBarry Smith 306c07bf074SBarry Smith PetscFunctionBegin; 307a06653b4SBarry Smith ierr = PCReset_MG(pc);CHKERRQ(ierr); 308a06653b4SBarry Smith if (mglevels) { 309a06653b4SBarry Smith n = mglevels[0]->levels; 310a06653b4SBarry Smith for (i=0; i<n; i++) { 311a06653b4SBarry Smith if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 3126bf464f9SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr); 313a06653b4SBarry Smith } 3146bf464f9SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr); 315a06653b4SBarry Smith ierr = PetscFree(mglevels[i]);CHKERRQ(ierr); 316a06653b4SBarry Smith } 317a06653b4SBarry Smith ierr = PetscFree(mg->levels);CHKERRQ(ierr); 318a06653b4SBarry Smith } 319c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 320f3fbd535SBarry Smith PetscFunctionReturn(0); 321f3fbd535SBarry Smith } 322f3fbd535SBarry Smith 323f3fbd535SBarry Smith 324f3fbd535SBarry Smith 32509573ac7SBarry Smith extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**); 32609573ac7SBarry Smith extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**); 32709573ac7SBarry Smith extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**); 328f3fbd535SBarry Smith 329f3fbd535SBarry Smith /* 330f3fbd535SBarry Smith PCApply_MG - Runs either an additive, multiplicative, Kaskadic 331f3fbd535SBarry Smith or full cycle of multigrid. 332f3fbd535SBarry Smith 333f3fbd535SBarry Smith Note: 334f3fbd535SBarry Smith A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle(). 335f3fbd535SBarry Smith */ 336f3fbd535SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x) 337f3fbd535SBarry Smith { 338f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 339f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 340f3fbd535SBarry Smith PetscErrorCode ierr; 341391689abSStefano Zampini PC tpc; 342f3fbd535SBarry Smith PetscInt levels = mglevels[0]->levels,i; 343391689abSStefano Zampini PetscBool changeu,changed; 344f3fbd535SBarry Smith 345f3fbd535SBarry Smith PetscFunctionBegin; 346a9db87a2SMatthew G Knepley if (mg->stageApply) {ierr = PetscLogStagePush(mg->stageApply);CHKERRQ(ierr);} 347e1d8e5deSBarry Smith /* When the DM is supplying the matrix then it will not exist until here */ 348298cc208SBarry Smith for (i=0; i<levels; i++) { 349e1d8e5deSBarry Smith if (!mglevels[i]->A) { 35023ee1639SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr); 351298cc208SBarry Smith ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr); 352e1d8e5deSBarry Smith } 353e1d8e5deSBarry Smith } 354e1d8e5deSBarry Smith 355391689abSStefano Zampini ierr = KSPGetPC(mglevels[levels-1]->smoothd,&tpc);CHKERRQ(ierr); 356391689abSStefano Zampini ierr = PCPreSolveChangeRHS(tpc,&changed);CHKERRQ(ierr); 357391689abSStefano Zampini ierr = KSPGetPC(mglevels[levels-1]->smoothu,&tpc);CHKERRQ(ierr); 358391689abSStefano Zampini ierr = PCPreSolveChangeRHS(tpc,&changeu);CHKERRQ(ierr); 359391689abSStefano Zampini if (!changeu && !changed) { 360391689abSStefano Zampini ierr = VecDestroy(&mglevels[levels-1]->b);CHKERRQ(ierr); 361f3fbd535SBarry Smith mglevels[levels-1]->b = b; 362391689abSStefano Zampini } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */ 363391689abSStefano Zampini if (!mglevels[levels-1]->b) { 364391689abSStefano Zampini Vec *vec; 365391689abSStefano Zampini 366391689abSStefano Zampini ierr = KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr); 367391689abSStefano Zampini mglevels[levels-1]->b = *vec; 368*7ae21d43SStefano Zampini ierr = PetscFree(vec);CHKERRQ(ierr); 369391689abSStefano Zampini } 370391689abSStefano Zampini ierr = VecCopy(b,mglevels[levels-1]->b);CHKERRQ(ierr); 371391689abSStefano Zampini } 372f3fbd535SBarry Smith mglevels[levels-1]->x = x; 373391689abSStefano Zampini 37431567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 375f3fbd535SBarry Smith ierr = VecSet(x,0.0);CHKERRQ(ierr); 37631567311SBarry Smith for (i=0; i<mg->cyclesperpcapply; i++) { 3770298fd71SBarry Smith ierr = PCMGMCycle_Private(pc,mglevels+levels-1,NULL);CHKERRQ(ierr); 378f3fbd535SBarry Smith } 3792fa5cd67SKarl Rupp } else if (mg->am == PC_MG_ADDITIVE) { 38031567311SBarry Smith ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr); 3812fa5cd67SKarl Rupp } else if (mg->am == PC_MG_KASKADE) { 38231567311SBarry Smith ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr); 3832fa5cd67SKarl Rupp } else { 384f3fbd535SBarry Smith ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr); 385f3fbd535SBarry Smith } 386a9db87a2SMatthew G Knepley if (mg->stageApply) {ierr = PetscLogStagePop();CHKERRQ(ierr);} 387391689abSStefano Zampini if (!changeu && !changed) mglevels[levels-1]->b = NULL; 388f3fbd535SBarry Smith PetscFunctionReturn(0); 389f3fbd535SBarry Smith } 390f3fbd535SBarry Smith 391f3fbd535SBarry Smith 3924416b707SBarry Smith PetscErrorCode PCSetFromOptions_MG(PetscOptionItems *PetscOptionsObject,PC pc) 393f3fbd535SBarry Smith { 394f3fbd535SBarry Smith PetscErrorCode ierr; 395710315b6SLawrence Mitchell PetscInt levels,cycles; 3962134b1e4SBarry Smith PetscBool flg; 397f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 3983d3eaba7SBarry Smith PC_MG_Levels **mglevels; 399f3fbd535SBarry Smith PCMGType mgtype; 400f3fbd535SBarry Smith PCMGCycleType mgctype; 4012134b1e4SBarry Smith PCMGGalerkinType gtype; 402f3fbd535SBarry Smith 403f3fbd535SBarry Smith PetscFunctionBegin; 4041d743356SStefano Zampini levels = PetscMax(mg->nlevels,1); 405e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Multigrid options");CHKERRQ(ierr); 406f3fbd535SBarry Smith ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr); 4071a97d7b6SStefano Zampini if (!flg && !mg->levels && pc->dm) { 408eb3f98d2SBarry Smith ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr); 409eb3f98d2SBarry Smith levels++; 4103aeaf226SBarry Smith mg->usedmfornumberoflevels = PETSC_TRUE; 411eb3f98d2SBarry Smith } 4120298fd71SBarry Smith ierr = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr); 4133aeaf226SBarry Smith mglevels = mg->levels; 4143aeaf226SBarry Smith 415f3fbd535SBarry Smith mgctype = (PCMGCycleType) mglevels[0]->cycles; 416f3fbd535SBarry Smith ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr); 417f3fbd535SBarry Smith if (flg) { 418f3fbd535SBarry Smith ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr); 4192fa5cd67SKarl Rupp } 4202134b1e4SBarry Smith gtype = mg->galerkin; 4212134b1e4SBarry Smith ierr = PetscOptionsEnum("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",PCMGGalerkinTypes,(PetscEnum)gtype,(PetscEnum*)>ype,&flg);CHKERRQ(ierr); 4222134b1e4SBarry Smith if (flg) { 4232134b1e4SBarry Smith ierr = PCMGSetGalerkin(pc,gtype);CHKERRQ(ierr); 424f3fbd535SBarry Smith } 425f56b1265SBarry Smith flg = PETSC_FALSE; 426f442ab6aSBarry Smith ierr = PetscOptionsBool("-pc_mg_distinct_smoothup","Create seperate smoothup KSP and append the prefix _up","PCMGSetDistinctSmoothUp",PETSC_FALSE,&flg,NULL);CHKERRQ(ierr); 427f442ab6aSBarry Smith if (flg) { 428f442ab6aSBarry Smith ierr = PCMGSetDistinctSmoothUp(pc);CHKERRQ(ierr); 429f442ab6aSBarry Smith } 43031567311SBarry Smith mgtype = mg->am; 431f3fbd535SBarry Smith ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr); 432f3fbd535SBarry Smith if (flg) { 433f3fbd535SBarry Smith ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr); 434f3fbd535SBarry Smith } 43531567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 4365363c1e0SLisandro Dalcin ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGMultiplicativeSetCycles",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr); 437f3fbd535SBarry Smith if (flg) { 438f3fbd535SBarry Smith ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr); 439f3fbd535SBarry Smith } 440f3fbd535SBarry Smith } 441f3fbd535SBarry Smith flg = PETSC_FALSE; 4420298fd71SBarry Smith ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,NULL);CHKERRQ(ierr); 443f3fbd535SBarry Smith if (flg) { 444f3fbd535SBarry Smith PetscInt i; 445f3fbd535SBarry Smith char eventname[128]; 4461a97d7b6SStefano Zampini 447f3fbd535SBarry Smith levels = mglevels[0]->levels; 448f3fbd535SBarry Smith for (i=0; i<levels; i++) { 449f3fbd535SBarry Smith sprintf(eventname,"MGSetup Level %d",(int)i); 45063e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr); 451f3fbd535SBarry Smith sprintf(eventname,"MGSmooth Level %d",(int)i); 45263e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr); 453f3fbd535SBarry Smith if (i) { 454f3fbd535SBarry Smith sprintf(eventname,"MGResid Level %d",(int)i); 45563e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr); 456f3fbd535SBarry Smith sprintf(eventname,"MGInterp Level %d",(int)i); 45763e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr); 458f3fbd535SBarry Smith } 459f3fbd535SBarry Smith } 460eec35531SMatthew G Knepley 461ec60ca73SSatish Balay #if defined(PETSC_USE_LOG) 462eec35531SMatthew G Knepley { 463eec35531SMatthew G Knepley const char *sname = "MG Apply"; 464eec35531SMatthew G Knepley PetscStageLog stageLog; 465eec35531SMatthew G Knepley PetscInt st; 466eec35531SMatthew G Knepley 467eec35531SMatthew G Knepley ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr); 468eec35531SMatthew G Knepley for (st = 0; st < stageLog->numStages; ++st) { 469eec35531SMatthew G Knepley PetscBool same; 470eec35531SMatthew G Knepley 471eec35531SMatthew G Knepley ierr = PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);CHKERRQ(ierr); 4722fa5cd67SKarl Rupp if (same) mg->stageApply = st; 473eec35531SMatthew G Knepley } 474eec35531SMatthew G Knepley if (!mg->stageApply) { 475eec35531SMatthew G Knepley ierr = PetscLogStageRegister(sname, &mg->stageApply);CHKERRQ(ierr); 476eec35531SMatthew G Knepley } 477eec35531SMatthew G Knepley } 478ec60ca73SSatish Balay #endif 479f3fbd535SBarry Smith } 480f3fbd535SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 481f3fbd535SBarry Smith PetscFunctionReturn(0); 482f3fbd535SBarry Smith } 483f3fbd535SBarry Smith 4846a6fc655SJed Brown const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0}; 4856a6fc655SJed Brown const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0}; 4862134b1e4SBarry Smith const char *const PCMGGalerkinTypes[] = {"both","pmat","mat","none","external","PCMGGalerkinType","PC_MG_GALERKIN",0}; 487f3fbd535SBarry Smith 4889804daf3SBarry Smith #include <petscdraw.h> 489f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer) 490f3fbd535SBarry Smith { 491f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 492f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 493f3fbd535SBarry Smith PetscErrorCode ierr; 494e3deeeafSJed Brown PetscInt levels = mglevels ? mglevels[0]->levels : 0,i; 495219b1045SBarry Smith PetscBool iascii,isbinary,isdraw; 496f3fbd535SBarry Smith 497f3fbd535SBarry Smith PetscFunctionBegin; 498251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 4995b0b0462SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 500219b1045SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 501f3fbd535SBarry Smith if (iascii) { 502e3deeeafSJed Brown const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown"; 503efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,cyclename);CHKERRQ(ierr); 50431567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 50531567311SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr); 506f3fbd535SBarry Smith } 5072134b1e4SBarry Smith if (mg->galerkin == PC_MG_GALERKIN_BOTH) { 508f3fbd535SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 5092134b1e4SBarry Smith } else if (mg->galerkin == PC_MG_GALERKIN_PMAT) { 5102134b1e4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices for pmat\n");CHKERRQ(ierr); 5112134b1e4SBarry Smith } else if (mg->galerkin == PC_MG_GALERKIN_MAT) { 5122134b1e4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices for mat\n");CHKERRQ(ierr); 5132134b1e4SBarry Smith } else if (mg->galerkin == PC_MG_GALERKIN_EXTERNAL) { 5142134b1e4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using externally compute Galerkin coarse grid matrices\n");CHKERRQ(ierr); 5154f66f45eSBarry Smith } else { 5164f66f45eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 517f3fbd535SBarry Smith } 5185adeb434SBarry Smith if (mg->view){ 5195adeb434SBarry Smith ierr = (*mg->view)(pc,viewer);CHKERRQ(ierr); 5205adeb434SBarry Smith } 521f3fbd535SBarry Smith for (i=0; i<levels; i++) { 522f3fbd535SBarry Smith if (!i) { 52389cce641SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);CHKERRQ(ierr); 524f3fbd535SBarry Smith } else { 52589cce641SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 526f3fbd535SBarry Smith } 527f3fbd535SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 528f3fbd535SBarry Smith ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 529f3fbd535SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 530f3fbd535SBarry Smith if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) { 531f3fbd535SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 532f3fbd535SBarry Smith } else if (i) { 53389cce641SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 534f3fbd535SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 535f3fbd535SBarry Smith ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 536f3fbd535SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 537f3fbd535SBarry Smith } 538f3fbd535SBarry Smith } 5395b0b0462SBarry Smith } else if (isbinary) { 5405b0b0462SBarry Smith for (i=levels-1; i>=0; i--) { 5415b0b0462SBarry Smith ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 5425b0b0462SBarry Smith if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) { 5435b0b0462SBarry Smith ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 5445b0b0462SBarry Smith } 5455b0b0462SBarry Smith } 546219b1045SBarry Smith } else if (isdraw) { 547219b1045SBarry Smith PetscDraw draw; 548b4375e8dSPeter Brune PetscReal x,w,y,bottom,th; 549219b1045SBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 550219b1045SBarry Smith ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr); 5510298fd71SBarry Smith ierr = PetscDrawStringGetSize(draw,NULL,&th);CHKERRQ(ierr); 552219b1045SBarry Smith bottom = y - th; 553219b1045SBarry Smith for (i=levels-1; i>=0; i--) { 554b4375e8dSPeter Brune if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) { 555219b1045SBarry Smith ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr); 556219b1045SBarry Smith ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 557219b1045SBarry Smith ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 558b4375e8dSPeter Brune } else { 559b4375e8dSPeter Brune w = 0.5*PetscMin(1.0-x,x); 560b4375e8dSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr); 561b4375e8dSPeter Brune ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 562b4375e8dSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 563b4375e8dSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr); 564b4375e8dSPeter Brune ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 565b4375e8dSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 566b4375e8dSPeter Brune } 5670298fd71SBarry Smith ierr = PetscDrawGetBoundingBox(draw,NULL,&bottom,NULL,NULL);CHKERRQ(ierr); 5681cd381d6SBarry Smith bottom -= th; 569219b1045SBarry Smith } 5705b0b0462SBarry Smith } 571f3fbd535SBarry Smith PetscFunctionReturn(0); 572f3fbd535SBarry Smith } 573f3fbd535SBarry Smith 574af0996ceSBarry Smith #include <petsc/private/dmimpl.h> 575af0996ceSBarry Smith #include <petsc/private/kspimpl.h> 576cab2e9ccSBarry Smith 577f3fbd535SBarry Smith /* 578f3fbd535SBarry Smith Calls setup for the KSP on each level 579f3fbd535SBarry Smith */ 580f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc) 581f3fbd535SBarry Smith { 582f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 583f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 584f3fbd535SBarry Smith PetscErrorCode ierr; 5851a97d7b6SStefano Zampini PetscInt i,n; 58698e04984SBarry Smith PC cpc; 5873db492dfSBarry Smith PetscBool dump = PETSC_FALSE,opsset,use_amat,missinginterpolate = PETSC_FALSE; 588f3fbd535SBarry Smith Mat dA,dB; 589f3fbd535SBarry Smith Vec tvec; 590218a07d4SBarry Smith DM *dms; 591649052a6SBarry Smith PetscViewer viewer = 0; 5922134b1e4SBarry Smith PetscBool dAeqdB = PETSC_FALSE, needRestricts = PETSC_FALSE; 593f3fbd535SBarry Smith 594f3fbd535SBarry Smith PetscFunctionBegin; 5951a97d7b6SStefano Zampini if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels with PCMGSetLevels() before setting up"); 5961a97d7b6SStefano Zampini n = mglevels[0]->levels; 59701bc414fSDmitry Karpeev /* FIX: Move this to PCSetFromOptions_MG? */ 5983aeaf226SBarry Smith if (mg->usedmfornumberoflevels) { 5993aeaf226SBarry Smith PetscInt levels; 6003aeaf226SBarry Smith ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr); 6013aeaf226SBarry Smith levels++; 6023aeaf226SBarry Smith if (levels > n) { /* the problem is now being solved on a finer grid */ 6030298fd71SBarry Smith ierr = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr); 6043aeaf226SBarry Smith n = levels; 6053aeaf226SBarry Smith ierr = PCSetFromOptions(pc);CHKERRQ(ierr); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */ 6063aeaf226SBarry Smith mglevels = mg->levels; 6073aeaf226SBarry Smith } 6083aeaf226SBarry Smith } 60998e04984SBarry Smith ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr); 6103aeaf226SBarry Smith 611f3fbd535SBarry Smith 612f3fbd535SBarry Smith /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */ 613f3fbd535SBarry Smith /* so use those from global PC */ 614f3fbd535SBarry Smith /* Is this what we always want? What if user wants to keep old one? */ 6150298fd71SBarry Smith ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset);CHKERRQ(ierr); 61698e04984SBarry Smith if (opsset) { 61798e04984SBarry Smith Mat mmat; 61823ee1639SBarry Smith ierr = KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat);CHKERRQ(ierr); 61998e04984SBarry Smith if (mmat == pc->pmat) opsset = PETSC_FALSE; 62098e04984SBarry Smith } 6214a5f07a5SMark F. Adams 6224a5f07a5SMark F. Adams if (!opsset) { 62371b23a65SMark F. Adams ierr = PCGetUseAmat(pc,&use_amat);CHKERRQ(ierr); 6244a5f07a5SMark F. Adams if(use_amat){ 625f3fbd535SBarry 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); 62623ee1639SBarry Smith ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat);CHKERRQ(ierr); 627f3fbd535SBarry Smith } 6284a5f07a5SMark F. Adams else { 6294a5f07a5SMark 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); 63023ee1639SBarry Smith ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat);CHKERRQ(ierr); 6314a5f07a5SMark F. Adams } 6324a5f07a5SMark F. Adams } 633f3fbd535SBarry Smith 6349c7ffb39SBarry Smith for (i=n-1; i>0; i--) { 6359c7ffb39SBarry Smith if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) { 6369c7ffb39SBarry Smith missinginterpolate = PETSC_TRUE; 6379c7ffb39SBarry Smith continue; 6389c7ffb39SBarry Smith } 6399c7ffb39SBarry Smith } 6402134b1e4SBarry Smith 6412134b1e4SBarry Smith ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB);CHKERRQ(ierr); 6422134b1e4SBarry Smith if (dA == dB) dAeqdB = PETSC_TRUE; 6432134b1e4SBarry Smith if ((mg->galerkin == PC_MG_GALERKIN_NONE) || (((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_MAT)) && !dAeqdB)) { 6442134b1e4SBarry Smith needRestricts = PETSC_TRUE; /* user must compute either mat, pmat, or both so must restrict x to coarser levels */ 6452134b1e4SBarry Smith } 6462134b1e4SBarry Smith 6472134b1e4SBarry Smith 6489c7ffb39SBarry Smith /* 6499c7ffb39SBarry Smith Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS) 6509c7ffb39SBarry Smith Skipping for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs? 6519c7ffb39SBarry Smith */ 6522134b1e4SBarry Smith if (missinginterpolate && pc->dm && mg->galerkin != PC_MG_GALERKIN_EXTERNAL && !pc->setupcalled) { 6532d2b81a6SBarry Smith /* construct the interpolation from the DMs */ 6542e499ae9SBarry Smith Mat p; 65573250ac0SBarry Smith Vec rscale; 656785e854fSJed Brown ierr = PetscMalloc1(n,&dms);CHKERRQ(ierr); 657218a07d4SBarry Smith dms[n-1] = pc->dm; 658ef1267afSMatthew G. Knepley /* Separately create them so we do not get DMKSP interference between levels */ 659ef1267afSMatthew G. Knepley for (i=n-2; i>-1; i--) {ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr);} 66041fce8fdSHong Zhang /* 66141fce8fdSHong Zhang Force the mat type of coarse level operator to be AIJ because usually we want to use LU for coarse level. 66241fce8fdSHong Zhang Notice that it can be overwritten by -mat_type because KSPSetUp() reads command line options. 66341fce8fdSHong Zhang But it is safe to use -dm_mat_type. 66441fce8fdSHong Zhang 66541fce8fdSHong Zhang The mat type should not be hardcoded like this, we need to find a better way. 66641fce8fdSHong Zhang ierr = DMSetMatType(dms[0],MATAIJ);CHKERRQ(ierr); 66741fce8fdSHong Zhang */ 668218a07d4SBarry Smith for (i=n-2; i>-1; i--) { 669942e3340SBarry Smith DMKSP kdm; 670eab52d2bSLawrence Mitchell PetscBool dmhasrestrict, dmhasinject; 6713c0c59f3SBarry Smith ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr); 6722134b1e4SBarry Smith if (!needRestricts) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);} 673942e3340SBarry Smith ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr); 674d1a3e677SJed Brown /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take 675d1a3e677SJed Brown * a bitwise OR of computing the matrix, RHS, and initial iterate. */ 6760298fd71SBarry Smith kdm->ops->computerhs = NULL; 6770298fd71SBarry Smith kdm->rhsctx = NULL; 67824c3aa18SBarry Smith if (!mglevels[i+1]->interpolate) { 679e727c939SJed Brown ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr); 6802d2b81a6SBarry Smith ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr); 68100ac8be1SBarry Smith if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);} 68273250ac0SBarry Smith ierr = VecDestroy(&rscale);CHKERRQ(ierr); 6836bf464f9SBarry Smith ierr = MatDestroy(&p);CHKERRQ(ierr); 684218a07d4SBarry Smith } 6853ad4599aSBarry Smith ierr = DMHasCreateRestriction(dms[i],&dmhasrestrict);CHKERRQ(ierr); 6863ad4599aSBarry Smith if (dmhasrestrict && !mglevels[i+1]->restrct){ 6873ad4599aSBarry Smith ierr = DMCreateRestriction(dms[i],dms[i+1],&p);CHKERRQ(ierr); 6883ad4599aSBarry Smith ierr = PCMGSetRestriction(pc,i+1,p);CHKERRQ(ierr); 6893ad4599aSBarry Smith ierr = MatDestroy(&p);CHKERRQ(ierr); 6903ad4599aSBarry Smith } 691eab52d2bSLawrence Mitchell ierr = DMHasCreateInjection(dms[i],&dmhasinject);CHKERRQ(ierr); 692eab52d2bSLawrence Mitchell if (dmhasinject && !mglevels[i+1]->inject){ 693eab52d2bSLawrence Mitchell ierr = DMCreateInjection(dms[i],dms[i+1],&p);CHKERRQ(ierr); 694eab52d2bSLawrence Mitchell ierr = PCMGSetInjection(pc,i+1,p);CHKERRQ(ierr); 695eab52d2bSLawrence Mitchell ierr = MatDestroy(&p);CHKERRQ(ierr); 696eab52d2bSLawrence Mitchell } 69724c3aa18SBarry Smith } 6982d2b81a6SBarry Smith 699ef1267afSMatthew G. Knepley for (i=n-2; i>-1; i--) {ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);} 7002d2b81a6SBarry Smith ierr = PetscFree(dms);CHKERRQ(ierr); 701ce4cda84SJed Brown } 702cab2e9ccSBarry Smith 703ce4cda84SJed Brown if (pc->dm && !pc->setupcalled) { 7042134b1e4SBarry Smith /* finest smoother also gets DM but it is not active, independent of whether galerkin==PC_MG_GALERKIN_EXTERNAL */ 705cab2e9ccSBarry Smith ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr); 706cab2e9ccSBarry Smith ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr); 707218a07d4SBarry Smith } 708218a07d4SBarry Smith 7092134b1e4SBarry Smith if (mg->galerkin < PC_MG_GALERKIN_NONE) { 7102134b1e4SBarry Smith Mat A,B; 7112134b1e4SBarry Smith PetscBool doA = PETSC_FALSE,doB = PETSC_FALSE; 7122134b1e4SBarry Smith MatReuse reuse = MAT_INITIAL_MATRIX; 7132134b1e4SBarry Smith 7142134b1e4SBarry Smith if ((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_BOTH)) doB = PETSC_TRUE; 7152134b1e4SBarry Smith if ((mg->galerkin == PC_MG_GALERKIN_MAT) || ((mg->galerkin == PC_MG_GALERKIN_BOTH) && (dA != dB))) doA = PETSC_TRUE; 7162134b1e4SBarry Smith if (pc->setupcalled) reuse = MAT_REUSE_MATRIX; 717f3fbd535SBarry Smith for (i=n-2; i>-1; i--) { 718b9367914SBarry 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"); 719b9367914SBarry Smith if (!mglevels[i+1]->interpolate) { 720b9367914SBarry Smith ierr = PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct);CHKERRQ(ierr); 721b9367914SBarry Smith } 722b9367914SBarry Smith if (!mglevels[i+1]->restrct) { 723b9367914SBarry Smith ierr = PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate);CHKERRQ(ierr); 724b9367914SBarry Smith } 7252134b1e4SBarry Smith if (reuse == MAT_REUSE_MATRIX) { 7262134b1e4SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,&A,&B);CHKERRQ(ierr); 7272134b1e4SBarry Smith } 7282134b1e4SBarry Smith if (doA) { 7292df6c5c3SPatrick Farrell ierr = MatGalerkin(mglevels[i+1]->restrct,dA,mglevels[i+1]->interpolate,reuse,1.0,&A);CHKERRQ(ierr); 7302134b1e4SBarry Smith } 7312134b1e4SBarry Smith if (doB) { 7322df6c5c3SPatrick Farrell ierr = MatGalerkin(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,reuse,1.0,&B);CHKERRQ(ierr); 733b9367914SBarry Smith } 7342134b1e4SBarry Smith /* the management of the PetscObjectReference() and PetscObjecDereference() below is rather delicate */ 7352134b1e4SBarry Smith if (!doA && dAeqdB) { 7362134b1e4SBarry Smith if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);} 7372134b1e4SBarry Smith A = B; 7382134b1e4SBarry Smith } else if (!doA && reuse == MAT_INITIAL_MATRIX ) { 7392134b1e4SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,&A,NULL);CHKERRQ(ierr); 7402134b1e4SBarry Smith ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 741b9367914SBarry Smith } 7422134b1e4SBarry Smith if (!doB && dAeqdB) { 7432134b1e4SBarry Smith if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);} 7442134b1e4SBarry Smith B = A; 7452134b1e4SBarry Smith } else if (!doB && reuse == MAT_INITIAL_MATRIX) { 74623ee1639SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&B);CHKERRQ(ierr); 7472134b1e4SBarry Smith ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 7487d537d92SJed Brown } 7492134b1e4SBarry Smith if (reuse == MAT_INITIAL_MATRIX) { 7502134b1e4SBarry Smith ierr = KSPSetOperators(mglevels[i]->smoothd,A,B);CHKERRQ(ierr); 7512134b1e4SBarry Smith ierr = PetscObjectDereference((PetscObject)A);CHKERRQ(ierr); 7522134b1e4SBarry Smith ierr = PetscObjectDereference((PetscObject)B);CHKERRQ(ierr); 7532134b1e4SBarry Smith } 7542134b1e4SBarry Smith dA = A; 755f3fbd535SBarry Smith dB = B; 756f3fbd535SBarry Smith } 757f3fbd535SBarry Smith } 7582134b1e4SBarry Smith if (needRestricts && pc->dm && pc->dm->x) { 759cab2e9ccSBarry Smith /* need to restrict Jacobian location to coarser meshes for evaluation */ 760cab2e9ccSBarry Smith for (i=n-2; i>-1; i--) { 761c88c5224SJed Brown Mat R; 762c88c5224SJed Brown Vec rscale; 763cab2e9ccSBarry Smith if (!mglevels[i]->smoothd->dm->x) { 764cab2e9ccSBarry Smith Vec *vecs; 7652a7a6963SBarry Smith ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vecs,0,NULL);CHKERRQ(ierr); 766cab2e9ccSBarry Smith mglevels[i]->smoothd->dm->x = vecs[0]; 767cab2e9ccSBarry Smith ierr = PetscFree(vecs);CHKERRQ(ierr); 768cab2e9ccSBarry Smith } 769c88c5224SJed Brown ierr = PCMGGetRestriction(pc,i+1,&R);CHKERRQ(ierr); 770c88c5224SJed Brown ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr); 771c88c5224SJed Brown ierr = MatRestrict(R,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);CHKERRQ(ierr); 772c88c5224SJed Brown ierr = VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,rscale);CHKERRQ(ierr); 773cab2e9ccSBarry Smith } 774f3fbd535SBarry Smith } 7752134b1e4SBarry Smith if (needRestricts && pc->dm) { 776caa4e7f2SJed Brown for (i=n-2; i>=0; i--) { 777ccceb50cSJed Brown DM dmfine,dmcoarse; 778caa4e7f2SJed Brown Mat Restrict,Inject; 779caa4e7f2SJed Brown Vec rscale; 780ccceb50cSJed Brown ierr = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr); 781ccceb50cSJed Brown ierr = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr); 782caa4e7f2SJed Brown ierr = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr); 783caa4e7f2SJed Brown ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr); 784eab52d2bSLawrence Mitchell ierr = PCMGGetInjection(pc,i+1,&Inject);CHKERRQ(ierr); 785caa4e7f2SJed Brown ierr = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr); 786caa4e7f2SJed Brown } 787caa4e7f2SJed Brown } 788f3fbd535SBarry Smith 789f3fbd535SBarry Smith if (!pc->setupcalled) { 790f3fbd535SBarry Smith for (i=0; i<n; i++) { 791f3fbd535SBarry Smith ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr); 792f3fbd535SBarry Smith } 793f3fbd535SBarry Smith for (i=1; i<n; i++) { 794f3fbd535SBarry Smith if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) { 795f3fbd535SBarry Smith ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr); 796f3fbd535SBarry Smith } 797f3fbd535SBarry Smith } 7983ad4599aSBarry Smith /* insure that if either interpolation or restriction is set the other other one is set */ 799f3fbd535SBarry Smith for (i=1; i<n; i++) { 8003ad4599aSBarry Smith ierr = PCMGGetInterpolation(pc,i,NULL);CHKERRQ(ierr); 8013ad4599aSBarry Smith ierr = PCMGGetRestriction(pc,i,NULL);CHKERRQ(ierr); 802f3fbd535SBarry Smith } 803f3fbd535SBarry Smith for (i=0; i<n-1; i++) { 804f3fbd535SBarry Smith if (!mglevels[i]->b) { 805f3fbd535SBarry Smith Vec *vec; 8062a7a6963SBarry Smith ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr); 807f3fbd535SBarry Smith ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr); 8086bf464f9SBarry Smith ierr = VecDestroy(vec);CHKERRQ(ierr); 809f3fbd535SBarry Smith ierr = PetscFree(vec);CHKERRQ(ierr); 810f3fbd535SBarry Smith } 811f3fbd535SBarry Smith if (!mglevels[i]->r && i) { 812f3fbd535SBarry Smith ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 813f3fbd535SBarry Smith ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr); 8146bf464f9SBarry Smith ierr = VecDestroy(&tvec);CHKERRQ(ierr); 815f3fbd535SBarry Smith } 816f3fbd535SBarry Smith if (!mglevels[i]->x) { 817f3fbd535SBarry Smith ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 818f3fbd535SBarry Smith ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr); 8196bf464f9SBarry Smith ierr = VecDestroy(&tvec);CHKERRQ(ierr); 820f3fbd535SBarry Smith } 821f3fbd535SBarry Smith } 822f3fbd535SBarry Smith if (n != 1 && !mglevels[n-1]->r) { 823f3fbd535SBarry Smith /* PCMGSetR() on the finest level if user did not supply it */ 824f3fbd535SBarry Smith Vec *vec; 8252a7a6963SBarry Smith ierr = KSPCreateVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr); 826f3fbd535SBarry Smith ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr); 8276bf464f9SBarry Smith ierr = VecDestroy(vec);CHKERRQ(ierr); 828f3fbd535SBarry Smith ierr = PetscFree(vec);CHKERRQ(ierr); 829f3fbd535SBarry Smith } 830f3fbd535SBarry Smith } 831f3fbd535SBarry Smith 83298e04984SBarry Smith if (pc->dm) { 83398e04984SBarry Smith /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */ 83498e04984SBarry Smith for (i=0; i<n-1; i++) { 83598e04984SBarry Smith if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX; 83698e04984SBarry Smith } 83798e04984SBarry Smith } 838f3fbd535SBarry Smith 839f3fbd535SBarry Smith for (i=1; i<n; i++) { 8402a21e185SBarry Smith if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1){ 841f3fbd535SBarry Smith /* if doing only down then initial guess is zero */ 842f3fbd535SBarry Smith ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 843f3fbd535SBarry Smith } 84463e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 845f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr); 846899639b0SHong Zhang if (mglevels[i]->smoothd->reason == KSP_DIVERGED_PCSETUP_FAILED) { 847899639b0SHong Zhang pc->failedreason = PC_SUBPC_ERROR; 848899639b0SHong Zhang } 84963e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 850d42688cbSBarry Smith if (!mglevels[i]->residual) { 851d42688cbSBarry Smith Mat mat; 85213842ffbSBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,&mat,NULL);CHKERRQ(ierr); 85354b2cd4bSJed Brown ierr = PCMGSetResidual(pc,i,PCMGResidualDefault,mat);CHKERRQ(ierr); 854d42688cbSBarry Smith } 855f3fbd535SBarry Smith } 856f3fbd535SBarry Smith for (i=1; i<n; i++) { 857f3fbd535SBarry Smith if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) { 858f3fbd535SBarry Smith Mat downmat,downpmat; 859f3fbd535SBarry Smith 860f3fbd535SBarry Smith /* check if operators have been set for up, if not use down operators to set them */ 8610298fd71SBarry Smith ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);CHKERRQ(ierr); 862f3fbd535SBarry Smith if (!opsset) { 86323ee1639SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);CHKERRQ(ierr); 86423ee1639SBarry Smith ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat);CHKERRQ(ierr); 865f3fbd535SBarry Smith } 866f3fbd535SBarry Smith 867f3fbd535SBarry Smith ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr); 86863e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 869f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr); 870899639b0SHong Zhang if (mglevels[i]->smoothu->reason == KSP_DIVERGED_PCSETUP_FAILED) { 871899639b0SHong Zhang pc->failedreason = PC_SUBPC_ERROR; 872899639b0SHong Zhang } 87363e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 874f3fbd535SBarry Smith } 875f3fbd535SBarry Smith } 876f3fbd535SBarry Smith 87763e6d426SJed Brown if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 878f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr); 879899639b0SHong Zhang if (mglevels[0]->smoothd->reason == KSP_DIVERGED_PCSETUP_FAILED) { 880899639b0SHong Zhang pc->failedreason = PC_SUBPC_ERROR; 881899639b0SHong Zhang } 88263e6d426SJed Brown if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 883f3fbd535SBarry Smith 884f3fbd535SBarry Smith /* 885f3fbd535SBarry Smith Dump the interpolation/restriction matrices plus the 886e3c5b3baSBarry Smith Jacobian/stiffness on each level. This allows MATLAB users to 887f3fbd535SBarry Smith easily check if the Galerkin condition A_c = R A_f R^T is satisfied. 888f3fbd535SBarry Smith 889f3fbd535SBarry Smith Only support one or the other at the same time. 890f3fbd535SBarry Smith */ 891f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER) 892c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);CHKERRQ(ierr); 893ce94432eSBarry Smith if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc)); 894f3fbd535SBarry Smith dump = PETSC_FALSE; 895f3fbd535SBarry Smith #endif 896c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);CHKERRQ(ierr); 897ce94432eSBarry Smith if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc)); 898f3fbd535SBarry Smith 899f3fbd535SBarry Smith if (viewer) { 900f3fbd535SBarry Smith for (i=1; i<n; i++) { 901f3fbd535SBarry Smith ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr); 902f3fbd535SBarry Smith } 903f3fbd535SBarry Smith for (i=0; i<n; i++) { 904f3fbd535SBarry Smith ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr); 905f3fbd535SBarry Smith ierr = MatView(pc->mat,viewer);CHKERRQ(ierr); 906f3fbd535SBarry Smith } 907f3fbd535SBarry Smith } 908f3fbd535SBarry Smith PetscFunctionReturn(0); 909f3fbd535SBarry Smith } 910f3fbd535SBarry Smith 911f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/ 912f3fbd535SBarry Smith 9134b9ad928SBarry Smith /*@ 91497177400SBarry Smith PCMGGetLevels - Gets the number of levels to use with MG. 9154b9ad928SBarry Smith 9164b9ad928SBarry Smith Not Collective 9174b9ad928SBarry Smith 9184b9ad928SBarry Smith Input Parameter: 9194b9ad928SBarry Smith . pc - the preconditioner context 9204b9ad928SBarry Smith 9214b9ad928SBarry Smith Output parameter: 9224b9ad928SBarry Smith . levels - the number of levels 9234b9ad928SBarry Smith 9244b9ad928SBarry Smith Level: advanced 9254b9ad928SBarry Smith 9264b9ad928SBarry Smith .keywords: MG, get, levels, multigrid 9274b9ad928SBarry Smith 92897177400SBarry Smith .seealso: PCMGSetLevels() 9294b9ad928SBarry Smith @*/ 9307087cfbeSBarry Smith PetscErrorCode PCMGGetLevels(PC pc,PetscInt *levels) 9314b9ad928SBarry Smith { 932f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 9334b9ad928SBarry Smith 9344b9ad928SBarry Smith PetscFunctionBegin; 9350700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 9364482741eSBarry Smith PetscValidIntPointer(levels,2); 937f3fbd535SBarry Smith *levels = mg->nlevels; 9384b9ad928SBarry Smith PetscFunctionReturn(0); 9394b9ad928SBarry Smith } 9404b9ad928SBarry Smith 9414b9ad928SBarry Smith /*@ 94297177400SBarry Smith PCMGSetType - Determines the form of multigrid to use: 9434b9ad928SBarry Smith multiplicative, additive, full, or the Kaskade algorithm. 9444b9ad928SBarry Smith 945ad4df100SBarry Smith Logically Collective on PC 9464b9ad928SBarry Smith 9474b9ad928SBarry Smith Input Parameters: 9484b9ad928SBarry Smith + pc - the preconditioner context 9499dcbbd2bSBarry Smith - form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE, 9509dcbbd2bSBarry Smith PC_MG_FULL, PC_MG_KASKADE 9514b9ad928SBarry Smith 9524b9ad928SBarry Smith Options Database Key: 9534b9ad928SBarry Smith . -pc_mg_type <form> - Sets <form>, one of multiplicative, 9544b9ad928SBarry Smith additive, full, kaskade 9554b9ad928SBarry Smith 9564b9ad928SBarry Smith Level: advanced 9574b9ad928SBarry Smith 9584b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid 9594b9ad928SBarry Smith 96097177400SBarry Smith .seealso: PCMGSetLevels() 9614b9ad928SBarry Smith @*/ 9627087cfbeSBarry Smith PetscErrorCode PCMGSetType(PC pc,PCMGType form) 9634b9ad928SBarry Smith { 964f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 9654b9ad928SBarry Smith 9664b9ad928SBarry Smith PetscFunctionBegin; 9670700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 968c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc,form,2); 96931567311SBarry Smith mg->am = form; 9709dcbbd2bSBarry Smith if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG; 971c60c7ad4SBarry Smith else pc->ops->applyrichardson = NULL; 972c60c7ad4SBarry Smith PetscFunctionReturn(0); 973c60c7ad4SBarry Smith } 974c60c7ad4SBarry Smith 975c60c7ad4SBarry Smith /*@ 976c60c7ad4SBarry Smith PCMGGetType - Determines the form of multigrid to use: 977c60c7ad4SBarry Smith multiplicative, additive, full, or the Kaskade algorithm. 978c60c7ad4SBarry Smith 979c60c7ad4SBarry Smith Logically Collective on PC 980c60c7ad4SBarry Smith 981c60c7ad4SBarry Smith Input Parameter: 982c60c7ad4SBarry Smith . pc - the preconditioner context 983c60c7ad4SBarry Smith 984c60c7ad4SBarry Smith Output Parameter: 985c60c7ad4SBarry Smith . type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE 986c60c7ad4SBarry Smith 987c60c7ad4SBarry Smith 988c60c7ad4SBarry Smith Level: advanced 989c60c7ad4SBarry Smith 990c60c7ad4SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid 991c60c7ad4SBarry Smith 992c60c7ad4SBarry Smith .seealso: PCMGSetLevels() 993c60c7ad4SBarry Smith @*/ 994c60c7ad4SBarry Smith PetscErrorCode PCMGGetType(PC pc,PCMGType *type) 995c60c7ad4SBarry Smith { 996c60c7ad4SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 997c60c7ad4SBarry Smith 998c60c7ad4SBarry Smith PetscFunctionBegin; 999c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1000c60c7ad4SBarry Smith *type = mg->am; 10014b9ad928SBarry Smith PetscFunctionReturn(0); 10024b9ad928SBarry Smith } 10034b9ad928SBarry Smith 10044b9ad928SBarry Smith /*@ 10050d353602SBarry Smith PCMGSetCycleType - Sets the type cycles to use. Use PCMGSetCycleTypeOnLevel() for more 10064b9ad928SBarry Smith complicated cycling. 10074b9ad928SBarry Smith 1008ad4df100SBarry Smith Logically Collective on PC 10094b9ad928SBarry Smith 10104b9ad928SBarry Smith Input Parameters: 1011c2be2410SBarry Smith + pc - the multigrid context 1012c1cbb1deSBarry Smith - n - either PC_MG_CYCLE_V or PC_MG_CYCLE_W 10134b9ad928SBarry Smith 10144b9ad928SBarry Smith Options Database Key: 1015c1cbb1deSBarry Smith . -pc_mg_cycle_type <v,w> - provide the cycle desired 10164b9ad928SBarry Smith 10174b9ad928SBarry Smith Level: advanced 10184b9ad928SBarry Smith 10194b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 10204b9ad928SBarry Smith 10210d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel() 10224b9ad928SBarry Smith @*/ 10237087cfbeSBarry Smith PetscErrorCode PCMGSetCycleType(PC pc,PCMGCycleType n) 10244b9ad928SBarry Smith { 1025f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1026f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 102779416396SBarry Smith PetscInt i,levels; 10284b9ad928SBarry Smith 10294b9ad928SBarry Smith PetscFunctionBegin; 10300700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 103169fbec6eSBarry Smith PetscValidLogicalCollectiveEnum(pc,n,2); 10321a97d7b6SStefano Zampini if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling"); 1033f3fbd535SBarry Smith levels = mglevels[0]->levels; 10342fa5cd67SKarl Rupp for (i=0; i<levels; i++) mglevels[i]->cycles = n; 10354b9ad928SBarry Smith PetscFunctionReturn(0); 10364b9ad928SBarry Smith } 10374b9ad928SBarry Smith 10388cc2d5dfSBarry Smith /*@ 10398cc2d5dfSBarry Smith PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step 10408cc2d5dfSBarry Smith of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used 10418cc2d5dfSBarry Smith 1042ad4df100SBarry Smith Logically Collective on PC 10438cc2d5dfSBarry Smith 10448cc2d5dfSBarry Smith Input Parameters: 10458cc2d5dfSBarry Smith + pc - the multigrid context 10468cc2d5dfSBarry Smith - n - number of cycles (default is 1) 10478cc2d5dfSBarry Smith 10488cc2d5dfSBarry Smith Options Database Key: 1049e1bc860dSBarry Smith . -pc_mg_multiplicative_cycles n 10508cc2d5dfSBarry Smith 10518cc2d5dfSBarry Smith Level: advanced 10528cc2d5dfSBarry Smith 105395452b02SPatrick Sanan Notes: 105495452b02SPatrick Sanan This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType() 10558cc2d5dfSBarry Smith 10568cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 10578cc2d5dfSBarry Smith 10588cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType() 10598cc2d5dfSBarry Smith @*/ 10607087cfbeSBarry Smith PetscErrorCode PCMGMultiplicativeSetCycles(PC pc,PetscInt n) 10618cc2d5dfSBarry Smith { 1062f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 10638cc2d5dfSBarry Smith 10648cc2d5dfSBarry Smith PetscFunctionBegin; 10650700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1066c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 10672a21e185SBarry Smith mg->cyclesperpcapply = n; 10688cc2d5dfSBarry Smith PetscFunctionReturn(0); 10698cc2d5dfSBarry Smith } 10708cc2d5dfSBarry Smith 10712134b1e4SBarry Smith PetscErrorCode PCMGSetGalerkin_MG(PC pc,PCMGGalerkinType use) 1072fb15c04eSBarry Smith { 1073fb15c04eSBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1074fb15c04eSBarry Smith 1075fb15c04eSBarry Smith PetscFunctionBegin; 10762134b1e4SBarry Smith mg->galerkin = use; 1077fb15c04eSBarry Smith PetscFunctionReturn(0); 1078fb15c04eSBarry Smith } 1079fb15c04eSBarry Smith 1080c2be2410SBarry Smith /*@ 108197177400SBarry Smith PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the 108282b87d87SMatthew G. Knepley finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i 1083c2be2410SBarry Smith 1084ad4df100SBarry Smith Logically Collective on PC 1085c2be2410SBarry Smith 1086c2be2410SBarry Smith Input Parameters: 1087c91913d3SJed Brown + pc - the multigrid context 10882134b1e4SBarry Smith - use - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, or PC_MG_GALERKIN_NONE 1089c2be2410SBarry Smith 1090c2be2410SBarry Smith Options Database Key: 10912134b1e4SBarry Smith . -pc_mg_galerkin <both,pmat,mat,none> 1092c2be2410SBarry Smith 1093c2be2410SBarry Smith Level: intermediate 1094c2be2410SBarry Smith 109595452b02SPatrick Sanan Notes: 109695452b02SPatrick Sanan Some codes that use PCMG such as PCGAMG use Galerkin internally while constructing the hierarchy and thus do not 10971c1aac46SBarry Smith use the PCMG construction of the coarser grids. 10981c1aac46SBarry Smith 1099c2be2410SBarry Smith .keywords: MG, set, Galerkin 1100c2be2410SBarry Smith 11012134b1e4SBarry Smith .seealso: PCMGGetGalerkin(), PCMGGalerkinType 11023fc8bf9cSBarry Smith 1103c2be2410SBarry Smith @*/ 11042134b1e4SBarry Smith PetscErrorCode PCMGSetGalerkin(PC pc,PCMGGalerkinType use) 1105c2be2410SBarry Smith { 1106fb15c04eSBarry Smith PetscErrorCode ierr; 1107c2be2410SBarry Smith 1108c2be2410SBarry Smith PetscFunctionBegin; 11090700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 11102134b1e4SBarry Smith ierr = PetscTryMethod(pc,"PCMGSetGalerkin_C",(PC,PCMGGalerkinType),(pc,use));CHKERRQ(ierr); 1111c2be2410SBarry Smith PetscFunctionReturn(0); 1112c2be2410SBarry Smith } 1113c2be2410SBarry Smith 11143fc8bf9cSBarry Smith /*@ 11153fc8bf9cSBarry Smith PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. 111682b87d87SMatthew G. Knepley A_i-1 = r_i * A_i * p_i 11173fc8bf9cSBarry Smith 11183fc8bf9cSBarry Smith Not Collective 11193fc8bf9cSBarry Smith 11203fc8bf9cSBarry Smith Input Parameter: 11213fc8bf9cSBarry Smith . pc - the multigrid context 11223fc8bf9cSBarry Smith 11233fc8bf9cSBarry Smith Output Parameter: 11242134b1e4SBarry Smith . galerkin - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, PC_MG_GALERKIN_NONE, or PC_MG_GALERKIN_EXTERNAL 11253fc8bf9cSBarry Smith 11263fc8bf9cSBarry Smith Level: intermediate 11273fc8bf9cSBarry Smith 11283fc8bf9cSBarry Smith .keywords: MG, set, Galerkin 11293fc8bf9cSBarry Smith 11302134b1e4SBarry Smith .seealso: PCMGSetGalerkin(), PCMGGalerkinType 11313fc8bf9cSBarry Smith 11323fc8bf9cSBarry Smith @*/ 11332134b1e4SBarry Smith PetscErrorCode PCMGGetGalerkin(PC pc,PCMGGalerkinType *galerkin) 11343fc8bf9cSBarry Smith { 1135f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 11363fc8bf9cSBarry Smith 11373fc8bf9cSBarry Smith PetscFunctionBegin; 11380700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 11392134b1e4SBarry Smith *galerkin = mg->galerkin; 11403fc8bf9cSBarry Smith PetscFunctionReturn(0); 11413fc8bf9cSBarry Smith } 11423fc8bf9cSBarry Smith 11434b9ad928SBarry Smith /*@ 114406792cafSBarry Smith PCMGSetNumberSmooth - Sets the number of pre and post-smoothing steps to use 1145710315b6SLawrence Mitchell on all levels. Use PCMGDistinctSmoothUp() to create separate up and down smoothers if you want different numbers of 1146710315b6SLawrence Mitchell pre- and post-smoothing steps. 114706792cafSBarry Smith 114806792cafSBarry Smith Logically Collective on PC 114906792cafSBarry Smith 115006792cafSBarry Smith Input Parameters: 115106792cafSBarry Smith + mg - the multigrid context 115206792cafSBarry Smith - n - the number of smoothing steps 115306792cafSBarry Smith 115406792cafSBarry Smith Options Database Key: 115506792cafSBarry Smith + -mg_levels_ksp_max_it <n> - Sets number of pre and post-smoothing steps 115606792cafSBarry Smith 115706792cafSBarry Smith Level: advanced 115806792cafSBarry Smith 115995452b02SPatrick Sanan Notes: 116095452b02SPatrick Sanan this does not set a value on the coarsest grid, since we assume that 116106792cafSBarry Smith there is no separate smooth up on the coarsest grid. 116206792cafSBarry Smith 116306792cafSBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid 116406792cafSBarry Smith 1165710315b6SLawrence Mitchell .seealso: PCMGSetDistinctSmoothUp() 116606792cafSBarry Smith @*/ 116706792cafSBarry Smith PetscErrorCode PCMGSetNumberSmooth(PC pc,PetscInt n) 116806792cafSBarry Smith { 116906792cafSBarry Smith PC_MG *mg = (PC_MG*)pc->data; 117006792cafSBarry Smith PC_MG_Levels **mglevels = mg->levels; 117106792cafSBarry Smith PetscErrorCode ierr; 117206792cafSBarry Smith PetscInt i,levels; 117306792cafSBarry Smith 117406792cafSBarry Smith PetscFunctionBegin; 117506792cafSBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 117606792cafSBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 11771a97d7b6SStefano Zampini if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling"); 117806792cafSBarry Smith levels = mglevels[0]->levels; 117906792cafSBarry Smith 118006792cafSBarry Smith for (i=1; i<levels; i++) { 118106792cafSBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 118206792cafSBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 118306792cafSBarry Smith mg->default_smoothu = n; 118406792cafSBarry Smith mg->default_smoothd = n; 118506792cafSBarry Smith } 118606792cafSBarry Smith PetscFunctionReturn(0); 118706792cafSBarry Smith } 118806792cafSBarry Smith 1189f442ab6aSBarry Smith /*@ 1190f442ab6aSBarry Smith PCMGSetDistinctSmoothUp - sets the up (post) smoother to be a seperate KSP from the down (pre) smoother on all levels 1191710315b6SLawrence Mitchell and adds the suffix _up to the options name 1192f442ab6aSBarry Smith 1193f442ab6aSBarry Smith Logically Collective on PC 1194f442ab6aSBarry Smith 1195f442ab6aSBarry Smith Input Parameters: 1196f442ab6aSBarry Smith . pc - the preconditioner context 1197f442ab6aSBarry Smith 1198f442ab6aSBarry Smith Options Database Key: 1199f442ab6aSBarry Smith . -pc_mg_distinct_smoothup 1200f442ab6aSBarry Smith 1201f442ab6aSBarry Smith Level: advanced 1202f442ab6aSBarry Smith 120395452b02SPatrick Sanan Notes: 120495452b02SPatrick Sanan this does not set a value on the coarsest grid, since we assume that 1205f442ab6aSBarry Smith there is no separate smooth up on the coarsest grid. 1206f442ab6aSBarry Smith 1207f442ab6aSBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid 1208f442ab6aSBarry Smith 1209710315b6SLawrence Mitchell .seealso: PCMGSetNumberSmooth() 1210f442ab6aSBarry Smith @*/ 1211f442ab6aSBarry Smith PetscErrorCode PCMGSetDistinctSmoothUp(PC pc) 1212f442ab6aSBarry Smith { 1213f442ab6aSBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1214f442ab6aSBarry Smith PC_MG_Levels **mglevels = mg->levels; 1215f442ab6aSBarry Smith PetscErrorCode ierr; 1216f442ab6aSBarry Smith PetscInt i,levels; 1217f442ab6aSBarry Smith KSP subksp; 1218f442ab6aSBarry Smith 1219f442ab6aSBarry Smith PetscFunctionBegin; 1220f442ab6aSBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1221f442ab6aSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling"); 1222f442ab6aSBarry Smith levels = mglevels[0]->levels; 1223f442ab6aSBarry Smith 1224f442ab6aSBarry Smith for (i=1; i<levels; i++) { 1225710315b6SLawrence Mitchell const char *prefix = NULL; 1226f442ab6aSBarry Smith /* make sure smoother up and down are different */ 1227f442ab6aSBarry Smith ierr = PCMGGetSmootherUp(pc,i,&subksp);CHKERRQ(ierr); 1228710315b6SLawrence Mitchell ierr = KSPGetOptionsPrefix(mglevels[i]->smoothd,&prefix);CHKERRQ(ierr); 1229710315b6SLawrence Mitchell ierr = KSPSetOptionsPrefix(subksp,prefix);CHKERRQ(ierr); 1230f442ab6aSBarry Smith ierr = KSPAppendOptionsPrefix(subksp,"up_");CHKERRQ(ierr); 1231f442ab6aSBarry Smith } 1232f442ab6aSBarry Smith PetscFunctionReturn(0); 1233f442ab6aSBarry Smith } 1234f442ab6aSBarry Smith 12354b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/ 12364b9ad928SBarry Smith 12373b09bd56SBarry Smith /*MC 1238ccb205f8SBarry Smith PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional 12393b09bd56SBarry Smith information about the coarser grid matrices and restriction/interpolation operators. 12403b09bd56SBarry Smith 12413b09bd56SBarry Smith Options Database Keys: 12423b09bd56SBarry Smith + -pc_mg_levels <nlevels> - number of levels including finest 1243391689abSStefano Zampini . -pc_mg_cycle_type <v,w> - provide the cycle desired 12448c1c2452SJed Brown . -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default 12453b09bd56SBarry Smith . -pc_mg_log - log information about time spent on each level of the solver 1246710315b6SLawrence Mitchell . -pc_mg_distinct_smoothup - configure up (after interpolation) and down (before restriction) smoothers separately (with different options prefixes) 12472134b1e4SBarry Smith . -pc_mg_galerkin <both,pmat,mat,none> - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R' 12488cf6037eSBarry Smith . -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1) 12498cf6037eSBarry Smith . -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 1250e3c5b3baSBarry Smith to the Socket viewer for reading from MATLAB. 12518cf6037eSBarry Smith - -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices 12528cf6037eSBarry Smith to the binary output file called binaryoutput 12533b09bd56SBarry Smith 125495452b02SPatrick Sanan Notes: 125595452b02SPatrick Sanan If one uses a Krylov method such GMRES or CG as the smoother than one must use KSPFGMRES, KSPGCG, or KSPRICHARDSON as the outer Krylov method 12563b09bd56SBarry Smith 12578cf6037eSBarry Smith When run with a single level the smoother options are used on that level NOT the coarse grid solver options 12588cf6037eSBarry Smith 125923067569SBarry Smith When run with KSPRICHARDSON the convergence test changes slightly if monitor is turned on. The iteration count may change slightly. This 126023067569SBarry Smith is because without monitoring the residual norm is computed WITHIN each multigrid cycle on the finest level after the pre-smoothing 126123067569SBarry Smith (because the residual has just been computed for the multigrid algorithm and is hence available for free) while with monitoring the 126223067569SBarry Smith residual is computed at the end of each cycle. 126323067569SBarry Smith 12643b09bd56SBarry Smith Level: intermediate 12653b09bd56SBarry Smith 12668f87f92bSBarry Smith Concepts: multigrid/multilevel 12673b09bd56SBarry Smith 12688cf6037eSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE 1269710315b6SLawrence Mitchell PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), 1270710315b6SLawrence Mitchell PCMGSetDistinctSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 127197177400SBarry Smith PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 12720d353602SBarry Smith PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 12733b09bd56SBarry Smith M*/ 12743b09bd56SBarry Smith 12758cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc) 12764b9ad928SBarry Smith { 1277f3fbd535SBarry Smith PC_MG *mg; 1278f3fbd535SBarry Smith PetscErrorCode ierr; 1279f3fbd535SBarry Smith 12804b9ad928SBarry Smith PetscFunctionBegin; 1281b00a9115SJed Brown ierr = PetscNewLog(pc,&mg);CHKERRQ(ierr); 1282f3fbd535SBarry Smith pc->data = (void*)mg; 1283f3fbd535SBarry Smith mg->nlevels = -1; 128410eca3edSLisandro Dalcin mg->am = PC_MG_MULTIPLICATIVE; 12852134b1e4SBarry Smith mg->galerkin = PC_MG_GALERKIN_NONE; 1286f3fbd535SBarry Smith 128737a44384SMark Adams pc->useAmat = PETSC_TRUE; 128837a44384SMark Adams 12894b9ad928SBarry Smith pc->ops->apply = PCApply_MG; 12904b9ad928SBarry Smith pc->ops->setup = PCSetUp_MG; 1291a06653b4SBarry Smith pc->ops->reset = PCReset_MG; 12924b9ad928SBarry Smith pc->ops->destroy = PCDestroy_MG; 12934b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_MG; 12944b9ad928SBarry Smith pc->ops->view = PCView_MG; 1295fb15c04eSBarry Smith 1296fb15c04eSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",PCMGSetGalerkin_MG);CHKERRQ(ierr); 12974b9ad928SBarry Smith PetscFunctionReturn(0); 12984b9ad928SBarry Smith } 1299