1dba47a55SKris Buschelman 24b9ad928SBarry Smith /* 34b9ad928SBarry Smith Defines the multigrid preconditioner interface. 44b9ad928SBarry Smith */ 5af0996ceSBarry Smith #include <petsc/private/pcmgimpl.h> /*I "petscksp.h" I*/ 641b6fd38SMatthew G. Knepley #include <petsc/private/kspimpl.h> 71e25c274SJed Brown #include <petscdm.h> 8391689abSStefano Zampini PETSC_INTERN PetscErrorCode PCPreSolveChangeRHS(PC,PetscBool*); 94b9ad928SBarry Smith 10f3b08a26SMatthew G. Knepley /* 11f3b08a26SMatthew G. Knepley Contains the list of registered coarse space construction routines 12f3b08a26SMatthew G. Knepley */ 13f3b08a26SMatthew G. Knepley PetscFunctionList PCMGCoarseList = NULL; 14f3b08a26SMatthew G. Knepley 1530b0564aSStefano Zampini PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG_Levels **mglevelsin,PetscBool transpose,PetscBool matapp,PCRichardsonConvergedReason *reason) 164b9ad928SBarry Smith { 1731567311SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1831567311SBarry Smith PC_MG_Levels *mgc,*mglevels = *mglevelsin; 1931567311SBarry Smith PetscInt cycles = (mglevels->level == 1) ? 1 : (PetscInt) mglevels->cycles; 204b9ad928SBarry Smith 214b9ad928SBarry Smith PetscFunctionBegin; 229566063dSJacob Faibussowitsch if (mglevels->eventsmoothsolve) PetscCall(PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0)); 23fcb023d4SJed Brown if (!transpose) { 2430b0564aSStefano Zampini if (matapp) { 259566063dSJacob Faibussowitsch PetscCall(KSPMatSolve(mglevels->smoothd,mglevels->B,mglevels->X)); /* pre-smooth */ 269566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(mglevels->smoothd,pc,NULL)); 2730b0564aSStefano Zampini } else { 289566063dSJacob Faibussowitsch PetscCall(KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x)); /* pre-smooth */ 299566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(mglevels->smoothd,pc,mglevels->x)); 3030b0564aSStefano Zampini } 31fcb023d4SJed Brown } else { 3228b400f6SJacob Faibussowitsch PetscCheck(!matapp,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Not supported"); 339566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(mglevels->smoothu,mglevels->b,mglevels->x)); /* transpose of post-smooth */ 349566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(mglevels->smoothu,pc,mglevels->x)); 35fcb023d4SJed Brown } 369566063dSJacob Faibussowitsch if (mglevels->eventsmoothsolve) PetscCall(PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0)); 3731567311SBarry Smith if (mglevels->level) { /* not the coarsest grid */ 389566063dSJacob Faibussowitsch if (mglevels->eventresidual) PetscCall(PetscLogEventBegin(mglevels->eventresidual,0,0,0,0)); 3930b0564aSStefano Zampini if (matapp && !mglevels->R) { 409566063dSJacob Faibussowitsch PetscCall(MatDuplicate(mglevels->B,MAT_DO_NOT_COPY_VALUES,&mglevels->R)); 4130b0564aSStefano Zampini } 42fcb023d4SJed Brown if (!transpose) { 439566063dSJacob Faibussowitsch if (matapp) PetscCall((*mglevels->matresidual)(mglevels->A,mglevels->B,mglevels->X,mglevels->R)); 449566063dSJacob Faibussowitsch else PetscCall((*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r)); 45fcb023d4SJed Brown } else { 469566063dSJacob Faibussowitsch if (matapp) PetscCall((*mglevels->matresidualtranspose)(mglevels->A,mglevels->B,mglevels->X,mglevels->R)); 479566063dSJacob Faibussowitsch else PetscCall((*mglevels->residualtranspose)(mglevels->A,mglevels->b,mglevels->x,mglevels->r)); 48fcb023d4SJed Brown } 499566063dSJacob Faibussowitsch if (mglevels->eventresidual) PetscCall(PetscLogEventEnd(mglevels->eventresidual,0,0,0,0)); 504b9ad928SBarry Smith 514b9ad928SBarry Smith /* if on finest level and have convergence criteria set */ 5231567311SBarry Smith if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) { 534b9ad928SBarry Smith PetscReal rnorm; 549566063dSJacob Faibussowitsch PetscCall(VecNorm(mglevels->r,NORM_2,&rnorm)); 554b9ad928SBarry Smith if (rnorm <= mg->ttol) { 5670441072SBarry Smith if (rnorm < mg->abstol) { 574d0a8057SBarry Smith *reason = PCRICHARDSON_CONVERGED_ATOL; 589566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc,"Linear solver has converged. Residual norm %g is less than absolute tolerance %g\n",(double)rnorm,(double)mg->abstol)); 594b9ad928SBarry Smith } else { 604d0a8057SBarry Smith *reason = PCRICHARDSON_CONVERGED_RTOL; 619566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc,"Linear solver has converged. Residual norm %g is less than relative tolerance times initial residual norm %g\n",(double)rnorm,(double)mg->ttol)); 624b9ad928SBarry Smith } 634b9ad928SBarry Smith PetscFunctionReturn(0); 644b9ad928SBarry Smith } 654b9ad928SBarry Smith } 664b9ad928SBarry Smith 6731567311SBarry Smith mgc = *(mglevelsin - 1); 689566063dSJacob Faibussowitsch if (mglevels->eventinterprestrict) PetscCall(PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0)); 69fcb023d4SJed Brown if (!transpose) { 709566063dSJacob Faibussowitsch if (matapp) PetscCall(MatMatRestrict(mglevels->restrct,mglevels->R,&mgc->B)); 719566063dSJacob Faibussowitsch else PetscCall(MatRestrict(mglevels->restrct,mglevels->r,mgc->b)); 72fcb023d4SJed Brown } else { 739566063dSJacob Faibussowitsch if (matapp) PetscCall(MatMatRestrict(mglevels->interpolate,mglevels->R,&mgc->B)); 749566063dSJacob Faibussowitsch else PetscCall(MatRestrict(mglevels->interpolate,mglevels->r,mgc->b)); 75fcb023d4SJed Brown } 769566063dSJacob Faibussowitsch if (mglevels->eventinterprestrict) PetscCall(PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0)); 7730b0564aSStefano Zampini if (matapp) { 7830b0564aSStefano Zampini if (!mgc->X) { 799566063dSJacob Faibussowitsch PetscCall(MatDuplicate(mgc->B,MAT_DO_NOT_COPY_VALUES,&mgc->X)); 8030b0564aSStefano Zampini } else { 819566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(mgc->X)); 8230b0564aSStefano Zampini } 8330b0564aSStefano Zampini } else { 849566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(mgc->x)); 8530b0564aSStefano Zampini } 864b9ad928SBarry Smith while (cycles--) { 879566063dSJacob Faibussowitsch PetscCall(PCMGMCycle_Private(pc,mglevelsin-1,transpose,matapp,reason)); 884b9ad928SBarry Smith } 899566063dSJacob Faibussowitsch if (mglevels->eventinterprestrict) PetscCall(PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0)); 90fcb023d4SJed Brown if (!transpose) { 919566063dSJacob Faibussowitsch if (matapp) PetscCall(MatMatInterpolateAdd(mglevels->interpolate,mgc->X,mglevels->X,&mglevels->X)); 929566063dSJacob Faibussowitsch else PetscCall(MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x)); 93fcb023d4SJed Brown } else { 949566063dSJacob Faibussowitsch PetscCall(MatInterpolateAdd(mglevels->restrct,mgc->x,mglevels->x,mglevels->x)); 95fcb023d4SJed Brown } 969566063dSJacob Faibussowitsch if (mglevels->eventinterprestrict) PetscCall(PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0)); 979566063dSJacob Faibussowitsch if (mglevels->eventsmoothsolve) PetscCall(PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0)); 98fcb023d4SJed Brown if (!transpose) { 9930b0564aSStefano Zampini if (matapp) { 1009566063dSJacob Faibussowitsch PetscCall(KSPMatSolve(mglevels->smoothu,mglevels->B,mglevels->X)); /* post smooth */ 1019566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(mglevels->smoothu,pc,NULL)); 10230b0564aSStefano Zampini } else { 1039566063dSJacob Faibussowitsch PetscCall(KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x)); /* post smooth */ 1049566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(mglevels->smoothu,pc,mglevels->x)); 10530b0564aSStefano Zampini } 106fcb023d4SJed Brown } else { 10728b400f6SJacob Faibussowitsch PetscCheck(!matapp,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Not supported"); 1089566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(mglevels->smoothd,mglevels->b,mglevels->x)); /* post smooth */ 1099566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(mglevels->smoothd,pc,mglevels->x)); 110fcb023d4SJed Brown } 11141b6fd38SMatthew G. Knepley if (mglevels->cr) { 11228b400f6SJacob Faibussowitsch PetscCheck(!matapp,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Not supported"); 11341b6fd38SMatthew G. Knepley /* TODO Turn on copy and turn off noisy if we have an exact solution 1149566063dSJacob Faibussowitsch PetscCall(VecCopy(mglevels->x, mglevels->crx)); 1159566063dSJacob Faibussowitsch PetscCall(VecCopy(mglevels->b, mglevels->crb)); */ 1169566063dSJacob Faibussowitsch PetscCall(KSPSetNoisy_Private(mglevels->crx)); 1179566063dSJacob Faibussowitsch PetscCall(KSPSolve(mglevels->cr,mglevels->crb,mglevels->crx)); /* compatible relaxation */ 1189566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(mglevels->cr,pc,mglevels->crx)); 11941b6fd38SMatthew G. Knepley } 1209566063dSJacob Faibussowitsch if (mglevels->eventsmoothsolve) PetscCall(PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0)); 1214b9ad928SBarry Smith } 1224b9ad928SBarry Smith PetscFunctionReturn(0); 1234b9ad928SBarry Smith } 1244b9ad928SBarry Smith 125ace3abfcSBarry 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) 1264b9ad928SBarry Smith { 127f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 128f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 129391689abSStefano Zampini PC tpc; 130391689abSStefano Zampini PetscBool changeu,changed; 131f3fbd535SBarry Smith PetscInt levels = mglevels[0]->levels,i; 1324b9ad928SBarry Smith 1334b9ad928SBarry Smith PetscFunctionBegin; 134a762d673SBarry Smith /* When the DM is supplying the matrix then it will not exist until here */ 135a762d673SBarry Smith for (i=0; i<levels; i++) { 136a762d673SBarry Smith if (!mglevels[i]->A) { 1379566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL)); 1389566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mglevels[i]->A)); 139a762d673SBarry Smith } 140a762d673SBarry Smith } 141391689abSStefano Zampini 1429566063dSJacob Faibussowitsch PetscCall(KSPGetPC(mglevels[levels-1]->smoothd,&tpc)); 1439566063dSJacob Faibussowitsch PetscCall(PCPreSolveChangeRHS(tpc,&changed)); 1449566063dSJacob Faibussowitsch PetscCall(KSPGetPC(mglevels[levels-1]->smoothu,&tpc)); 1459566063dSJacob Faibussowitsch PetscCall(PCPreSolveChangeRHS(tpc,&changeu)); 146391689abSStefano Zampini if (!changed && !changeu) { 1479566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[levels-1]->b)); 148f3fbd535SBarry Smith mglevels[levels-1]->b = b; 149391689abSStefano Zampini } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */ 150391689abSStefano Zampini if (!mglevels[levels-1]->b) { 151391689abSStefano Zampini Vec *vec; 152391689abSStefano Zampini 1539566063dSJacob Faibussowitsch PetscCall(KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL)); 154391689abSStefano Zampini mglevels[levels-1]->b = *vec; 1559566063dSJacob Faibussowitsch PetscCall(PetscFree(vec)); 156391689abSStefano Zampini } 1579566063dSJacob Faibussowitsch PetscCall(VecCopy(b,mglevels[levels-1]->b)); 158391689abSStefano Zampini } 159f3fbd535SBarry Smith mglevels[levels-1]->x = x; 1604b9ad928SBarry Smith 16131567311SBarry Smith mg->rtol = rtol; 16231567311SBarry Smith mg->abstol = abstol; 16331567311SBarry Smith mg->dtol = dtol; 1644b9ad928SBarry Smith if (rtol) { 1654b9ad928SBarry Smith /* compute initial residual norm for relative convergence test */ 1664b9ad928SBarry Smith PetscReal rnorm; 1677319c654SBarry Smith if (zeroguess) { 1689566063dSJacob Faibussowitsch PetscCall(VecNorm(b,NORM_2,&rnorm)); 1697319c654SBarry Smith } else { 1709566063dSJacob Faibussowitsch PetscCall((*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w)); 1719566063dSJacob Faibussowitsch PetscCall(VecNorm(w,NORM_2,&rnorm)); 1727319c654SBarry Smith } 17331567311SBarry Smith mg->ttol = PetscMax(rtol*rnorm,abstol); 1742fa5cd67SKarl Rupp } else if (abstol) mg->ttol = abstol; 1752fa5cd67SKarl Rupp else mg->ttol = 0.0; 1764b9ad928SBarry Smith 1774d0a8057SBarry Smith /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't 178336babb1SJed Brown stop prematurely due to small residual */ 1794d0a8057SBarry Smith for (i=1; i<levels; i++) { 1809566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT)); 181f3fbd535SBarry Smith if (mglevels[i]->smoothu != mglevels[i]->smoothd) { 18223067569SBarry Smith /* For Richardson the initial guess is nonzero since it is solving in each cycle the original system not just applying as a preconditioner */ 1839566063dSJacob Faibussowitsch PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE)); 1849566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT)); 1854b9ad928SBarry Smith } 1864d0a8057SBarry Smith } 1874d0a8057SBarry Smith 1884d0a8057SBarry Smith *reason = (PCRichardsonConvergedReason)0; 1894d0a8057SBarry Smith for (i=0; i<its; i++) { 1909566063dSJacob Faibussowitsch PetscCall(PCMGMCycle_Private(pc,mglevels+levels-1,PETSC_FALSE,PETSC_FALSE,reason)); 1914d0a8057SBarry Smith if (*reason) break; 1924d0a8057SBarry Smith } 1934d0a8057SBarry Smith if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS; 1944d0a8057SBarry Smith *outits = i; 195391689abSStefano Zampini if (!changed && !changeu) mglevels[levels-1]->b = NULL; 1964b9ad928SBarry Smith PetscFunctionReturn(0); 1974b9ad928SBarry Smith } 1984b9ad928SBarry Smith 1993aeaf226SBarry Smith PetscErrorCode PCReset_MG(PC pc) 2003aeaf226SBarry Smith { 2013aeaf226SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 2023aeaf226SBarry Smith PC_MG_Levels **mglevels = mg->levels; 203*2b3cbbdaSStefano Zampini PetscInt i,n; 2043aeaf226SBarry Smith 2053aeaf226SBarry Smith PetscFunctionBegin; 2063aeaf226SBarry Smith if (mglevels) { 2073aeaf226SBarry Smith n = mglevels[0]->levels; 2083aeaf226SBarry Smith for (i=0; i<n-1; i++) { 2099566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[i+1]->r)); 2109566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[i]->b)); 2119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[i]->x)); 2129566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[i+1]->R)); 2139566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[i]->B)); 2149566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[i]->X)); 2159566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[i]->crx)); 2169566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[i]->crb)); 2179566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[i+1]->restrct)); 2189566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[i+1]->interpolate)); 2199566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[i+1]->inject)); 2209566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[i+1]->rscale)); 2213aeaf226SBarry Smith } 2229566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[n-1]->crx)); 2239566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[n-1]->crb)); 224391689abSStefano Zampini /* this is not null only if the smoother on the finest level 225391689abSStefano Zampini changes the rhs during PreSolve */ 2269566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[n-1]->b)); 2279566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[n-1]->B)); 2283aeaf226SBarry Smith 2293aeaf226SBarry Smith for (i=0; i<n; i++) { 230*2b3cbbdaSStefano Zampini PetscCall(MatDestroy(&mglevels[i]->coarseSpace)); 2319566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[i]->A)); 2323aeaf226SBarry Smith if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 2339566063dSJacob Faibussowitsch PetscCall(KSPReset(mglevels[i]->smoothd)); 2343aeaf226SBarry Smith } 2359566063dSJacob Faibussowitsch PetscCall(KSPReset(mglevels[i]->smoothu)); 2369566063dSJacob Faibussowitsch if (mglevels[i]->cr) PetscCall(KSPReset(mglevels[i]->cr)); 2373aeaf226SBarry Smith } 238f3b08a26SMatthew G. Knepley mg->Nc = 0; 2393aeaf226SBarry Smith } 2403aeaf226SBarry Smith PetscFunctionReturn(0); 2413aeaf226SBarry Smith } 2423aeaf226SBarry Smith 24341b6fd38SMatthew G. Knepley /* Implementing CR 24441b6fd38SMatthew G. Knepley 24541b6fd38SMatthew G. Knepley We only want to make corrections that ``do not change'' the coarse solution. What we mean by not changing is that if I prolong my coarse solution to the fine grid and then inject that fine solution back to the coarse grid, I get the same answer. Injection is what Brannick calls R. We want the complementary projector to Inj, which we will call S, after Brannick, so that Inj S = 0. Now the orthogonal projector onto the range of Inj^T is 24641b6fd38SMatthew G. Knepley 24741b6fd38SMatthew G. Knepley Inj^T (Inj Inj^T)^{-1} Inj 24841b6fd38SMatthew G. Knepley 24941b6fd38SMatthew G. Knepley and if Inj is a VecScatter, as it is now in PETSc, we have 25041b6fd38SMatthew G. Knepley 25141b6fd38SMatthew G. Knepley Inj^T Inj 25241b6fd38SMatthew G. Knepley 25341b6fd38SMatthew G. Knepley and 25441b6fd38SMatthew G. Knepley 25541b6fd38SMatthew G. Knepley S = I - Inj^T Inj 25641b6fd38SMatthew G. Knepley 25741b6fd38SMatthew G. Knepley since 25841b6fd38SMatthew G. Knepley 25941b6fd38SMatthew G. Knepley Inj S = Inj - (Inj Inj^T) Inj = 0. 26041b6fd38SMatthew G. Knepley 26141b6fd38SMatthew G. Knepley Brannick suggests 26241b6fd38SMatthew G. Knepley 26341b6fd38SMatthew G. Knepley A \to S^T A S \qquad\mathrm{and}\qquad M \to S^T M S 26441b6fd38SMatthew G. Knepley 26541b6fd38SMatthew G. Knepley but I do not think his :math:`S^T S = I` is correct. Our S is an orthogonal projector, so :math:`S^T S = S^2 = S`. We will use 26641b6fd38SMatthew G. Knepley 26741b6fd38SMatthew G. Knepley M^{-1} A \to S M^{-1} A S 26841b6fd38SMatthew G. Knepley 26941b6fd38SMatthew G. Knepley In fact, since it is somewhat hard in PETSc to do the symmetric application, we will just apply S on the left. 27041b6fd38SMatthew G. Knepley 27141b6fd38SMatthew G. Knepley Check: || Inj P - I ||_F < tol 27241b6fd38SMatthew G. Knepley Check: In general, Inj Inj^T = I 27341b6fd38SMatthew G. Knepley */ 27441b6fd38SMatthew G. Knepley 27541b6fd38SMatthew G. Knepley typedef struct { 27641b6fd38SMatthew G. Knepley PC mg; /* The PCMG object */ 27741b6fd38SMatthew G. Knepley PetscInt l; /* The multigrid level for this solver */ 27841b6fd38SMatthew G. Knepley Mat Inj; /* The injection matrix */ 27941b6fd38SMatthew G. Knepley Mat S; /* I - Inj^T Inj */ 28041b6fd38SMatthew G. Knepley } CRContext; 28141b6fd38SMatthew G. Knepley 28241b6fd38SMatthew G. Knepley static PetscErrorCode CRSetup_Private(PC pc) 28341b6fd38SMatthew G. Knepley { 28441b6fd38SMatthew G. Knepley CRContext *ctx; 28541b6fd38SMatthew G. Knepley Mat It; 28641b6fd38SMatthew G. Knepley 28741b6fd38SMatthew G. Knepley PetscFunctionBeginUser; 2889566063dSJacob Faibussowitsch PetscCall(PCShellGetContext(pc, &ctx)); 2899566063dSJacob Faibussowitsch PetscCall(PCMGGetInjection(ctx->mg, ctx->l, &It)); 29028b400f6SJacob Faibussowitsch PetscCheck(It,PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_WRONGSTATE, "CR requires that injection be defined for this PCMG"); 2919566063dSJacob Faibussowitsch PetscCall(MatCreateTranspose(It, &ctx->Inj)); 2929566063dSJacob Faibussowitsch PetscCall(MatCreateNormal(ctx->Inj, &ctx->S)); 2939566063dSJacob Faibussowitsch PetscCall(MatScale(ctx->S, -1.0)); 2949566063dSJacob Faibussowitsch PetscCall(MatShift(ctx->S, 1.0)); 29541b6fd38SMatthew G. Knepley PetscFunctionReturn(0); 29641b6fd38SMatthew G. Knepley } 29741b6fd38SMatthew G. Knepley 29841b6fd38SMatthew G. Knepley static PetscErrorCode CRApply_Private(PC pc, Vec x, Vec y) 29941b6fd38SMatthew G. Knepley { 30041b6fd38SMatthew G. Knepley CRContext *ctx; 30141b6fd38SMatthew G. Knepley 30241b6fd38SMatthew G. Knepley PetscFunctionBeginUser; 3039566063dSJacob Faibussowitsch PetscCall(PCShellGetContext(pc, &ctx)); 3049566063dSJacob Faibussowitsch PetscCall(MatMult(ctx->S, x, y)); 30541b6fd38SMatthew G. Knepley PetscFunctionReturn(0); 30641b6fd38SMatthew G. Knepley } 30741b6fd38SMatthew G. Knepley 30841b6fd38SMatthew G. Knepley static PetscErrorCode CRDestroy_Private(PC pc) 30941b6fd38SMatthew G. Knepley { 31041b6fd38SMatthew G. Knepley CRContext *ctx; 31141b6fd38SMatthew G. Knepley 31241b6fd38SMatthew G. Knepley PetscFunctionBeginUser; 3139566063dSJacob Faibussowitsch PetscCall(PCShellGetContext(pc, &ctx)); 3149566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx->Inj)); 3159566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx->S)); 3169566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 3179566063dSJacob Faibussowitsch PetscCall(PCShellSetContext(pc, NULL)); 31841b6fd38SMatthew G. Knepley PetscFunctionReturn(0); 31941b6fd38SMatthew G. Knepley } 32041b6fd38SMatthew G. Knepley 32141b6fd38SMatthew G. Knepley static PetscErrorCode CreateCR_Private(PC pc, PetscInt l, PC *cr) 32241b6fd38SMatthew G. Knepley { 32341b6fd38SMatthew G. Knepley CRContext *ctx; 32441b6fd38SMatthew G. Knepley 32541b6fd38SMatthew G. Knepley PetscFunctionBeginUser; 3269566063dSJacob Faibussowitsch PetscCall(PCCreate(PetscObjectComm((PetscObject) pc), cr)); 3279566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) *cr, "S (complementary projector to injection)")); 3289566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(1, &ctx)); 32941b6fd38SMatthew G. Knepley ctx->mg = pc; 33041b6fd38SMatthew G. Knepley ctx->l = l; 3319566063dSJacob Faibussowitsch PetscCall(PCSetType(*cr, PCSHELL)); 3329566063dSJacob Faibussowitsch PetscCall(PCShellSetContext(*cr, ctx)); 3339566063dSJacob Faibussowitsch PetscCall(PCShellSetApply(*cr, CRApply_Private)); 3349566063dSJacob Faibussowitsch PetscCall(PCShellSetSetUp(*cr, CRSetup_Private)); 3359566063dSJacob Faibussowitsch PetscCall(PCShellSetDestroy(*cr, CRDestroy_Private)); 33641b6fd38SMatthew G. Knepley PetscFunctionReturn(0); 33741b6fd38SMatthew G. Knepley } 33841b6fd38SMatthew G. Knepley 33997d33e41SMatthew G. Knepley PetscErrorCode PCMGSetLevels_MG(PC pc,PetscInt levels,MPI_Comm *comms) 3404b9ad928SBarry Smith { 341f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 342ce94432eSBarry Smith MPI_Comm comm; 3433aeaf226SBarry Smith PC_MG_Levels **mglevels = mg->levels; 34410eca3edSLisandro Dalcin PCMGType mgtype = mg->am; 34510167fecSBarry Smith PetscInt mgctype = (PetscInt) PC_MG_CYCLE_V; 346f3fbd535SBarry Smith PetscInt i; 347f3fbd535SBarry Smith PetscMPIInt size; 348f3fbd535SBarry Smith const char *prefix; 349f3fbd535SBarry Smith PC ipc; 3503aeaf226SBarry Smith PetscInt n; 3514b9ad928SBarry Smith 3524b9ad928SBarry Smith PetscFunctionBegin; 3530700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 354548767e0SJed Brown PetscValidLogicalCollectiveInt(pc,levels,2); 355548767e0SJed Brown if (mg->nlevels == levels) PetscFunctionReturn(0); 3569566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)pc,&comm)); 3573aeaf226SBarry Smith if (mglevels) { 35810eca3edSLisandro Dalcin mgctype = mglevels[0]->cycles; 3593aeaf226SBarry Smith /* changing the number of levels so free up the previous stuff */ 3609566063dSJacob Faibussowitsch PetscCall(PCReset_MG(pc)); 3613aeaf226SBarry Smith n = mglevels[0]->levels; 3623aeaf226SBarry Smith for (i=0; i<n; i++) { 3633aeaf226SBarry Smith if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 3649566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&mglevels[i]->smoothd)); 3653aeaf226SBarry Smith } 3669566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&mglevels[i]->smoothu)); 3679566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&mglevels[i]->cr)); 3689566063dSJacob Faibussowitsch PetscCall(PetscFree(mglevels[i])); 3693aeaf226SBarry Smith } 3709566063dSJacob Faibussowitsch PetscCall(PetscFree(mg->levels)); 3713aeaf226SBarry Smith } 372f3fbd535SBarry Smith 373f3fbd535SBarry Smith mg->nlevels = levels; 374f3fbd535SBarry Smith 3759566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(levels,&mglevels)); 3769566063dSJacob Faibussowitsch PetscCall(PetscLogObjectMemory((PetscObject)pc,levels*(sizeof(PC_MG*)))); 377f3fbd535SBarry Smith 3789566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc,&prefix)); 379f3fbd535SBarry Smith 380a9db87a2SMatthew G Knepley mg->stageApply = 0; 381f3fbd535SBarry Smith for (i=0; i<levels; i++) { 3829566063dSJacob Faibussowitsch PetscCall(PetscNewLog(pc,&mglevels[i])); 3832fa5cd67SKarl Rupp 384f3fbd535SBarry Smith mglevels[i]->level = i; 385f3fbd535SBarry Smith mglevels[i]->levels = levels; 38610eca3edSLisandro Dalcin mglevels[i]->cycles = mgctype; 387336babb1SJed Brown mg->default_smoothu = 2; 388336babb1SJed Brown mg->default_smoothd = 2; 38963e6d426SJed Brown mglevels[i]->eventsmoothsetup = 0; 39063e6d426SJed Brown mglevels[i]->eventsmoothsolve = 0; 39163e6d426SJed Brown mglevels[i]->eventresidual = 0; 39263e6d426SJed Brown mglevels[i]->eventinterprestrict = 0; 393f3fbd535SBarry Smith 394f3fbd535SBarry Smith if (comms) comm = comms[i]; 395d5a8f7e6SBarry Smith if (comm != MPI_COMM_NULL) { 3969566063dSJacob Faibussowitsch PetscCall(KSPCreate(comm,&mglevels[i]->smoothd)); 3979566063dSJacob Faibussowitsch PetscCall(KSPSetErrorIfNotConverged(mglevels[i]->smoothd,pc->erroriffailure)); 3989566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i)); 3999566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix)); 4009566063dSJacob Faibussowitsch PetscCall(PetscObjectComposedDataSetInt((PetscObject) mglevels[i]->smoothd, PetscMGLevelId, mglevels[i]->level)); 4010ee9a668SBarry Smith if (i || levels == 1) { 4020ee9a668SBarry Smith char tprefix[128]; 4030ee9a668SBarry Smith 4049566063dSJacob Faibussowitsch PetscCall(KSPSetType(mglevels[i]->smoothd,KSPCHEBYSHEV)); 4059566063dSJacob Faibussowitsch PetscCall(KSPSetConvergenceTest(mglevels[i]->smoothd,KSPConvergedSkip,NULL,NULL)); 4069566063dSJacob Faibussowitsch PetscCall(KSPSetNormType(mglevels[i]->smoothd,KSP_NORM_NONE)); 4079566063dSJacob Faibussowitsch PetscCall(KSPGetPC(mglevels[i]->smoothd,&ipc)); 4089566063dSJacob Faibussowitsch PetscCall(PCSetType(ipc,PCSOR)); 4099566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg->default_smoothd)); 410f3fbd535SBarry Smith 4119566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(tprefix,128,"mg_levels_%d_",(int)i)); 4129566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix)); 4130ee9a668SBarry Smith } else { 4149566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_")); 415f3fbd535SBarry Smith 4169dbfc187SHong Zhang /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */ 4179566063dSJacob Faibussowitsch PetscCall(KSPSetType(mglevels[0]->smoothd,KSPPREONLY)); 4189566063dSJacob Faibussowitsch PetscCall(KSPGetPC(mglevels[0]->smoothd,&ipc)); 4199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 420f3fbd535SBarry Smith if (size > 1) { 4219566063dSJacob Faibussowitsch PetscCall(PCSetType(ipc,PCREDUNDANT)); 422f3fbd535SBarry Smith } else { 4239566063dSJacob Faibussowitsch PetscCall(PCSetType(ipc,PCLU)); 424f3fbd535SBarry Smith } 4259566063dSJacob Faibussowitsch PetscCall(PCFactorSetShiftType(ipc,MAT_SHIFT_INBLOCKS)); 426f3fbd535SBarry Smith } 4279566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[i]->smoothd)); 428d5a8f7e6SBarry Smith } 429f3fbd535SBarry Smith mglevels[i]->smoothu = mglevels[i]->smoothd; 43031567311SBarry Smith mg->rtol = 0.0; 43131567311SBarry Smith mg->abstol = 0.0; 43231567311SBarry Smith mg->dtol = 0.0; 43331567311SBarry Smith mg->ttol = 0.0; 43431567311SBarry Smith mg->cyclesperpcapply = 1; 435f3fbd535SBarry Smith } 436f3fbd535SBarry Smith mg->levels = mglevels; 4379566063dSJacob Faibussowitsch PetscCall(PCMGSetType(pc,mgtype)); 4384b9ad928SBarry Smith PetscFunctionReturn(0); 4394b9ad928SBarry Smith } 4404b9ad928SBarry Smith 44197d33e41SMatthew G. Knepley /*@C 44297d33e41SMatthew G. Knepley PCMGSetLevels - Sets the number of levels to use with MG. 44397d33e41SMatthew G. Knepley Must be called before any other MG routine. 44497d33e41SMatthew G. Knepley 44597d33e41SMatthew G. Knepley Logically Collective on PC 44697d33e41SMatthew G. Knepley 44797d33e41SMatthew G. Knepley Input Parameters: 44897d33e41SMatthew G. Knepley + pc - the preconditioner context 44997d33e41SMatthew G. Knepley . levels - the number of levels 45097d33e41SMatthew G. Knepley - comms - optional communicators for each level; this is to allow solving the coarser problems 451d5a8f7e6SBarry Smith on smaller sets of processes. For processes that are not included in the computation 45205552d4bSJunchao Zhang you must pass MPI_COMM_NULL. Use comms = NULL to specify that all processes 45305552d4bSJunchao Zhang should participate in each level of problem. 45497d33e41SMatthew G. Knepley 45597d33e41SMatthew G. Knepley Level: intermediate 45697d33e41SMatthew G. Knepley 45797d33e41SMatthew G. Knepley Notes: 45897d33e41SMatthew G. Knepley If the number of levels is one then the multigrid uses the -mg_levels prefix 45997d33e41SMatthew G. Knepley for setting the level options rather than the -mg_coarse prefix. 46097d33e41SMatthew G. Knepley 461d5a8f7e6SBarry Smith You can free the information in comms after this routine is called. 462d5a8f7e6SBarry Smith 463d5a8f7e6SBarry Smith The array of MPI communicators must contain MPI_COMM_NULL for those ranks that at each level 464d5a8f7e6SBarry Smith are not participating in the coarser solve. For example, with 2 levels and 1 and 2 ranks on 465d5a8f7e6SBarry Smith the two levels, rank 0 in the original communicator will pass in an array of 2 communicators 466d5a8f7e6SBarry Smith of size 2 and 1, while rank 1 in the original communicator will pass in array of 2 communicators 467d5a8f7e6SBarry Smith the first of size 2 and the second of value MPI_COMM_NULL since the rank 1 does not participate 468d5a8f7e6SBarry Smith in the coarse grid solve. 469d5a8f7e6SBarry Smith 470d5a8f7e6SBarry Smith Since each coarser level may have a new MPI_Comm with fewer ranks than the previous, one 471d5a8f7e6SBarry Smith must take special care in providing the restriction and interpolation operation. We recommend 472d5a8f7e6SBarry Smith providing these as two step operations; first perform a standard restriction or interpolation on 473d5a8f7e6SBarry Smith the full number of ranks for that level and then use an MPI call to copy the resulting vector 47405552d4bSJunchao Zhang array entries (after calls to VecGetArray()) to the smaller or larger number of ranks, note in both 475d5a8f7e6SBarry Smith cases the MPI calls must be made on the larger of the two communicators. Traditional MPI send and 476d5a8f7e6SBarry Smith recieves or MPI_AlltoAllv() could be used to do the reshuffling of the vector entries. 477d5a8f7e6SBarry Smith 47805552d4bSJunchao Zhang Fortran Notes: 47905552d4bSJunchao Zhang Use comms = PETSC_NULL_MPI_COMM as the equivalent of NULL in the C interface. Note PETSC_NULL_MPI_COMM 48005552d4bSJunchao Zhang is not MPI_COMM_NULL. It is more like PETSC_NULL_INTEGER, PETSC_NULL_REAL etc. 481d5a8f7e6SBarry Smith 482db781477SPatrick Sanan .seealso: `PCMGSetType()`, `PCMGGetLevels()` 48397d33e41SMatthew G. Knepley @*/ 48497d33e41SMatthew G. Knepley PetscErrorCode PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms) 48597d33e41SMatthew G. Knepley { 48697d33e41SMatthew G. Knepley PetscFunctionBegin; 48797d33e41SMatthew G. Knepley PetscValidHeaderSpecific(pc,PC_CLASSID,1); 48897d33e41SMatthew G. Knepley if (comms) PetscValidPointer(comms,3); 489cac4c232SBarry Smith PetscTryMethod(pc,"PCMGSetLevels_C",(PC,PetscInt,MPI_Comm*),(pc,levels,comms)); 49097d33e41SMatthew G. Knepley PetscFunctionReturn(0); 49197d33e41SMatthew G. Knepley } 49297d33e41SMatthew G. Knepley 493c07bf074SBarry Smith PetscErrorCode PCDestroy_MG(PC pc) 494c07bf074SBarry Smith { 495a06653b4SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 496a06653b4SBarry Smith PC_MG_Levels **mglevels = mg->levels; 497a06653b4SBarry Smith PetscInt i,n; 498c07bf074SBarry Smith 499c07bf074SBarry Smith PetscFunctionBegin; 5009566063dSJacob Faibussowitsch PetscCall(PCReset_MG(pc)); 501a06653b4SBarry Smith if (mglevels) { 502a06653b4SBarry Smith n = mglevels[0]->levels; 503a06653b4SBarry Smith for (i=0; i<n; i++) { 504a06653b4SBarry Smith if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 5059566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&mglevels[i]->smoothd)); 506a06653b4SBarry Smith } 5079566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&mglevels[i]->smoothu)); 5089566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&mglevels[i]->cr)); 5099566063dSJacob Faibussowitsch PetscCall(PetscFree(mglevels[i])); 510a06653b4SBarry Smith } 5119566063dSJacob Faibussowitsch PetscCall(PetscFree(mg->levels)); 512a06653b4SBarry Smith } 5139566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 5149566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGetInterpolations_C",NULL)); 5159566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGetCoarseOperators_C",NULL)); 516*2b3cbbdaSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",NULL)); 517*2b3cbbdaSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGGetLevels_C",NULL)); 518*2b3cbbdaSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGSetLevels_C",NULL)); 519*2b3cbbdaSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGetInterpolations_C",NULL)); 520*2b3cbbdaSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGetCoarseOperators_C",NULL)); 521*2b3cbbdaSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGSetAdaptInterpolation_C",NULL)); 522*2b3cbbdaSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGGetAdaptInterpolation_C",NULL)); 523*2b3cbbdaSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGSetAdaptCR_C",NULL)); 524*2b3cbbdaSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGGetAdaptCR_C",NULL)); 525*2b3cbbdaSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGSetAdaptCoarseSpaceType_C",NULL)); 526*2b3cbbdaSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGGetAdaptCoarseSpaceType_C",NULL)); 527f3fbd535SBarry Smith PetscFunctionReturn(0); 528f3fbd535SBarry Smith } 529f3fbd535SBarry Smith 530f3fbd535SBarry Smith /* 531f3fbd535SBarry Smith PCApply_MG - Runs either an additive, multiplicative, Kaskadic 532f3fbd535SBarry Smith or full cycle of multigrid. 533f3fbd535SBarry Smith 534f3fbd535SBarry Smith Note: 535f3fbd535SBarry Smith A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle(). 536f3fbd535SBarry Smith */ 53730b0564aSStefano Zampini static PetscErrorCode PCApply_MG_Internal(PC pc,Vec b,Vec x,Mat B,Mat X,PetscBool transpose) 538f3fbd535SBarry Smith { 539f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 540f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 541391689abSStefano Zampini PC tpc; 542f3fbd535SBarry Smith PetscInt levels = mglevels[0]->levels,i; 54330b0564aSStefano Zampini PetscBool changeu,changed,matapp; 544f3fbd535SBarry Smith 545f3fbd535SBarry Smith PetscFunctionBegin; 54630b0564aSStefano Zampini matapp = (PetscBool)(B && X); 5479566063dSJacob Faibussowitsch if (mg->stageApply) PetscCall(PetscLogStagePush(mg->stageApply)); 548e1d8e5deSBarry Smith /* When the DM is supplying the matrix then it will not exist until here */ 549298cc208SBarry Smith for (i=0; i<levels; i++) { 550e1d8e5deSBarry Smith if (!mglevels[i]->A) { 5519566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL)); 5529566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mglevels[i]->A)); 553e1d8e5deSBarry Smith } 554e1d8e5deSBarry Smith } 555e1d8e5deSBarry Smith 5569566063dSJacob Faibussowitsch PetscCall(KSPGetPC(mglevels[levels-1]->smoothd,&tpc)); 5579566063dSJacob Faibussowitsch PetscCall(PCPreSolveChangeRHS(tpc,&changed)); 5589566063dSJacob Faibussowitsch PetscCall(KSPGetPC(mglevels[levels-1]->smoothu,&tpc)); 5599566063dSJacob Faibussowitsch PetscCall(PCPreSolveChangeRHS(tpc,&changeu)); 560391689abSStefano Zampini if (!changeu && !changed) { 56130b0564aSStefano Zampini if (matapp) { 5629566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[levels-1]->B)); 56330b0564aSStefano Zampini mglevels[levels-1]->B = B; 56430b0564aSStefano Zampini } else { 5659566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[levels-1]->b)); 566f3fbd535SBarry Smith mglevels[levels-1]->b = b; 56730b0564aSStefano Zampini } 568391689abSStefano Zampini } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */ 56930b0564aSStefano Zampini if (matapp) { 57030b0564aSStefano Zampini if (mglevels[levels-1]->B) { 57130b0564aSStefano Zampini PetscInt N1,N2; 57230b0564aSStefano Zampini PetscBool flg; 57330b0564aSStefano Zampini 5749566063dSJacob Faibussowitsch PetscCall(MatGetSize(mglevels[levels-1]->B,NULL,&N1)); 5759566063dSJacob Faibussowitsch PetscCall(MatGetSize(B,NULL,&N2)); 5769566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)mglevels[levels-1]->B,((PetscObject)B)->type_name,&flg)); 57730b0564aSStefano Zampini if (N1 != N2 || !flg) { 5789566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[levels-1]->B)); 57930b0564aSStefano Zampini } 58030b0564aSStefano Zampini } 58130b0564aSStefano Zampini if (!mglevels[levels-1]->B) { 5829566063dSJacob Faibussowitsch PetscCall(MatDuplicate(B,MAT_COPY_VALUES,&mglevels[levels-1]->B)); 58330b0564aSStefano Zampini } else { 5849566063dSJacob Faibussowitsch PetscCall(MatCopy(B,mglevels[levels-1]->B,SAME_NONZERO_PATTERN)); 58530b0564aSStefano Zampini } 58630b0564aSStefano Zampini } else { 587391689abSStefano Zampini if (!mglevels[levels-1]->b) { 588391689abSStefano Zampini Vec *vec; 589391689abSStefano Zampini 5909566063dSJacob Faibussowitsch PetscCall(KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL)); 591391689abSStefano Zampini mglevels[levels-1]->b = *vec; 5929566063dSJacob Faibussowitsch PetscCall(PetscFree(vec)); 593391689abSStefano Zampini } 5949566063dSJacob Faibussowitsch PetscCall(VecCopy(b,mglevels[levels-1]->b)); 595391689abSStefano Zampini } 59630b0564aSStefano Zampini } 59730b0564aSStefano Zampini if (matapp) { mglevels[levels-1]->X = X; } 59830b0564aSStefano Zampini else { mglevels[levels-1]->x = x; } 59930b0564aSStefano Zampini 60030b0564aSStefano Zampini /* If coarser Xs are present, it means we have already block applied the PC at least once 60130b0564aSStefano Zampini Reset operators if sizes/type do no match */ 60230b0564aSStefano Zampini if (matapp && levels > 1 && mglevels[levels-2]->X) { 60330b0564aSStefano Zampini PetscInt Xc,Bc; 60430b0564aSStefano Zampini PetscBool flg; 60530b0564aSStefano Zampini 6069566063dSJacob Faibussowitsch PetscCall(MatGetSize(mglevels[levels-2]->X,NULL,&Xc)); 6079566063dSJacob Faibussowitsch PetscCall(MatGetSize(mglevels[levels-1]->B,NULL,&Bc)); 6089566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)mglevels[levels-2]->X,((PetscObject)mglevels[levels-1]->X)->type_name,&flg)); 60930b0564aSStefano Zampini if (Xc != Bc || !flg) { 6109566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[levels-1]->R)); 61130b0564aSStefano Zampini for (i=0;i<levels-1;i++) { 6129566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[i]->R)); 6139566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[i]->B)); 6149566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[i]->X)); 61530b0564aSStefano Zampini } 61630b0564aSStefano Zampini } 61730b0564aSStefano Zampini } 618391689abSStefano Zampini 61931567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 6209566063dSJacob Faibussowitsch if (matapp) PetscCall(MatZeroEntries(X)); 6219566063dSJacob Faibussowitsch else PetscCall(VecZeroEntries(x)); 62231567311SBarry Smith for (i=0; i<mg->cyclesperpcapply; i++) { 6239566063dSJacob Faibussowitsch PetscCall(PCMGMCycle_Private(pc,mglevels+levels-1,transpose,matapp,NULL)); 624f3fbd535SBarry Smith } 6252fa5cd67SKarl Rupp } else if (mg->am == PC_MG_ADDITIVE) { 6269566063dSJacob Faibussowitsch PetscCall(PCMGACycle_Private(pc,mglevels,transpose,matapp)); 6272fa5cd67SKarl Rupp } else if (mg->am == PC_MG_KASKADE) { 6289566063dSJacob Faibussowitsch PetscCall(PCMGKCycle_Private(pc,mglevels,transpose,matapp)); 6292fa5cd67SKarl Rupp } else { 6309566063dSJacob Faibussowitsch PetscCall(PCMGFCycle_Private(pc,mglevels,transpose,matapp)); 631f3fbd535SBarry Smith } 6329566063dSJacob Faibussowitsch if (mg->stageApply) PetscCall(PetscLogStagePop()); 63330b0564aSStefano Zampini if (!changeu && !changed) { 63430b0564aSStefano Zampini if (matapp) { mglevels[levels-1]->B = NULL; } 63530b0564aSStefano Zampini else { mglevels[levels-1]->b = NULL; } 63630b0564aSStefano Zampini } 637f3fbd535SBarry Smith PetscFunctionReturn(0); 638f3fbd535SBarry Smith } 639f3fbd535SBarry Smith 640fcb023d4SJed Brown static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x) 641fcb023d4SJed Brown { 642fcb023d4SJed Brown PetscFunctionBegin; 6439566063dSJacob Faibussowitsch PetscCall(PCApply_MG_Internal(pc,b,x,NULL,NULL,PETSC_FALSE)); 644fcb023d4SJed Brown PetscFunctionReturn(0); 645fcb023d4SJed Brown } 646fcb023d4SJed Brown 647fcb023d4SJed Brown static PetscErrorCode PCApplyTranspose_MG(PC pc,Vec b,Vec x) 648fcb023d4SJed Brown { 649fcb023d4SJed Brown PetscFunctionBegin; 6509566063dSJacob Faibussowitsch PetscCall(PCApply_MG_Internal(pc,b,x,NULL,NULL,PETSC_TRUE)); 65130b0564aSStefano Zampini PetscFunctionReturn(0); 65230b0564aSStefano Zampini } 65330b0564aSStefano Zampini 65430b0564aSStefano Zampini static PetscErrorCode PCMatApply_MG(PC pc,Mat b,Mat x) 65530b0564aSStefano Zampini { 65630b0564aSStefano Zampini PetscFunctionBegin; 6579566063dSJacob Faibussowitsch PetscCall(PCApply_MG_Internal(pc,NULL,NULL,b,x,PETSC_FALSE)); 658fcb023d4SJed Brown PetscFunctionReturn(0); 659fcb023d4SJed Brown } 660f3fbd535SBarry Smith 6614416b707SBarry Smith PetscErrorCode PCSetFromOptions_MG(PetscOptionItems *PetscOptionsObject,PC pc) 662f3fbd535SBarry Smith { 663710315b6SLawrence Mitchell PetscInt levels,cycles; 664f3b08a26SMatthew G. Knepley PetscBool flg, flg2; 665f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 6663d3eaba7SBarry Smith PC_MG_Levels **mglevels; 667f3fbd535SBarry Smith PCMGType mgtype; 668f3fbd535SBarry Smith PCMGCycleType mgctype; 6692134b1e4SBarry Smith PCMGGalerkinType gtype; 670*2b3cbbdaSStefano Zampini PCMGCoarseSpaceType coarseSpaceType; 671f3fbd535SBarry Smith 672f3fbd535SBarry Smith PetscFunctionBegin; 6731d743356SStefano Zampini levels = PetscMax(mg->nlevels,1); 674d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"Multigrid options"); 6759566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg)); 6761a97d7b6SStefano Zampini if (!flg && !mg->levels && pc->dm) { 6779566063dSJacob Faibussowitsch PetscCall(DMGetRefineLevel(pc->dm,&levels)); 678eb3f98d2SBarry Smith levels++; 6793aeaf226SBarry Smith mg->usedmfornumberoflevels = PETSC_TRUE; 680eb3f98d2SBarry Smith } 6819566063dSJacob Faibussowitsch PetscCall(PCMGSetLevels(pc,levels,NULL)); 6823aeaf226SBarry Smith mglevels = mg->levels; 6833aeaf226SBarry Smith 684f3fbd535SBarry Smith mgctype = (PCMGCycleType) mglevels[0]->cycles; 6859566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg)); 686f3fbd535SBarry Smith if (flg) { 6879566063dSJacob Faibussowitsch PetscCall(PCMGSetCycleType(pc,mgctype)); 6882fa5cd67SKarl Rupp } 6892134b1e4SBarry Smith gtype = mg->galerkin; 6909566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",PCMGGalerkinTypes,(PetscEnum)gtype,(PetscEnum*)>ype,&flg)); 6912134b1e4SBarry Smith if (flg) { 6929566063dSJacob Faibussowitsch PetscCall(PCMGSetGalerkin(pc,gtype)); 693f3fbd535SBarry Smith } 694*2b3cbbdaSStefano Zampini coarseSpaceType = mg->coarseSpaceType; 695*2b3cbbdaSStefano Zampini PetscCall(PetscOptionsEnum("-pc_mg_adapt_interp_coarse_space","Type of adaptive coarse space: none, polynomial, harmonic, eigenvector, generalized_eigenvector, gdsw","PCMGSetAdaptCoarseSpaceType",PCMGCoarseSpaceTypes,(PetscEnum)coarseSpaceType,(PetscEnum*)&coarseSpaceType,&flg)); 696*2b3cbbdaSStefano Zampini if (flg) PetscCall(PCMGSetAdaptCoarseSpaceType(pc,coarseSpaceType)); 6979566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_mg_adapt_interp_n","Size of the coarse space for adaptive interpolation","PCMGSetCoarseSpace",mg->Nc,&mg->Nc,&flg)); 6989566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_mg_mesp_monitor","Monitor the multilevel eigensolver","PCMGSetAdaptInterpolation",PETSC_FALSE,&mg->mespMonitor,&flg)); 69941b6fd38SMatthew G. Knepley flg2 = PETSC_FALSE; 7009566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_mg_adapt_cr","Monitor coarse space quality using Compatible Relaxation (CR)","PCMGSetAdaptCR",PETSC_FALSE,&flg2,&flg)); 7019566063dSJacob Faibussowitsch if (flg) PetscCall(PCMGSetAdaptCR(pc, flg2)); 702f56b1265SBarry Smith flg = PETSC_FALSE; 7039566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_mg_distinct_smoothup","Create separate smoothup KSP and append the prefix _up","PCMGSetDistinctSmoothUp",PETSC_FALSE,&flg,NULL)); 704f442ab6aSBarry Smith if (flg) { 7059566063dSJacob Faibussowitsch PetscCall(PCMGSetDistinctSmoothUp(pc)); 706f442ab6aSBarry Smith } 70731567311SBarry Smith mgtype = mg->am; 7089566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg)); 709f3fbd535SBarry Smith if (flg) { 7109566063dSJacob Faibussowitsch PetscCall(PCMGSetType(pc,mgtype)); 711f3fbd535SBarry Smith } 71231567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 7139566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGMultiplicativeSetCycles",mg->cyclesperpcapply,&cycles,&flg)); 714f3fbd535SBarry Smith if (flg) { 7159566063dSJacob Faibussowitsch PetscCall(PCMGMultiplicativeSetCycles(pc,cycles)); 716f3fbd535SBarry Smith } 717f3fbd535SBarry Smith } 718f3fbd535SBarry Smith flg = PETSC_FALSE; 7199566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,NULL)); 720f3fbd535SBarry Smith if (flg) { 721f3fbd535SBarry Smith PetscInt i; 722f3fbd535SBarry Smith char eventname[128]; 7231a97d7b6SStefano Zampini 724f3fbd535SBarry Smith levels = mglevels[0]->levels; 725f3fbd535SBarry Smith for (i=0; i<levels; i++) { 726f3fbd535SBarry Smith sprintf(eventname,"MGSetup Level %d",(int)i); 7279566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup)); 728f3fbd535SBarry Smith sprintf(eventname,"MGSmooth Level %d",(int)i); 7299566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve)); 730f3fbd535SBarry Smith if (i) { 731f3fbd535SBarry Smith sprintf(eventname,"MGResid Level %d",(int)i); 7329566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual)); 733f3fbd535SBarry Smith sprintf(eventname,"MGInterp Level %d",(int)i); 7349566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict)); 735f3fbd535SBarry Smith } 736f3fbd535SBarry Smith } 737eec35531SMatthew G Knepley 738ec60ca73SSatish Balay #if defined(PETSC_USE_LOG) 739eec35531SMatthew G Knepley { 740eec35531SMatthew G Knepley const char *sname = "MG Apply"; 741eec35531SMatthew G Knepley PetscStageLog stageLog; 742eec35531SMatthew G Knepley PetscInt st; 743eec35531SMatthew G Knepley 7449566063dSJacob Faibussowitsch PetscCall(PetscLogGetStageLog(&stageLog)); 745eec35531SMatthew G Knepley for (st = 0; st < stageLog->numStages; ++st) { 746eec35531SMatthew G Knepley PetscBool same; 747eec35531SMatthew G Knepley 7489566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(stageLog->stageInfo[st].name, sname, &same)); 7492fa5cd67SKarl Rupp if (same) mg->stageApply = st; 750eec35531SMatthew G Knepley } 751eec35531SMatthew G Knepley if (!mg->stageApply) { 7529566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister(sname, &mg->stageApply)); 753eec35531SMatthew G Knepley } 754eec35531SMatthew G Knepley } 755ec60ca73SSatish Balay #endif 756f3fbd535SBarry Smith } 757d0609cedSBarry Smith PetscOptionsHeadEnd(); 758f3b08a26SMatthew G. Knepley /* Check option consistency */ 7599566063dSJacob Faibussowitsch PetscCall(PCMGGetGalerkin(pc, >ype)); 7609566063dSJacob Faibussowitsch PetscCall(PCMGGetAdaptInterpolation(pc, &flg)); 76108401ef6SPierre Jolivet PetscCheck(!flg || !(gtype >= PC_MG_GALERKIN_NONE),PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_INCOMP, "Must use Galerkin coarse operators when adapting the interpolator"); 762f3fbd535SBarry Smith PetscFunctionReturn(0); 763f3fbd535SBarry Smith } 764f3fbd535SBarry Smith 7650a545947SLisandro Dalcin const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",NULL}; 7660a545947SLisandro Dalcin const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",NULL}; 7670a545947SLisandro Dalcin const char *const PCMGGalerkinTypes[] = {"both","pmat","mat","none","external","PCMGGalerkinType","PC_MG_GALERKIN",NULL}; 768*2b3cbbdaSStefano Zampini const char *const PCMGCoarseSpaceTypes[] = {"none","polynomial","harmonic","eigenvector","generalized_eigenvector","gdsw","PCMGCoarseSpaceType","PCMG_ADAPT_NONE",NULL}; 769f3fbd535SBarry Smith 7709804daf3SBarry Smith #include <petscdraw.h> 771f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer) 772f3fbd535SBarry Smith { 773f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 774f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 775e3deeeafSJed Brown PetscInt levels = mglevels ? mglevels[0]->levels : 0,i; 776219b1045SBarry Smith PetscBool iascii,isbinary,isdraw; 777f3fbd535SBarry Smith 778f3fbd535SBarry Smith PetscFunctionBegin; 7799566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 7809566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary)); 7819566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw)); 782f3fbd535SBarry Smith if (iascii) { 783e3deeeafSJed Brown const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown"; 78463a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," type is %s, levels=%" PetscInt_FMT " cycles=%s\n", PCMGTypes[mg->am],levels,cyclename)); 78531567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 78663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Cycles per PCApply=%" PetscInt_FMT "\n",mg->cyclesperpcapply)); 787f3fbd535SBarry Smith } 7882134b1e4SBarry Smith if (mg->galerkin == PC_MG_GALERKIN_BOTH) { 7899566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices\n")); 7902134b1e4SBarry Smith } else if (mg->galerkin == PC_MG_GALERKIN_PMAT) { 7919566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices for pmat\n")); 7922134b1e4SBarry Smith } else if (mg->galerkin == PC_MG_GALERKIN_MAT) { 7939566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices for mat\n")); 7942134b1e4SBarry Smith } else if (mg->galerkin == PC_MG_GALERKIN_EXTERNAL) { 7959566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Using externally compute Galerkin coarse grid matrices\n")); 7964f66f45eSBarry Smith } else { 7979566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Not using Galerkin computed coarse grid matrices\n")); 798f3fbd535SBarry Smith } 7995adeb434SBarry Smith if (mg->view) { 8009566063dSJacob Faibussowitsch PetscCall((*mg->view)(pc,viewer)); 8015adeb434SBarry Smith } 802f3fbd535SBarry Smith for (i=0; i<levels; i++) { 80363a3b9bcSJacob Faibussowitsch if (i) { 80463a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %" PetscInt_FMT " -------------------------------\n",i)); 805f3fbd535SBarry Smith } else { 80663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %" PetscInt_FMT " -------------------------------\n",i)); 807f3fbd535SBarry Smith } 8089566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 8099566063dSJacob Faibussowitsch PetscCall(KSPView(mglevels[i]->smoothd,viewer)); 8109566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 811f3fbd535SBarry Smith if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) { 8129566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n")); 813f3fbd535SBarry Smith } else if (i) { 81463a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %" PetscInt_FMT " -------------------------------\n",i)); 8159566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 8169566063dSJacob Faibussowitsch PetscCall(KSPView(mglevels[i]->smoothu,viewer)); 8179566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 818f3fbd535SBarry Smith } 81941b6fd38SMatthew G. Knepley if (i && mglevels[i]->cr) { 82063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"CR solver on level %" PetscInt_FMT " -------------------------------\n",i)); 8219566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 8229566063dSJacob Faibussowitsch PetscCall(KSPView(mglevels[i]->cr,viewer)); 8239566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 82441b6fd38SMatthew G. Knepley } 825f3fbd535SBarry Smith } 8265b0b0462SBarry Smith } else if (isbinary) { 8275b0b0462SBarry Smith for (i=levels-1; i>=0; i--) { 8289566063dSJacob Faibussowitsch PetscCall(KSPView(mglevels[i]->smoothd,viewer)); 8295b0b0462SBarry Smith if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) { 8309566063dSJacob Faibussowitsch PetscCall(KSPView(mglevels[i]->smoothu,viewer)); 8315b0b0462SBarry Smith } 8325b0b0462SBarry Smith } 833219b1045SBarry Smith } else if (isdraw) { 834219b1045SBarry Smith PetscDraw draw; 835b4375e8dSPeter Brune PetscReal x,w,y,bottom,th; 8369566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(viewer,0,&draw)); 8379566063dSJacob Faibussowitsch PetscCall(PetscDrawGetCurrentPoint(draw,&x,&y)); 8389566063dSJacob Faibussowitsch PetscCall(PetscDrawStringGetSize(draw,NULL,&th)); 839219b1045SBarry Smith bottom = y - th; 840219b1045SBarry Smith for (i=levels-1; i>=0; i--) { 841b4375e8dSPeter Brune if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) { 8429566063dSJacob Faibussowitsch PetscCall(PetscDrawPushCurrentPoint(draw,x,bottom)); 8439566063dSJacob Faibussowitsch PetscCall(KSPView(mglevels[i]->smoothd,viewer)); 8449566063dSJacob Faibussowitsch PetscCall(PetscDrawPopCurrentPoint(draw)); 845b4375e8dSPeter Brune } else { 846b4375e8dSPeter Brune w = 0.5*PetscMin(1.0-x,x); 8479566063dSJacob Faibussowitsch PetscCall(PetscDrawPushCurrentPoint(draw,x+w,bottom)); 8489566063dSJacob Faibussowitsch PetscCall(KSPView(mglevels[i]->smoothd,viewer)); 8499566063dSJacob Faibussowitsch PetscCall(PetscDrawPopCurrentPoint(draw)); 8509566063dSJacob Faibussowitsch PetscCall(PetscDrawPushCurrentPoint(draw,x-w,bottom)); 8519566063dSJacob Faibussowitsch PetscCall(KSPView(mglevels[i]->smoothu,viewer)); 8529566063dSJacob Faibussowitsch PetscCall(PetscDrawPopCurrentPoint(draw)); 853b4375e8dSPeter Brune } 8549566063dSJacob Faibussowitsch PetscCall(PetscDrawGetBoundingBox(draw,NULL,&bottom,NULL,NULL)); 8551cd381d6SBarry Smith bottom -= th; 856219b1045SBarry Smith } 8575b0b0462SBarry Smith } 858f3fbd535SBarry Smith PetscFunctionReturn(0); 859f3fbd535SBarry Smith } 860f3fbd535SBarry Smith 861af0996ceSBarry Smith #include <petsc/private/kspimpl.h> 862cab2e9ccSBarry Smith 863f3fbd535SBarry Smith /* 864f3fbd535SBarry Smith Calls setup for the KSP on each level 865f3fbd535SBarry Smith */ 866f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc) 867f3fbd535SBarry Smith { 868f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 869f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 8701a97d7b6SStefano Zampini PetscInt i,n; 87198e04984SBarry Smith PC cpc; 8723db492dfSBarry Smith PetscBool dump = PETSC_FALSE,opsset,use_amat,missinginterpolate = PETSC_FALSE; 873f3fbd535SBarry Smith Mat dA,dB; 874f3fbd535SBarry Smith Vec tvec; 875218a07d4SBarry Smith DM *dms; 8760a545947SLisandro Dalcin PetscViewer viewer = NULL; 87741b6fd38SMatthew G. Knepley PetscBool dAeqdB = PETSC_FALSE, needRestricts = PETSC_FALSE, doCR = PETSC_FALSE; 878*2b3cbbdaSStefano Zampini PetscBool adaptInterpolation = mg->adaptInterpolation; 879f3fbd535SBarry Smith 880f3fbd535SBarry Smith PetscFunctionBegin; 88128b400f6SJacob Faibussowitsch PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels with PCMGSetLevels() before setting up"); 8821a97d7b6SStefano Zampini n = mglevels[0]->levels; 88301bc414fSDmitry Karpeev /* FIX: Move this to PCSetFromOptions_MG? */ 8843aeaf226SBarry Smith if (mg->usedmfornumberoflevels) { 8853aeaf226SBarry Smith PetscInt levels; 8869566063dSJacob Faibussowitsch PetscCall(DMGetRefineLevel(pc->dm,&levels)); 8873aeaf226SBarry Smith levels++; 8883aeaf226SBarry Smith if (levels > n) { /* the problem is now being solved on a finer grid */ 8899566063dSJacob Faibussowitsch PetscCall(PCMGSetLevels(pc,levels,NULL)); 8903aeaf226SBarry Smith n = levels; 8919566063dSJacob Faibussowitsch PetscCall(PCSetFromOptions(pc)); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */ 8923aeaf226SBarry Smith mglevels = mg->levels; 8933aeaf226SBarry Smith } 8943aeaf226SBarry Smith } 8959566063dSJacob Faibussowitsch PetscCall(KSPGetPC(mglevels[0]->smoothd,&cpc)); 8963aeaf226SBarry Smith 897f3fbd535SBarry Smith /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */ 898f3fbd535SBarry Smith /* so use those from global PC */ 899f3fbd535SBarry Smith /* Is this what we always want? What if user wants to keep old one? */ 9009566063dSJacob Faibussowitsch PetscCall(KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset)); 90198e04984SBarry Smith if (opsset) { 90298e04984SBarry Smith Mat mmat; 9039566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat)); 90498e04984SBarry Smith if (mmat == pc->pmat) opsset = PETSC_FALSE; 90598e04984SBarry Smith } 9064a5f07a5SMark F. Adams 90741b6fd38SMatthew G. Knepley /* Create CR solvers */ 9089566063dSJacob Faibussowitsch PetscCall(PCMGGetAdaptCR(pc, &doCR)); 90941b6fd38SMatthew G. Knepley if (doCR) { 91041b6fd38SMatthew G. Knepley const char *prefix; 91141b6fd38SMatthew G. Knepley 9129566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &prefix)); 91341b6fd38SMatthew G. Knepley for (i = 1; i < n; ++i) { 91441b6fd38SMatthew G. Knepley PC ipc, cr; 91541b6fd38SMatthew G. Knepley char crprefix[128]; 91641b6fd38SMatthew G. Knepley 9179566063dSJacob Faibussowitsch PetscCall(KSPCreate(PetscObjectComm((PetscObject) pc), &mglevels[i]->cr)); 9189566063dSJacob Faibussowitsch PetscCall(KSPSetErrorIfNotConverged(mglevels[i]->cr, PETSC_FALSE)); 9199566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject) mglevels[i]->cr, (PetscObject) pc, n-i)); 9209566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(mglevels[i]->cr, prefix)); 9219566063dSJacob Faibussowitsch PetscCall(PetscObjectComposedDataSetInt((PetscObject) mglevels[i]->cr, PetscMGLevelId, mglevels[i]->level)); 9229566063dSJacob Faibussowitsch PetscCall(KSPSetType(mglevels[i]->cr, KSPCHEBYSHEV)); 9239566063dSJacob Faibussowitsch PetscCall(KSPSetConvergenceTest(mglevels[i]->cr, KSPConvergedSkip, NULL, NULL)); 9249566063dSJacob Faibussowitsch PetscCall(KSPSetNormType(mglevels[i]->cr, KSP_NORM_PRECONDITIONED)); 9259566063dSJacob Faibussowitsch PetscCall(KSPGetPC(mglevels[i]->cr, &ipc)); 92641b6fd38SMatthew G. Knepley 9279566063dSJacob Faibussowitsch PetscCall(PCSetType(ipc, PCCOMPOSITE)); 9289566063dSJacob Faibussowitsch PetscCall(PCCompositeSetType(ipc, PC_COMPOSITE_MULTIPLICATIVE)); 9299566063dSJacob Faibussowitsch PetscCall(PCCompositeAddPCType(ipc, PCSOR)); 9309566063dSJacob Faibussowitsch PetscCall(CreateCR_Private(pc, i, &cr)); 9319566063dSJacob Faibussowitsch PetscCall(PCCompositeAddPC(ipc, cr)); 9329566063dSJacob Faibussowitsch PetscCall(PCDestroy(&cr)); 93341b6fd38SMatthew G. Knepley 9349566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(mglevels[i]->cr, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, mg->default_smoothd)); 9359566063dSJacob Faibussowitsch PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->cr, PETSC_TRUE)); 9369566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(crprefix, 128, "mg_levels_%d_cr_", (int) i)); 9379566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(mglevels[i]->cr, crprefix)); 9389566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject) pc, (PetscObject) mglevels[i]->cr)); 93941b6fd38SMatthew G. Knepley } 94041b6fd38SMatthew G. Knepley } 94141b6fd38SMatthew G. Knepley 9424a5f07a5SMark F. Adams if (!opsset) { 9439566063dSJacob Faibussowitsch PetscCall(PCGetUseAmat(pc,&use_amat)); 9444a5f07a5SMark F. Adams if (use_amat) { 9459566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc,"Using outer operators to define finest grid operator \n because PCMGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n")); 9469566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat)); 94769858f1bSStefano Zampini } else { 9489566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc,"Using matrix (pmat) operators to define finest grid operator \n because PCMGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n")); 9499566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat)); 9504a5f07a5SMark F. Adams } 9514a5f07a5SMark F. Adams } 952f3fbd535SBarry Smith 9539c7ffb39SBarry Smith for (i=n-1; i>0; i--) { 9549c7ffb39SBarry Smith if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) { 9559c7ffb39SBarry Smith missinginterpolate = PETSC_TRUE; 956*2b3cbbdaSStefano Zampini break; 9579c7ffb39SBarry Smith } 9589c7ffb39SBarry Smith } 9592134b1e4SBarry Smith 9609566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB)); 9612134b1e4SBarry Smith if (dA == dB) dAeqdB = PETSC_TRUE; 962*2b3cbbdaSStefano Zampini if (mg->galerkin == PC_MG_GALERKIN_NONE || ((mg->galerkin == PC_MG_GALERKIN_PMAT || mg->galerkin == PC_MG_GALERKIN_MAT) && !dAeqdB)) { 9632134b1e4SBarry Smith needRestricts = PETSC_TRUE; /* user must compute either mat, pmat, or both so must restrict x to coarser levels */ 9642134b1e4SBarry Smith } 9652134b1e4SBarry Smith 966*2b3cbbdaSStefano Zampini if (pc->dm && !pc->setupcalled) { 967*2b3cbbdaSStefano Zampini /* finest smoother also gets DM but it is not active, independent of whether galerkin==PC_MG_GALERKIN_EXTERNAL */ 968*2b3cbbdaSStefano Zampini PetscCall(KSPSetDM(mglevels[n-1]->smoothd,pc->dm)); 969*2b3cbbdaSStefano Zampini PetscCall(KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE)); 970*2b3cbbdaSStefano Zampini if (mglevels[n-1]->smoothd != mglevels[n-1]->smoothu) { 971*2b3cbbdaSStefano Zampini PetscCall(KSPSetDM(mglevels[n-1]->smoothu,pc->dm)); 972*2b3cbbdaSStefano Zampini PetscCall(KSPSetDMActive(mglevels[n-1]->smoothu,PETSC_FALSE)); 973*2b3cbbdaSStefano Zampini } 974*2b3cbbdaSStefano Zampini if (mglevels[n-1]->cr) { 975*2b3cbbdaSStefano Zampini PetscCall(KSPSetDM(mglevels[n-1]->cr,pc->dm)); 976*2b3cbbdaSStefano Zampini PetscCall(KSPSetDMActive(mglevels[n-1]->cr,PETSC_FALSE)); 977*2b3cbbdaSStefano Zampini } 978*2b3cbbdaSStefano Zampini } 979*2b3cbbdaSStefano Zampini 9809c7ffb39SBarry Smith /* 9819c7ffb39SBarry Smith Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS) 982*2b3cbbdaSStefano Zampini Skipping for externally managed hierarchy (such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs? 9839c7ffb39SBarry Smith */ 9842134b1e4SBarry Smith if (missinginterpolate && pc->dm && mg->galerkin != PC_MG_GALERKIN_EXTERNAL && !pc->setupcalled) { 9852d2b81a6SBarry Smith /* construct the interpolation from the DMs */ 9862e499ae9SBarry Smith Mat p; 98773250ac0SBarry Smith Vec rscale; 9889566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n,&dms)); 989218a07d4SBarry Smith dms[n-1] = pc->dm; 990ef1267afSMatthew G. Knepley /* Separately create them so we do not get DMKSP interference between levels */ 9919566063dSJacob Faibussowitsch for (i=n-2; i>-1; i--) PetscCall(DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i])); 992218a07d4SBarry Smith for (i=n-2; i>-1; i--) { 993942e3340SBarry Smith DMKSP kdm; 994eab52d2bSLawrence Mitchell PetscBool dmhasrestrict, dmhasinject; 9959566063dSJacob Faibussowitsch PetscCall(KSPSetDM(mglevels[i]->smoothd,dms[i])); 9969566063dSJacob Faibussowitsch if (!needRestricts) PetscCall(KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE)); 997c27ee7a3SPatrick Farrell if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 9989566063dSJacob Faibussowitsch PetscCall(KSPSetDM(mglevels[i]->smoothu,dms[i])); 9999566063dSJacob Faibussowitsch if (!needRestricts) PetscCall(KSPSetDMActive(mglevels[i]->smoothu,PETSC_FALSE)); 1000c27ee7a3SPatrick Farrell } 100141b6fd38SMatthew G. Knepley if (mglevels[i]->cr) { 10029566063dSJacob Faibussowitsch PetscCall(KSPSetDM(mglevels[i]->cr,dms[i])); 10039566063dSJacob Faibussowitsch if (!needRestricts) PetscCall(KSPSetDMActive(mglevels[i]->cr,PETSC_FALSE)); 100441b6fd38SMatthew G. Knepley } 10059566063dSJacob Faibussowitsch PetscCall(DMGetDMKSPWrite(dms[i],&kdm)); 1006d1a3e677SJed Brown /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take 1007d1a3e677SJed Brown * a bitwise OR of computing the matrix, RHS, and initial iterate. */ 10080298fd71SBarry Smith kdm->ops->computerhs = NULL; 10090298fd71SBarry Smith kdm->rhsctx = NULL; 101024c3aa18SBarry Smith if (!mglevels[i+1]->interpolate) { 10119566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale)); 10129566063dSJacob Faibussowitsch PetscCall(PCMGSetInterpolation(pc,i+1,p)); 10139566063dSJacob Faibussowitsch if (rscale) PetscCall(PCMGSetRScale(pc,i+1,rscale)); 10149566063dSJacob Faibussowitsch PetscCall(VecDestroy(&rscale)); 10159566063dSJacob Faibussowitsch PetscCall(MatDestroy(&p)); 1016218a07d4SBarry Smith } 10179566063dSJacob Faibussowitsch PetscCall(DMHasCreateRestriction(dms[i],&dmhasrestrict)); 10183ad4599aSBarry Smith if (dmhasrestrict && !mglevels[i+1]->restrct) { 10199566063dSJacob Faibussowitsch PetscCall(DMCreateRestriction(dms[i],dms[i+1],&p)); 10209566063dSJacob Faibussowitsch PetscCall(PCMGSetRestriction(pc,i+1,p)); 10219566063dSJacob Faibussowitsch PetscCall(MatDestroy(&p)); 10223ad4599aSBarry Smith } 10239566063dSJacob Faibussowitsch PetscCall(DMHasCreateInjection(dms[i],&dmhasinject)); 1024eab52d2bSLawrence Mitchell if (dmhasinject && !mglevels[i+1]->inject) { 10259566063dSJacob Faibussowitsch PetscCall(DMCreateInjection(dms[i],dms[i+1],&p)); 10269566063dSJacob Faibussowitsch PetscCall(PCMGSetInjection(pc,i+1,p)); 10279566063dSJacob Faibussowitsch PetscCall(MatDestroy(&p)); 1028eab52d2bSLawrence Mitchell } 102924c3aa18SBarry Smith } 10302d2b81a6SBarry Smith 10319566063dSJacob Faibussowitsch for (i=n-2; i>-1; i--) PetscCall(DMDestroy(&dms[i])); 10329566063dSJacob Faibussowitsch PetscCall(PetscFree(dms)); 1033ce4cda84SJed Brown } 1034cab2e9ccSBarry Smith 10352134b1e4SBarry Smith if (mg->galerkin < PC_MG_GALERKIN_NONE) { 10362134b1e4SBarry Smith Mat A,B; 10372134b1e4SBarry Smith PetscBool doA = PETSC_FALSE,doB = PETSC_FALSE; 10382134b1e4SBarry Smith MatReuse reuse = MAT_INITIAL_MATRIX; 10392134b1e4SBarry Smith 1040*2b3cbbdaSStefano Zampini if (mg->galerkin == PC_MG_GALERKIN_PMAT || mg->galerkin == PC_MG_GALERKIN_BOTH) doB = PETSC_TRUE; 1041*2b3cbbdaSStefano Zampini if (mg->galerkin == PC_MG_GALERKIN_MAT || (mg->galerkin == PC_MG_GALERKIN_BOTH && dA != dB)) doA = PETSC_TRUE; 10422134b1e4SBarry Smith if (pc->setupcalled) reuse = MAT_REUSE_MATRIX; 1043f3fbd535SBarry Smith for (i=n-2; i>-1; i--) { 1044*2b3cbbdaSStefano Zampini if (!mglevels[i+1]->restrct && !mglevels[i+1]->interpolate) { /* see if we can compute a coarse space */ 1045*2b3cbbdaSStefano Zampini PetscCall(PCMGComputeCoarseSpace_Internal(pc, i+1, mg->coarseSpaceType, mg->Nc, NULL, &mglevels[i+1]->coarseSpace)); 1046*2b3cbbdaSStefano Zampini PetscCall(PCMGSetInterpolation(pc,i+1,mglevels[i+1]->coarseSpace)); 1047*2b3cbbdaSStefano Zampini } 10487827d75bSBarry Smith PetscCheck(mglevels[i+1]->restrct || mglevels[i+1]->interpolate,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must provide interpolation or restriction for each MG level except level 0"); 1049b9367914SBarry Smith if (!mglevels[i+1]->interpolate) { 10509566063dSJacob Faibussowitsch PetscCall(PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct)); 1051b9367914SBarry Smith } 1052b9367914SBarry Smith if (!mglevels[i+1]->restrct) { 10539566063dSJacob Faibussowitsch PetscCall(PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate)); 1054b9367914SBarry Smith } 10552134b1e4SBarry Smith if (reuse == MAT_REUSE_MATRIX) { 10569566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[i]->smoothd,&A,&B)); 10572134b1e4SBarry Smith } 10582134b1e4SBarry Smith if (doA) { 10599566063dSJacob Faibussowitsch PetscCall(MatGalerkin(mglevels[i+1]->restrct,dA,mglevels[i+1]->interpolate,reuse,1.0,&A)); 10602134b1e4SBarry Smith } 10612134b1e4SBarry Smith if (doB) { 10629566063dSJacob Faibussowitsch PetscCall(MatGalerkin(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,reuse,1.0,&B)); 1063b9367914SBarry Smith } 10642134b1e4SBarry Smith /* the management of the PetscObjectReference() and PetscObjecDereference() below is rather delicate */ 10652134b1e4SBarry Smith if (!doA && dAeqdB) { 10669566063dSJacob Faibussowitsch if (reuse == MAT_INITIAL_MATRIX) PetscCall(PetscObjectReference((PetscObject)B)); 10672134b1e4SBarry Smith A = B; 10682134b1e4SBarry Smith } else if (!doA && reuse == MAT_INITIAL_MATRIX) { 10699566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[i]->smoothd,&A,NULL)); 10709566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)A)); 1071b9367914SBarry Smith } 10722134b1e4SBarry Smith if (!doB && dAeqdB) { 10739566063dSJacob Faibussowitsch if (reuse == MAT_INITIAL_MATRIX) PetscCall(PetscObjectReference((PetscObject)A)); 10742134b1e4SBarry Smith B = A; 10752134b1e4SBarry Smith } else if (!doB && reuse == MAT_INITIAL_MATRIX) { 10769566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[i]->smoothd,NULL,&B)); 10779566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)B)); 10787d537d92SJed Brown } 10792134b1e4SBarry Smith if (reuse == MAT_INITIAL_MATRIX) { 10809566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(mglevels[i]->smoothd,A,B)); 10819566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)A)); 10829566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)B)); 10832134b1e4SBarry Smith } 10842134b1e4SBarry Smith dA = A; 1085f3fbd535SBarry Smith dB = B; 1086f3fbd535SBarry Smith } 1087f3fbd535SBarry Smith } 1088f3b08a26SMatthew G. Knepley 1089f3b08a26SMatthew G. Knepley /* Adapt interpolation matrices */ 1090*2b3cbbdaSStefano Zampini if (adaptInterpolation) { 1091f3b08a26SMatthew G. Knepley for (i = 0; i < n; ++i) { 1092*2b3cbbdaSStefano Zampini if (!mglevels[i]->coarseSpace) { 10939566063dSJacob Faibussowitsch PetscCall(PCMGComputeCoarseSpace_Internal(pc, i, mg->coarseSpaceType, mg->Nc, !i ? NULL : mglevels[i-1]->coarseSpace, &mglevels[i]->coarseSpace)); 1094*2b3cbbdaSStefano Zampini } 1095*2b3cbbdaSStefano Zampini if (i) PetscCall(PCMGAdaptInterpolator_Internal(pc, i, mglevels[i-1]->smoothu, mglevels[i]->smoothu, mglevels[i-1]->coarseSpace, mglevels[i]->coarseSpace)); 1096f3b08a26SMatthew G. Knepley } 1097f3b08a26SMatthew G. Knepley for (i = n-2; i > -1; --i) { 10989566063dSJacob Faibussowitsch PetscCall(PCMGRecomputeLevelOperators_Internal(pc, i)); 1099f3b08a26SMatthew G. Knepley } 1100f3b08a26SMatthew G. Knepley } 1101f3b08a26SMatthew G. Knepley 11022134b1e4SBarry Smith if (needRestricts && pc->dm) { 1103caa4e7f2SJed Brown for (i=n-2; i>=0; i--) { 1104ccceb50cSJed Brown DM dmfine,dmcoarse; 1105caa4e7f2SJed Brown Mat Restrict,Inject; 1106caa4e7f2SJed Brown Vec rscale; 11079566063dSJacob Faibussowitsch PetscCall(KSPGetDM(mglevels[i+1]->smoothd,&dmfine)); 11089566063dSJacob Faibussowitsch PetscCall(KSPGetDM(mglevels[i]->smoothd,&dmcoarse)); 11099566063dSJacob Faibussowitsch PetscCall(PCMGGetRestriction(pc,i+1,&Restrict)); 11109566063dSJacob Faibussowitsch PetscCall(PCMGGetRScale(pc,i+1,&rscale)); 11119566063dSJacob Faibussowitsch PetscCall(PCMGGetInjection(pc,i+1,&Inject)); 11129566063dSJacob Faibussowitsch PetscCall(DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse)); 1113caa4e7f2SJed Brown } 1114caa4e7f2SJed Brown } 1115f3fbd535SBarry Smith 1116f3fbd535SBarry Smith if (!pc->setupcalled) { 1117f3fbd535SBarry Smith for (i=0; i<n; i++) { 11189566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(mglevels[i]->smoothd)); 1119f3fbd535SBarry Smith } 1120f3fbd535SBarry Smith for (i=1; i<n; i++) { 1121f3fbd535SBarry Smith if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) { 11229566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(mglevels[i]->smoothu)); 1123f3fbd535SBarry Smith } 112441b6fd38SMatthew G. Knepley if (mglevels[i]->cr) { 11259566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(mglevels[i]->cr)); 112641b6fd38SMatthew G. Knepley } 1127f3fbd535SBarry Smith } 11283ad4599aSBarry Smith /* insure that if either interpolation or restriction is set the other other one is set */ 1129f3fbd535SBarry Smith for (i=1; i<n; i++) { 11309566063dSJacob Faibussowitsch PetscCall(PCMGGetInterpolation(pc,i,NULL)); 11319566063dSJacob Faibussowitsch PetscCall(PCMGGetRestriction(pc,i,NULL)); 1132f3fbd535SBarry Smith } 1133f3fbd535SBarry Smith for (i=0; i<n-1; i++) { 1134f3fbd535SBarry Smith if (!mglevels[i]->b) { 1135f3fbd535SBarry Smith Vec *vec; 11369566063dSJacob Faibussowitsch PetscCall(KSPCreateVecs(mglevels[i]->smoothd,1,&vec,0,NULL)); 11379566063dSJacob Faibussowitsch PetscCall(PCMGSetRhs(pc,i,*vec)); 11389566063dSJacob Faibussowitsch PetscCall(VecDestroy(vec)); 11399566063dSJacob Faibussowitsch PetscCall(PetscFree(vec)); 1140f3fbd535SBarry Smith } 1141f3fbd535SBarry Smith if (!mglevels[i]->r && i) { 11429566063dSJacob Faibussowitsch PetscCall(VecDuplicate(mglevels[i]->b,&tvec)); 11439566063dSJacob Faibussowitsch PetscCall(PCMGSetR(pc,i,tvec)); 11449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tvec)); 1145f3fbd535SBarry Smith } 1146f3fbd535SBarry Smith if (!mglevels[i]->x) { 11479566063dSJacob Faibussowitsch PetscCall(VecDuplicate(mglevels[i]->b,&tvec)); 11489566063dSJacob Faibussowitsch PetscCall(PCMGSetX(pc,i,tvec)); 11499566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tvec)); 1150f3fbd535SBarry Smith } 115141b6fd38SMatthew G. Knepley if (doCR) { 11529566063dSJacob Faibussowitsch PetscCall(VecDuplicate(mglevels[i]->b,&mglevels[i]->crx)); 11539566063dSJacob Faibussowitsch PetscCall(VecDuplicate(mglevels[i]->b,&mglevels[i]->crb)); 11549566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(mglevels[i]->crb)); 115541b6fd38SMatthew G. Knepley } 1156f3fbd535SBarry Smith } 1157f3fbd535SBarry Smith if (n != 1 && !mglevels[n-1]->r) { 1158f3fbd535SBarry Smith /* PCMGSetR() on the finest level if user did not supply it */ 1159f3fbd535SBarry Smith Vec *vec; 11609566063dSJacob Faibussowitsch PetscCall(KSPCreateVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL)); 11619566063dSJacob Faibussowitsch PetscCall(PCMGSetR(pc,n-1,*vec)); 11629566063dSJacob Faibussowitsch PetscCall(VecDestroy(vec)); 11639566063dSJacob Faibussowitsch PetscCall(PetscFree(vec)); 1164f3fbd535SBarry Smith } 116541b6fd38SMatthew G. Knepley if (doCR) { 11669566063dSJacob Faibussowitsch PetscCall(VecDuplicate(mglevels[n-1]->r, &mglevels[n-1]->crx)); 11679566063dSJacob Faibussowitsch PetscCall(VecDuplicate(mglevels[n-1]->r, &mglevels[n-1]->crb)); 11689566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(mglevels[n-1]->crb)); 116941b6fd38SMatthew G. Knepley } 1170f3fbd535SBarry Smith } 1171f3fbd535SBarry Smith 117298e04984SBarry Smith if (pc->dm) { 117398e04984SBarry Smith /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */ 117498e04984SBarry Smith for (i=0; i<n-1; i++) { 117598e04984SBarry Smith if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX; 117698e04984SBarry Smith } 117798e04984SBarry Smith } 117808549ed4SJed Brown // We got here (PCSetUp_MG) because the matrix has changed, which means the smoother needs to be set up again (e.g., 117908549ed4SJed Brown // new diagonal for Jacobi). Setting it here allows it to be logged under PCSetUp rather than deep inside a PCApply. 118008549ed4SJed Brown if (mglevels[n-1]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[n-1]->smoothd->setupstage = KSP_SETUP_NEWMATRIX; 1181f3fbd535SBarry Smith 1182f3fbd535SBarry Smith for (i=1; i<n; i++) { 11832a21e185SBarry Smith if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1) { 1184f3fbd535SBarry Smith /* if doing only down then initial guess is zero */ 11859566063dSJacob Faibussowitsch PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE)); 1186f3fbd535SBarry Smith } 11879566063dSJacob Faibussowitsch if (mglevels[i]->cr) PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->cr,PETSC_TRUE)); 11889566063dSJacob Faibussowitsch if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0)); 11899566063dSJacob Faibussowitsch PetscCall(KSPSetUp(mglevels[i]->smoothd)); 1190c0decd05SBarry Smith if (mglevels[i]->smoothd->reason == KSP_DIVERGED_PC_FAILED) { 1191899639b0SHong Zhang pc->failedreason = PC_SUBPC_ERROR; 1192899639b0SHong Zhang } 11939566063dSJacob Faibussowitsch if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0)); 1194d42688cbSBarry Smith if (!mglevels[i]->residual) { 1195d42688cbSBarry Smith Mat mat; 11969566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[i]->smoothd,&mat,NULL)); 11979566063dSJacob Faibussowitsch PetscCall(PCMGSetResidual(pc,i,PCMGResidualDefault,mat)); 1198d42688cbSBarry Smith } 1199fcb023d4SJed Brown if (!mglevels[i]->residualtranspose) { 1200fcb023d4SJed Brown Mat mat; 12019566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[i]->smoothd,&mat,NULL)); 12029566063dSJacob Faibussowitsch PetscCall(PCMGSetResidualTranspose(pc,i,PCMGResidualTransposeDefault,mat)); 1203fcb023d4SJed Brown } 1204f3fbd535SBarry Smith } 1205f3fbd535SBarry Smith for (i=1; i<n; i++) { 1206f3fbd535SBarry Smith if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) { 1207f3fbd535SBarry Smith Mat downmat,downpmat; 1208f3fbd535SBarry Smith 1209f3fbd535SBarry Smith /* check if operators have been set for up, if not use down operators to set them */ 12109566063dSJacob Faibussowitsch PetscCall(KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL)); 1211f3fbd535SBarry Smith if (!opsset) { 12129566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat)); 12139566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat)); 1214f3fbd535SBarry Smith } 1215f3fbd535SBarry Smith 12169566063dSJacob Faibussowitsch PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE)); 12179566063dSJacob Faibussowitsch if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0)); 12189566063dSJacob Faibussowitsch PetscCall(KSPSetUp(mglevels[i]->smoothu)); 1219c0decd05SBarry Smith if (mglevels[i]->smoothu->reason == KSP_DIVERGED_PC_FAILED) { 1220899639b0SHong Zhang pc->failedreason = PC_SUBPC_ERROR; 1221899639b0SHong Zhang } 12229566063dSJacob Faibussowitsch if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0)); 1223f3fbd535SBarry Smith } 122441b6fd38SMatthew G. Knepley if (mglevels[i]->cr) { 122541b6fd38SMatthew G. Knepley Mat downmat,downpmat; 122641b6fd38SMatthew G. Knepley 122741b6fd38SMatthew G. Knepley /* check if operators have been set for up, if not use down operators to set them */ 12289566063dSJacob Faibussowitsch PetscCall(KSPGetOperatorsSet(mglevels[i]->cr,&opsset,NULL)); 122941b6fd38SMatthew G. Knepley if (!opsset) { 12309566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat)); 12319566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(mglevels[i]->cr,downmat,downpmat)); 123241b6fd38SMatthew G. Knepley } 123341b6fd38SMatthew G. Knepley 12349566063dSJacob Faibussowitsch PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->cr,PETSC_TRUE)); 12359566063dSJacob Faibussowitsch if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0)); 12369566063dSJacob Faibussowitsch PetscCall(KSPSetUp(mglevels[i]->cr)); 123741b6fd38SMatthew G. Knepley if (mglevels[i]->cr->reason == KSP_DIVERGED_PC_FAILED) { 123841b6fd38SMatthew G. Knepley pc->failedreason = PC_SUBPC_ERROR; 123941b6fd38SMatthew G. Knepley } 12409566063dSJacob Faibussowitsch if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0)); 124141b6fd38SMatthew G. Knepley } 1242f3fbd535SBarry Smith } 1243f3fbd535SBarry Smith 12449566063dSJacob Faibussowitsch if (mglevels[0]->eventsmoothsetup) PetscCall(PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0)); 12459566063dSJacob Faibussowitsch PetscCall(KSPSetUp(mglevels[0]->smoothd)); 1246c0decd05SBarry Smith if (mglevels[0]->smoothd->reason == KSP_DIVERGED_PC_FAILED) { 1247899639b0SHong Zhang pc->failedreason = PC_SUBPC_ERROR; 1248899639b0SHong Zhang } 12499566063dSJacob Faibussowitsch if (mglevels[0]->eventsmoothsetup) PetscCall(PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0)); 1250f3fbd535SBarry Smith 1251f3fbd535SBarry Smith /* 1252f3fbd535SBarry Smith Dump the interpolation/restriction matrices plus the 1253e3c5b3baSBarry Smith Jacobian/stiffness on each level. This allows MATLAB users to 1254f3fbd535SBarry Smith easily check if the Galerkin condition A_c = R A_f R^T is satisfied. 1255f3fbd535SBarry Smith 1256f3fbd535SBarry Smith Only support one or the other at the same time. 1257f3fbd535SBarry Smith */ 1258f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER) 12599566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL)); 1260ce94432eSBarry Smith if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc)); 1261f3fbd535SBarry Smith dump = PETSC_FALSE; 1262f3fbd535SBarry Smith #endif 12639566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL)); 1264ce94432eSBarry Smith if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc)); 1265f3fbd535SBarry Smith 1266f3fbd535SBarry Smith if (viewer) { 1267f3fbd535SBarry Smith for (i=1; i<n; i++) { 12689566063dSJacob Faibussowitsch PetscCall(MatView(mglevels[i]->restrct,viewer)); 1269f3fbd535SBarry Smith } 1270f3fbd535SBarry Smith for (i=0; i<n; i++) { 12719566063dSJacob Faibussowitsch PetscCall(KSPGetPC(mglevels[i]->smoothd,&pc)); 12729566063dSJacob Faibussowitsch PetscCall(MatView(pc->mat,viewer)); 1273f3fbd535SBarry Smith } 1274f3fbd535SBarry Smith } 1275f3fbd535SBarry Smith PetscFunctionReturn(0); 1276f3fbd535SBarry Smith } 1277f3fbd535SBarry Smith 1278f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/ 1279f3fbd535SBarry Smith 128097d33e41SMatthew G. Knepley PetscErrorCode PCMGGetLevels_MG(PC pc, PetscInt *levels) 128197d33e41SMatthew G. Knepley { 128297d33e41SMatthew G. Knepley PC_MG *mg = (PC_MG *) pc->data; 128397d33e41SMatthew G. Knepley 128497d33e41SMatthew G. Knepley PetscFunctionBegin; 128597d33e41SMatthew G. Knepley *levels = mg->nlevels; 128697d33e41SMatthew G. Knepley PetscFunctionReturn(0); 128797d33e41SMatthew G. Knepley } 128897d33e41SMatthew G. Knepley 12894b9ad928SBarry Smith /*@ 129097177400SBarry Smith PCMGGetLevels - Gets the number of levels to use with MG. 12914b9ad928SBarry Smith 12924b9ad928SBarry Smith Not Collective 12934b9ad928SBarry Smith 12944b9ad928SBarry Smith Input Parameter: 12954b9ad928SBarry Smith . pc - the preconditioner context 12964b9ad928SBarry Smith 12974b9ad928SBarry Smith Output parameter: 12984b9ad928SBarry Smith . levels - the number of levels 12994b9ad928SBarry Smith 13004b9ad928SBarry Smith Level: advanced 13014b9ad928SBarry Smith 1302db781477SPatrick Sanan .seealso: `PCMGSetLevels()` 13034b9ad928SBarry Smith @*/ 13047087cfbeSBarry Smith PetscErrorCode PCMGGetLevels(PC pc,PetscInt *levels) 13054b9ad928SBarry Smith { 13064b9ad928SBarry Smith PetscFunctionBegin; 13070700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 13084482741eSBarry Smith PetscValidIntPointer(levels,2); 130997d33e41SMatthew G. Knepley *levels = 0; 1310cac4c232SBarry Smith PetscTryMethod(pc,"PCMGGetLevels_C",(PC,PetscInt*),(pc,levels)); 13114b9ad928SBarry Smith PetscFunctionReturn(0); 13124b9ad928SBarry Smith } 13134b9ad928SBarry Smith 13144b9ad928SBarry Smith /*@ 1315e7d4b4cbSMark Adams PCMGGetGridComplexity - compute operator and grid complexity of MG hierarchy 1316e7d4b4cbSMark Adams 1317e7d4b4cbSMark Adams Input Parameter: 1318e7d4b4cbSMark Adams . pc - the preconditioner context 1319e7d4b4cbSMark Adams 132090db8557SMark Adams Output Parameters: 132190db8557SMark Adams + gc - grid complexity = sum_i(n_i) / n_0 132290db8557SMark Adams - oc - operator complexity = sum_i(nnz_i) / nnz_0 1323e7d4b4cbSMark Adams 1324e7d4b4cbSMark Adams Level: advanced 1325e7d4b4cbSMark Adams 1326db781477SPatrick Sanan .seealso: `PCMGGetLevels()` 1327e7d4b4cbSMark Adams @*/ 1328e7d4b4cbSMark Adams PetscErrorCode PCMGGetGridComplexity(PC pc, PetscReal *gc, PetscReal *oc) 1329e7d4b4cbSMark Adams { 1330e7d4b4cbSMark Adams PC_MG *mg = (PC_MG*)pc->data; 1331e7d4b4cbSMark Adams PC_MG_Levels **mglevels = mg->levels; 1332e7d4b4cbSMark Adams PetscInt lev,N; 1333e7d4b4cbSMark Adams PetscLogDouble nnz0 = 0, sgc = 0, soc = 0, n0 = 0; 1334e7d4b4cbSMark Adams MatInfo info; 1335e7d4b4cbSMark Adams 1336e7d4b4cbSMark Adams PetscFunctionBegin; 133790db8557SMark Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 133890db8557SMark Adams if (gc) PetscValidRealPointer(gc,2); 133990db8557SMark Adams if (oc) PetscValidRealPointer(oc,3); 1340e7d4b4cbSMark Adams if (!pc->setupcalled) { 134190db8557SMark Adams if (gc) *gc = 0; 134290db8557SMark Adams if (oc) *oc = 0; 1343e7d4b4cbSMark Adams PetscFunctionReturn(0); 1344e7d4b4cbSMark Adams } 1345e7d4b4cbSMark Adams PetscCheck(mg->nlevels > 0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MG has no levels"); 1346e7d4b4cbSMark Adams for (lev=0; lev<mg->nlevels; lev++) { 1347e7d4b4cbSMark Adams Mat dB; 13489566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[lev]->smoothd,NULL,&dB)); 13499566063dSJacob Faibussowitsch PetscCall(MatGetInfo(dB,MAT_GLOBAL_SUM,&info)); /* global reduction */ 13509566063dSJacob Faibussowitsch PetscCall(MatGetSize(dB,&N,NULL)); 1351e7d4b4cbSMark Adams sgc += N; 1352e7d4b4cbSMark Adams soc += info.nz_used; 1353e7d4b4cbSMark Adams if (lev==mg->nlevels-1) {nnz0 = info.nz_used; n0 = N;} 1354e7d4b4cbSMark Adams } 135590db8557SMark Adams if (n0 > 0 && gc) *gc = (PetscReal)(sgc/n0); 1356e7d4b4cbSMark Adams else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Number for grid points on finest level is not available"); 135790db8557SMark Adams if (nnz0 > 0 && oc) *oc = (PetscReal)(soc/nnz0); 1358e7d4b4cbSMark Adams PetscFunctionReturn(0); 1359e7d4b4cbSMark Adams } 1360e7d4b4cbSMark Adams 1361e7d4b4cbSMark Adams /*@ 136297177400SBarry Smith PCMGSetType - Determines the form of multigrid to use: 13634b9ad928SBarry Smith multiplicative, additive, full, or the Kaskade algorithm. 13644b9ad928SBarry Smith 1365ad4df100SBarry Smith Logically Collective on PC 13664b9ad928SBarry Smith 13674b9ad928SBarry Smith Input Parameters: 13684b9ad928SBarry Smith + pc - the preconditioner context 13699dcbbd2bSBarry Smith - form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE, 13709dcbbd2bSBarry Smith PC_MG_FULL, PC_MG_KASKADE 13714b9ad928SBarry Smith 13724b9ad928SBarry Smith Options Database Key: 13734b9ad928SBarry Smith . -pc_mg_type <form> - Sets <form>, one of multiplicative, 13744b9ad928SBarry Smith additive, full, kaskade 13754b9ad928SBarry Smith 13764b9ad928SBarry Smith Level: advanced 13774b9ad928SBarry Smith 1378db781477SPatrick Sanan .seealso: `PCMGSetLevels()` 13794b9ad928SBarry Smith @*/ 13807087cfbeSBarry Smith PetscErrorCode PCMGSetType(PC pc,PCMGType form) 13814b9ad928SBarry Smith { 1382f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 13834b9ad928SBarry Smith 13844b9ad928SBarry Smith PetscFunctionBegin; 13850700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1386c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc,form,2); 138731567311SBarry Smith mg->am = form; 13889dcbbd2bSBarry Smith if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG; 1389c60c7ad4SBarry Smith else pc->ops->applyrichardson = NULL; 1390c60c7ad4SBarry Smith PetscFunctionReturn(0); 1391c60c7ad4SBarry Smith } 1392c60c7ad4SBarry Smith 1393c60c7ad4SBarry Smith /*@ 1394c60c7ad4SBarry Smith PCMGGetType - Determines the form of multigrid to use: 1395c60c7ad4SBarry Smith multiplicative, additive, full, or the Kaskade algorithm. 1396c60c7ad4SBarry Smith 1397c60c7ad4SBarry Smith Logically Collective on PC 1398c60c7ad4SBarry Smith 1399c60c7ad4SBarry Smith Input Parameter: 1400c60c7ad4SBarry Smith . pc - the preconditioner context 1401c60c7ad4SBarry Smith 1402c60c7ad4SBarry Smith Output Parameter: 1403c60c7ad4SBarry Smith . type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE 1404c60c7ad4SBarry Smith 1405c60c7ad4SBarry Smith Level: advanced 1406c60c7ad4SBarry Smith 1407db781477SPatrick Sanan .seealso: `PCMGSetLevels()` 1408c60c7ad4SBarry Smith @*/ 1409c60c7ad4SBarry Smith PetscErrorCode PCMGGetType(PC pc,PCMGType *type) 1410c60c7ad4SBarry Smith { 1411c60c7ad4SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1412c60c7ad4SBarry Smith 1413c60c7ad4SBarry Smith PetscFunctionBegin; 1414c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1415c60c7ad4SBarry Smith *type = mg->am; 14164b9ad928SBarry Smith PetscFunctionReturn(0); 14174b9ad928SBarry Smith } 14184b9ad928SBarry Smith 14194b9ad928SBarry Smith /*@ 14200d353602SBarry Smith PCMGSetCycleType - Sets the type cycles to use. Use PCMGSetCycleTypeOnLevel() for more 14214b9ad928SBarry Smith complicated cycling. 14224b9ad928SBarry Smith 1423ad4df100SBarry Smith Logically Collective on PC 14244b9ad928SBarry Smith 14254b9ad928SBarry Smith Input Parameters: 1426c2be2410SBarry Smith + pc - the multigrid context 1427c1cbb1deSBarry Smith - n - either PC_MG_CYCLE_V or PC_MG_CYCLE_W 14284b9ad928SBarry Smith 14294b9ad928SBarry Smith Options Database Key: 1430c1cbb1deSBarry Smith . -pc_mg_cycle_type <v,w> - provide the cycle desired 14314b9ad928SBarry Smith 14324b9ad928SBarry Smith Level: advanced 14334b9ad928SBarry Smith 1434db781477SPatrick Sanan .seealso: `PCMGSetCycleTypeOnLevel()` 14354b9ad928SBarry Smith @*/ 14367087cfbeSBarry Smith PetscErrorCode PCMGSetCycleType(PC pc,PCMGCycleType n) 14374b9ad928SBarry Smith { 1438f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1439f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 144079416396SBarry Smith PetscInt i,levels; 14414b9ad928SBarry Smith 14424b9ad928SBarry Smith PetscFunctionBegin; 14430700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 144469fbec6eSBarry Smith PetscValidLogicalCollectiveEnum(pc,n,2); 144528b400f6SJacob Faibussowitsch PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling"); 1446f3fbd535SBarry Smith levels = mglevels[0]->levels; 14472fa5cd67SKarl Rupp for (i=0; i<levels; i++) mglevels[i]->cycles = n; 14484b9ad928SBarry Smith PetscFunctionReturn(0); 14494b9ad928SBarry Smith } 14504b9ad928SBarry Smith 14518cc2d5dfSBarry Smith /*@ 14528cc2d5dfSBarry Smith PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step 14538cc2d5dfSBarry Smith of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used 14548cc2d5dfSBarry Smith 1455ad4df100SBarry Smith Logically Collective on PC 14568cc2d5dfSBarry Smith 14578cc2d5dfSBarry Smith Input Parameters: 14588cc2d5dfSBarry Smith + pc - the multigrid context 14598cc2d5dfSBarry Smith - n - number of cycles (default is 1) 14608cc2d5dfSBarry Smith 14618cc2d5dfSBarry Smith Options Database Key: 1462147403d9SBarry Smith . -pc_mg_multiplicative_cycles n - set the number of cycles 14638cc2d5dfSBarry Smith 14648cc2d5dfSBarry Smith Level: advanced 14658cc2d5dfSBarry Smith 146695452b02SPatrick Sanan Notes: 146795452b02SPatrick Sanan This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType() 14688cc2d5dfSBarry Smith 1469db781477SPatrick Sanan .seealso: `PCMGSetCycleTypeOnLevel()`, `PCMGSetCycleType()` 14708cc2d5dfSBarry Smith @*/ 14717087cfbeSBarry Smith PetscErrorCode PCMGMultiplicativeSetCycles(PC pc,PetscInt n) 14728cc2d5dfSBarry Smith { 1473f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 14748cc2d5dfSBarry Smith 14758cc2d5dfSBarry Smith PetscFunctionBegin; 14760700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1477c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 14782a21e185SBarry Smith mg->cyclesperpcapply = n; 14798cc2d5dfSBarry Smith PetscFunctionReturn(0); 14808cc2d5dfSBarry Smith } 14818cc2d5dfSBarry Smith 14822134b1e4SBarry Smith PetscErrorCode PCMGSetGalerkin_MG(PC pc,PCMGGalerkinType use) 1483fb15c04eSBarry Smith { 1484fb15c04eSBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1485fb15c04eSBarry Smith 1486fb15c04eSBarry Smith PetscFunctionBegin; 14872134b1e4SBarry Smith mg->galerkin = use; 1488fb15c04eSBarry Smith PetscFunctionReturn(0); 1489fb15c04eSBarry Smith } 1490fb15c04eSBarry Smith 1491c2be2410SBarry Smith /*@ 149297177400SBarry Smith PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the 149382b87d87SMatthew G. Knepley finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i 1494c2be2410SBarry Smith 1495ad4df100SBarry Smith Logically Collective on PC 1496c2be2410SBarry Smith 1497c2be2410SBarry Smith Input Parameters: 1498c91913d3SJed Brown + pc - the multigrid context 14992134b1e4SBarry Smith - use - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, or PC_MG_GALERKIN_NONE 1500c2be2410SBarry Smith 1501c2be2410SBarry Smith Options Database Key: 1502147403d9SBarry Smith . -pc_mg_galerkin <both,pmat,mat,none> - set the matrices to form via the Galerkin process 1503c2be2410SBarry Smith 1504c2be2410SBarry Smith Level: intermediate 1505c2be2410SBarry Smith 150695452b02SPatrick Sanan Notes: 150795452b02SPatrick Sanan Some codes that use PCMG such as PCGAMG use Galerkin internally while constructing the hierarchy and thus do not 15081c1aac46SBarry Smith use the PCMG construction of the coarser grids. 15091c1aac46SBarry Smith 1510db781477SPatrick Sanan .seealso: `PCMGGetGalerkin()`, `PCMGGalerkinType` 15113fc8bf9cSBarry Smith 1512c2be2410SBarry Smith @*/ 15132134b1e4SBarry Smith PetscErrorCode PCMGSetGalerkin(PC pc,PCMGGalerkinType use) 1514c2be2410SBarry Smith { 1515c2be2410SBarry Smith PetscFunctionBegin; 15160700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1517cac4c232SBarry Smith PetscTryMethod(pc,"PCMGSetGalerkin_C",(PC,PCMGGalerkinType),(pc,use)); 1518c2be2410SBarry Smith PetscFunctionReturn(0); 1519c2be2410SBarry Smith } 1520c2be2410SBarry Smith 15213fc8bf9cSBarry Smith /*@ 15223fc8bf9cSBarry Smith PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. 152382b87d87SMatthew G. Knepley A_i-1 = r_i * A_i * p_i 15243fc8bf9cSBarry Smith 15253fc8bf9cSBarry Smith Not Collective 15263fc8bf9cSBarry Smith 15273fc8bf9cSBarry Smith Input Parameter: 15283fc8bf9cSBarry Smith . pc - the multigrid context 15293fc8bf9cSBarry Smith 15303fc8bf9cSBarry Smith Output Parameter: 15312134b1e4SBarry 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 15323fc8bf9cSBarry Smith 15333fc8bf9cSBarry Smith Level: intermediate 15343fc8bf9cSBarry Smith 1535db781477SPatrick Sanan .seealso: `PCMGSetGalerkin()`, `PCMGGalerkinType` 15363fc8bf9cSBarry Smith 15373fc8bf9cSBarry Smith @*/ 15382134b1e4SBarry Smith PetscErrorCode PCMGGetGalerkin(PC pc,PCMGGalerkinType *galerkin) 15393fc8bf9cSBarry Smith { 1540f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 15413fc8bf9cSBarry Smith 15423fc8bf9cSBarry Smith PetscFunctionBegin; 15430700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 15442134b1e4SBarry Smith *galerkin = mg->galerkin; 15453fc8bf9cSBarry Smith PetscFunctionReturn(0); 15463fc8bf9cSBarry Smith } 15473fc8bf9cSBarry Smith 1548f3b08a26SMatthew G. Knepley PetscErrorCode PCMGSetAdaptInterpolation_MG(PC pc, PetscBool adapt) 1549f3b08a26SMatthew G. Knepley { 1550f3b08a26SMatthew G. Knepley PC_MG *mg = (PC_MG *) pc->data; 1551f3b08a26SMatthew G. Knepley 1552f3b08a26SMatthew G. Knepley PetscFunctionBegin; 1553f3b08a26SMatthew G. Knepley mg->adaptInterpolation = adapt; 1554f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 1555f3b08a26SMatthew G. Knepley } 1556f3b08a26SMatthew G. Knepley 1557f3b08a26SMatthew G. Knepley PetscErrorCode PCMGGetAdaptInterpolation_MG(PC pc, PetscBool *adapt) 1558f3b08a26SMatthew G. Knepley { 1559f3b08a26SMatthew G. Knepley PC_MG *mg = (PC_MG *) pc->data; 1560f3b08a26SMatthew G. Knepley 1561f3b08a26SMatthew G. Knepley PetscFunctionBegin; 1562f3b08a26SMatthew G. Knepley *adapt = mg->adaptInterpolation; 1563f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 1564f3b08a26SMatthew G. Knepley } 1565f3b08a26SMatthew G. Knepley 1566*2b3cbbdaSStefano Zampini PetscErrorCode PCMGSetAdaptCoarseSpaceType_MG(PC pc, PCMGCoarseSpaceType ctype) 1567*2b3cbbdaSStefano Zampini { 1568*2b3cbbdaSStefano Zampini PC_MG *mg = (PC_MG *) pc->data; 1569*2b3cbbdaSStefano Zampini 1570*2b3cbbdaSStefano Zampini PetscFunctionBegin; 1571*2b3cbbdaSStefano Zampini mg->adaptInterpolation = ctype != PCMG_ADAPT_NONE ? PETSC_TRUE : PETSC_FALSE; 1572*2b3cbbdaSStefano Zampini mg->coarseSpaceType = ctype; 1573*2b3cbbdaSStefano Zampini PetscFunctionReturn(0); 1574*2b3cbbdaSStefano Zampini } 1575*2b3cbbdaSStefano Zampini 1576*2b3cbbdaSStefano Zampini PetscErrorCode PCMGGetAdaptCoarseSpaceType_MG(PC pc, PCMGCoarseSpaceType *ctype) 1577*2b3cbbdaSStefano Zampini { 1578*2b3cbbdaSStefano Zampini PC_MG *mg = (PC_MG *) pc->data; 1579*2b3cbbdaSStefano Zampini 1580*2b3cbbdaSStefano Zampini PetscFunctionBegin; 1581*2b3cbbdaSStefano Zampini *ctype = mg->coarseSpaceType; 1582*2b3cbbdaSStefano Zampini PetscFunctionReturn(0); 1583*2b3cbbdaSStefano Zampini } 1584*2b3cbbdaSStefano Zampini 158541b6fd38SMatthew G. Knepley PetscErrorCode PCMGSetAdaptCR_MG(PC pc, PetscBool cr) 158641b6fd38SMatthew G. Knepley { 158741b6fd38SMatthew G. Knepley PC_MG *mg = (PC_MG *) pc->data; 158841b6fd38SMatthew G. Knepley 158941b6fd38SMatthew G. Knepley PetscFunctionBegin; 159041b6fd38SMatthew G. Knepley mg->compatibleRelaxation = cr; 159141b6fd38SMatthew G. Knepley PetscFunctionReturn(0); 159241b6fd38SMatthew G. Knepley } 159341b6fd38SMatthew G. Knepley 159441b6fd38SMatthew G. Knepley PetscErrorCode PCMGGetAdaptCR_MG(PC pc, PetscBool *cr) 159541b6fd38SMatthew G. Knepley { 159641b6fd38SMatthew G. Knepley PC_MG *mg = (PC_MG *) pc->data; 159741b6fd38SMatthew G. Knepley 159841b6fd38SMatthew G. Knepley PetscFunctionBegin; 159941b6fd38SMatthew G. Knepley *cr = mg->compatibleRelaxation; 160041b6fd38SMatthew G. Knepley PetscFunctionReturn(0); 160141b6fd38SMatthew G. Knepley } 160241b6fd38SMatthew G. Knepley 1603*2b3cbbdaSStefano Zampini /*@C 1604*2b3cbbdaSStefano Zampini PCMGSetAdaptCoarseSpaceType - Set the type of adaptive coarse space. 1605*2b3cbbdaSStefano Zampini 1606*2b3cbbdaSStefano Zampini Adapts or creates the interpolator based upon a vector space which should be accurately captured by the next coarser mesh, and thus accurately interpolated. 1607*2b3cbbdaSStefano Zampini 1608*2b3cbbdaSStefano Zampini Logically Collective on PC 1609*2b3cbbdaSStefano Zampini 1610*2b3cbbdaSStefano Zampini Input Parameters: 1611*2b3cbbdaSStefano Zampini + pc - the multigrid context 1612*2b3cbbdaSStefano Zampini - ctype - the type of coarse space 1613*2b3cbbdaSStefano Zampini 1614*2b3cbbdaSStefano Zampini Options Database Keys: 1615*2b3cbbdaSStefano Zampini + -pc_mg_adapt_interp_n <int> - The number of modes to use 1616*2b3cbbdaSStefano Zampini - -pc_mg_adapt_interp_coarse_space <type> - The type of coarse space: none, polynomial, harmonic, eigenvector, generalized_eigenvector, gdsw 1617*2b3cbbdaSStefano Zampini 1618*2b3cbbdaSStefano Zampini Level: intermediate 1619*2b3cbbdaSStefano Zampini 1620*2b3cbbdaSStefano Zampini .keywords: MG, set, Galerkin 1621*2b3cbbdaSStefano Zampini .seealso: `PCMGCoarseSpaceType`, `PCMGGetAdaptCoarseSpaceType()`, `PCMGSetGalerkin()` 1622*2b3cbbdaSStefano Zampini @*/ 1623*2b3cbbdaSStefano Zampini PetscErrorCode PCMGSetAdaptCoarseSpaceType(PC pc, PCMGCoarseSpaceType ctype) 1624*2b3cbbdaSStefano Zampini { 1625*2b3cbbdaSStefano Zampini PetscFunctionBegin; 1626*2b3cbbdaSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1627*2b3cbbdaSStefano Zampini PetscValidLogicalCollectiveEnum(pc,ctype,2); 1628*2b3cbbdaSStefano Zampini PetscTryMethod(pc,"PCMGSetAdaptCoarseSpaceType_C",(PC,PCMGCoarseSpaceType),(pc,ctype)); 1629*2b3cbbdaSStefano Zampini PetscFunctionReturn(0); 1630*2b3cbbdaSStefano Zampini } 1631*2b3cbbdaSStefano Zampini 1632*2b3cbbdaSStefano Zampini /*@C 1633*2b3cbbdaSStefano Zampini PCMGGetAdaptCoarseSpaceType - Get the type of adaptive coarse space. 1634*2b3cbbdaSStefano Zampini 1635*2b3cbbdaSStefano Zampini Not Collective 1636*2b3cbbdaSStefano Zampini 1637*2b3cbbdaSStefano Zampini Input Parameter: 1638*2b3cbbdaSStefano Zampini . pc - the multigrid context 1639*2b3cbbdaSStefano Zampini 1640*2b3cbbdaSStefano Zampini Output Parameter: 1641*2b3cbbdaSStefano Zampini . ctype - the type of coarse space 1642*2b3cbbdaSStefano Zampini 1643*2b3cbbdaSStefano Zampini Level: intermediate 1644*2b3cbbdaSStefano Zampini 1645*2b3cbbdaSStefano Zampini .keywords: MG, Get, Galerkin 1646*2b3cbbdaSStefano Zampini .seealso: `PCMGCoarseSpaceType`, `PCMGSetAdaptCoarseSpaceType()`, `PCMGSetGalerkin()` 1647*2b3cbbdaSStefano Zampini @*/ 1648*2b3cbbdaSStefano Zampini PetscErrorCode PCMGGetAdaptCoarseSpaceType(PC pc, PCMGCoarseSpaceType *ctype) 1649*2b3cbbdaSStefano Zampini { 1650*2b3cbbdaSStefano Zampini PetscFunctionBegin; 1651*2b3cbbdaSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1652*2b3cbbdaSStefano Zampini PetscValidPointer(ctype,2); 1653*2b3cbbdaSStefano Zampini PetscUseMethod(pc,"PCMGGetAdaptCoarseSpaceType_C",(PC,PCMGCoarseSpaceType*),(pc,ctype)); 1654*2b3cbbdaSStefano Zampini PetscFunctionReturn(0); 1655*2b3cbbdaSStefano Zampini } 1656*2b3cbbdaSStefano Zampini 1657*2b3cbbdaSStefano Zampini /* MATT: REMOVE? */ 1658f3b08a26SMatthew G. Knepley /*@ 1659f3b08a26SMatthew G. Knepley PCMGSetAdaptInterpolation - Adapt the interpolator based upon a vector space which should be accurately captured by the next coarser mesh, and thus accurately interpolated. 1660f3b08a26SMatthew G. Knepley 1661f3b08a26SMatthew G. Knepley Logically Collective on PC 1662f3b08a26SMatthew G. Knepley 1663f3b08a26SMatthew G. Knepley Input Parameters: 1664f3b08a26SMatthew G. Knepley + pc - the multigrid context 1665f3b08a26SMatthew G. Knepley - adapt - flag for adaptation of the interpolator 1666f3b08a26SMatthew G. Knepley 1667f3b08a26SMatthew G. Knepley Options Database Keys: 1668f3b08a26SMatthew G. Knepley + -pc_mg_adapt_interp - Turn on adaptation 1669f3b08a26SMatthew G. Knepley . -pc_mg_adapt_interp_n <int> - The number of modes to use, should be divisible by dimension 1670f3b08a26SMatthew G. Knepley - -pc_mg_adapt_interp_coarse_space <type> - The type of coarse space: polynomial, harmonic, eigenvector, generalized_eigenvector 1671f3b08a26SMatthew G. Knepley 1672f3b08a26SMatthew G. Knepley Level: intermediate 1673f3b08a26SMatthew G. Knepley 1674f3b08a26SMatthew G. Knepley .keywords: MG, set, Galerkin 1675db781477SPatrick Sanan .seealso: `PCMGGetAdaptInterpolation()`, `PCMGSetGalerkin()` 1676f3b08a26SMatthew G. Knepley @*/ 1677f3b08a26SMatthew G. Knepley PetscErrorCode PCMGSetAdaptInterpolation(PC pc, PetscBool adapt) 1678f3b08a26SMatthew G. Knepley { 1679f3b08a26SMatthew G. Knepley PetscFunctionBegin; 1680f3b08a26SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1681cac4c232SBarry Smith PetscTryMethod(pc,"PCMGSetAdaptInterpolation_C",(PC,PetscBool),(pc,adapt)); 1682f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 1683f3b08a26SMatthew G. Knepley } 1684f3b08a26SMatthew G. Knepley 1685f3b08a26SMatthew G. Knepley /*@ 1686f3b08a26SMatthew G. Knepley PCMGGetAdaptInterpolation - Get the flag to adapt the interpolator based upon a vector space which should be accurately captured by the next coarser mesh, and thus accurately interpolated. 1687f3b08a26SMatthew G. Knepley 1688f3b08a26SMatthew G. Knepley Not collective 1689f3b08a26SMatthew G. Knepley 1690f3b08a26SMatthew G. Knepley Input Parameter: 1691f3b08a26SMatthew G. Knepley . pc - the multigrid context 1692f3b08a26SMatthew G. Knepley 1693f3b08a26SMatthew G. Knepley Output Parameter: 1694f3b08a26SMatthew G. Knepley . adapt - flag for adaptation of the interpolator 1695f3b08a26SMatthew G. Knepley 1696f3b08a26SMatthew G. Knepley Level: intermediate 1697f3b08a26SMatthew G. Knepley 1698f3b08a26SMatthew G. Knepley .keywords: MG, set, Galerkin 1699db781477SPatrick Sanan .seealso: `PCMGSetAdaptInterpolation()`, `PCMGSetGalerkin()` 1700f3b08a26SMatthew G. Knepley @*/ 1701f3b08a26SMatthew G. Knepley PetscErrorCode PCMGGetAdaptInterpolation(PC pc, PetscBool *adapt) 1702f3b08a26SMatthew G. Knepley { 1703f3b08a26SMatthew G. Knepley PetscFunctionBegin; 1704f3b08a26SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1705dadcf809SJacob Faibussowitsch PetscValidBoolPointer(adapt, 2); 1706cac4c232SBarry Smith PetscUseMethod(pc,"PCMGGetAdaptInterpolation_C",(PC,PetscBool*),(pc,adapt)); 1707f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 1708f3b08a26SMatthew G. Knepley } 1709f3b08a26SMatthew G. Knepley 17104b9ad928SBarry Smith /*@ 171141b6fd38SMatthew G. Knepley PCMGSetAdaptCR - Monitor the coarse space quality using an auxiliary solve with compatible relaxation. 171241b6fd38SMatthew G. Knepley 171341b6fd38SMatthew G. Knepley Logically Collective on PC 171441b6fd38SMatthew G. Knepley 171541b6fd38SMatthew G. Knepley Input Parameters: 171641b6fd38SMatthew G. Knepley + pc - the multigrid context 171741b6fd38SMatthew G. Knepley - cr - flag for compatible relaxation 171841b6fd38SMatthew G. Knepley 171941b6fd38SMatthew G. Knepley Options Database Keys: 172041b6fd38SMatthew G. Knepley . -pc_mg_adapt_cr - Turn on compatible relaxation 172141b6fd38SMatthew G. Knepley 172241b6fd38SMatthew G. Knepley Level: intermediate 172341b6fd38SMatthew G. Knepley 172441b6fd38SMatthew G. Knepley .keywords: MG, set, Galerkin 1725db781477SPatrick Sanan .seealso: `PCMGGetAdaptCR()`, `PCMGSetAdaptInterpolation()`, `PCMGSetGalerkin()` 172641b6fd38SMatthew G. Knepley @*/ 172741b6fd38SMatthew G. Knepley PetscErrorCode PCMGSetAdaptCR(PC pc, PetscBool cr) 172841b6fd38SMatthew G. Knepley { 172941b6fd38SMatthew G. Knepley PetscFunctionBegin; 173041b6fd38SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1731cac4c232SBarry Smith PetscTryMethod(pc,"PCMGSetAdaptCR_C",(PC,PetscBool),(pc,cr)); 173241b6fd38SMatthew G. Knepley PetscFunctionReturn(0); 173341b6fd38SMatthew G. Knepley } 173441b6fd38SMatthew G. Knepley 173541b6fd38SMatthew G. Knepley /*@ 173641b6fd38SMatthew G. Knepley PCMGGetAdaptCR - Get the flag to monitor coarse space quality using an auxiliary solve with compatible relaxation. 173741b6fd38SMatthew G. Knepley 173841b6fd38SMatthew G. Knepley Not collective 173941b6fd38SMatthew G. Knepley 174041b6fd38SMatthew G. Knepley Input Parameter: 174141b6fd38SMatthew G. Knepley . pc - the multigrid context 174241b6fd38SMatthew G. Knepley 174341b6fd38SMatthew G. Knepley Output Parameter: 174441b6fd38SMatthew G. Knepley . cr - flag for compatible relaxaion 174541b6fd38SMatthew G. Knepley 174641b6fd38SMatthew G. Knepley Level: intermediate 174741b6fd38SMatthew G. Knepley 174841b6fd38SMatthew G. Knepley .keywords: MG, set, Galerkin 1749db781477SPatrick Sanan .seealso: `PCMGSetAdaptCR()`, `PCMGGetAdaptInterpolation()`, `PCMGSetGalerkin()` 175041b6fd38SMatthew G. Knepley @*/ 175141b6fd38SMatthew G. Knepley PetscErrorCode PCMGGetAdaptCR(PC pc, PetscBool *cr) 175241b6fd38SMatthew G. Knepley { 175341b6fd38SMatthew G. Knepley PetscFunctionBegin; 175441b6fd38SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1755dadcf809SJacob Faibussowitsch PetscValidBoolPointer(cr, 2); 1756cac4c232SBarry Smith PetscUseMethod(pc,"PCMGGetAdaptCR_C",(PC,PetscBool*),(pc,cr)); 175741b6fd38SMatthew G. Knepley PetscFunctionReturn(0); 175841b6fd38SMatthew G. Knepley } 175941b6fd38SMatthew G. Knepley 176041b6fd38SMatthew G. Knepley /*@ 176106792cafSBarry Smith PCMGSetNumberSmooth - Sets the number of pre and post-smoothing steps to use 1762710315b6SLawrence Mitchell on all levels. Use PCMGDistinctSmoothUp() to create separate up and down smoothers if you want different numbers of 1763710315b6SLawrence Mitchell pre- and post-smoothing steps. 176406792cafSBarry Smith 176506792cafSBarry Smith Logically Collective on PC 176606792cafSBarry Smith 176706792cafSBarry Smith Input Parameters: 176806792cafSBarry Smith + mg - the multigrid context 176906792cafSBarry Smith - n - the number of smoothing steps 177006792cafSBarry Smith 177106792cafSBarry Smith Options Database Key: 1772a2b725a8SWilliam Gropp . -mg_levels_ksp_max_it <n> - Sets number of pre and post-smoothing steps 177306792cafSBarry Smith 177406792cafSBarry Smith Level: advanced 177506792cafSBarry Smith 177695452b02SPatrick Sanan Notes: 177795452b02SPatrick Sanan this does not set a value on the coarsest grid, since we assume that 177806792cafSBarry Smith there is no separate smooth up on the coarsest grid. 177906792cafSBarry Smith 1780db781477SPatrick Sanan .seealso: `PCMGSetDistinctSmoothUp()` 178106792cafSBarry Smith @*/ 178206792cafSBarry Smith PetscErrorCode PCMGSetNumberSmooth(PC pc,PetscInt n) 178306792cafSBarry Smith { 178406792cafSBarry Smith PC_MG *mg = (PC_MG*)pc->data; 178506792cafSBarry Smith PC_MG_Levels **mglevels = mg->levels; 178606792cafSBarry Smith PetscInt i,levels; 178706792cafSBarry Smith 178806792cafSBarry Smith PetscFunctionBegin; 178906792cafSBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 179006792cafSBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 179128b400f6SJacob Faibussowitsch PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling"); 179206792cafSBarry Smith levels = mglevels[0]->levels; 179306792cafSBarry Smith 179406792cafSBarry Smith for (i=1; i<levels; i++) { 17959566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n)); 17969566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n)); 179706792cafSBarry Smith mg->default_smoothu = n; 179806792cafSBarry Smith mg->default_smoothd = n; 179906792cafSBarry Smith } 180006792cafSBarry Smith PetscFunctionReturn(0); 180106792cafSBarry Smith } 180206792cafSBarry Smith 1803f442ab6aSBarry Smith /*@ 18048e5aa403SBarry Smith PCMGSetDistinctSmoothUp - sets the up (post) smoother to be a separate KSP from the down (pre) smoother on all levels 1805710315b6SLawrence Mitchell and adds the suffix _up to the options name 1806f442ab6aSBarry Smith 1807f442ab6aSBarry Smith Logically Collective on PC 1808f442ab6aSBarry Smith 1809f442ab6aSBarry Smith Input Parameters: 1810f442ab6aSBarry Smith . pc - the preconditioner context 1811f442ab6aSBarry Smith 1812f442ab6aSBarry Smith Options Database Key: 1813147403d9SBarry Smith . -pc_mg_distinct_smoothup <bool> - use distinct smoothing objects 1814f442ab6aSBarry Smith 1815f442ab6aSBarry Smith Level: advanced 1816f442ab6aSBarry Smith 181795452b02SPatrick Sanan Notes: 181895452b02SPatrick Sanan this does not set a value on the coarsest grid, since we assume that 1819f442ab6aSBarry Smith there is no separate smooth up on the coarsest grid. 1820f442ab6aSBarry Smith 1821db781477SPatrick Sanan .seealso: `PCMGSetNumberSmooth()` 1822f442ab6aSBarry Smith @*/ 1823f442ab6aSBarry Smith PetscErrorCode PCMGSetDistinctSmoothUp(PC pc) 1824f442ab6aSBarry Smith { 1825f442ab6aSBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1826f442ab6aSBarry Smith PC_MG_Levels **mglevels = mg->levels; 1827f442ab6aSBarry Smith PetscInt i,levels; 1828f442ab6aSBarry Smith KSP subksp; 1829f442ab6aSBarry Smith 1830f442ab6aSBarry Smith PetscFunctionBegin; 1831f442ab6aSBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 183228b400f6SJacob Faibussowitsch PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling"); 1833f442ab6aSBarry Smith levels = mglevels[0]->levels; 1834f442ab6aSBarry Smith 1835f442ab6aSBarry Smith for (i=1; i<levels; i++) { 1836710315b6SLawrence Mitchell const char *prefix = NULL; 1837f442ab6aSBarry Smith /* make sure smoother up and down are different */ 18389566063dSJacob Faibussowitsch PetscCall(PCMGGetSmootherUp(pc,i,&subksp)); 18399566063dSJacob Faibussowitsch PetscCall(KSPGetOptionsPrefix(mglevels[i]->smoothd,&prefix)); 18409566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(subksp,prefix)); 18419566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(subksp,"up_")); 1842f442ab6aSBarry Smith } 1843f442ab6aSBarry Smith PetscFunctionReturn(0); 1844f442ab6aSBarry Smith } 1845f442ab6aSBarry Smith 184607a4832bSFande Kong /* No new matrices are created, and the coarse operator matrices are the references to the original ones */ 184707a4832bSFande Kong PetscErrorCode PCGetInterpolations_MG(PC pc,PetscInt *num_levels,Mat *interpolations[]) 184807a4832bSFande Kong { 184907a4832bSFande Kong PC_MG *mg = (PC_MG*)pc->data; 185007a4832bSFande Kong PC_MG_Levels **mglevels = mg->levels; 185107a4832bSFande Kong Mat *mat; 185207a4832bSFande Kong PetscInt l; 185307a4832bSFande Kong 185407a4832bSFande Kong PetscFunctionBegin; 185528b400f6SJacob Faibussowitsch PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 18569566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(mg->nlevels,&mat)); 185707a4832bSFande Kong for (l=1; l< mg->nlevels; l++) { 185807a4832bSFande Kong mat[l-1] = mglevels[l]->interpolate; 18599566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat[l-1])); 186007a4832bSFande Kong } 186107a4832bSFande Kong *num_levels = mg->nlevels; 186207a4832bSFande Kong *interpolations = mat; 186307a4832bSFande Kong PetscFunctionReturn(0); 186407a4832bSFande Kong } 186507a4832bSFande Kong 186607a4832bSFande Kong /* No new matrices are created, and the coarse operator matrices are the references to the original ones */ 186707a4832bSFande Kong PetscErrorCode PCGetCoarseOperators_MG(PC pc,PetscInt *num_levels,Mat *coarseOperators[]) 186807a4832bSFande Kong { 186907a4832bSFande Kong PC_MG *mg = (PC_MG*)pc->data; 187007a4832bSFande Kong PC_MG_Levels **mglevels = mg->levels; 187107a4832bSFande Kong PetscInt l; 187207a4832bSFande Kong Mat *mat; 187307a4832bSFande Kong 187407a4832bSFande Kong PetscFunctionBegin; 187528b400f6SJacob Faibussowitsch PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 18769566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(mg->nlevels,&mat)); 187707a4832bSFande Kong for (l=0; l<mg->nlevels-1; l++) { 18789566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[l]->smoothd,NULL,&(mat[l]))); 18799566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat[l])); 188007a4832bSFande Kong } 188107a4832bSFande Kong *num_levels = mg->nlevels; 188207a4832bSFande Kong *coarseOperators = mat; 188307a4832bSFande Kong PetscFunctionReturn(0); 188407a4832bSFande Kong } 188507a4832bSFande Kong 1886f3b08a26SMatthew G. Knepley /*@C 1887f3b08a26SMatthew G. Knepley PCMGRegisterCoarseSpaceConstructor - Adds a method to the PCMG package for coarse space construction. 1888f3b08a26SMatthew G. Knepley 1889f3b08a26SMatthew G. Knepley Not collective 1890f3b08a26SMatthew G. Knepley 1891f3b08a26SMatthew G. Knepley Input Parameters: 1892f3b08a26SMatthew G. Knepley + name - name of the constructor 1893f3b08a26SMatthew G. Knepley - function - constructor routine 1894f3b08a26SMatthew G. Knepley 1895f3b08a26SMatthew G. Knepley Notes: 1896f3b08a26SMatthew G. Knepley Calling sequence for the routine: 1897*2b3cbbdaSStefano Zampini $ my_csp(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, Mat initGuess, Mat *coarseSp) 1898f3b08a26SMatthew G. Knepley $ pc - The PC object 1899f3b08a26SMatthew G. Knepley $ l - The multigrid level, 0 is the coarse level 1900f3b08a26SMatthew G. Knepley $ dm - The DM for this level 1901f3b08a26SMatthew G. Knepley $ smooth - The level smoother 1902f3b08a26SMatthew G. Knepley $ Nc - The size of the coarse space 1903f3b08a26SMatthew G. Knepley $ initGuess - Basis for an initial guess for the space 1904f3b08a26SMatthew G. Knepley $ coarseSp - A basis for the computed coarse space 1905f3b08a26SMatthew G. Knepley 1906f3b08a26SMatthew G. Knepley Level: advanced 1907f3b08a26SMatthew G. Knepley 1908db781477SPatrick Sanan .seealso: `PCMGGetCoarseSpaceConstructor()`, `PCRegister()` 1909f3b08a26SMatthew G. Knepley @*/ 1910*2b3cbbdaSStefano Zampini PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char name[], PetscErrorCode (*function)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat*)) 1911f3b08a26SMatthew G. Knepley { 1912f3b08a26SMatthew G. Knepley PetscFunctionBegin; 19139566063dSJacob Faibussowitsch PetscCall(PCInitializePackage()); 19149566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&PCMGCoarseList,name,function)); 1915f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 1916f3b08a26SMatthew G. Knepley } 1917f3b08a26SMatthew G. Knepley 1918f3b08a26SMatthew G. Knepley /*@C 1919f3b08a26SMatthew G. Knepley PCMGGetCoarseSpaceConstructor - Returns the given coarse space construction method. 1920f3b08a26SMatthew G. Knepley 1921f3b08a26SMatthew G. Knepley Not collective 1922f3b08a26SMatthew G. Knepley 1923f3b08a26SMatthew G. Knepley Input Parameter: 1924f3b08a26SMatthew G. Knepley . name - name of the constructor 1925f3b08a26SMatthew G. Knepley 1926f3b08a26SMatthew G. Knepley Output Parameter: 1927f3b08a26SMatthew G. Knepley . function - constructor routine 1928f3b08a26SMatthew G. Knepley 1929f3b08a26SMatthew G. Knepley Notes: 1930f3b08a26SMatthew G. Knepley Calling sequence for the routine: 1931*2b3cbbdaSStefano Zampini $ my_csp(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, Mat initGuess, Mat *coarseSp) 1932f3b08a26SMatthew G. Knepley $ pc - The PC object 1933f3b08a26SMatthew G. Knepley $ l - The multigrid level, 0 is the coarse level 1934f3b08a26SMatthew G. Knepley $ dm - The DM for this level 1935f3b08a26SMatthew G. Knepley $ smooth - The level smoother 1936f3b08a26SMatthew G. Knepley $ Nc - The size of the coarse space 1937f3b08a26SMatthew G. Knepley $ initGuess - Basis for an initial guess for the space 1938f3b08a26SMatthew G. Knepley $ coarseSp - A basis for the computed coarse space 1939f3b08a26SMatthew G. Knepley 1940f3b08a26SMatthew G. Knepley Level: advanced 1941f3b08a26SMatthew G. Knepley 1942db781477SPatrick Sanan .seealso: `PCMGRegisterCoarseSpaceConstructor()`, `PCRegister()` 1943f3b08a26SMatthew G. Knepley @*/ 1944*2b3cbbdaSStefano Zampini PetscErrorCode PCMGGetCoarseSpaceConstructor(const char name[], PetscErrorCode (**function)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat*)) 1945f3b08a26SMatthew G. Knepley { 1946f3b08a26SMatthew G. Knepley PetscFunctionBegin; 19479566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(PCMGCoarseList,name,function)); 1948f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 1949f3b08a26SMatthew G. Knepley } 1950f3b08a26SMatthew G. Knepley 19514b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/ 19524b9ad928SBarry Smith 19533b09bd56SBarry Smith /*MC 1954ccb205f8SBarry Smith PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional 19553b09bd56SBarry Smith information about the coarser grid matrices and restriction/interpolation operators. 19563b09bd56SBarry Smith 19573b09bd56SBarry Smith Options Database Keys: 19583b09bd56SBarry Smith + -pc_mg_levels <nlevels> - number of levels including finest 1959391689abSStefano Zampini . -pc_mg_cycle_type <v,w> - provide the cycle desired 19608c1c2452SJed Brown . -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default 19613b09bd56SBarry Smith . -pc_mg_log - log information about time spent on each level of the solver 1962710315b6SLawrence Mitchell . -pc_mg_distinct_smoothup - configure up (after interpolation) and down (before restriction) smoothers separately (with different options prefixes) 19632134b1e4SBarry Smith . -pc_mg_galerkin <both,pmat,mat,none> - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R' 19648cf6037eSBarry Smith . -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1) 19658cf6037eSBarry Smith . -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 1966e3c5b3baSBarry Smith to the Socket viewer for reading from MATLAB. 19678cf6037eSBarry Smith - -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices 19688cf6037eSBarry Smith to the binary output file called binaryoutput 19693b09bd56SBarry Smith 197095452b02SPatrick Sanan Notes: 19710b3c753eSRichard Tran Mills If one uses a Krylov method such GMRES or CG as the smoother then one must use KSPFGMRES, KSPGCR, or KSPRICHARDSON as the outer Krylov method 19723b09bd56SBarry Smith 19738cf6037eSBarry Smith When run with a single level the smoother options are used on that level NOT the coarse grid solver options 19748cf6037eSBarry Smith 197523067569SBarry Smith When run with KSPRICHARDSON the convergence test changes slightly if monitor is turned on. The iteration count may change slightly. This 197623067569SBarry Smith is because without monitoring the residual norm is computed WITHIN each multigrid cycle on the finest level after the pre-smoothing 197723067569SBarry Smith (because the residual has just been computed for the multigrid algorithm and is hence available for free) while with monitoring the 197823067569SBarry Smith residual is computed at the end of each cycle. 197923067569SBarry Smith 19803b09bd56SBarry Smith Level: intermediate 19813b09bd56SBarry Smith 1982db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCMGType`, `PCEXOTIC`, `PCGAMG`, `PCML`, `PCHYPRE` 1983db781477SPatrick Sanan `PCMGSetLevels()`, `PCMGGetLevels()`, `PCMGSetType()`, `PCMGSetCycleType()`, 1984db781477SPatrick Sanan `PCMGSetDistinctSmoothUp()`, `PCMGGetCoarseSolve()`, `PCMGSetResidual()`, `PCMGSetInterpolation()`, 1985db781477SPatrick Sanan `PCMGSetRestriction()`, `PCMGGetSmoother()`, `PCMGGetSmootherUp()`, `PCMGGetSmootherDown()`, 1986db781477SPatrick Sanan `PCMGSetCycleTypeOnLevel()`, `PCMGSetRhs()`, `PCMGSetX()`, `PCMGSetR()` 19873b09bd56SBarry Smith M*/ 19883b09bd56SBarry Smith 19898cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc) 19904b9ad928SBarry Smith { 1991f3fbd535SBarry Smith PC_MG *mg; 1992f3fbd535SBarry Smith 19934b9ad928SBarry Smith PetscFunctionBegin; 19949566063dSJacob Faibussowitsch PetscCall(PetscNewLog(pc,&mg)); 19953ec1f749SStefano Zampini pc->data = mg; 1996f3fbd535SBarry Smith mg->nlevels = -1; 199710eca3edSLisandro Dalcin mg->am = PC_MG_MULTIPLICATIVE; 19982134b1e4SBarry Smith mg->galerkin = PC_MG_GALERKIN_NONE; 1999f3b08a26SMatthew G. Knepley mg->adaptInterpolation = PETSC_FALSE; 2000f3b08a26SMatthew G. Knepley mg->Nc = -1; 2001f3b08a26SMatthew G. Knepley mg->eigenvalue = -1; 2002f3fbd535SBarry Smith 200337a44384SMark Adams pc->useAmat = PETSC_TRUE; 200437a44384SMark Adams 20054b9ad928SBarry Smith pc->ops->apply = PCApply_MG; 2006fcb023d4SJed Brown pc->ops->applytranspose = PCApplyTranspose_MG; 200730b0564aSStefano Zampini pc->ops->matapply = PCMatApply_MG; 20084b9ad928SBarry Smith pc->ops->setup = PCSetUp_MG; 2009a06653b4SBarry Smith pc->ops->reset = PCReset_MG; 20104b9ad928SBarry Smith pc->ops->destroy = PCDestroy_MG; 20114b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_MG; 20124b9ad928SBarry Smith pc->ops->view = PCView_MG; 2013fb15c04eSBarry Smith 20149566063dSJacob Faibussowitsch PetscCall(PetscObjectComposedDataRegister(&mg->eigenvalue)); 20159566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",PCMGSetGalerkin_MG)); 20169566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGGetLevels_C",PCMGGetLevels_MG)); 20179566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGSetLevels_C",PCMGSetLevels_MG)); 20189566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGetInterpolations_C",PCGetInterpolations_MG)); 20199566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGetCoarseOperators_C",PCGetCoarseOperators_MG)); 20209566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGSetAdaptInterpolation_C",PCMGSetAdaptInterpolation_MG)); 20219566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGGetAdaptInterpolation_C",PCMGGetAdaptInterpolation_MG)); 20229566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGSetAdaptCR_C",PCMGSetAdaptCR_MG)); 20239566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGGetAdaptCR_C",PCMGGetAdaptCR_MG)); 2024*2b3cbbdaSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGSetAdaptCoarseSpaceType_C",PCMGSetAdaptCoarseSpaceType_MG)); 2025*2b3cbbdaSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGGetAdaptCoarseSpaceType_C",PCMGGetAdaptCoarseSpaceType_MG)); 20264b9ad928SBarry Smith PetscFunctionReturn(0); 20274b9ad928SBarry Smith } 2028