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 159371c9d4SSatish Balay PetscErrorCode PCMGMCycle_Private(PC pc, PC_MG_Levels **mglevelsin, PetscBool transpose, PetscBool matapp, PCRichardsonConvergedReason *reason) { 1631567311SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 1731567311SBarry Smith PC_MG_Levels *mgc, *mglevels = *mglevelsin; 1831567311SBarry Smith PetscInt cycles = (mglevels->level == 1) ? 1 : (PetscInt)mglevels->cycles; 194b9ad928SBarry Smith 204b9ad928SBarry Smith PetscFunctionBegin; 219566063dSJacob Faibussowitsch if (mglevels->eventsmoothsolve) PetscCall(PetscLogEventBegin(mglevels->eventsmoothsolve, 0, 0, 0, 0)); 22fcb023d4SJed Brown if (!transpose) { 2330b0564aSStefano Zampini if (matapp) { 249566063dSJacob Faibussowitsch PetscCall(KSPMatSolve(mglevels->smoothd, mglevels->B, mglevels->X)); /* pre-smooth */ 259566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(mglevels->smoothd, pc, NULL)); 2630b0564aSStefano Zampini } else { 279566063dSJacob Faibussowitsch PetscCall(KSPSolve(mglevels->smoothd, mglevels->b, mglevels->x)); /* pre-smooth */ 289566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(mglevels->smoothd, pc, mglevels->x)); 2930b0564aSStefano Zampini } 30fcb023d4SJed Brown } else { 3128b400f6SJacob Faibussowitsch PetscCheck(!matapp, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Not supported"); 329566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(mglevels->smoothu, mglevels->b, mglevels->x)); /* transpose of post-smooth */ 339566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(mglevels->smoothu, pc, mglevels->x)); 34fcb023d4SJed Brown } 359566063dSJacob Faibussowitsch if (mglevels->eventsmoothsolve) PetscCall(PetscLogEventEnd(mglevels->eventsmoothsolve, 0, 0, 0, 0)); 3631567311SBarry Smith if (mglevels->level) { /* not the coarsest grid */ 379566063dSJacob Faibussowitsch if (mglevels->eventresidual) PetscCall(PetscLogEventBegin(mglevels->eventresidual, 0, 0, 0, 0)); 3848a46eb9SPierre Jolivet if (matapp && !mglevels->R) PetscCall(MatDuplicate(mglevels->B, MAT_DO_NOT_COPY_VALUES, &mglevels->R)); 39fcb023d4SJed Brown if (!transpose) { 409566063dSJacob Faibussowitsch if (matapp) PetscCall((*mglevels->matresidual)(mglevels->A, mglevels->B, mglevels->X, mglevels->R)); 419566063dSJacob Faibussowitsch else PetscCall((*mglevels->residual)(mglevels->A, mglevels->b, mglevels->x, mglevels->r)); 42fcb023d4SJed Brown } else { 439566063dSJacob Faibussowitsch if (matapp) PetscCall((*mglevels->matresidualtranspose)(mglevels->A, mglevels->B, mglevels->X, mglevels->R)); 449566063dSJacob Faibussowitsch else PetscCall((*mglevels->residualtranspose)(mglevels->A, mglevels->b, mglevels->x, mglevels->r)); 45fcb023d4SJed Brown } 469566063dSJacob Faibussowitsch if (mglevels->eventresidual) PetscCall(PetscLogEventEnd(mglevels->eventresidual, 0, 0, 0, 0)); 474b9ad928SBarry Smith 484b9ad928SBarry Smith /* if on finest level and have convergence criteria set */ 4931567311SBarry Smith if (mglevels->level == mglevels->levels - 1 && mg->ttol && reason) { 504b9ad928SBarry Smith PetscReal rnorm; 519566063dSJacob Faibussowitsch PetscCall(VecNorm(mglevels->r, NORM_2, &rnorm)); 524b9ad928SBarry Smith if (rnorm <= mg->ttol) { 5370441072SBarry Smith if (rnorm < mg->abstol) { 544d0a8057SBarry Smith *reason = PCRICHARDSON_CONVERGED_ATOL; 559566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "Linear solver has converged. Residual norm %g is less than absolute tolerance %g\n", (double)rnorm, (double)mg->abstol)); 564b9ad928SBarry Smith } else { 574d0a8057SBarry Smith *reason = PCRICHARDSON_CONVERGED_RTOL; 589566063dSJacob 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)); 594b9ad928SBarry Smith } 604b9ad928SBarry Smith PetscFunctionReturn(0); 614b9ad928SBarry Smith } 624b9ad928SBarry Smith } 634b9ad928SBarry Smith 6431567311SBarry Smith mgc = *(mglevelsin - 1); 659566063dSJacob Faibussowitsch if (mglevels->eventinterprestrict) PetscCall(PetscLogEventBegin(mglevels->eventinterprestrict, 0, 0, 0, 0)); 66fcb023d4SJed Brown if (!transpose) { 679566063dSJacob Faibussowitsch if (matapp) PetscCall(MatMatRestrict(mglevels->restrct, mglevels->R, &mgc->B)); 689566063dSJacob Faibussowitsch else PetscCall(MatRestrict(mglevels->restrct, mglevels->r, mgc->b)); 69fcb023d4SJed Brown } else { 709566063dSJacob Faibussowitsch if (matapp) PetscCall(MatMatRestrict(mglevels->interpolate, mglevels->R, &mgc->B)); 719566063dSJacob Faibussowitsch else PetscCall(MatRestrict(mglevels->interpolate, mglevels->r, mgc->b)); 72fcb023d4SJed Brown } 739566063dSJacob Faibussowitsch if (mglevels->eventinterprestrict) PetscCall(PetscLogEventEnd(mglevels->eventinterprestrict, 0, 0, 0, 0)); 7430b0564aSStefano Zampini if (matapp) { 7530b0564aSStefano Zampini if (!mgc->X) { 769566063dSJacob Faibussowitsch PetscCall(MatDuplicate(mgc->B, MAT_DO_NOT_COPY_VALUES, &mgc->X)); 7730b0564aSStefano Zampini } else { 789566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(mgc->X)); 7930b0564aSStefano Zampini } 8030b0564aSStefano Zampini } else { 819566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(mgc->x)); 8230b0564aSStefano Zampini } 8348a46eb9SPierre Jolivet while (cycles--) PetscCall(PCMGMCycle_Private(pc, mglevelsin - 1, transpose, matapp, reason)); 849566063dSJacob Faibussowitsch if (mglevels->eventinterprestrict) PetscCall(PetscLogEventBegin(mglevels->eventinterprestrict, 0, 0, 0, 0)); 85fcb023d4SJed Brown if (!transpose) { 869566063dSJacob Faibussowitsch if (matapp) PetscCall(MatMatInterpolateAdd(mglevels->interpolate, mgc->X, mglevels->X, &mglevels->X)); 879566063dSJacob Faibussowitsch else PetscCall(MatInterpolateAdd(mglevels->interpolate, mgc->x, mglevels->x, mglevels->x)); 88fcb023d4SJed Brown } else { 899566063dSJacob Faibussowitsch PetscCall(MatInterpolateAdd(mglevels->restrct, mgc->x, mglevels->x, mglevels->x)); 90fcb023d4SJed Brown } 919566063dSJacob Faibussowitsch if (mglevels->eventinterprestrict) PetscCall(PetscLogEventEnd(mglevels->eventinterprestrict, 0, 0, 0, 0)); 929566063dSJacob Faibussowitsch if (mglevels->eventsmoothsolve) PetscCall(PetscLogEventBegin(mglevels->eventsmoothsolve, 0, 0, 0, 0)); 93fcb023d4SJed Brown if (!transpose) { 9430b0564aSStefano Zampini if (matapp) { 959566063dSJacob Faibussowitsch PetscCall(KSPMatSolve(mglevels->smoothu, mglevels->B, mglevels->X)); /* post smooth */ 969566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(mglevels->smoothu, pc, NULL)); 9730b0564aSStefano Zampini } else { 989566063dSJacob Faibussowitsch PetscCall(KSPSolve(mglevels->smoothu, mglevels->b, mglevels->x)); /* post smooth */ 999566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(mglevels->smoothu, pc, mglevels->x)); 10030b0564aSStefano Zampini } 101fcb023d4SJed Brown } else { 10228b400f6SJacob Faibussowitsch PetscCheck(!matapp, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Not supported"); 1039566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(mglevels->smoothd, mglevels->b, mglevels->x)); /* post smooth */ 1049566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(mglevels->smoothd, pc, mglevels->x)); 105fcb023d4SJed Brown } 10641b6fd38SMatthew G. Knepley if (mglevels->cr) { 10728b400f6SJacob Faibussowitsch PetscCheck(!matapp, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Not supported"); 10841b6fd38SMatthew G. Knepley /* TODO Turn on copy and turn off noisy if we have an exact solution 1099566063dSJacob Faibussowitsch PetscCall(VecCopy(mglevels->x, mglevels->crx)); 1109566063dSJacob Faibussowitsch PetscCall(VecCopy(mglevels->b, mglevels->crb)); */ 1119566063dSJacob Faibussowitsch PetscCall(KSPSetNoisy_Private(mglevels->crx)); 1129566063dSJacob Faibussowitsch PetscCall(KSPSolve(mglevels->cr, mglevels->crb, mglevels->crx)); /* compatible relaxation */ 1139566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(mglevels->cr, pc, mglevels->crx)); 11441b6fd38SMatthew G. Knepley } 1159566063dSJacob Faibussowitsch if (mglevels->eventsmoothsolve) PetscCall(PetscLogEventEnd(mglevels->eventsmoothsolve, 0, 0, 0, 0)); 1164b9ad928SBarry Smith } 1174b9ad928SBarry Smith PetscFunctionReturn(0); 1184b9ad928SBarry Smith } 1194b9ad928SBarry Smith 1209371c9d4SSatish Balay 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) { 121f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 122f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 123391689abSStefano Zampini PC tpc; 124391689abSStefano Zampini PetscBool changeu, changed; 125f3fbd535SBarry Smith PetscInt levels = mglevels[0]->levels, i; 1264b9ad928SBarry Smith 1274b9ad928SBarry Smith PetscFunctionBegin; 128a762d673SBarry Smith /* When the DM is supplying the matrix then it will not exist until here */ 129a762d673SBarry Smith for (i = 0; i < levels; i++) { 130a762d673SBarry Smith if (!mglevels[i]->A) { 1319566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[i]->smoothu, &mglevels[i]->A, NULL)); 1329566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mglevels[i]->A)); 133a762d673SBarry Smith } 134a762d673SBarry Smith } 135391689abSStefano Zampini 1369566063dSJacob Faibussowitsch PetscCall(KSPGetPC(mglevels[levels - 1]->smoothd, &tpc)); 1379566063dSJacob Faibussowitsch PetscCall(PCPreSolveChangeRHS(tpc, &changed)); 1389566063dSJacob Faibussowitsch PetscCall(KSPGetPC(mglevels[levels - 1]->smoothu, &tpc)); 1399566063dSJacob Faibussowitsch PetscCall(PCPreSolveChangeRHS(tpc, &changeu)); 140391689abSStefano Zampini if (!changed && !changeu) { 1419566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[levels - 1]->b)); 142f3fbd535SBarry Smith mglevels[levels - 1]->b = b; 143391689abSStefano Zampini } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */ 144391689abSStefano Zampini if (!mglevels[levels - 1]->b) { 145391689abSStefano Zampini Vec *vec; 146391689abSStefano Zampini 1479566063dSJacob Faibussowitsch PetscCall(KSPCreateVecs(mglevels[levels - 1]->smoothd, 1, &vec, 0, NULL)); 148391689abSStefano Zampini mglevels[levels - 1]->b = *vec; 1499566063dSJacob Faibussowitsch PetscCall(PetscFree(vec)); 150391689abSStefano Zampini } 1519566063dSJacob Faibussowitsch PetscCall(VecCopy(b, mglevels[levels - 1]->b)); 152391689abSStefano Zampini } 153f3fbd535SBarry Smith mglevels[levels - 1]->x = x; 1544b9ad928SBarry Smith 15531567311SBarry Smith mg->rtol = rtol; 15631567311SBarry Smith mg->abstol = abstol; 15731567311SBarry Smith mg->dtol = dtol; 1584b9ad928SBarry Smith if (rtol) { 1594b9ad928SBarry Smith /* compute initial residual norm for relative convergence test */ 1604b9ad928SBarry Smith PetscReal rnorm; 1617319c654SBarry Smith if (zeroguess) { 1629566063dSJacob Faibussowitsch PetscCall(VecNorm(b, NORM_2, &rnorm)); 1637319c654SBarry Smith } else { 1649566063dSJacob Faibussowitsch PetscCall((*mglevels[levels - 1]->residual)(mglevels[levels - 1]->A, b, x, w)); 1659566063dSJacob Faibussowitsch PetscCall(VecNorm(w, NORM_2, &rnorm)); 1667319c654SBarry Smith } 16731567311SBarry Smith mg->ttol = PetscMax(rtol * rnorm, abstol); 1682fa5cd67SKarl Rupp } else if (abstol) mg->ttol = abstol; 1692fa5cd67SKarl Rupp else mg->ttol = 0.0; 1704b9ad928SBarry Smith 1714d0a8057SBarry Smith /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't 172336babb1SJed Brown stop prematurely due to small residual */ 1734d0a8057SBarry Smith for (i = 1; i < levels; i++) { 1749566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(mglevels[i]->smoothu, 0, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 175f3fbd535SBarry Smith if (mglevels[i]->smoothu != mglevels[i]->smoothd) { 17623067569SBarry Smith /* For Richardson the initial guess is nonzero since it is solving in each cycle the original system not just applying as a preconditioner */ 1779566063dSJacob Faibussowitsch PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->smoothd, PETSC_TRUE)); 1789566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(mglevels[i]->smoothd, 0, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 1794b9ad928SBarry Smith } 1804d0a8057SBarry Smith } 1814d0a8057SBarry Smith 1824d0a8057SBarry Smith *reason = (PCRichardsonConvergedReason)0; 1834d0a8057SBarry Smith for (i = 0; i < its; i++) { 1849566063dSJacob Faibussowitsch PetscCall(PCMGMCycle_Private(pc, mglevels + levels - 1, PETSC_FALSE, PETSC_FALSE, reason)); 1854d0a8057SBarry Smith if (*reason) break; 1864d0a8057SBarry Smith } 1874d0a8057SBarry Smith if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS; 1884d0a8057SBarry Smith *outits = i; 189391689abSStefano Zampini if (!changed && !changeu) mglevels[levels - 1]->b = NULL; 1904b9ad928SBarry Smith PetscFunctionReturn(0); 1914b9ad928SBarry Smith } 1924b9ad928SBarry Smith 1939371c9d4SSatish Balay PetscErrorCode PCReset_MG(PC pc) { 1943aeaf226SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 1953aeaf226SBarry Smith PC_MG_Levels **mglevels = mg->levels; 1962b3cbbdaSStefano Zampini PetscInt i, n; 1973aeaf226SBarry Smith 1983aeaf226SBarry Smith PetscFunctionBegin; 1993aeaf226SBarry Smith if (mglevels) { 2003aeaf226SBarry Smith n = mglevels[0]->levels; 2013aeaf226SBarry Smith for (i = 0; i < n - 1; i++) { 2029566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[i + 1]->r)); 2039566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[i]->b)); 2049566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[i]->x)); 2059566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[i + 1]->R)); 2069566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[i]->B)); 2079566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[i]->X)); 2089566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[i]->crx)); 2099566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[i]->crb)); 2109566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[i + 1]->restrct)); 2119566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[i + 1]->interpolate)); 2129566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[i + 1]->inject)); 2139566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[i + 1]->rscale)); 2143aeaf226SBarry Smith } 2159566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[n - 1]->crx)); 2169566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[n - 1]->crb)); 217391689abSStefano Zampini /* this is not null only if the smoother on the finest level 218391689abSStefano Zampini changes the rhs during PreSolve */ 2199566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[n - 1]->b)); 2209566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[n - 1]->B)); 2213aeaf226SBarry Smith 2223aeaf226SBarry Smith for (i = 0; i < n; i++) { 2232b3cbbdaSStefano Zampini PetscCall(MatDestroy(&mglevels[i]->coarseSpace)); 2249566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[i]->A)); 22548a46eb9SPierre Jolivet if (mglevels[i]->smoothd != mglevels[i]->smoothu) PetscCall(KSPReset(mglevels[i]->smoothd)); 2269566063dSJacob Faibussowitsch PetscCall(KSPReset(mglevels[i]->smoothu)); 2279566063dSJacob Faibussowitsch if (mglevels[i]->cr) PetscCall(KSPReset(mglevels[i]->cr)); 2283aeaf226SBarry Smith } 229f3b08a26SMatthew G. Knepley mg->Nc = 0; 2303aeaf226SBarry Smith } 2313aeaf226SBarry Smith PetscFunctionReturn(0); 2323aeaf226SBarry Smith } 2333aeaf226SBarry Smith 23441b6fd38SMatthew G. Knepley /* Implementing CR 23541b6fd38SMatthew G. Knepley 23641b6fd38SMatthew 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 23741b6fd38SMatthew G. Knepley 23841b6fd38SMatthew G. Knepley Inj^T (Inj Inj^T)^{-1} Inj 23941b6fd38SMatthew G. Knepley 24041b6fd38SMatthew G. Knepley and if Inj is a VecScatter, as it is now in PETSc, we have 24141b6fd38SMatthew G. Knepley 24241b6fd38SMatthew G. Knepley Inj^T Inj 24341b6fd38SMatthew G. Knepley 24441b6fd38SMatthew G. Knepley and 24541b6fd38SMatthew G. Knepley 24641b6fd38SMatthew G. Knepley S = I - Inj^T Inj 24741b6fd38SMatthew G. Knepley 24841b6fd38SMatthew G. Knepley since 24941b6fd38SMatthew G. Knepley 25041b6fd38SMatthew G. Knepley Inj S = Inj - (Inj Inj^T) Inj = 0. 25141b6fd38SMatthew G. Knepley 25241b6fd38SMatthew G. Knepley Brannick suggests 25341b6fd38SMatthew G. Knepley 25441b6fd38SMatthew G. Knepley A \to S^T A S \qquad\mathrm{and}\qquad M \to S^T M S 25541b6fd38SMatthew G. Knepley 25641b6fd38SMatthew 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 25741b6fd38SMatthew G. Knepley 25841b6fd38SMatthew G. Knepley M^{-1} A \to S M^{-1} A S 25941b6fd38SMatthew G. Knepley 26041b6fd38SMatthew G. Knepley In fact, since it is somewhat hard in PETSc to do the symmetric application, we will just apply S on the left. 26141b6fd38SMatthew G. Knepley 26241b6fd38SMatthew G. Knepley Check: || Inj P - I ||_F < tol 26341b6fd38SMatthew G. Knepley Check: In general, Inj Inj^T = I 26441b6fd38SMatthew G. Knepley */ 26541b6fd38SMatthew G. Knepley 26641b6fd38SMatthew G. Knepley typedef struct { 26741b6fd38SMatthew G. Knepley PC mg; /* The PCMG object */ 26841b6fd38SMatthew G. Knepley PetscInt l; /* The multigrid level for this solver */ 26941b6fd38SMatthew G. Knepley Mat Inj; /* The injection matrix */ 27041b6fd38SMatthew G. Knepley Mat S; /* I - Inj^T Inj */ 27141b6fd38SMatthew G. Knepley } CRContext; 27241b6fd38SMatthew G. Knepley 2739371c9d4SSatish Balay static PetscErrorCode CRSetup_Private(PC pc) { 27441b6fd38SMatthew G. Knepley CRContext *ctx; 27541b6fd38SMatthew G. Knepley Mat It; 27641b6fd38SMatthew G. Knepley 27741b6fd38SMatthew G. Knepley PetscFunctionBeginUser; 2789566063dSJacob Faibussowitsch PetscCall(PCShellGetContext(pc, &ctx)); 2799566063dSJacob Faibussowitsch PetscCall(PCMGGetInjection(ctx->mg, ctx->l, &It)); 28028b400f6SJacob Faibussowitsch PetscCheck(It, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "CR requires that injection be defined for this PCMG"); 2819566063dSJacob Faibussowitsch PetscCall(MatCreateTranspose(It, &ctx->Inj)); 2829566063dSJacob Faibussowitsch PetscCall(MatCreateNormal(ctx->Inj, &ctx->S)); 2839566063dSJacob Faibussowitsch PetscCall(MatScale(ctx->S, -1.0)); 2849566063dSJacob Faibussowitsch PetscCall(MatShift(ctx->S, 1.0)); 28541b6fd38SMatthew G. Knepley PetscFunctionReturn(0); 28641b6fd38SMatthew G. Knepley } 28741b6fd38SMatthew G. Knepley 2889371c9d4SSatish Balay static PetscErrorCode CRApply_Private(PC pc, Vec x, Vec y) { 28941b6fd38SMatthew G. Knepley CRContext *ctx; 29041b6fd38SMatthew G. Knepley 29141b6fd38SMatthew G. Knepley PetscFunctionBeginUser; 2929566063dSJacob Faibussowitsch PetscCall(PCShellGetContext(pc, &ctx)); 2939566063dSJacob Faibussowitsch PetscCall(MatMult(ctx->S, x, y)); 29441b6fd38SMatthew G. Knepley PetscFunctionReturn(0); 29541b6fd38SMatthew G. Knepley } 29641b6fd38SMatthew G. Knepley 2979371c9d4SSatish Balay static PetscErrorCode CRDestroy_Private(PC pc) { 29841b6fd38SMatthew G. Knepley CRContext *ctx; 29941b6fd38SMatthew G. Knepley 30041b6fd38SMatthew G. Knepley PetscFunctionBeginUser; 3019566063dSJacob Faibussowitsch PetscCall(PCShellGetContext(pc, &ctx)); 3029566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx->Inj)); 3039566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx->S)); 3049566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 3059566063dSJacob Faibussowitsch PetscCall(PCShellSetContext(pc, NULL)); 30641b6fd38SMatthew G. Knepley PetscFunctionReturn(0); 30741b6fd38SMatthew G. Knepley } 30841b6fd38SMatthew G. Knepley 3099371c9d4SSatish Balay static PetscErrorCode CreateCR_Private(PC pc, PetscInt l, PC *cr) { 31041b6fd38SMatthew G. Knepley CRContext *ctx; 31141b6fd38SMatthew G. Knepley 31241b6fd38SMatthew G. Knepley PetscFunctionBeginUser; 3139566063dSJacob Faibussowitsch PetscCall(PCCreate(PetscObjectComm((PetscObject)pc), cr)); 3149566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*cr, "S (complementary projector to injection)")); 3159566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(1, &ctx)); 31641b6fd38SMatthew G. Knepley ctx->mg = pc; 31741b6fd38SMatthew G. Knepley ctx->l = l; 3189566063dSJacob Faibussowitsch PetscCall(PCSetType(*cr, PCSHELL)); 3199566063dSJacob Faibussowitsch PetscCall(PCShellSetContext(*cr, ctx)); 3209566063dSJacob Faibussowitsch PetscCall(PCShellSetApply(*cr, CRApply_Private)); 3219566063dSJacob Faibussowitsch PetscCall(PCShellSetSetUp(*cr, CRSetup_Private)); 3229566063dSJacob Faibussowitsch PetscCall(PCShellSetDestroy(*cr, CRDestroy_Private)); 32341b6fd38SMatthew G. Knepley PetscFunctionReturn(0); 32441b6fd38SMatthew G. Knepley } 32541b6fd38SMatthew G. Knepley 3269371c9d4SSatish Balay PetscErrorCode PCMGSetLevels_MG(PC pc, PetscInt levels, MPI_Comm *comms) { 327f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 328ce94432eSBarry Smith MPI_Comm comm; 3293aeaf226SBarry Smith PC_MG_Levels **mglevels = mg->levels; 33010eca3edSLisandro Dalcin PCMGType mgtype = mg->am; 33110167fecSBarry Smith PetscInt mgctype = (PetscInt)PC_MG_CYCLE_V; 332f3fbd535SBarry Smith PetscInt i; 333f3fbd535SBarry Smith PetscMPIInt size; 334f3fbd535SBarry Smith const char *prefix; 335f3fbd535SBarry Smith PC ipc; 3363aeaf226SBarry Smith PetscInt n; 3374b9ad928SBarry Smith 3384b9ad928SBarry Smith PetscFunctionBegin; 3390700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 340548767e0SJed Brown PetscValidLogicalCollectiveInt(pc, levels, 2); 341548767e0SJed Brown if (mg->nlevels == levels) PetscFunctionReturn(0); 3429566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)pc, &comm)); 3433aeaf226SBarry Smith if (mglevels) { 34410eca3edSLisandro Dalcin mgctype = mglevels[0]->cycles; 3453aeaf226SBarry Smith /* changing the number of levels so free up the previous stuff */ 3469566063dSJacob Faibussowitsch PetscCall(PCReset_MG(pc)); 3473aeaf226SBarry Smith n = mglevels[0]->levels; 3483aeaf226SBarry Smith for (i = 0; i < n; i++) { 34948a46eb9SPierre Jolivet if (mglevels[i]->smoothd != mglevels[i]->smoothu) PetscCall(KSPDestroy(&mglevels[i]->smoothd)); 3509566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&mglevels[i]->smoothu)); 3519566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&mglevels[i]->cr)); 3529566063dSJacob Faibussowitsch PetscCall(PetscFree(mglevels[i])); 3533aeaf226SBarry Smith } 3549566063dSJacob Faibussowitsch PetscCall(PetscFree(mg->levels)); 3553aeaf226SBarry Smith } 356f3fbd535SBarry Smith 357f3fbd535SBarry Smith mg->nlevels = levels; 358f3fbd535SBarry Smith 3599566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(levels, &mglevels)); 3609566063dSJacob Faibussowitsch PetscCall(PetscLogObjectMemory((PetscObject)pc, levels * (sizeof(PC_MG *)))); 361f3fbd535SBarry Smith 3629566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &prefix)); 363f3fbd535SBarry Smith 364a9db87a2SMatthew G Knepley mg->stageApply = 0; 365f3fbd535SBarry Smith for (i = 0; i < levels; i++) { 3669566063dSJacob Faibussowitsch PetscCall(PetscNewLog(pc, &mglevels[i])); 3672fa5cd67SKarl Rupp 368f3fbd535SBarry Smith mglevels[i]->level = i; 369f3fbd535SBarry Smith mglevels[i]->levels = levels; 37010eca3edSLisandro Dalcin mglevels[i]->cycles = mgctype; 371336babb1SJed Brown mg->default_smoothu = 2; 372336babb1SJed Brown mg->default_smoothd = 2; 37363e6d426SJed Brown mglevels[i]->eventsmoothsetup = 0; 37463e6d426SJed Brown mglevels[i]->eventsmoothsolve = 0; 37563e6d426SJed Brown mglevels[i]->eventresidual = 0; 37663e6d426SJed Brown mglevels[i]->eventinterprestrict = 0; 377f3fbd535SBarry Smith 378f3fbd535SBarry Smith if (comms) comm = comms[i]; 379d5a8f7e6SBarry Smith if (comm != MPI_COMM_NULL) { 3809566063dSJacob Faibussowitsch PetscCall(KSPCreate(comm, &mglevels[i]->smoothd)); 3819566063dSJacob Faibussowitsch PetscCall(KSPSetErrorIfNotConverged(mglevels[i]->smoothd, pc->erroriffailure)); 3829566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd, (PetscObject)pc, levels - i)); 3839566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(mglevels[i]->smoothd, prefix)); 3849566063dSJacob Faibussowitsch PetscCall(PetscObjectComposedDataSetInt((PetscObject)mglevels[i]->smoothd, PetscMGLevelId, mglevels[i]->level)); 3850ee9a668SBarry Smith if (i || levels == 1) { 3860ee9a668SBarry Smith char tprefix[128]; 3870ee9a668SBarry Smith 3889566063dSJacob Faibussowitsch PetscCall(KSPSetType(mglevels[i]->smoothd, KSPCHEBYSHEV)); 3899566063dSJacob Faibussowitsch PetscCall(KSPSetConvergenceTest(mglevels[i]->smoothd, KSPConvergedSkip, NULL, NULL)); 3909566063dSJacob Faibussowitsch PetscCall(KSPSetNormType(mglevels[i]->smoothd, KSP_NORM_NONE)); 3919566063dSJacob Faibussowitsch PetscCall(KSPGetPC(mglevels[i]->smoothd, &ipc)); 3929566063dSJacob Faibussowitsch PetscCall(PCSetType(ipc, PCSOR)); 3939566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(mglevels[i]->smoothd, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, mg->default_smoothd)); 394f3fbd535SBarry Smith 3959566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(tprefix, 128, "mg_levels_%d_", (int)i)); 3969566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(mglevels[i]->smoothd, tprefix)); 3970ee9a668SBarry Smith } else { 3989566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(mglevels[0]->smoothd, "mg_coarse_")); 399f3fbd535SBarry Smith 4009dbfc187SHong Zhang /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */ 4019566063dSJacob Faibussowitsch PetscCall(KSPSetType(mglevels[0]->smoothd, KSPPREONLY)); 4029566063dSJacob Faibussowitsch PetscCall(KSPGetPC(mglevels[0]->smoothd, &ipc)); 4039566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 404f3fbd535SBarry Smith if (size > 1) { 4059566063dSJacob Faibussowitsch PetscCall(PCSetType(ipc, PCREDUNDANT)); 406f3fbd535SBarry Smith } else { 4079566063dSJacob Faibussowitsch PetscCall(PCSetType(ipc, PCLU)); 408f3fbd535SBarry Smith } 4099566063dSJacob Faibussowitsch PetscCall(PCFactorSetShiftType(ipc, MAT_SHIFT_INBLOCKS)); 410f3fbd535SBarry Smith } 4119566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)mglevels[i]->smoothd)); 412d5a8f7e6SBarry Smith } 413f3fbd535SBarry Smith mglevels[i]->smoothu = mglevels[i]->smoothd; 41431567311SBarry Smith mg->rtol = 0.0; 41531567311SBarry Smith mg->abstol = 0.0; 41631567311SBarry Smith mg->dtol = 0.0; 41731567311SBarry Smith mg->ttol = 0.0; 41831567311SBarry Smith mg->cyclesperpcapply = 1; 419f3fbd535SBarry Smith } 420f3fbd535SBarry Smith mg->levels = mglevels; 4219566063dSJacob Faibussowitsch PetscCall(PCMGSetType(pc, mgtype)); 4224b9ad928SBarry Smith PetscFunctionReturn(0); 4234b9ad928SBarry Smith } 4244b9ad928SBarry Smith 42597d33e41SMatthew G. Knepley /*@C 426*f1580f4eSBarry Smith PCMGSetLevels - Sets the number of levels to use with `PCMG`. 427*f1580f4eSBarry Smith Must be called before any other `PCMG` routine. 42897d33e41SMatthew G. Knepley 429*f1580f4eSBarry Smith Logically Collective on pc 43097d33e41SMatthew G. Knepley 43197d33e41SMatthew G. Knepley Input Parameters: 43297d33e41SMatthew G. Knepley + pc - the preconditioner context 43397d33e41SMatthew G. Knepley . levels - the number of levels 43497d33e41SMatthew G. Knepley - comms - optional communicators for each level; this is to allow solving the coarser problems 435d5a8f7e6SBarry Smith on smaller sets of processes. For processes that are not included in the computation 436*f1580f4eSBarry Smith you must pass `MPI_COMM_NULL`. Use comms = NULL to specify that all processes 43705552d4bSJunchao Zhang should participate in each level of problem. 43897d33e41SMatthew G. Knepley 43997d33e41SMatthew G. Knepley Level: intermediate 44097d33e41SMatthew G. Knepley 44197d33e41SMatthew G. Knepley Notes: 44297d33e41SMatthew G. Knepley If the number of levels is one then the multigrid uses the -mg_levels prefix 44397d33e41SMatthew G. Knepley for setting the level options rather than the -mg_coarse prefix. 44497d33e41SMatthew G. Knepley 445d5a8f7e6SBarry Smith You can free the information in comms after this routine is called. 446d5a8f7e6SBarry Smith 447*f1580f4eSBarry Smith The array of MPI communicators must contain `MPI_COMM_NULL` for those ranks that at each level 448d5a8f7e6SBarry Smith are not participating in the coarser solve. For example, with 2 levels and 1 and 2 ranks on 449d5a8f7e6SBarry Smith the two levels, rank 0 in the original communicator will pass in an array of 2 communicators 450d5a8f7e6SBarry Smith of size 2 and 1, while rank 1 in the original communicator will pass in array of 2 communicators 451*f1580f4eSBarry Smith the first of size 2 and the second of value `MPI_COMM_NULL` since the rank 1 does not participate 452d5a8f7e6SBarry Smith in the coarse grid solve. 453d5a8f7e6SBarry Smith 454*f1580f4eSBarry Smith Since each coarser level may have a new `MPI_Comm` with fewer ranks than the previous, one 455d5a8f7e6SBarry Smith must take special care in providing the restriction and interpolation operation. We recommend 456d5a8f7e6SBarry Smith providing these as two step operations; first perform a standard restriction or interpolation on 457d5a8f7e6SBarry Smith the full number of ranks for that level and then use an MPI call to copy the resulting vector 45805552d4bSJunchao Zhang array entries (after calls to VecGetArray()) to the smaller or larger number of ranks, note in both 459d5a8f7e6SBarry Smith cases the MPI calls must be made on the larger of the two communicators. Traditional MPI send and 460d5a8f7e6SBarry Smith recieves or MPI_AlltoAllv() could be used to do the reshuffling of the vector entries. 461d5a8f7e6SBarry Smith 462*f1580f4eSBarry Smith Fortran Note: 463*f1580f4eSBarry Smith Use comms = `PETSC_NULL_MPI_COMM` as the equivalent of NULL in the C interface. Note `PETSC_NULL_MPI_COMM` 464*f1580f4eSBarry Smith is not `MPI_COMM_NULL`. It is more like `PETSC_NULL_INTEGER`, `PETSC_NULL_REAL` etc. 465d5a8f7e6SBarry Smith 466db781477SPatrick Sanan .seealso: `PCMGSetType()`, `PCMGGetLevels()` 46797d33e41SMatthew G. Knepley @*/ 4689371c9d4SSatish Balay PetscErrorCode PCMGSetLevels(PC pc, PetscInt levels, MPI_Comm *comms) { 46997d33e41SMatthew G. Knepley PetscFunctionBegin; 47097d33e41SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 47197d33e41SMatthew G. Knepley if (comms) PetscValidPointer(comms, 3); 472cac4c232SBarry Smith PetscTryMethod(pc, "PCMGSetLevels_C", (PC, PetscInt, MPI_Comm *), (pc, levels, comms)); 47397d33e41SMatthew G. Knepley PetscFunctionReturn(0); 47497d33e41SMatthew G. Knepley } 47597d33e41SMatthew G. Knepley 4769371c9d4SSatish Balay PetscErrorCode PCDestroy_MG(PC pc) { 477a06653b4SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 478a06653b4SBarry Smith PC_MG_Levels **mglevels = mg->levels; 479a06653b4SBarry Smith PetscInt i, n; 480c07bf074SBarry Smith 481c07bf074SBarry Smith PetscFunctionBegin; 4829566063dSJacob Faibussowitsch PetscCall(PCReset_MG(pc)); 483a06653b4SBarry Smith if (mglevels) { 484a06653b4SBarry Smith n = mglevels[0]->levels; 485a06653b4SBarry Smith for (i = 0; i < n; i++) { 48648a46eb9SPierre Jolivet if (mglevels[i]->smoothd != mglevels[i]->smoothu) PetscCall(KSPDestroy(&mglevels[i]->smoothd)); 4879566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&mglevels[i]->smoothu)); 4889566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&mglevels[i]->cr)); 4899566063dSJacob Faibussowitsch PetscCall(PetscFree(mglevels[i])); 490a06653b4SBarry Smith } 4919566063dSJacob Faibussowitsch PetscCall(PetscFree(mg->levels)); 492a06653b4SBarry Smith } 4939566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 4949566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGetInterpolations_C", NULL)); 4959566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGetCoarseOperators_C", NULL)); 4962b3cbbdaSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetGalerkin_C", NULL)); 4972b3cbbdaSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetLevels_C", NULL)); 4982b3cbbdaSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetLevels_C", NULL)); 4992b3cbbdaSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGetInterpolations_C", NULL)); 5002b3cbbdaSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGetCoarseOperators_C", NULL)); 5012b3cbbdaSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetAdaptInterpolation_C", NULL)); 5022b3cbbdaSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetAdaptInterpolation_C", NULL)); 5032b3cbbdaSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetAdaptCR_C", NULL)); 5042b3cbbdaSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetAdaptCR_C", NULL)); 5052b3cbbdaSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetAdaptCoarseSpaceType_C", NULL)); 5062b3cbbdaSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetAdaptCoarseSpaceType_C", NULL)); 507f3fbd535SBarry Smith PetscFunctionReturn(0); 508f3fbd535SBarry Smith } 509f3fbd535SBarry Smith 510f3fbd535SBarry Smith /* 511f3fbd535SBarry Smith PCApply_MG - Runs either an additive, multiplicative, Kaskadic 512f3fbd535SBarry Smith or full cycle of multigrid. 513f3fbd535SBarry Smith 514f3fbd535SBarry Smith Note: 515f3fbd535SBarry Smith A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle(). 516f3fbd535SBarry Smith */ 5179371c9d4SSatish Balay static PetscErrorCode PCApply_MG_Internal(PC pc, Vec b, Vec x, Mat B, Mat X, PetscBool transpose) { 518f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 519f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 520391689abSStefano Zampini PC tpc; 521f3fbd535SBarry Smith PetscInt levels = mglevels[0]->levels, i; 52230b0564aSStefano Zampini PetscBool changeu, changed, matapp; 523f3fbd535SBarry Smith 524f3fbd535SBarry Smith PetscFunctionBegin; 52530b0564aSStefano Zampini matapp = (PetscBool)(B && X); 5269566063dSJacob Faibussowitsch if (mg->stageApply) PetscCall(PetscLogStagePush(mg->stageApply)); 527e1d8e5deSBarry Smith /* When the DM is supplying the matrix then it will not exist until here */ 528298cc208SBarry Smith for (i = 0; i < levels; i++) { 529e1d8e5deSBarry Smith if (!mglevels[i]->A) { 5309566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[i]->smoothu, &mglevels[i]->A, NULL)); 5319566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mglevels[i]->A)); 532e1d8e5deSBarry Smith } 533e1d8e5deSBarry Smith } 534e1d8e5deSBarry Smith 5359566063dSJacob Faibussowitsch PetscCall(KSPGetPC(mglevels[levels - 1]->smoothd, &tpc)); 5369566063dSJacob Faibussowitsch PetscCall(PCPreSolveChangeRHS(tpc, &changed)); 5379566063dSJacob Faibussowitsch PetscCall(KSPGetPC(mglevels[levels - 1]->smoothu, &tpc)); 5389566063dSJacob Faibussowitsch PetscCall(PCPreSolveChangeRHS(tpc, &changeu)); 539391689abSStefano Zampini if (!changeu && !changed) { 54030b0564aSStefano Zampini if (matapp) { 5419566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[levels - 1]->B)); 54230b0564aSStefano Zampini mglevels[levels - 1]->B = B; 54330b0564aSStefano Zampini } else { 5449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[levels - 1]->b)); 545f3fbd535SBarry Smith mglevels[levels - 1]->b = b; 54630b0564aSStefano Zampini } 547391689abSStefano Zampini } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */ 54830b0564aSStefano Zampini if (matapp) { 54930b0564aSStefano Zampini if (mglevels[levels - 1]->B) { 55030b0564aSStefano Zampini PetscInt N1, N2; 55130b0564aSStefano Zampini PetscBool flg; 55230b0564aSStefano Zampini 5539566063dSJacob Faibussowitsch PetscCall(MatGetSize(mglevels[levels - 1]->B, NULL, &N1)); 5549566063dSJacob Faibussowitsch PetscCall(MatGetSize(B, NULL, &N2)); 5559566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)mglevels[levels - 1]->B, ((PetscObject)B)->type_name, &flg)); 55648a46eb9SPierre Jolivet if (N1 != N2 || !flg) PetscCall(MatDestroy(&mglevels[levels - 1]->B)); 55730b0564aSStefano Zampini } 55830b0564aSStefano Zampini if (!mglevels[levels - 1]->B) { 5599566063dSJacob Faibussowitsch PetscCall(MatDuplicate(B, MAT_COPY_VALUES, &mglevels[levels - 1]->B)); 56030b0564aSStefano Zampini } else { 5619566063dSJacob Faibussowitsch PetscCall(MatCopy(B, mglevels[levels - 1]->B, SAME_NONZERO_PATTERN)); 56230b0564aSStefano Zampini } 56330b0564aSStefano Zampini } else { 564391689abSStefano Zampini if (!mglevels[levels - 1]->b) { 565391689abSStefano Zampini Vec *vec; 566391689abSStefano Zampini 5679566063dSJacob Faibussowitsch PetscCall(KSPCreateVecs(mglevels[levels - 1]->smoothd, 1, &vec, 0, NULL)); 568391689abSStefano Zampini mglevels[levels - 1]->b = *vec; 5699566063dSJacob Faibussowitsch PetscCall(PetscFree(vec)); 570391689abSStefano Zampini } 5719566063dSJacob Faibussowitsch PetscCall(VecCopy(b, mglevels[levels - 1]->b)); 572391689abSStefano Zampini } 57330b0564aSStefano Zampini } 5749371c9d4SSatish Balay if (matapp) { 5759371c9d4SSatish Balay mglevels[levels - 1]->X = X; 5769371c9d4SSatish Balay } else { 5779371c9d4SSatish Balay mglevels[levels - 1]->x = x; 5789371c9d4SSatish Balay } 57930b0564aSStefano Zampini 58030b0564aSStefano Zampini /* If coarser Xs are present, it means we have already block applied the PC at least once 58130b0564aSStefano Zampini Reset operators if sizes/type do no match */ 58230b0564aSStefano Zampini if (matapp && levels > 1 && mglevels[levels - 2]->X) { 58330b0564aSStefano Zampini PetscInt Xc, Bc; 58430b0564aSStefano Zampini PetscBool flg; 58530b0564aSStefano Zampini 5869566063dSJacob Faibussowitsch PetscCall(MatGetSize(mglevels[levels - 2]->X, NULL, &Xc)); 5879566063dSJacob Faibussowitsch PetscCall(MatGetSize(mglevels[levels - 1]->B, NULL, &Bc)); 5889566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)mglevels[levels - 2]->X, ((PetscObject)mglevels[levels - 1]->X)->type_name, &flg)); 58930b0564aSStefano Zampini if (Xc != Bc || !flg) { 5909566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[levels - 1]->R)); 59130b0564aSStefano Zampini for (i = 0; i < levels - 1; i++) { 5929566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[i]->R)); 5939566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[i]->B)); 5949566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[i]->X)); 59530b0564aSStefano Zampini } 59630b0564aSStefano Zampini } 59730b0564aSStefano Zampini } 598391689abSStefano Zampini 59931567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 6009566063dSJacob Faibussowitsch if (matapp) PetscCall(MatZeroEntries(X)); 6019566063dSJacob Faibussowitsch else PetscCall(VecZeroEntries(x)); 60248a46eb9SPierre Jolivet for (i = 0; i < mg->cyclesperpcapply; i++) PetscCall(PCMGMCycle_Private(pc, mglevels + levels - 1, transpose, matapp, NULL)); 6032fa5cd67SKarl Rupp } else if (mg->am == PC_MG_ADDITIVE) { 6049566063dSJacob Faibussowitsch PetscCall(PCMGACycle_Private(pc, mglevels, transpose, matapp)); 6052fa5cd67SKarl Rupp } else if (mg->am == PC_MG_KASKADE) { 6069566063dSJacob Faibussowitsch PetscCall(PCMGKCycle_Private(pc, mglevels, transpose, matapp)); 6072fa5cd67SKarl Rupp } else { 6089566063dSJacob Faibussowitsch PetscCall(PCMGFCycle_Private(pc, mglevels, transpose, matapp)); 609f3fbd535SBarry Smith } 6109566063dSJacob Faibussowitsch if (mg->stageApply) PetscCall(PetscLogStagePop()); 61130b0564aSStefano Zampini if (!changeu && !changed) { 6129371c9d4SSatish Balay if (matapp) { 6139371c9d4SSatish Balay mglevels[levels - 1]->B = NULL; 6149371c9d4SSatish Balay } else { 6159371c9d4SSatish Balay mglevels[levels - 1]->b = NULL; 6169371c9d4SSatish Balay } 61730b0564aSStefano Zampini } 618f3fbd535SBarry Smith PetscFunctionReturn(0); 619f3fbd535SBarry Smith } 620f3fbd535SBarry Smith 6219371c9d4SSatish Balay static PetscErrorCode PCApply_MG(PC pc, Vec b, Vec x) { 622fcb023d4SJed Brown PetscFunctionBegin; 6239566063dSJacob Faibussowitsch PetscCall(PCApply_MG_Internal(pc, b, x, NULL, NULL, PETSC_FALSE)); 624fcb023d4SJed Brown PetscFunctionReturn(0); 625fcb023d4SJed Brown } 626fcb023d4SJed Brown 6279371c9d4SSatish Balay static PetscErrorCode PCApplyTranspose_MG(PC pc, Vec b, Vec x) { 628fcb023d4SJed Brown PetscFunctionBegin; 6299566063dSJacob Faibussowitsch PetscCall(PCApply_MG_Internal(pc, b, x, NULL, NULL, PETSC_TRUE)); 63030b0564aSStefano Zampini PetscFunctionReturn(0); 63130b0564aSStefano Zampini } 63230b0564aSStefano Zampini 6339371c9d4SSatish Balay static PetscErrorCode PCMatApply_MG(PC pc, Mat b, Mat x) { 63430b0564aSStefano Zampini PetscFunctionBegin; 6359566063dSJacob Faibussowitsch PetscCall(PCApply_MG_Internal(pc, NULL, NULL, b, x, PETSC_FALSE)); 636fcb023d4SJed Brown PetscFunctionReturn(0); 637fcb023d4SJed Brown } 638f3fbd535SBarry Smith 6399371c9d4SSatish Balay PetscErrorCode PCSetFromOptions_MG(PC pc, PetscOptionItems *PetscOptionsObject) { 640710315b6SLawrence Mitchell PetscInt levels, cycles; 641f3b08a26SMatthew G. Knepley PetscBool flg, flg2; 642f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 6433d3eaba7SBarry Smith PC_MG_Levels **mglevels; 644f3fbd535SBarry Smith PCMGType mgtype; 645f3fbd535SBarry Smith PCMGCycleType mgctype; 6462134b1e4SBarry Smith PCMGGalerkinType gtype; 6472b3cbbdaSStefano Zampini PCMGCoarseSpaceType coarseSpaceType; 648f3fbd535SBarry Smith 649f3fbd535SBarry Smith PetscFunctionBegin; 6501d743356SStefano Zampini levels = PetscMax(mg->nlevels, 1); 651d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Multigrid options"); 6529566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_mg_levels", "Number of Levels", "PCMGSetLevels", levels, &levels, &flg)); 6531a97d7b6SStefano Zampini if (!flg && !mg->levels && pc->dm) { 6549566063dSJacob Faibussowitsch PetscCall(DMGetRefineLevel(pc->dm, &levels)); 655eb3f98d2SBarry Smith levels++; 6563aeaf226SBarry Smith mg->usedmfornumberoflevels = PETSC_TRUE; 657eb3f98d2SBarry Smith } 6589566063dSJacob Faibussowitsch PetscCall(PCMGSetLevels(pc, levels, NULL)); 6593aeaf226SBarry Smith mglevels = mg->levels; 6603aeaf226SBarry Smith 661f3fbd535SBarry Smith mgctype = (PCMGCycleType)mglevels[0]->cycles; 6629566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_mg_cycle_type", "V cycle or for W-cycle", "PCMGSetCycleType", PCMGCycleTypes, (PetscEnum)mgctype, (PetscEnum *)&mgctype, &flg)); 6631baa6e33SBarry Smith if (flg) PetscCall(PCMGSetCycleType(pc, mgctype)); 6642134b1e4SBarry Smith gtype = mg->galerkin; 6659566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_mg_galerkin", "Use Galerkin process to compute coarser operators", "PCMGSetGalerkin", PCMGGalerkinTypes, (PetscEnum)gtype, (PetscEnum *)>ype, &flg)); 6661baa6e33SBarry Smith if (flg) PetscCall(PCMGSetGalerkin(pc, gtype)); 6672b3cbbdaSStefano Zampini coarseSpaceType = mg->coarseSpaceType; 6682b3cbbdaSStefano 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)); 6692b3cbbdaSStefano Zampini if (flg) PetscCall(PCMGSetAdaptCoarseSpaceType(pc, coarseSpaceType)); 6709566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_mg_adapt_interp_n", "Size of the coarse space for adaptive interpolation", "PCMGSetCoarseSpace", mg->Nc, &mg->Nc, &flg)); 6719566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_mg_mesp_monitor", "Monitor the multilevel eigensolver", "PCMGSetAdaptInterpolation", PETSC_FALSE, &mg->mespMonitor, &flg)); 67241b6fd38SMatthew G. Knepley flg2 = PETSC_FALSE; 6739566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_mg_adapt_cr", "Monitor coarse space quality using Compatible Relaxation (CR)", "PCMGSetAdaptCR", PETSC_FALSE, &flg2, &flg)); 6749566063dSJacob Faibussowitsch if (flg) PetscCall(PCMGSetAdaptCR(pc, flg2)); 675f56b1265SBarry Smith flg = PETSC_FALSE; 6769566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_mg_distinct_smoothup", "Create separate smoothup KSP and append the prefix _up", "PCMGSetDistinctSmoothUp", PETSC_FALSE, &flg, NULL)); 6771baa6e33SBarry Smith if (flg) PetscCall(PCMGSetDistinctSmoothUp(pc)); 67831567311SBarry Smith mgtype = mg->am; 6799566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_mg_type", "Multigrid type", "PCMGSetType", PCMGTypes, (PetscEnum)mgtype, (PetscEnum *)&mgtype, &flg)); 6801baa6e33SBarry Smith if (flg) PetscCall(PCMGSetType(pc, mgtype)); 68131567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 6829566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_mg_multiplicative_cycles", "Number of cycles for each preconditioner step", "PCMGMultiplicativeSetCycles", mg->cyclesperpcapply, &cycles, &flg)); 6831baa6e33SBarry Smith if (flg) PetscCall(PCMGMultiplicativeSetCycles(pc, cycles)); 684f3fbd535SBarry Smith } 685f3fbd535SBarry Smith flg = PETSC_FALSE; 6869566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_mg_log", "Log times for each multigrid level", "None", flg, &flg, NULL)); 687f3fbd535SBarry Smith if (flg) { 688f3fbd535SBarry Smith PetscInt i; 689f3fbd535SBarry Smith char eventname[128]; 6901a97d7b6SStefano Zampini 691f3fbd535SBarry Smith levels = mglevels[0]->levels; 692f3fbd535SBarry Smith for (i = 0; i < levels; i++) { 693f3fbd535SBarry Smith sprintf(eventname, "MGSetup Level %d", (int)i); 6949566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(eventname, ((PetscObject)pc)->classid, &mglevels[i]->eventsmoothsetup)); 695f3fbd535SBarry Smith sprintf(eventname, "MGSmooth Level %d", (int)i); 6969566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(eventname, ((PetscObject)pc)->classid, &mglevels[i]->eventsmoothsolve)); 697f3fbd535SBarry Smith if (i) { 698f3fbd535SBarry Smith sprintf(eventname, "MGResid Level %d", (int)i); 6999566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(eventname, ((PetscObject)pc)->classid, &mglevels[i]->eventresidual)); 700f3fbd535SBarry Smith sprintf(eventname, "MGInterp Level %d", (int)i); 7019566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(eventname, ((PetscObject)pc)->classid, &mglevels[i]->eventinterprestrict)); 702f3fbd535SBarry Smith } 703f3fbd535SBarry Smith } 704eec35531SMatthew G Knepley 705ec60ca73SSatish Balay #if defined(PETSC_USE_LOG) 706eec35531SMatthew G Knepley { 707eec35531SMatthew G Knepley const char *sname = "MG Apply"; 708eec35531SMatthew G Knepley PetscStageLog stageLog; 709eec35531SMatthew G Knepley PetscInt st; 710eec35531SMatthew G Knepley 7119566063dSJacob Faibussowitsch PetscCall(PetscLogGetStageLog(&stageLog)); 712eec35531SMatthew G Knepley for (st = 0; st < stageLog->numStages; ++st) { 713eec35531SMatthew G Knepley PetscBool same; 714eec35531SMatthew G Knepley 7159566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(stageLog->stageInfo[st].name, sname, &same)); 7162fa5cd67SKarl Rupp if (same) mg->stageApply = st; 717eec35531SMatthew G Knepley } 71848a46eb9SPierre Jolivet if (!mg->stageApply) PetscCall(PetscLogStageRegister(sname, &mg->stageApply)); 719eec35531SMatthew G Knepley } 720ec60ca73SSatish Balay #endif 721f3fbd535SBarry Smith } 722d0609cedSBarry Smith PetscOptionsHeadEnd(); 723f3b08a26SMatthew G. Knepley /* Check option consistency */ 7249566063dSJacob Faibussowitsch PetscCall(PCMGGetGalerkin(pc, >ype)); 7259566063dSJacob Faibussowitsch PetscCall(PCMGGetAdaptInterpolation(pc, &flg)); 72608401ef6SPierre Jolivet PetscCheck(!flg || !(gtype >= PC_MG_GALERKIN_NONE), PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Must use Galerkin coarse operators when adapting the interpolator"); 727f3fbd535SBarry Smith PetscFunctionReturn(0); 728f3fbd535SBarry Smith } 729f3fbd535SBarry Smith 7300a545947SLisandro Dalcin const char *const PCMGTypes[] = {"MULTIPLICATIVE", "ADDITIVE", "FULL", "KASKADE", "PCMGType", "PC_MG", NULL}; 7310a545947SLisandro Dalcin const char *const PCMGCycleTypes[] = {"invalid", "v", "w", "PCMGCycleType", "PC_MG_CYCLE", NULL}; 7320a545947SLisandro Dalcin const char *const PCMGGalerkinTypes[] = {"both", "pmat", "mat", "none", "external", "PCMGGalerkinType", "PC_MG_GALERKIN", NULL}; 7332b3cbbdaSStefano Zampini const char *const PCMGCoarseSpaceTypes[] = {"none", "polynomial", "harmonic", "eigenvector", "generalized_eigenvector", "gdsw", "PCMGCoarseSpaceType", "PCMG_ADAPT_NONE", NULL}; 734f3fbd535SBarry Smith 7359804daf3SBarry Smith #include <petscdraw.h> 7369371c9d4SSatish Balay PetscErrorCode PCView_MG(PC pc, PetscViewer viewer) { 737f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 738f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 739e3deeeafSJed Brown PetscInt levels = mglevels ? mglevels[0]->levels : 0, i; 740219b1045SBarry Smith PetscBool iascii, isbinary, isdraw; 741f3fbd535SBarry Smith 742f3fbd535SBarry Smith PetscFunctionBegin; 7439566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 7449566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary)); 7459566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw)); 746f3fbd535SBarry Smith if (iascii) { 747e3deeeafSJed Brown const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown"; 74863a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " type is %s, levels=%" PetscInt_FMT " cycles=%s\n", PCMGTypes[mg->am], levels, cyclename)); 74948a46eb9SPierre Jolivet if (mg->am == PC_MG_MULTIPLICATIVE) PetscCall(PetscViewerASCIIPrintf(viewer, " Cycles per PCApply=%" PetscInt_FMT "\n", mg->cyclesperpcapply)); 7502134b1e4SBarry Smith if (mg->galerkin == PC_MG_GALERKIN_BOTH) { 7519566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Using Galerkin computed coarse grid matrices\n")); 7522134b1e4SBarry Smith } else if (mg->galerkin == PC_MG_GALERKIN_PMAT) { 7539566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Using Galerkin computed coarse grid matrices for pmat\n")); 7542134b1e4SBarry Smith } else if (mg->galerkin == PC_MG_GALERKIN_MAT) { 7559566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Using Galerkin computed coarse grid matrices for mat\n")); 7562134b1e4SBarry Smith } else if (mg->galerkin == PC_MG_GALERKIN_EXTERNAL) { 7579566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Using externally compute Galerkin coarse grid matrices\n")); 7584f66f45eSBarry Smith } else { 7599566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Not using Galerkin computed coarse grid matrices\n")); 760f3fbd535SBarry Smith } 7611baa6e33SBarry Smith if (mg->view) PetscCall((*mg->view)(pc, viewer)); 762f3fbd535SBarry Smith for (i = 0; i < levels; i++) { 76363a3b9bcSJacob Faibussowitsch if (i) { 76463a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Down solver (pre-smoother) on level %" PetscInt_FMT " -------------------------------\n", i)); 765f3fbd535SBarry Smith } else { 76663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Coarse grid solver -- level %" PetscInt_FMT " -------------------------------\n", i)); 767f3fbd535SBarry Smith } 7689566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 7699566063dSJacob Faibussowitsch PetscCall(KSPView(mglevels[i]->smoothd, viewer)); 7709566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 771f3fbd535SBarry Smith if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) { 7729566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Up solver (post-smoother) same as down solver (pre-smoother)\n")); 773f3fbd535SBarry Smith } else if (i) { 77463a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Up solver (post-smoother) on level %" PetscInt_FMT " -------------------------------\n", i)); 7759566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 7769566063dSJacob Faibussowitsch PetscCall(KSPView(mglevels[i]->smoothu, viewer)); 7779566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 778f3fbd535SBarry Smith } 77941b6fd38SMatthew G. Knepley if (i && mglevels[i]->cr) { 78063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "CR solver on level %" PetscInt_FMT " -------------------------------\n", i)); 7819566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 7829566063dSJacob Faibussowitsch PetscCall(KSPView(mglevels[i]->cr, viewer)); 7839566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 78441b6fd38SMatthew G. Knepley } 785f3fbd535SBarry Smith } 7865b0b0462SBarry Smith } else if (isbinary) { 7875b0b0462SBarry Smith for (i = levels - 1; i >= 0; i--) { 7889566063dSJacob Faibussowitsch PetscCall(KSPView(mglevels[i]->smoothd, viewer)); 78948a46eb9SPierre Jolivet if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) PetscCall(KSPView(mglevels[i]->smoothu, viewer)); 7905b0b0462SBarry Smith } 791219b1045SBarry Smith } else if (isdraw) { 792219b1045SBarry Smith PetscDraw draw; 793b4375e8dSPeter Brune PetscReal x, w, y, bottom, th; 7949566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw)); 7959566063dSJacob Faibussowitsch PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y)); 7969566063dSJacob Faibussowitsch PetscCall(PetscDrawStringGetSize(draw, NULL, &th)); 797219b1045SBarry Smith bottom = y - th; 798219b1045SBarry Smith for (i = levels - 1; i >= 0; i--) { 799b4375e8dSPeter Brune if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) { 8009566063dSJacob Faibussowitsch PetscCall(PetscDrawPushCurrentPoint(draw, x, bottom)); 8019566063dSJacob Faibussowitsch PetscCall(KSPView(mglevels[i]->smoothd, viewer)); 8029566063dSJacob Faibussowitsch PetscCall(PetscDrawPopCurrentPoint(draw)); 803b4375e8dSPeter Brune } else { 804b4375e8dSPeter Brune w = 0.5 * PetscMin(1.0 - x, x); 8059566063dSJacob Faibussowitsch PetscCall(PetscDrawPushCurrentPoint(draw, x + w, bottom)); 8069566063dSJacob Faibussowitsch PetscCall(KSPView(mglevels[i]->smoothd, viewer)); 8079566063dSJacob Faibussowitsch PetscCall(PetscDrawPopCurrentPoint(draw)); 8089566063dSJacob Faibussowitsch PetscCall(PetscDrawPushCurrentPoint(draw, x - w, bottom)); 8099566063dSJacob Faibussowitsch PetscCall(KSPView(mglevels[i]->smoothu, viewer)); 8109566063dSJacob Faibussowitsch PetscCall(PetscDrawPopCurrentPoint(draw)); 811b4375e8dSPeter Brune } 8129566063dSJacob Faibussowitsch PetscCall(PetscDrawGetBoundingBox(draw, NULL, &bottom, NULL, NULL)); 8131cd381d6SBarry Smith bottom -= th; 814219b1045SBarry Smith } 8155b0b0462SBarry Smith } 816f3fbd535SBarry Smith PetscFunctionReturn(0); 817f3fbd535SBarry Smith } 818f3fbd535SBarry Smith 819af0996ceSBarry Smith #include <petsc/private/kspimpl.h> 820cab2e9ccSBarry Smith 821f3fbd535SBarry Smith /* 822f3fbd535SBarry Smith Calls setup for the KSP on each level 823f3fbd535SBarry Smith */ 8249371c9d4SSatish Balay PetscErrorCode PCSetUp_MG(PC pc) { 825f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 826f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 8271a97d7b6SStefano Zampini PetscInt i, n; 82898e04984SBarry Smith PC cpc; 8293db492dfSBarry Smith PetscBool dump = PETSC_FALSE, opsset, use_amat, missinginterpolate = PETSC_FALSE; 830f3fbd535SBarry Smith Mat dA, dB; 831f3fbd535SBarry Smith Vec tvec; 832218a07d4SBarry Smith DM *dms; 8330a545947SLisandro Dalcin PetscViewer viewer = NULL; 83441b6fd38SMatthew G. Knepley PetscBool dAeqdB = PETSC_FALSE, needRestricts = PETSC_FALSE, doCR = PETSC_FALSE; 8352b3cbbdaSStefano Zampini PetscBool adaptInterpolation = mg->adaptInterpolation; 836f3fbd535SBarry Smith 837f3fbd535SBarry Smith PetscFunctionBegin; 83828b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels with PCMGSetLevels() before setting up"); 8391a97d7b6SStefano Zampini n = mglevels[0]->levels; 84001bc414fSDmitry Karpeev /* FIX: Move this to PCSetFromOptions_MG? */ 8413aeaf226SBarry Smith if (mg->usedmfornumberoflevels) { 8423aeaf226SBarry Smith PetscInt levels; 8439566063dSJacob Faibussowitsch PetscCall(DMGetRefineLevel(pc->dm, &levels)); 8443aeaf226SBarry Smith levels++; 8453aeaf226SBarry Smith if (levels > n) { /* the problem is now being solved on a finer grid */ 8469566063dSJacob Faibussowitsch PetscCall(PCMGSetLevels(pc, levels, NULL)); 8473aeaf226SBarry Smith n = levels; 8489566063dSJacob Faibussowitsch PetscCall(PCSetFromOptions(pc)); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */ 8493aeaf226SBarry Smith mglevels = mg->levels; 8503aeaf226SBarry Smith } 8513aeaf226SBarry Smith } 8529566063dSJacob Faibussowitsch PetscCall(KSPGetPC(mglevels[0]->smoothd, &cpc)); 8533aeaf226SBarry Smith 854f3fbd535SBarry Smith /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */ 855f3fbd535SBarry Smith /* so use those from global PC */ 856f3fbd535SBarry Smith /* Is this what we always want? What if user wants to keep old one? */ 8579566063dSJacob Faibussowitsch PetscCall(KSPGetOperatorsSet(mglevels[n - 1]->smoothd, NULL, &opsset)); 85898e04984SBarry Smith if (opsset) { 85998e04984SBarry Smith Mat mmat; 8609566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[n - 1]->smoothd, NULL, &mmat)); 86198e04984SBarry Smith if (mmat == pc->pmat) opsset = PETSC_FALSE; 86298e04984SBarry Smith } 8634a5f07a5SMark F. Adams 86441b6fd38SMatthew G. Knepley /* Create CR solvers */ 8659566063dSJacob Faibussowitsch PetscCall(PCMGGetAdaptCR(pc, &doCR)); 86641b6fd38SMatthew G. Knepley if (doCR) { 86741b6fd38SMatthew G. Knepley const char *prefix; 86841b6fd38SMatthew G. Knepley 8699566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &prefix)); 87041b6fd38SMatthew G. Knepley for (i = 1; i < n; ++i) { 87141b6fd38SMatthew G. Knepley PC ipc, cr; 87241b6fd38SMatthew G. Knepley char crprefix[128]; 87341b6fd38SMatthew G. Knepley 8749566063dSJacob Faibussowitsch PetscCall(KSPCreate(PetscObjectComm((PetscObject)pc), &mglevels[i]->cr)); 8759566063dSJacob Faibussowitsch PetscCall(KSPSetErrorIfNotConverged(mglevels[i]->cr, PETSC_FALSE)); 8769566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->cr, (PetscObject)pc, n - i)); 8779566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(mglevels[i]->cr, prefix)); 8789566063dSJacob Faibussowitsch PetscCall(PetscObjectComposedDataSetInt((PetscObject)mglevels[i]->cr, PetscMGLevelId, mglevels[i]->level)); 8799566063dSJacob Faibussowitsch PetscCall(KSPSetType(mglevels[i]->cr, KSPCHEBYSHEV)); 8809566063dSJacob Faibussowitsch PetscCall(KSPSetConvergenceTest(mglevels[i]->cr, KSPConvergedSkip, NULL, NULL)); 8819566063dSJacob Faibussowitsch PetscCall(KSPSetNormType(mglevels[i]->cr, KSP_NORM_PRECONDITIONED)); 8829566063dSJacob Faibussowitsch PetscCall(KSPGetPC(mglevels[i]->cr, &ipc)); 88341b6fd38SMatthew G. Knepley 8849566063dSJacob Faibussowitsch PetscCall(PCSetType(ipc, PCCOMPOSITE)); 8859566063dSJacob Faibussowitsch PetscCall(PCCompositeSetType(ipc, PC_COMPOSITE_MULTIPLICATIVE)); 8869566063dSJacob Faibussowitsch PetscCall(PCCompositeAddPCType(ipc, PCSOR)); 8879566063dSJacob Faibussowitsch PetscCall(CreateCR_Private(pc, i, &cr)); 8889566063dSJacob Faibussowitsch PetscCall(PCCompositeAddPC(ipc, cr)); 8899566063dSJacob Faibussowitsch PetscCall(PCDestroy(&cr)); 89041b6fd38SMatthew G. Knepley 8919566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(mglevels[i]->cr, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, mg->default_smoothd)); 8929566063dSJacob Faibussowitsch PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->cr, PETSC_TRUE)); 8939566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(crprefix, 128, "mg_levels_%d_cr_", (int)i)); 8949566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(mglevels[i]->cr, crprefix)); 8959566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)mglevels[i]->cr)); 89641b6fd38SMatthew G. Knepley } 89741b6fd38SMatthew G. Knepley } 89841b6fd38SMatthew G. Knepley 8994a5f07a5SMark F. Adams if (!opsset) { 9009566063dSJacob Faibussowitsch PetscCall(PCGetUseAmat(pc, &use_amat)); 9014a5f07a5SMark F. Adams if (use_amat) { 9029566063dSJacob 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")); 9039566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(mglevels[n - 1]->smoothd, pc->mat, pc->pmat)); 90469858f1bSStefano Zampini } else { 9059566063dSJacob 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")); 9069566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(mglevels[n - 1]->smoothd, pc->pmat, pc->pmat)); 9074a5f07a5SMark F. Adams } 9084a5f07a5SMark F. Adams } 909f3fbd535SBarry Smith 9109c7ffb39SBarry Smith for (i = n - 1; i > 0; i--) { 9119c7ffb39SBarry Smith if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) { 9129c7ffb39SBarry Smith missinginterpolate = PETSC_TRUE; 9132b3cbbdaSStefano Zampini break; 9149c7ffb39SBarry Smith } 9159c7ffb39SBarry Smith } 9162134b1e4SBarry Smith 9179566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[n - 1]->smoothd, &dA, &dB)); 9182134b1e4SBarry Smith if (dA == dB) dAeqdB = PETSC_TRUE; 9192b3cbbdaSStefano Zampini if (mg->galerkin == PC_MG_GALERKIN_NONE || ((mg->galerkin == PC_MG_GALERKIN_PMAT || mg->galerkin == PC_MG_GALERKIN_MAT) && !dAeqdB)) { 9202134b1e4SBarry Smith needRestricts = PETSC_TRUE; /* user must compute either mat, pmat, or both so must restrict x to coarser levels */ 9212134b1e4SBarry Smith } 9222134b1e4SBarry Smith 9232b3cbbdaSStefano Zampini if (pc->dm && !pc->setupcalled) { 9242b3cbbdaSStefano Zampini /* finest smoother also gets DM but it is not active, independent of whether galerkin==PC_MG_GALERKIN_EXTERNAL */ 9252b3cbbdaSStefano Zampini PetscCall(KSPSetDM(mglevels[n - 1]->smoothd, pc->dm)); 9262b3cbbdaSStefano Zampini PetscCall(KSPSetDMActive(mglevels[n - 1]->smoothd, PETSC_FALSE)); 9272b3cbbdaSStefano Zampini if (mglevels[n - 1]->smoothd != mglevels[n - 1]->smoothu) { 9282b3cbbdaSStefano Zampini PetscCall(KSPSetDM(mglevels[n - 1]->smoothu, pc->dm)); 9292b3cbbdaSStefano Zampini PetscCall(KSPSetDMActive(mglevels[n - 1]->smoothu, PETSC_FALSE)); 9302b3cbbdaSStefano Zampini } 9312b3cbbdaSStefano Zampini if (mglevels[n - 1]->cr) { 9322b3cbbdaSStefano Zampini PetscCall(KSPSetDM(mglevels[n - 1]->cr, pc->dm)); 9332b3cbbdaSStefano Zampini PetscCall(KSPSetDMActive(mglevels[n - 1]->cr, PETSC_FALSE)); 9342b3cbbdaSStefano Zampini } 9352b3cbbdaSStefano Zampini } 9362b3cbbdaSStefano Zampini 9379c7ffb39SBarry Smith /* 9389c7ffb39SBarry Smith Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS) 9392b3cbbdaSStefano Zampini Skipping for externally managed hierarchy (such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs? 9409c7ffb39SBarry Smith */ 94132fe681dSStefano Zampini if (missinginterpolate && mg->galerkin != PC_MG_GALERKIN_EXTERNAL && !pc->setupcalled) { 94232fe681dSStefano Zampini /* first see if we can compute a coarse space */ 94332fe681dSStefano Zampini if (mg->coarseSpaceType == PCMG_ADAPT_GDSW) { 94432fe681dSStefano Zampini for (i = n - 2; i > -1; i--) { 94532fe681dSStefano Zampini if (!mglevels[i + 1]->restrct && !mglevels[i + 1]->interpolate) { 94632fe681dSStefano Zampini PetscCall(PCMGComputeCoarseSpace_Internal(pc, i + 1, mg->coarseSpaceType, mg->Nc, NULL, &mglevels[i + 1]->coarseSpace)); 94732fe681dSStefano Zampini PetscCall(PCMGSetInterpolation(pc, i + 1, mglevels[i + 1]->coarseSpace)); 94832fe681dSStefano Zampini } 94932fe681dSStefano Zampini } 95032fe681dSStefano Zampini } else { /* construct the interpolation from the DMs */ 9512e499ae9SBarry Smith Mat p; 95273250ac0SBarry Smith Vec rscale; 9539566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &dms)); 954218a07d4SBarry Smith dms[n - 1] = pc->dm; 955ef1267afSMatthew G. Knepley /* Separately create them so we do not get DMKSP interference between levels */ 9569566063dSJacob Faibussowitsch for (i = n - 2; i > -1; i--) PetscCall(DMCoarsen(dms[i + 1], MPI_COMM_NULL, &dms[i])); 957218a07d4SBarry Smith for (i = n - 2; i > -1; i--) { 958942e3340SBarry Smith DMKSP kdm; 959eab52d2bSLawrence Mitchell PetscBool dmhasrestrict, dmhasinject; 9609566063dSJacob Faibussowitsch PetscCall(KSPSetDM(mglevels[i]->smoothd, dms[i])); 9619566063dSJacob Faibussowitsch if (!needRestricts) PetscCall(KSPSetDMActive(mglevels[i]->smoothd, PETSC_FALSE)); 962c27ee7a3SPatrick Farrell if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 9639566063dSJacob Faibussowitsch PetscCall(KSPSetDM(mglevels[i]->smoothu, dms[i])); 9649566063dSJacob Faibussowitsch if (!needRestricts) PetscCall(KSPSetDMActive(mglevels[i]->smoothu, PETSC_FALSE)); 965c27ee7a3SPatrick Farrell } 96641b6fd38SMatthew G. Knepley if (mglevels[i]->cr) { 9679566063dSJacob Faibussowitsch PetscCall(KSPSetDM(mglevels[i]->cr, dms[i])); 9689566063dSJacob Faibussowitsch if (!needRestricts) PetscCall(KSPSetDMActive(mglevels[i]->cr, PETSC_FALSE)); 96941b6fd38SMatthew G. Knepley } 9709566063dSJacob Faibussowitsch PetscCall(DMGetDMKSPWrite(dms[i], &kdm)); 971d1a3e677SJed Brown /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take 972d1a3e677SJed Brown * a bitwise OR of computing the matrix, RHS, and initial iterate. */ 9730298fd71SBarry Smith kdm->ops->computerhs = NULL; 9740298fd71SBarry Smith kdm->rhsctx = NULL; 97524c3aa18SBarry Smith if (!mglevels[i + 1]->interpolate) { 9769566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(dms[i], dms[i + 1], &p, &rscale)); 9779566063dSJacob Faibussowitsch PetscCall(PCMGSetInterpolation(pc, i + 1, p)); 9789566063dSJacob Faibussowitsch if (rscale) PetscCall(PCMGSetRScale(pc, i + 1, rscale)); 9799566063dSJacob Faibussowitsch PetscCall(VecDestroy(&rscale)); 9809566063dSJacob Faibussowitsch PetscCall(MatDestroy(&p)); 981218a07d4SBarry Smith } 9829566063dSJacob Faibussowitsch PetscCall(DMHasCreateRestriction(dms[i], &dmhasrestrict)); 9833ad4599aSBarry Smith if (dmhasrestrict && !mglevels[i + 1]->restrct) { 9849566063dSJacob Faibussowitsch PetscCall(DMCreateRestriction(dms[i], dms[i + 1], &p)); 9859566063dSJacob Faibussowitsch PetscCall(PCMGSetRestriction(pc, i + 1, p)); 9869566063dSJacob Faibussowitsch PetscCall(MatDestroy(&p)); 9873ad4599aSBarry Smith } 9889566063dSJacob Faibussowitsch PetscCall(DMHasCreateInjection(dms[i], &dmhasinject)); 989eab52d2bSLawrence Mitchell if (dmhasinject && !mglevels[i + 1]->inject) { 9909566063dSJacob Faibussowitsch PetscCall(DMCreateInjection(dms[i], dms[i + 1], &p)); 9919566063dSJacob Faibussowitsch PetscCall(PCMGSetInjection(pc, i + 1, p)); 9929566063dSJacob Faibussowitsch PetscCall(MatDestroy(&p)); 993eab52d2bSLawrence Mitchell } 99424c3aa18SBarry Smith } 9952d2b81a6SBarry Smith 9969566063dSJacob Faibussowitsch for (i = n - 2; i > -1; i--) PetscCall(DMDestroy(&dms[i])); 9979566063dSJacob Faibussowitsch PetscCall(PetscFree(dms)); 998ce4cda84SJed Brown } 99932fe681dSStefano Zampini } 1000cab2e9ccSBarry Smith 10012134b1e4SBarry Smith if (mg->galerkin < PC_MG_GALERKIN_NONE) { 10022134b1e4SBarry Smith Mat A, B; 10032134b1e4SBarry Smith PetscBool doA = PETSC_FALSE, doB = PETSC_FALSE; 10042134b1e4SBarry Smith MatReuse reuse = MAT_INITIAL_MATRIX; 10052134b1e4SBarry Smith 10062b3cbbdaSStefano Zampini if (mg->galerkin == PC_MG_GALERKIN_PMAT || mg->galerkin == PC_MG_GALERKIN_BOTH) doB = PETSC_TRUE; 10072b3cbbdaSStefano Zampini if (mg->galerkin == PC_MG_GALERKIN_MAT || (mg->galerkin == PC_MG_GALERKIN_BOTH && dA != dB)) doA = PETSC_TRUE; 10082134b1e4SBarry Smith if (pc->setupcalled) reuse = MAT_REUSE_MATRIX; 1009f3fbd535SBarry Smith for (i = n - 2; i > -1; i--) { 10107827d75bSBarry 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"); 101148a46eb9SPierre Jolivet if (!mglevels[i + 1]->interpolate) PetscCall(PCMGSetInterpolation(pc, i + 1, mglevels[i + 1]->restrct)); 101248a46eb9SPierre Jolivet if (!mglevels[i + 1]->restrct) PetscCall(PCMGSetRestriction(pc, i + 1, mglevels[i + 1]->interpolate)); 101348a46eb9SPierre Jolivet if (reuse == MAT_REUSE_MATRIX) PetscCall(KSPGetOperators(mglevels[i]->smoothd, &A, &B)); 101448a46eb9SPierre Jolivet if (doA) PetscCall(MatGalerkin(mglevels[i + 1]->restrct, dA, mglevels[i + 1]->interpolate, reuse, 1.0, &A)); 101548a46eb9SPierre Jolivet if (doB) PetscCall(MatGalerkin(mglevels[i + 1]->restrct, dB, mglevels[i + 1]->interpolate, reuse, 1.0, &B)); 10162134b1e4SBarry Smith /* the management of the PetscObjectReference() and PetscObjecDereference() below is rather delicate */ 10172134b1e4SBarry Smith if (!doA && dAeqdB) { 10189566063dSJacob Faibussowitsch if (reuse == MAT_INITIAL_MATRIX) PetscCall(PetscObjectReference((PetscObject)B)); 10192134b1e4SBarry Smith A = B; 10202134b1e4SBarry Smith } else if (!doA && reuse == MAT_INITIAL_MATRIX) { 10219566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[i]->smoothd, &A, NULL)); 10229566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)A)); 1023b9367914SBarry Smith } 10242134b1e4SBarry Smith if (!doB && dAeqdB) { 10259566063dSJacob Faibussowitsch if (reuse == MAT_INITIAL_MATRIX) PetscCall(PetscObjectReference((PetscObject)A)); 10262134b1e4SBarry Smith B = A; 10272134b1e4SBarry Smith } else if (!doB && reuse == MAT_INITIAL_MATRIX) { 10289566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[i]->smoothd, NULL, &B)); 10299566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)B)); 10307d537d92SJed Brown } 10312134b1e4SBarry Smith if (reuse == MAT_INITIAL_MATRIX) { 10329566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(mglevels[i]->smoothd, A, B)); 10339566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)A)); 10349566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)B)); 10352134b1e4SBarry Smith } 10362134b1e4SBarry Smith dA = A; 1037f3fbd535SBarry Smith dB = B; 1038f3fbd535SBarry Smith } 1039f3fbd535SBarry Smith } 1040f3b08a26SMatthew G. Knepley 1041f3b08a26SMatthew G. Knepley /* Adapt interpolation matrices */ 10422b3cbbdaSStefano Zampini if (adaptInterpolation) { 1043f3b08a26SMatthew G. Knepley for (i = 0; i < n; ++i) { 104448a46eb9SPierre Jolivet if (!mglevels[i]->coarseSpace) PetscCall(PCMGComputeCoarseSpace_Internal(pc, i, mg->coarseSpaceType, mg->Nc, !i ? NULL : mglevels[i - 1]->coarseSpace, &mglevels[i]->coarseSpace)); 10452b3cbbdaSStefano Zampini if (i) PetscCall(PCMGAdaptInterpolator_Internal(pc, i, mglevels[i - 1]->smoothu, mglevels[i]->smoothu, mglevels[i - 1]->coarseSpace, mglevels[i]->coarseSpace)); 1046f3b08a26SMatthew G. Knepley } 104748a46eb9SPierre Jolivet for (i = n - 2; i > -1; --i) PetscCall(PCMGRecomputeLevelOperators_Internal(pc, i)); 1048f3b08a26SMatthew G. Knepley } 1049f3b08a26SMatthew G. Knepley 10502134b1e4SBarry Smith if (needRestricts && pc->dm) { 1051caa4e7f2SJed Brown for (i = n - 2; i >= 0; i--) { 1052ccceb50cSJed Brown DM dmfine, dmcoarse; 1053caa4e7f2SJed Brown Mat Restrict, Inject; 1054caa4e7f2SJed Brown Vec rscale; 10559566063dSJacob Faibussowitsch PetscCall(KSPGetDM(mglevels[i + 1]->smoothd, &dmfine)); 10569566063dSJacob Faibussowitsch PetscCall(KSPGetDM(mglevels[i]->smoothd, &dmcoarse)); 10579566063dSJacob Faibussowitsch PetscCall(PCMGGetRestriction(pc, i + 1, &Restrict)); 10589566063dSJacob Faibussowitsch PetscCall(PCMGGetRScale(pc, i + 1, &rscale)); 10599566063dSJacob Faibussowitsch PetscCall(PCMGGetInjection(pc, i + 1, &Inject)); 10609566063dSJacob Faibussowitsch PetscCall(DMRestrict(dmfine, Restrict, rscale, Inject, dmcoarse)); 1061caa4e7f2SJed Brown } 1062caa4e7f2SJed Brown } 1063f3fbd535SBarry Smith 1064f3fbd535SBarry Smith if (!pc->setupcalled) { 106548a46eb9SPierre Jolivet for (i = 0; i < n; i++) PetscCall(KSPSetFromOptions(mglevels[i]->smoothd)); 1066f3fbd535SBarry Smith for (i = 1; i < n; i++) { 106748a46eb9SPierre Jolivet if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) PetscCall(KSPSetFromOptions(mglevels[i]->smoothu)); 106848a46eb9SPierre Jolivet if (mglevels[i]->cr) PetscCall(KSPSetFromOptions(mglevels[i]->cr)); 1069f3fbd535SBarry Smith } 10703ad4599aSBarry Smith /* insure that if either interpolation or restriction is set the other other one is set */ 1071f3fbd535SBarry Smith for (i = 1; i < n; i++) { 10729566063dSJacob Faibussowitsch PetscCall(PCMGGetInterpolation(pc, i, NULL)); 10739566063dSJacob Faibussowitsch PetscCall(PCMGGetRestriction(pc, i, NULL)); 1074f3fbd535SBarry Smith } 1075f3fbd535SBarry Smith for (i = 0; i < n - 1; i++) { 1076f3fbd535SBarry Smith if (!mglevels[i]->b) { 1077f3fbd535SBarry Smith Vec *vec; 10789566063dSJacob Faibussowitsch PetscCall(KSPCreateVecs(mglevels[i]->smoothd, 1, &vec, 0, NULL)); 10799566063dSJacob Faibussowitsch PetscCall(PCMGSetRhs(pc, i, *vec)); 10809566063dSJacob Faibussowitsch PetscCall(VecDestroy(vec)); 10819566063dSJacob Faibussowitsch PetscCall(PetscFree(vec)); 1082f3fbd535SBarry Smith } 1083f3fbd535SBarry Smith if (!mglevels[i]->r && i) { 10849566063dSJacob Faibussowitsch PetscCall(VecDuplicate(mglevels[i]->b, &tvec)); 10859566063dSJacob Faibussowitsch PetscCall(PCMGSetR(pc, i, tvec)); 10869566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tvec)); 1087f3fbd535SBarry Smith } 1088f3fbd535SBarry Smith if (!mglevels[i]->x) { 10899566063dSJacob Faibussowitsch PetscCall(VecDuplicate(mglevels[i]->b, &tvec)); 10909566063dSJacob Faibussowitsch PetscCall(PCMGSetX(pc, i, tvec)); 10919566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tvec)); 1092f3fbd535SBarry Smith } 109341b6fd38SMatthew G. Knepley if (doCR) { 10949566063dSJacob Faibussowitsch PetscCall(VecDuplicate(mglevels[i]->b, &mglevels[i]->crx)); 10959566063dSJacob Faibussowitsch PetscCall(VecDuplicate(mglevels[i]->b, &mglevels[i]->crb)); 10969566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(mglevels[i]->crb)); 109741b6fd38SMatthew G. Knepley } 1098f3fbd535SBarry Smith } 1099f3fbd535SBarry Smith if (n != 1 && !mglevels[n - 1]->r) { 1100f3fbd535SBarry Smith /* PCMGSetR() on the finest level if user did not supply it */ 1101f3fbd535SBarry Smith Vec *vec; 11029566063dSJacob Faibussowitsch PetscCall(KSPCreateVecs(mglevels[n - 1]->smoothd, 1, &vec, 0, NULL)); 11039566063dSJacob Faibussowitsch PetscCall(PCMGSetR(pc, n - 1, *vec)); 11049566063dSJacob Faibussowitsch PetscCall(VecDestroy(vec)); 11059566063dSJacob Faibussowitsch PetscCall(PetscFree(vec)); 1106f3fbd535SBarry Smith } 110741b6fd38SMatthew G. Knepley if (doCR) { 11089566063dSJacob Faibussowitsch PetscCall(VecDuplicate(mglevels[n - 1]->r, &mglevels[n - 1]->crx)); 11099566063dSJacob Faibussowitsch PetscCall(VecDuplicate(mglevels[n - 1]->r, &mglevels[n - 1]->crb)); 11109566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(mglevels[n - 1]->crb)); 111141b6fd38SMatthew G. Knepley } 1112f3fbd535SBarry Smith } 1113f3fbd535SBarry Smith 111498e04984SBarry Smith if (pc->dm) { 111598e04984SBarry Smith /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */ 111698e04984SBarry Smith for (i = 0; i < n - 1; i++) { 111798e04984SBarry Smith if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX; 111898e04984SBarry Smith } 111998e04984SBarry Smith } 112008549ed4SJed Brown // We got here (PCSetUp_MG) because the matrix has changed, which means the smoother needs to be set up again (e.g., 112108549ed4SJed Brown // new diagonal for Jacobi). Setting it here allows it to be logged under PCSetUp rather than deep inside a PCApply. 112208549ed4SJed Brown if (mglevels[n - 1]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[n - 1]->smoothd->setupstage = KSP_SETUP_NEWMATRIX; 1123f3fbd535SBarry Smith 1124f3fbd535SBarry Smith for (i = 1; i < n; i++) { 11252a21e185SBarry Smith if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1) { 1126f3fbd535SBarry Smith /* if doing only down then initial guess is zero */ 11279566063dSJacob Faibussowitsch PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->smoothd, PETSC_TRUE)); 1128f3fbd535SBarry Smith } 11299566063dSJacob Faibussowitsch if (mglevels[i]->cr) PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->cr, PETSC_TRUE)); 11309566063dSJacob Faibussowitsch if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventBegin(mglevels[i]->eventsmoothsetup, 0, 0, 0, 0)); 11319566063dSJacob Faibussowitsch PetscCall(KSPSetUp(mglevels[i]->smoothd)); 1132ad540459SPierre Jolivet if (mglevels[i]->smoothd->reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR; 11339566063dSJacob Faibussowitsch if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventEnd(mglevels[i]->eventsmoothsetup, 0, 0, 0, 0)); 1134d42688cbSBarry Smith if (!mglevels[i]->residual) { 1135d42688cbSBarry Smith Mat mat; 11369566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[i]->smoothd, &mat, NULL)); 11379566063dSJacob Faibussowitsch PetscCall(PCMGSetResidual(pc, i, PCMGResidualDefault, mat)); 1138d42688cbSBarry Smith } 1139fcb023d4SJed Brown if (!mglevels[i]->residualtranspose) { 1140fcb023d4SJed Brown Mat mat; 11419566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[i]->smoothd, &mat, NULL)); 11429566063dSJacob Faibussowitsch PetscCall(PCMGSetResidualTranspose(pc, i, PCMGResidualTransposeDefault, mat)); 1143fcb023d4SJed Brown } 1144f3fbd535SBarry Smith } 1145f3fbd535SBarry Smith for (i = 1; i < n; i++) { 1146f3fbd535SBarry Smith if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) { 1147f3fbd535SBarry Smith Mat downmat, downpmat; 1148f3fbd535SBarry Smith 1149f3fbd535SBarry Smith /* check if operators have been set for up, if not use down operators to set them */ 11509566063dSJacob Faibussowitsch PetscCall(KSPGetOperatorsSet(mglevels[i]->smoothu, &opsset, NULL)); 1151f3fbd535SBarry Smith if (!opsset) { 11529566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[i]->smoothd, &downmat, &downpmat)); 11539566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(mglevels[i]->smoothu, downmat, downpmat)); 1154f3fbd535SBarry Smith } 1155f3fbd535SBarry Smith 11569566063dSJacob Faibussowitsch PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->smoothu, PETSC_TRUE)); 11579566063dSJacob Faibussowitsch if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventBegin(mglevels[i]->eventsmoothsetup, 0, 0, 0, 0)); 11589566063dSJacob Faibussowitsch PetscCall(KSPSetUp(mglevels[i]->smoothu)); 1159ad540459SPierre Jolivet if (mglevels[i]->smoothu->reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR; 11609566063dSJacob Faibussowitsch if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventEnd(mglevels[i]->eventsmoothsetup, 0, 0, 0, 0)); 1161f3fbd535SBarry Smith } 116241b6fd38SMatthew G. Knepley if (mglevels[i]->cr) { 116341b6fd38SMatthew G. Knepley Mat downmat, downpmat; 116441b6fd38SMatthew G. Knepley 116541b6fd38SMatthew G. Knepley /* check if operators have been set for up, if not use down operators to set them */ 11669566063dSJacob Faibussowitsch PetscCall(KSPGetOperatorsSet(mglevels[i]->cr, &opsset, NULL)); 116741b6fd38SMatthew G. Knepley if (!opsset) { 11689566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[i]->smoothd, &downmat, &downpmat)); 11699566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(mglevels[i]->cr, downmat, downpmat)); 117041b6fd38SMatthew G. Knepley } 117141b6fd38SMatthew G. Knepley 11729566063dSJacob Faibussowitsch PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->cr, PETSC_TRUE)); 11739566063dSJacob Faibussowitsch if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventBegin(mglevels[i]->eventsmoothsetup, 0, 0, 0, 0)); 11749566063dSJacob Faibussowitsch PetscCall(KSPSetUp(mglevels[i]->cr)); 1175ad540459SPierre Jolivet if (mglevels[i]->cr->reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR; 11769566063dSJacob Faibussowitsch if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventEnd(mglevels[i]->eventsmoothsetup, 0, 0, 0, 0)); 117741b6fd38SMatthew G. Knepley } 1178f3fbd535SBarry Smith } 1179f3fbd535SBarry Smith 11809566063dSJacob Faibussowitsch if (mglevels[0]->eventsmoothsetup) PetscCall(PetscLogEventBegin(mglevels[0]->eventsmoothsetup, 0, 0, 0, 0)); 11819566063dSJacob Faibussowitsch PetscCall(KSPSetUp(mglevels[0]->smoothd)); 1182ad540459SPierre Jolivet if (mglevels[0]->smoothd->reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR; 11839566063dSJacob Faibussowitsch if (mglevels[0]->eventsmoothsetup) PetscCall(PetscLogEventEnd(mglevels[0]->eventsmoothsetup, 0, 0, 0, 0)); 1184f3fbd535SBarry Smith 1185f3fbd535SBarry Smith /* 1186f3fbd535SBarry Smith Dump the interpolation/restriction matrices plus the 1187e3c5b3baSBarry Smith Jacobian/stiffness on each level. This allows MATLAB users to 1188f3fbd535SBarry Smith easily check if the Galerkin condition A_c = R A_f R^T is satisfied. 1189f3fbd535SBarry Smith 1190f3fbd535SBarry Smith Only support one or the other at the same time. 1191f3fbd535SBarry Smith */ 1192f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER) 11939566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_mg_dump_matlab", &dump, NULL)); 1194ce94432eSBarry Smith if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc)); 1195f3fbd535SBarry Smith dump = PETSC_FALSE; 1196f3fbd535SBarry Smith #endif 11979566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_mg_dump_binary", &dump, NULL)); 1198ce94432eSBarry Smith if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc)); 1199f3fbd535SBarry Smith 1200f3fbd535SBarry Smith if (viewer) { 120148a46eb9SPierre Jolivet for (i = 1; i < n; i++) PetscCall(MatView(mglevels[i]->restrct, viewer)); 1202f3fbd535SBarry Smith for (i = 0; i < n; i++) { 12039566063dSJacob Faibussowitsch PetscCall(KSPGetPC(mglevels[i]->smoothd, &pc)); 12049566063dSJacob Faibussowitsch PetscCall(MatView(pc->mat, viewer)); 1205f3fbd535SBarry Smith } 1206f3fbd535SBarry Smith } 1207f3fbd535SBarry Smith PetscFunctionReturn(0); 1208f3fbd535SBarry Smith } 1209f3fbd535SBarry Smith 1210f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/ 1211f3fbd535SBarry Smith 12129371c9d4SSatish Balay PetscErrorCode PCMGGetLevels_MG(PC pc, PetscInt *levels) { 121397d33e41SMatthew G. Knepley PC_MG *mg = (PC_MG *)pc->data; 121497d33e41SMatthew G. Knepley 121597d33e41SMatthew G. Knepley PetscFunctionBegin; 121697d33e41SMatthew G. Knepley *levels = mg->nlevels; 121797d33e41SMatthew G. Knepley PetscFunctionReturn(0); 121897d33e41SMatthew G. Knepley } 121997d33e41SMatthew G. Knepley 12204b9ad928SBarry Smith /*@ 1221*f1580f4eSBarry Smith PCMGGetLevels - Gets the number of levels to use with `PCMG`. 12224b9ad928SBarry Smith 12234b9ad928SBarry Smith Not Collective 12244b9ad928SBarry Smith 12254b9ad928SBarry Smith Input Parameter: 12264b9ad928SBarry Smith . pc - the preconditioner context 12274b9ad928SBarry Smith 12284b9ad928SBarry Smith Output parameter: 12294b9ad928SBarry Smith . levels - the number of levels 12304b9ad928SBarry Smith 12314b9ad928SBarry Smith Level: advanced 12324b9ad928SBarry Smith 1233*f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetLevels()` 12344b9ad928SBarry Smith @*/ 12359371c9d4SSatish Balay PetscErrorCode PCMGGetLevels(PC pc, PetscInt *levels) { 12364b9ad928SBarry Smith PetscFunctionBegin; 12370700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 12384482741eSBarry Smith PetscValidIntPointer(levels, 2); 123997d33e41SMatthew G. Knepley *levels = 0; 1240cac4c232SBarry Smith PetscTryMethod(pc, "PCMGGetLevels_C", (PC, PetscInt *), (pc, levels)); 12414b9ad928SBarry Smith PetscFunctionReturn(0); 12424b9ad928SBarry Smith } 12434b9ad928SBarry Smith 12444b9ad928SBarry Smith /*@ 1245*f1580f4eSBarry Smith PCMGGetGridComplexity - compute operator and grid complexity of the `PCMG` hierarchy 1246e7d4b4cbSMark Adams 1247e7d4b4cbSMark Adams Input Parameter: 1248e7d4b4cbSMark Adams . pc - the preconditioner context 1249e7d4b4cbSMark Adams 125090db8557SMark Adams Output Parameters: 125190db8557SMark Adams + gc - grid complexity = sum_i(n_i) / n_0 125290db8557SMark Adams - oc - operator complexity = sum_i(nnz_i) / nnz_0 1253e7d4b4cbSMark Adams 1254e7d4b4cbSMark Adams Level: advanced 1255e7d4b4cbSMark Adams 1256*f1580f4eSBarry Smith Note: 1257*f1580f4eSBarry Smith This is often call the operator complexity in multigrid literature 1258*f1580f4eSBarry Smith 1259*f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGGetLevels()`, `PCMGSetLevels()` 1260e7d4b4cbSMark Adams @*/ 12619371c9d4SSatish Balay PetscErrorCode PCMGGetGridComplexity(PC pc, PetscReal *gc, PetscReal *oc) { 1262e7d4b4cbSMark Adams PC_MG *mg = (PC_MG *)pc->data; 1263e7d4b4cbSMark Adams PC_MG_Levels **mglevels = mg->levels; 1264e7d4b4cbSMark Adams PetscInt lev, N; 1265e7d4b4cbSMark Adams PetscLogDouble nnz0 = 0, sgc = 0, soc = 0, n0 = 0; 1266e7d4b4cbSMark Adams MatInfo info; 1267e7d4b4cbSMark Adams 1268e7d4b4cbSMark Adams PetscFunctionBegin; 126990db8557SMark Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 127090db8557SMark Adams if (gc) PetscValidRealPointer(gc, 2); 127190db8557SMark Adams if (oc) PetscValidRealPointer(oc, 3); 1272e7d4b4cbSMark Adams if (!pc->setupcalled) { 127390db8557SMark Adams if (gc) *gc = 0; 127490db8557SMark Adams if (oc) *oc = 0; 1275e7d4b4cbSMark Adams PetscFunctionReturn(0); 1276e7d4b4cbSMark Adams } 1277e7d4b4cbSMark Adams PetscCheck(mg->nlevels > 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MG has no levels"); 1278e7d4b4cbSMark Adams for (lev = 0; lev < mg->nlevels; lev++) { 1279e7d4b4cbSMark Adams Mat dB; 12809566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[lev]->smoothd, NULL, &dB)); 12819566063dSJacob Faibussowitsch PetscCall(MatGetInfo(dB, MAT_GLOBAL_SUM, &info)); /* global reduction */ 12829566063dSJacob Faibussowitsch PetscCall(MatGetSize(dB, &N, NULL)); 1283e7d4b4cbSMark Adams sgc += N; 1284e7d4b4cbSMark Adams soc += info.nz_used; 12859371c9d4SSatish Balay if (lev == mg->nlevels - 1) { 12869371c9d4SSatish Balay nnz0 = info.nz_used; 12879371c9d4SSatish Balay n0 = N; 12889371c9d4SSatish Balay } 1289e7d4b4cbSMark Adams } 129090db8557SMark Adams if (n0 > 0 && gc) *gc = (PetscReal)(sgc / n0); 1291e7d4b4cbSMark Adams else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number for grid points on finest level is not available"); 129290db8557SMark Adams if (nnz0 > 0 && oc) *oc = (PetscReal)(soc / nnz0); 1293e7d4b4cbSMark Adams PetscFunctionReturn(0); 1294e7d4b4cbSMark Adams } 1295e7d4b4cbSMark Adams 1296e7d4b4cbSMark Adams /*@ 129797177400SBarry Smith PCMGSetType - Determines the form of multigrid to use: 12984b9ad928SBarry Smith multiplicative, additive, full, or the Kaskade algorithm. 12994b9ad928SBarry Smith 1300*f1580f4eSBarry Smith Logically Collective on pc 13014b9ad928SBarry Smith 13024b9ad928SBarry Smith Input Parameters: 13034b9ad928SBarry Smith + pc - the preconditioner context 1304*f1580f4eSBarry Smith - form - multigrid form, one of `PC_MG_MULTIPLICATIVE`, `PC_MG_ADDITIVE`, `PC_MG_FULL`, `PC_MG_KASKADE` 13054b9ad928SBarry Smith 13064b9ad928SBarry Smith Options Database Key: 13074b9ad928SBarry Smith . -pc_mg_type <form> - Sets <form>, one of multiplicative, 13084b9ad928SBarry Smith additive, full, kaskade 13094b9ad928SBarry Smith 13104b9ad928SBarry Smith Level: advanced 13114b9ad928SBarry Smith 1312*f1580f4eSBarry Smith .seealso: `PCMGType`, `PCMG`, `PCMGGetLevels()`, `PCMGSetLevels()`, `PCMGGetType()`, `PCMGCycleType` 13134b9ad928SBarry Smith @*/ 13149371c9d4SSatish Balay PetscErrorCode PCMGSetType(PC pc, PCMGType form) { 1315f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 13164b9ad928SBarry Smith 13174b9ad928SBarry Smith PetscFunctionBegin; 13180700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1319c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc, form, 2); 132031567311SBarry Smith mg->am = form; 13219dcbbd2bSBarry Smith if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG; 1322c60c7ad4SBarry Smith else pc->ops->applyrichardson = NULL; 1323c60c7ad4SBarry Smith PetscFunctionReturn(0); 1324c60c7ad4SBarry Smith } 1325c60c7ad4SBarry Smith 1326c60c7ad4SBarry Smith /*@ 1327*f1580f4eSBarry Smith PCMGGetType - Finds the form of multigrid the `PCMG` is using multiplicative, additive, full, or the Kaskade algorithm. 1328c60c7ad4SBarry Smith 1329*f1580f4eSBarry Smith Logically Collective on pc 1330c60c7ad4SBarry Smith 1331c60c7ad4SBarry Smith Input Parameter: 1332c60c7ad4SBarry Smith . pc - the preconditioner context 1333c60c7ad4SBarry Smith 1334c60c7ad4SBarry Smith Output Parameter: 1335*f1580f4eSBarry Smith . type - one of `PC_MG_MULTIPLICATIVE`, `PC_MG_ADDITIVE`, `PC_MG_FULL`, `PC_MG_KASKADE`, `PCMGCycleType` 1336c60c7ad4SBarry Smith 1337c60c7ad4SBarry Smith Level: advanced 1338c60c7ad4SBarry Smith 1339*f1580f4eSBarry Smith .seealso: `PCMGType`, `PCMG`, `PCMGGetLevels()`, `PCMGSetLevels()`, `PCMGSetType()` 1340c60c7ad4SBarry Smith @*/ 13419371c9d4SSatish Balay PetscErrorCode PCMGGetType(PC pc, PCMGType *type) { 1342c60c7ad4SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 1343c60c7ad4SBarry Smith 1344c60c7ad4SBarry Smith PetscFunctionBegin; 1345c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1346c60c7ad4SBarry Smith *type = mg->am; 13474b9ad928SBarry Smith PetscFunctionReturn(0); 13484b9ad928SBarry Smith } 13494b9ad928SBarry Smith 13504b9ad928SBarry Smith /*@ 1351*f1580f4eSBarry Smith PCMGSetCycleType - Sets the type cycles to use. Use `PCMGSetCycleTypeOnLevel()` for more 13524b9ad928SBarry Smith complicated cycling. 13534b9ad928SBarry Smith 1354*f1580f4eSBarry Smith Logically Collective on pc 13554b9ad928SBarry Smith 13564b9ad928SBarry Smith Input Parameters: 1357c2be2410SBarry Smith + pc - the multigrid context 1358*f1580f4eSBarry Smith - n - either `PC_MG_CYCLE_V` or `PC_MG_CYCLE_W` 13594b9ad928SBarry Smith 13604b9ad928SBarry Smith Options Database Key: 1361c1cbb1deSBarry Smith . -pc_mg_cycle_type <v,w> - provide the cycle desired 13624b9ad928SBarry Smith 13634b9ad928SBarry Smith Level: advanced 13644b9ad928SBarry Smith 1365*f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetCycleTypeOnLevel()`, `PCMGType`, `PCMGCycleType` 13664b9ad928SBarry Smith @*/ 13679371c9d4SSatish Balay PetscErrorCode PCMGSetCycleType(PC pc, PCMGCycleType n) { 1368f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 1369f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 137079416396SBarry Smith PetscInt i, levels; 13714b9ad928SBarry Smith 13724b9ad928SBarry Smith PetscFunctionBegin; 13730700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 137469fbec6eSBarry Smith PetscValidLogicalCollectiveEnum(pc, n, 2); 137528b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ORDER, "Must set MG levels with PCMGSetLevels() before calling"); 1376f3fbd535SBarry Smith levels = mglevels[0]->levels; 13772fa5cd67SKarl Rupp for (i = 0; i < levels; i++) mglevels[i]->cycles = n; 13784b9ad928SBarry Smith PetscFunctionReturn(0); 13794b9ad928SBarry Smith } 13804b9ad928SBarry Smith 13818cc2d5dfSBarry Smith /*@ 13828cc2d5dfSBarry Smith PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step 1383*f1580f4eSBarry Smith of multigrid when `PCMGType` is `PC_MG_MULTIPLICATIVE` 13848cc2d5dfSBarry Smith 1385*f1580f4eSBarry Smith Logically Collective on pc 13868cc2d5dfSBarry Smith 13878cc2d5dfSBarry Smith Input Parameters: 13888cc2d5dfSBarry Smith + pc - the multigrid context 13898cc2d5dfSBarry Smith - n - number of cycles (default is 1) 13908cc2d5dfSBarry Smith 13918cc2d5dfSBarry Smith Options Database Key: 1392147403d9SBarry Smith . -pc_mg_multiplicative_cycles n - set the number of cycles 13938cc2d5dfSBarry Smith 13948cc2d5dfSBarry Smith Level: advanced 13958cc2d5dfSBarry Smith 1396*f1580f4eSBarry Smith Note: 1397*f1580f4eSBarry Smith This is not associated with setting a v or w cycle, that is set with `PCMGSetCycleType()` 13988cc2d5dfSBarry Smith 1399*f1580f4eSBarry Smith .seealso: `PCMGSetCycleTypeOnLevel()`, `PCMGSetCycleType()`, `PCMGCycleType`, `PCMGType` 14008cc2d5dfSBarry Smith @*/ 14019371c9d4SSatish Balay PetscErrorCode PCMGMultiplicativeSetCycles(PC pc, PetscInt n) { 1402f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 14038cc2d5dfSBarry Smith 14048cc2d5dfSBarry Smith PetscFunctionBegin; 14050700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1406c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc, n, 2); 14072a21e185SBarry Smith mg->cyclesperpcapply = n; 14088cc2d5dfSBarry Smith PetscFunctionReturn(0); 14098cc2d5dfSBarry Smith } 14108cc2d5dfSBarry Smith 14119371c9d4SSatish Balay PetscErrorCode PCMGSetGalerkin_MG(PC pc, PCMGGalerkinType use) { 1412fb15c04eSBarry Smith PC_MG *mg = (PC_MG *)pc->data; 1413fb15c04eSBarry Smith 1414fb15c04eSBarry Smith PetscFunctionBegin; 14152134b1e4SBarry Smith mg->galerkin = use; 1416fb15c04eSBarry Smith PetscFunctionReturn(0); 1417fb15c04eSBarry Smith } 1418fb15c04eSBarry Smith 1419c2be2410SBarry Smith /*@ 142097177400SBarry Smith PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the 142182b87d87SMatthew G. Knepley finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i 1422c2be2410SBarry Smith 1423*f1580f4eSBarry Smith Logically Collective on pc 1424c2be2410SBarry Smith 1425c2be2410SBarry Smith Input Parameters: 1426c91913d3SJed Brown + pc - the multigrid context 1427*f1580f4eSBarry Smith - use - one of `PC_MG_GALERKIN_BOTH`, `PC_MG_GALERKIN_PMAT`, `PC_MG_GALERKIN_MAT`, or `PC_MG_GALERKIN_NONE` 1428c2be2410SBarry Smith 1429c2be2410SBarry Smith Options Database Key: 1430147403d9SBarry Smith . -pc_mg_galerkin <both,pmat,mat,none> - set the matrices to form via the Galerkin process 1431c2be2410SBarry Smith 1432c2be2410SBarry Smith Level: intermediate 1433c2be2410SBarry Smith 1434*f1580f4eSBarry Smith Note: 1435*f1580f4eSBarry Smith Some codes that use `PCMG` such as `PCGAMG` use Galerkin internally while constructing the hierarchy and thus do not 1436*f1580f4eSBarry Smith use the `PCMG` construction of the coarser grids. 14371c1aac46SBarry Smith 1438*f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGGetGalerkin()`, `PCMGGalerkinType` 1439c2be2410SBarry Smith @*/ 14409371c9d4SSatish Balay PetscErrorCode PCMGSetGalerkin(PC pc, PCMGGalerkinType use) { 1441c2be2410SBarry Smith PetscFunctionBegin; 14420700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1443cac4c232SBarry Smith PetscTryMethod(pc, "PCMGSetGalerkin_C", (PC, PCMGGalerkinType), (pc, use)); 1444c2be2410SBarry Smith PetscFunctionReturn(0); 1445c2be2410SBarry Smith } 1446c2be2410SBarry Smith 14473fc8bf9cSBarry Smith /*@ 1448*f1580f4eSBarry Smith PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. A_i-1 = r_i * A_i * p_i 14493fc8bf9cSBarry Smith 14503fc8bf9cSBarry Smith Not Collective 14513fc8bf9cSBarry Smith 14523fc8bf9cSBarry Smith Input Parameter: 14533fc8bf9cSBarry Smith . pc - the multigrid context 14543fc8bf9cSBarry Smith 14553fc8bf9cSBarry Smith Output Parameter: 1456*f1580f4eSBarry 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` 14573fc8bf9cSBarry Smith 14583fc8bf9cSBarry Smith Level: intermediate 14593fc8bf9cSBarry Smith 1460*f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetGalerkin()`, `PCMGGalerkinType` 14613fc8bf9cSBarry Smith @*/ 14629371c9d4SSatish Balay PetscErrorCode PCMGGetGalerkin(PC pc, PCMGGalerkinType *galerkin) { 1463f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 14643fc8bf9cSBarry Smith 14653fc8bf9cSBarry Smith PetscFunctionBegin; 14660700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 14672134b1e4SBarry Smith *galerkin = mg->galerkin; 14683fc8bf9cSBarry Smith PetscFunctionReturn(0); 14693fc8bf9cSBarry Smith } 14703fc8bf9cSBarry Smith 14719371c9d4SSatish Balay PetscErrorCode PCMGSetAdaptInterpolation_MG(PC pc, PetscBool adapt) { 1472f3b08a26SMatthew G. Knepley PC_MG *mg = (PC_MG *)pc->data; 1473f3b08a26SMatthew G. Knepley 1474f3b08a26SMatthew G. Knepley PetscFunctionBegin; 1475f3b08a26SMatthew G. Knepley mg->adaptInterpolation = adapt; 1476f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 1477f3b08a26SMatthew G. Knepley } 1478f3b08a26SMatthew G. Knepley 14799371c9d4SSatish Balay PetscErrorCode PCMGGetAdaptInterpolation_MG(PC pc, PetscBool *adapt) { 1480f3b08a26SMatthew G. Knepley PC_MG *mg = (PC_MG *)pc->data; 1481f3b08a26SMatthew G. Knepley 1482f3b08a26SMatthew G. Knepley PetscFunctionBegin; 1483f3b08a26SMatthew G. Knepley *adapt = mg->adaptInterpolation; 1484f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 1485f3b08a26SMatthew G. Knepley } 1486f3b08a26SMatthew G. Knepley 14879371c9d4SSatish Balay PetscErrorCode PCMGSetAdaptCoarseSpaceType_MG(PC pc, PCMGCoarseSpaceType ctype) { 14882b3cbbdaSStefano Zampini PC_MG *mg = (PC_MG *)pc->data; 14892b3cbbdaSStefano Zampini 14902b3cbbdaSStefano Zampini PetscFunctionBegin; 14912b3cbbdaSStefano Zampini mg->adaptInterpolation = ctype != PCMG_ADAPT_NONE ? PETSC_TRUE : PETSC_FALSE; 14922b3cbbdaSStefano Zampini mg->coarseSpaceType = ctype; 14932b3cbbdaSStefano Zampini PetscFunctionReturn(0); 14942b3cbbdaSStefano Zampini } 14952b3cbbdaSStefano Zampini 14969371c9d4SSatish Balay PetscErrorCode PCMGGetAdaptCoarseSpaceType_MG(PC pc, PCMGCoarseSpaceType *ctype) { 14972b3cbbdaSStefano Zampini PC_MG *mg = (PC_MG *)pc->data; 14982b3cbbdaSStefano Zampini 14992b3cbbdaSStefano Zampini PetscFunctionBegin; 15002b3cbbdaSStefano Zampini *ctype = mg->coarseSpaceType; 15012b3cbbdaSStefano Zampini PetscFunctionReturn(0); 15022b3cbbdaSStefano Zampini } 15032b3cbbdaSStefano Zampini 15049371c9d4SSatish Balay PetscErrorCode PCMGSetAdaptCR_MG(PC pc, PetscBool cr) { 150541b6fd38SMatthew G. Knepley PC_MG *mg = (PC_MG *)pc->data; 150641b6fd38SMatthew G. Knepley 150741b6fd38SMatthew G. Knepley PetscFunctionBegin; 150841b6fd38SMatthew G. Knepley mg->compatibleRelaxation = cr; 150941b6fd38SMatthew G. Knepley PetscFunctionReturn(0); 151041b6fd38SMatthew G. Knepley } 151141b6fd38SMatthew G. Knepley 15129371c9d4SSatish Balay PetscErrorCode PCMGGetAdaptCR_MG(PC pc, PetscBool *cr) { 151341b6fd38SMatthew G. Knepley PC_MG *mg = (PC_MG *)pc->data; 151441b6fd38SMatthew G. Knepley 151541b6fd38SMatthew G. Knepley PetscFunctionBegin; 151641b6fd38SMatthew G. Knepley *cr = mg->compatibleRelaxation; 151741b6fd38SMatthew G. Knepley PetscFunctionReturn(0); 151841b6fd38SMatthew G. Knepley } 151941b6fd38SMatthew G. Knepley 15202b3cbbdaSStefano Zampini /*@C 15212b3cbbdaSStefano Zampini PCMGSetAdaptCoarseSpaceType - Set the type of adaptive coarse space. 15222b3cbbdaSStefano Zampini 15232b3cbbdaSStefano 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. 15242b3cbbdaSStefano Zampini 1525*f1580f4eSBarry Smith Logically Collective on pc 15262b3cbbdaSStefano Zampini 15272b3cbbdaSStefano Zampini Input Parameters: 15282b3cbbdaSStefano Zampini + pc - the multigrid context 15292b3cbbdaSStefano Zampini - ctype - the type of coarse space 15302b3cbbdaSStefano Zampini 15312b3cbbdaSStefano Zampini Options Database Keys: 15322b3cbbdaSStefano Zampini + -pc_mg_adapt_interp_n <int> - The number of modes to use 15332b3cbbdaSStefano Zampini - -pc_mg_adapt_interp_coarse_space <type> - The type of coarse space: none, polynomial, harmonic, eigenvector, generalized_eigenvector, gdsw 15342b3cbbdaSStefano Zampini 15352b3cbbdaSStefano Zampini Level: intermediate 15362b3cbbdaSStefano Zampini 1537*f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGCoarseSpaceType`, `PCMGGetAdaptCoarseSpaceType()`, `PCMGSetGalerkin()`, `PCMGSetAdaptInterpolation()` 15382b3cbbdaSStefano Zampini @*/ 15399371c9d4SSatish Balay PetscErrorCode PCMGSetAdaptCoarseSpaceType(PC pc, PCMGCoarseSpaceType ctype) { 15402b3cbbdaSStefano Zampini PetscFunctionBegin; 15412b3cbbdaSStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 15422b3cbbdaSStefano Zampini PetscValidLogicalCollectiveEnum(pc, ctype, 2); 15432b3cbbdaSStefano Zampini PetscTryMethod(pc, "PCMGSetAdaptCoarseSpaceType_C", (PC, PCMGCoarseSpaceType), (pc, ctype)); 15442b3cbbdaSStefano Zampini PetscFunctionReturn(0); 15452b3cbbdaSStefano Zampini } 15462b3cbbdaSStefano Zampini 15472b3cbbdaSStefano Zampini /*@C 15482b3cbbdaSStefano Zampini PCMGGetAdaptCoarseSpaceType - Get the type of adaptive coarse space. 15492b3cbbdaSStefano Zampini 15502b3cbbdaSStefano Zampini Not Collective 15512b3cbbdaSStefano Zampini 15522b3cbbdaSStefano Zampini Input Parameter: 15532b3cbbdaSStefano Zampini . pc - the multigrid context 15542b3cbbdaSStefano Zampini 15552b3cbbdaSStefano Zampini Output Parameter: 15562b3cbbdaSStefano Zampini . ctype - the type of coarse space 15572b3cbbdaSStefano Zampini 15582b3cbbdaSStefano Zampini Level: intermediate 15592b3cbbdaSStefano Zampini 1560*f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGCoarseSpaceType`, `PCMGSetAdaptCoarseSpaceType()`, `PCMGSetGalerkin()`, `PCMGSetAdaptInterpolation()` 15612b3cbbdaSStefano Zampini @*/ 15629371c9d4SSatish Balay PetscErrorCode PCMGGetAdaptCoarseSpaceType(PC pc, PCMGCoarseSpaceType *ctype) { 15632b3cbbdaSStefano Zampini PetscFunctionBegin; 15642b3cbbdaSStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 15652b3cbbdaSStefano Zampini PetscValidPointer(ctype, 2); 15662b3cbbdaSStefano Zampini PetscUseMethod(pc, "PCMGGetAdaptCoarseSpaceType_C", (PC, PCMGCoarseSpaceType *), (pc, ctype)); 15672b3cbbdaSStefano Zampini PetscFunctionReturn(0); 15682b3cbbdaSStefano Zampini } 15692b3cbbdaSStefano Zampini 15702b3cbbdaSStefano Zampini /* MATT: REMOVE? */ 1571f3b08a26SMatthew G. Knepley /*@ 1572f3b08a26SMatthew 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. 1573f3b08a26SMatthew G. Knepley 1574*f1580f4eSBarry Smith Logically Collective on pc 1575f3b08a26SMatthew G. Knepley 1576f3b08a26SMatthew G. Knepley Input Parameters: 1577f3b08a26SMatthew G. Knepley + pc - the multigrid context 1578f3b08a26SMatthew G. Knepley - adapt - flag for adaptation of the interpolator 1579f3b08a26SMatthew G. Knepley 1580f3b08a26SMatthew G. Knepley Options Database Keys: 1581f3b08a26SMatthew G. Knepley + -pc_mg_adapt_interp - Turn on adaptation 1582f3b08a26SMatthew G. Knepley . -pc_mg_adapt_interp_n <int> - The number of modes to use, should be divisible by dimension 1583f3b08a26SMatthew G. Knepley - -pc_mg_adapt_interp_coarse_space <type> - The type of coarse space: polynomial, harmonic, eigenvector, generalized_eigenvector 1584f3b08a26SMatthew G. Knepley 1585f3b08a26SMatthew G. Knepley Level: intermediate 1586f3b08a26SMatthew G. Knepley 1587*f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGGetAdaptInterpolation()`, `PCMGSetGalerkin()`, `PCMGGetAdaptCoarseSpaceType()`, `PCMGSetAdaptCoarseSpaceType()` 1588f3b08a26SMatthew G. Knepley @*/ 15899371c9d4SSatish Balay PetscErrorCode PCMGSetAdaptInterpolation(PC pc, PetscBool adapt) { 1590f3b08a26SMatthew G. Knepley PetscFunctionBegin; 1591f3b08a26SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1592cac4c232SBarry Smith PetscTryMethod(pc, "PCMGSetAdaptInterpolation_C", (PC, PetscBool), (pc, adapt)); 1593f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 1594f3b08a26SMatthew G. Knepley } 1595f3b08a26SMatthew G. Knepley 1596f3b08a26SMatthew G. Knepley /*@ 1597*f1580f4eSBarry Smith PCMGGetAdaptInterpolation - Get the flag to adapt the interpolator based upon a vector space which should be accurately captured by the next coarser mesh, 1598*f1580f4eSBarry Smith and thus accurately interpolated. 1599f3b08a26SMatthew G. Knepley 1600f3b08a26SMatthew G. Knepley Not collective 1601f3b08a26SMatthew G. Knepley 1602f3b08a26SMatthew G. Knepley Input Parameter: 1603f3b08a26SMatthew G. Knepley . pc - the multigrid context 1604f3b08a26SMatthew G. Knepley 1605f3b08a26SMatthew G. Knepley Output Parameter: 1606f3b08a26SMatthew G. Knepley . adapt - flag for adaptation of the interpolator 1607f3b08a26SMatthew G. Knepley 1608f3b08a26SMatthew G. Knepley Level: intermediate 1609f3b08a26SMatthew G. Knepley 1610*f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetAdaptInterpolation()`, `PCMGSetGalerkin()`, `PCMGGetAdaptCoarseSpaceType()`, `PCMGSetAdaptCoarseSpaceType()` 1611f3b08a26SMatthew G. Knepley @*/ 16129371c9d4SSatish Balay PetscErrorCode PCMGGetAdaptInterpolation(PC pc, PetscBool *adapt) { 1613f3b08a26SMatthew G. Knepley PetscFunctionBegin; 1614f3b08a26SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1615dadcf809SJacob Faibussowitsch PetscValidBoolPointer(adapt, 2); 1616cac4c232SBarry Smith PetscUseMethod(pc, "PCMGGetAdaptInterpolation_C", (PC, PetscBool *), (pc, adapt)); 1617f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 1618f3b08a26SMatthew G. Knepley } 1619f3b08a26SMatthew G. Knepley 16204b9ad928SBarry Smith /*@ 162141b6fd38SMatthew G. Knepley PCMGSetAdaptCR - Monitor the coarse space quality using an auxiliary solve with compatible relaxation. 162241b6fd38SMatthew G. Knepley 1623*f1580f4eSBarry Smith Logically Collective on pc 162441b6fd38SMatthew G. Knepley 162541b6fd38SMatthew G. Knepley Input Parameters: 162641b6fd38SMatthew G. Knepley + pc - the multigrid context 162741b6fd38SMatthew G. Knepley - cr - flag for compatible relaxation 162841b6fd38SMatthew G. Knepley 1629*f1580f4eSBarry Smith Options Database Key: 163041b6fd38SMatthew G. Knepley . -pc_mg_adapt_cr - Turn on compatible relaxation 163141b6fd38SMatthew G. Knepley 163241b6fd38SMatthew G. Knepley Level: intermediate 163341b6fd38SMatthew G. Knepley 1634*f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGGetAdaptCR()`, `PCMGSetAdaptInterpolation()`, `PCMGSetGalerkin()`, `PCMGGetAdaptCoarseSpaceType()`, `PCMGSetAdaptCoarseSpaceType()` 163541b6fd38SMatthew G. Knepley @*/ 16369371c9d4SSatish Balay PetscErrorCode PCMGSetAdaptCR(PC pc, PetscBool cr) { 163741b6fd38SMatthew G. Knepley PetscFunctionBegin; 163841b6fd38SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1639cac4c232SBarry Smith PetscTryMethod(pc, "PCMGSetAdaptCR_C", (PC, PetscBool), (pc, cr)); 164041b6fd38SMatthew G. Knepley PetscFunctionReturn(0); 164141b6fd38SMatthew G. Knepley } 164241b6fd38SMatthew G. Knepley 164341b6fd38SMatthew G. Knepley /*@ 164441b6fd38SMatthew G. Knepley PCMGGetAdaptCR - Get the flag to monitor coarse space quality using an auxiliary solve with compatible relaxation. 164541b6fd38SMatthew G. Knepley 164641b6fd38SMatthew G. Knepley Not collective 164741b6fd38SMatthew G. Knepley 164841b6fd38SMatthew G. Knepley Input Parameter: 164941b6fd38SMatthew G. Knepley . pc - the multigrid context 165041b6fd38SMatthew G. Knepley 165141b6fd38SMatthew G. Knepley Output Parameter: 165241b6fd38SMatthew G. Knepley . cr - flag for compatible relaxaion 165341b6fd38SMatthew G. Knepley 165441b6fd38SMatthew G. Knepley Level: intermediate 165541b6fd38SMatthew G. Knepley 1656*f1580f4eSBarry Smith .seealso: `PCMGSetAdaptCR()`, `PCMGGetAdaptInterpolation()`, `PCMGSetGalerkin()`, `PCMGGetAdaptCoarseSpaceType()`, `PCMGSetAdaptCoarseSpaceType()` 165741b6fd38SMatthew G. Knepley @*/ 16589371c9d4SSatish Balay PetscErrorCode PCMGGetAdaptCR(PC pc, PetscBool *cr) { 165941b6fd38SMatthew G. Knepley PetscFunctionBegin; 166041b6fd38SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1661dadcf809SJacob Faibussowitsch PetscValidBoolPointer(cr, 2); 1662cac4c232SBarry Smith PetscUseMethod(pc, "PCMGGetAdaptCR_C", (PC, PetscBool *), (pc, cr)); 166341b6fd38SMatthew G. Knepley PetscFunctionReturn(0); 166441b6fd38SMatthew G. Knepley } 166541b6fd38SMatthew G. Knepley 166641b6fd38SMatthew G. Knepley /*@ 166706792cafSBarry Smith PCMGSetNumberSmooth - Sets the number of pre and post-smoothing steps to use 1668*f1580f4eSBarry Smith on all levels. Use `PCMGDistinctSmoothUp()` to create separate up and down smoothers if you want different numbers of 1669710315b6SLawrence Mitchell pre- and post-smoothing steps. 167006792cafSBarry Smith 1671*f1580f4eSBarry Smith Logically Collective on pc 167206792cafSBarry Smith 167306792cafSBarry Smith Input Parameters: 167406792cafSBarry Smith + mg - the multigrid context 167506792cafSBarry Smith - n - the number of smoothing steps 167606792cafSBarry Smith 167706792cafSBarry Smith Options Database Key: 1678a2b725a8SWilliam Gropp . -mg_levels_ksp_max_it <n> - Sets number of pre and post-smoothing steps 167906792cafSBarry Smith 168006792cafSBarry Smith Level: advanced 168106792cafSBarry Smith 1682*f1580f4eSBarry Smith Note: 1683*f1580f4eSBarry Smith This does not set a value on the coarsest grid, since we assume that there is no separate smooth up on the coarsest grid. 168406792cafSBarry Smith 1685*f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetDistinctSmoothUp()` 168606792cafSBarry Smith @*/ 16879371c9d4SSatish Balay PetscErrorCode PCMGSetNumberSmooth(PC pc, PetscInt n) { 168806792cafSBarry Smith PC_MG *mg = (PC_MG *)pc->data; 168906792cafSBarry Smith PC_MG_Levels **mglevels = mg->levels; 169006792cafSBarry Smith PetscInt i, levels; 169106792cafSBarry Smith 169206792cafSBarry Smith PetscFunctionBegin; 169306792cafSBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 169406792cafSBarry Smith PetscValidLogicalCollectiveInt(pc, n, 2); 169528b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ORDER, "Must set MG levels with PCMGSetLevels() before calling"); 169606792cafSBarry Smith levels = mglevels[0]->levels; 169706792cafSBarry Smith 169806792cafSBarry Smith for (i = 1; i < levels; i++) { 16999566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(mglevels[i]->smoothu, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, n)); 17009566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(mglevels[i]->smoothd, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, n)); 170106792cafSBarry Smith mg->default_smoothu = n; 170206792cafSBarry Smith mg->default_smoothd = n; 170306792cafSBarry Smith } 170406792cafSBarry Smith PetscFunctionReturn(0); 170506792cafSBarry Smith } 170606792cafSBarry Smith 1707f442ab6aSBarry Smith /*@ 1708*f1580f4eSBarry Smith PCMGSetDistinctSmoothUp - sets the up (post) smoother to be a separate `KSP` from the down (pre) smoother on all levels 1709710315b6SLawrence Mitchell and adds the suffix _up to the options name 1710f442ab6aSBarry Smith 1711*f1580f4eSBarry Smith Logically Collective on pc 1712f442ab6aSBarry Smith 1713*f1580f4eSBarry Smith Input Parameter: 1714f442ab6aSBarry Smith . pc - the preconditioner context 1715f442ab6aSBarry Smith 1716f442ab6aSBarry Smith Options Database Key: 1717147403d9SBarry Smith . -pc_mg_distinct_smoothup <bool> - use distinct smoothing objects 1718f442ab6aSBarry Smith 1719f442ab6aSBarry Smith Level: advanced 1720f442ab6aSBarry Smith 1721*f1580f4eSBarry Smith Note: 1722*f1580f4eSBarry Smith This does not set a value on the coarsest grid, since we assume that there is no separate smooth up on the coarsest grid. 1723f442ab6aSBarry Smith 1724*f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetNumberSmooth()` 1725f442ab6aSBarry Smith @*/ 17269371c9d4SSatish Balay PetscErrorCode PCMGSetDistinctSmoothUp(PC pc) { 1727f442ab6aSBarry Smith PC_MG *mg = (PC_MG *)pc->data; 1728f442ab6aSBarry Smith PC_MG_Levels **mglevels = mg->levels; 1729f442ab6aSBarry Smith PetscInt i, levels; 1730f442ab6aSBarry Smith KSP subksp; 1731f442ab6aSBarry Smith 1732f442ab6aSBarry Smith PetscFunctionBegin; 1733f442ab6aSBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 173428b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ORDER, "Must set MG levels with PCMGSetLevels() before calling"); 1735f442ab6aSBarry Smith levels = mglevels[0]->levels; 1736f442ab6aSBarry Smith 1737f442ab6aSBarry Smith for (i = 1; i < levels; i++) { 1738710315b6SLawrence Mitchell const char *prefix = NULL; 1739f442ab6aSBarry Smith /* make sure smoother up and down are different */ 17409566063dSJacob Faibussowitsch PetscCall(PCMGGetSmootherUp(pc, i, &subksp)); 17419566063dSJacob Faibussowitsch PetscCall(KSPGetOptionsPrefix(mglevels[i]->smoothd, &prefix)); 17429566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(subksp, prefix)); 17439566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(subksp, "up_")); 1744f442ab6aSBarry Smith } 1745f442ab6aSBarry Smith PetscFunctionReturn(0); 1746f442ab6aSBarry Smith } 1747f442ab6aSBarry Smith 174807a4832bSFande Kong /* No new matrices are created, and the coarse operator matrices are the references to the original ones */ 17499371c9d4SSatish Balay PetscErrorCode PCGetInterpolations_MG(PC pc, PetscInt *num_levels, Mat *interpolations[]) { 175007a4832bSFande Kong PC_MG *mg = (PC_MG *)pc->data; 175107a4832bSFande Kong PC_MG_Levels **mglevels = mg->levels; 175207a4832bSFande Kong Mat *mat; 175307a4832bSFande Kong PetscInt l; 175407a4832bSFande Kong 175507a4832bSFande Kong PetscFunctionBegin; 175628b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling"); 17579566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(mg->nlevels, &mat)); 175807a4832bSFande Kong for (l = 1; l < mg->nlevels; l++) { 175907a4832bSFande Kong mat[l - 1] = mglevels[l]->interpolate; 17609566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat[l - 1])); 176107a4832bSFande Kong } 176207a4832bSFande Kong *num_levels = mg->nlevels; 176307a4832bSFande Kong *interpolations = mat; 176407a4832bSFande Kong PetscFunctionReturn(0); 176507a4832bSFande Kong } 176607a4832bSFande Kong 176707a4832bSFande Kong /* No new matrices are created, and the coarse operator matrices are the references to the original ones */ 17689371c9d4SSatish Balay PetscErrorCode PCGetCoarseOperators_MG(PC pc, PetscInt *num_levels, Mat *coarseOperators[]) { 176907a4832bSFande Kong PC_MG *mg = (PC_MG *)pc->data; 177007a4832bSFande Kong PC_MG_Levels **mglevels = mg->levels; 177107a4832bSFande Kong PetscInt l; 177207a4832bSFande Kong Mat *mat; 177307a4832bSFande Kong 177407a4832bSFande Kong PetscFunctionBegin; 177528b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling"); 17769566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(mg->nlevels, &mat)); 177707a4832bSFande Kong for (l = 0; l < mg->nlevels - 1; l++) { 17789566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[l]->smoothd, NULL, &(mat[l]))); 17799566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat[l])); 178007a4832bSFande Kong } 178107a4832bSFande Kong *num_levels = mg->nlevels; 178207a4832bSFande Kong *coarseOperators = mat; 178307a4832bSFande Kong PetscFunctionReturn(0); 178407a4832bSFande Kong } 178507a4832bSFande Kong 1786f3b08a26SMatthew G. Knepley /*@C 1787*f1580f4eSBarry Smith PCMGRegisterCoarseSpaceConstructor - Adds a method to the `PCMG` package for coarse space construction. 1788f3b08a26SMatthew G. Knepley 1789f3b08a26SMatthew G. Knepley Not collective 1790f3b08a26SMatthew G. Knepley 1791f3b08a26SMatthew G. Knepley Input Parameters: 1792f3b08a26SMatthew G. Knepley + name - name of the constructor 1793f3b08a26SMatthew G. Knepley - function - constructor routine 1794f3b08a26SMatthew G. Knepley 1795f3b08a26SMatthew G. Knepley Calling sequence for the routine: 17962b3cbbdaSStefano Zampini $ my_csp(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, Mat initGuess, Mat *coarseSp) 1797*f1580f4eSBarry Smith + pc - The PC object 1798*f1580f4eSBarry Smith . l - The multigrid level, 0 is the coarse level 1799*f1580f4eSBarry Smith . dm - The DM for this level 1800*f1580f4eSBarry Smith . smooth - The level smoother 1801*f1580f4eSBarry Smith . Nc - The size of the coarse space 1802*f1580f4eSBarry Smith . initGuess - Basis for an initial guess for the space 1803*f1580f4eSBarry Smith - coarseSp - A basis for the computed coarse space 1804f3b08a26SMatthew G. Knepley 1805f3b08a26SMatthew G. Knepley Level: advanced 1806f3b08a26SMatthew G. Knepley 1807*f1580f4eSBarry Smith Developer Note: 1808*f1580f4eSBarry Smith How come this is not used by `PCGAMG`? 1809*f1580f4eSBarry Smith 1810*f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGGetCoarseSpaceConstructor()`, `PCRegister()` 1811f3b08a26SMatthew G. Knepley @*/ 18129371c9d4SSatish Balay PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char name[], PetscErrorCode (*function)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat *)) { 1813f3b08a26SMatthew G. Knepley PetscFunctionBegin; 18149566063dSJacob Faibussowitsch PetscCall(PCInitializePackage()); 18159566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&PCMGCoarseList, name, function)); 1816f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 1817f3b08a26SMatthew G. Knepley } 1818f3b08a26SMatthew G. Knepley 1819f3b08a26SMatthew G. Knepley /*@C 1820f3b08a26SMatthew G. Knepley PCMGGetCoarseSpaceConstructor - Returns the given coarse space construction method. 1821f3b08a26SMatthew G. Knepley 1822f3b08a26SMatthew G. Knepley Not collective 1823f3b08a26SMatthew G. Knepley 1824f3b08a26SMatthew G. Knepley Input Parameter: 1825f3b08a26SMatthew G. Knepley . name - name of the constructor 1826f3b08a26SMatthew G. Knepley 1827f3b08a26SMatthew G. Knepley Output Parameter: 1828f3b08a26SMatthew G. Knepley . function - constructor routine 1829f3b08a26SMatthew G. Knepley 1830f3b08a26SMatthew G. Knepley Level: advanced 1831f3b08a26SMatthew G. Knepley 1832*f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGRegisterCoarseSpaceConstructor()`, `PCRegister()` 1833f3b08a26SMatthew G. Knepley @*/ 18349371c9d4SSatish Balay PetscErrorCode PCMGGetCoarseSpaceConstructor(const char name[], PetscErrorCode (**function)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat *)) { 1835f3b08a26SMatthew G. Knepley PetscFunctionBegin; 18369566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(PCMGCoarseList, name, function)); 1837f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 1838f3b08a26SMatthew G. Knepley } 1839f3b08a26SMatthew G. Knepley 18403b09bd56SBarry Smith /*MC 1841ccb205f8SBarry Smith PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional 18423b09bd56SBarry Smith information about the coarser grid matrices and restriction/interpolation operators. 18433b09bd56SBarry Smith 18443b09bd56SBarry Smith Options Database Keys: 18453b09bd56SBarry Smith + -pc_mg_levels <nlevels> - number of levels including finest 1846391689abSStefano Zampini . -pc_mg_cycle_type <v,w> - provide the cycle desired 18478c1c2452SJed Brown . -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default 18483b09bd56SBarry Smith . -pc_mg_log - log information about time spent on each level of the solver 1849710315b6SLawrence Mitchell . -pc_mg_distinct_smoothup - configure up (after interpolation) and down (before restriction) smoothers separately (with different options prefixes) 18502134b1e4SBarry Smith . -pc_mg_galerkin <both,pmat,mat,none> - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R' 18518cf6037eSBarry Smith . -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1) 18528cf6037eSBarry Smith . -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 1853e3c5b3baSBarry Smith to the Socket viewer for reading from MATLAB. 18548cf6037eSBarry Smith - -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices 18558cf6037eSBarry Smith to the binary output file called binaryoutput 18563b09bd56SBarry Smith 185795452b02SPatrick Sanan Notes: 1858*f1580f4eSBarry Smith If one uses a Krylov method such `KSPGMRES` or `KSPCG` as the smoother then one must use `KSPFGMRES`, `KSPGCR`, or `KSPRICHARDSON` as the outer Krylov method 18593b09bd56SBarry Smith 18608cf6037eSBarry Smith When run with a single level the smoother options are used on that level NOT the coarse grid solver options 18618cf6037eSBarry Smith 1862*f1580f4eSBarry Smith When run with `KSPRICHARDSON` the convergence test changes slightly if monitor is turned on. The iteration count may change slightly. This 186323067569SBarry Smith is because without monitoring the residual norm is computed WITHIN each multigrid cycle on the finest level after the pre-smoothing 186423067569SBarry Smith (because the residual has just been computed for the multigrid algorithm and is hence available for free) while with monitoring the 186523067569SBarry Smith residual is computed at the end of each cycle. 186623067569SBarry Smith 18673b09bd56SBarry Smith Level: intermediate 18683b09bd56SBarry Smith 1869db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCMGType`, `PCEXOTIC`, `PCGAMG`, `PCML`, `PCHYPRE` 1870db781477SPatrick Sanan `PCMGSetLevels()`, `PCMGGetLevels()`, `PCMGSetType()`, `PCMGSetCycleType()`, 1871db781477SPatrick Sanan `PCMGSetDistinctSmoothUp()`, `PCMGGetCoarseSolve()`, `PCMGSetResidual()`, `PCMGSetInterpolation()`, 1872db781477SPatrick Sanan `PCMGSetRestriction()`, `PCMGGetSmoother()`, `PCMGGetSmootherUp()`, `PCMGGetSmootherDown()`, 1873*f1580f4eSBarry Smith `PCMGSetCycleTypeOnLevel()`, `PCMGSetRhs()`, `PCMGSetX()`, `PCMGSetR()`, 1874*f1580f4eSBarry Smith `PCMGSetAdaptCR()`, `PCMGGetAdaptInterpolation()`, `PCMGSetGalerkin()`, `PCMGGetAdaptCoarseSpaceType()`, `PCMGSetAdaptCoarseSpaceType()` 18753b09bd56SBarry Smith M*/ 18763b09bd56SBarry Smith 18779371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc) { 1878f3fbd535SBarry Smith PC_MG *mg; 1879f3fbd535SBarry Smith 18804b9ad928SBarry Smith PetscFunctionBegin; 18819566063dSJacob Faibussowitsch PetscCall(PetscNewLog(pc, &mg)); 18823ec1f749SStefano Zampini pc->data = mg; 1883f3fbd535SBarry Smith mg->nlevels = -1; 188410eca3edSLisandro Dalcin mg->am = PC_MG_MULTIPLICATIVE; 18852134b1e4SBarry Smith mg->galerkin = PC_MG_GALERKIN_NONE; 1886f3b08a26SMatthew G. Knepley mg->adaptInterpolation = PETSC_FALSE; 1887f3b08a26SMatthew G. Knepley mg->Nc = -1; 1888f3b08a26SMatthew G. Knepley mg->eigenvalue = -1; 1889f3fbd535SBarry Smith 189037a44384SMark Adams pc->useAmat = PETSC_TRUE; 189137a44384SMark Adams 18924b9ad928SBarry Smith pc->ops->apply = PCApply_MG; 1893fcb023d4SJed Brown pc->ops->applytranspose = PCApplyTranspose_MG; 189430b0564aSStefano Zampini pc->ops->matapply = PCMatApply_MG; 18954b9ad928SBarry Smith pc->ops->setup = PCSetUp_MG; 1896a06653b4SBarry Smith pc->ops->reset = PCReset_MG; 18974b9ad928SBarry Smith pc->ops->destroy = PCDestroy_MG; 18984b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_MG; 18994b9ad928SBarry Smith pc->ops->view = PCView_MG; 1900fb15c04eSBarry Smith 19019566063dSJacob Faibussowitsch PetscCall(PetscObjectComposedDataRegister(&mg->eigenvalue)); 19029566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetGalerkin_C", PCMGSetGalerkin_MG)); 19039566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetLevels_C", PCMGGetLevels_MG)); 19049566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetLevels_C", PCMGSetLevels_MG)); 19059566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGetInterpolations_C", PCGetInterpolations_MG)); 19069566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGetCoarseOperators_C", PCGetCoarseOperators_MG)); 19079566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetAdaptInterpolation_C", PCMGSetAdaptInterpolation_MG)); 19089566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetAdaptInterpolation_C", PCMGGetAdaptInterpolation_MG)); 19099566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetAdaptCR_C", PCMGSetAdaptCR_MG)); 19109566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetAdaptCR_C", PCMGGetAdaptCR_MG)); 19112b3cbbdaSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetAdaptCoarseSpaceType_C", PCMGSetAdaptCoarseSpaceType_MG)); 19122b3cbbdaSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetAdaptCoarseSpaceType_C", PCMGGetAdaptCoarseSpaceType_MG)); 19134b9ad928SBarry Smith PetscFunctionReturn(0); 19144b9ad928SBarry Smith } 1915