14b9ad928SBarry Smith /* 24b9ad928SBarry Smith Defines the multigrid preconditioner interface. 34b9ad928SBarry Smith */ 4af0996ceSBarry Smith #include <petsc/private/pcmgimpl.h> /*I "petscksp.h" I*/ 541b6fd38SMatthew G. Knepley #include <petsc/private/kspimpl.h> 61e25c274SJed Brown #include <petscdm.h> 7391689abSStefano Zampini PETSC_INTERN PetscErrorCode PCPreSolveChangeRHS(PC, PetscBool *); 84b9ad928SBarry Smith 9f3b08a26SMatthew G. Knepley /* 10f3b08a26SMatthew G. Knepley Contains the list of registered coarse space construction routines 11f3b08a26SMatthew G. Knepley */ 12f3b08a26SMatthew G. Knepley PetscFunctionList PCMGCoarseList = NULL; 13f3b08a26SMatthew G. Knepley 14d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGMCycle_Private(PC pc, PC_MG_Levels **mglevelsin, PetscBool transpose, PetscBool matapp, PCRichardsonConvergedReason *reason) 15d71ae5a4SJacob Faibussowitsch { 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 } 603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 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 } 1173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1184b9ad928SBarry Smith } 1194b9ad928SBarry Smith 120d71ae5a4SJacob Faibussowitsch 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) 121d71ae5a4SJacob Faibussowitsch { 122f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 123f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 124391689abSStefano Zampini PC tpc; 125391689abSStefano Zampini PetscBool changeu, changed; 126f3fbd535SBarry Smith PetscInt levels = mglevels[0]->levels, i; 1274b9ad928SBarry Smith 1284b9ad928SBarry Smith PetscFunctionBegin; 129a762d673SBarry Smith /* When the DM is supplying the matrix then it will not exist until here */ 130a762d673SBarry Smith for (i = 0; i < levels; i++) { 131a762d673SBarry Smith if (!mglevels[i]->A) { 1329566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[i]->smoothu, &mglevels[i]->A, NULL)); 1339566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mglevels[i]->A)); 134a762d673SBarry Smith } 135a762d673SBarry Smith } 136391689abSStefano Zampini 1379566063dSJacob Faibussowitsch PetscCall(KSPGetPC(mglevels[levels - 1]->smoothd, &tpc)); 1389566063dSJacob Faibussowitsch PetscCall(PCPreSolveChangeRHS(tpc, &changed)); 1399566063dSJacob Faibussowitsch PetscCall(KSPGetPC(mglevels[levels - 1]->smoothu, &tpc)); 1409566063dSJacob Faibussowitsch PetscCall(PCPreSolveChangeRHS(tpc, &changeu)); 141391689abSStefano Zampini if (!changed && !changeu) { 1429566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[levels - 1]->b)); 143f3fbd535SBarry Smith mglevels[levels - 1]->b = b; 144391689abSStefano Zampini } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */ 145391689abSStefano Zampini if (!mglevels[levels - 1]->b) { 146391689abSStefano Zampini Vec *vec; 147391689abSStefano Zampini 1489566063dSJacob Faibussowitsch PetscCall(KSPCreateVecs(mglevels[levels - 1]->smoothd, 1, &vec, 0, NULL)); 149391689abSStefano Zampini mglevels[levels - 1]->b = *vec; 1509566063dSJacob Faibussowitsch PetscCall(PetscFree(vec)); 151391689abSStefano Zampini } 1529566063dSJacob Faibussowitsch PetscCall(VecCopy(b, mglevels[levels - 1]->b)); 153391689abSStefano Zampini } 154f3fbd535SBarry Smith mglevels[levels - 1]->x = x; 1554b9ad928SBarry Smith 15631567311SBarry Smith mg->rtol = rtol; 15731567311SBarry Smith mg->abstol = abstol; 15831567311SBarry Smith mg->dtol = dtol; 1594b9ad928SBarry Smith if (rtol) { 1604b9ad928SBarry Smith /* compute initial residual norm for relative convergence test */ 1614b9ad928SBarry Smith PetscReal rnorm; 1627319c654SBarry Smith if (zeroguess) { 1639566063dSJacob Faibussowitsch PetscCall(VecNorm(b, NORM_2, &rnorm)); 1647319c654SBarry Smith } else { 1659566063dSJacob Faibussowitsch PetscCall((*mglevels[levels - 1]->residual)(mglevels[levels - 1]->A, b, x, w)); 1669566063dSJacob Faibussowitsch PetscCall(VecNorm(w, NORM_2, &rnorm)); 1677319c654SBarry Smith } 16831567311SBarry Smith mg->ttol = PetscMax(rtol * rnorm, abstol); 1692fa5cd67SKarl Rupp } else if (abstol) mg->ttol = abstol; 1702fa5cd67SKarl Rupp else mg->ttol = 0.0; 1714b9ad928SBarry Smith 1724d0a8057SBarry Smith /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't 173336babb1SJed Brown stop prematurely due to small residual */ 1744d0a8057SBarry Smith for (i = 1; i < levels; i++) { 1759566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(mglevels[i]->smoothu, 0, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 176f3fbd535SBarry Smith if (mglevels[i]->smoothu != mglevels[i]->smoothd) { 17723067569SBarry Smith /* For Richardson the initial guess is nonzero since it is solving in each cycle the original system not just applying as a preconditioner */ 1789566063dSJacob Faibussowitsch PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->smoothd, PETSC_TRUE)); 1799566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(mglevels[i]->smoothd, 0, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 1804b9ad928SBarry Smith } 1814d0a8057SBarry Smith } 1824d0a8057SBarry Smith 1834d0a8057SBarry Smith *reason = (PCRichardsonConvergedReason)0; 1844d0a8057SBarry Smith for (i = 0; i < its; i++) { 1859566063dSJacob Faibussowitsch PetscCall(PCMGMCycle_Private(pc, mglevels + levels - 1, PETSC_FALSE, PETSC_FALSE, reason)); 1864d0a8057SBarry Smith if (*reason) break; 1874d0a8057SBarry Smith } 1884d0a8057SBarry Smith if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS; 1894d0a8057SBarry Smith *outits = i; 190391689abSStefano Zampini if (!changed && !changeu) mglevels[levels - 1]->b = NULL; 1913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1924b9ad928SBarry Smith } 1934b9ad928SBarry Smith 194d71ae5a4SJacob Faibussowitsch PetscErrorCode PCReset_MG(PC pc) 195d71ae5a4SJacob Faibussowitsch { 1963aeaf226SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 1973aeaf226SBarry Smith PC_MG_Levels **mglevels = mg->levels; 1982b3cbbdaSStefano Zampini PetscInt i, n; 1993aeaf226SBarry Smith 2003aeaf226SBarry Smith PetscFunctionBegin; 2013aeaf226SBarry Smith if (mglevels) { 2023aeaf226SBarry Smith n = mglevels[0]->levels; 2033aeaf226SBarry Smith for (i = 0; i < n - 1; i++) { 2049566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[i + 1]->r)); 2059566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[i]->b)); 2069566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[i]->x)); 2079566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[i + 1]->R)); 2089566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[i]->B)); 2099566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[i]->X)); 2109566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[i]->crx)); 2119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[i]->crb)); 2129566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[i + 1]->restrct)); 2139566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[i + 1]->interpolate)); 2149566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[i + 1]->inject)); 2159566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[i + 1]->rscale)); 2163aeaf226SBarry Smith } 2179566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[n - 1]->crx)); 2189566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[n - 1]->crb)); 219391689abSStefano Zampini /* this is not null only if the smoother on the finest level 220391689abSStefano Zampini changes the rhs during PreSolve */ 2219566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[n - 1]->b)); 2229566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[n - 1]->B)); 2233aeaf226SBarry Smith 2243aeaf226SBarry Smith for (i = 0; i < n; i++) { 2252b3cbbdaSStefano Zampini PetscCall(MatDestroy(&mglevels[i]->coarseSpace)); 2269566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[i]->A)); 22748a46eb9SPierre Jolivet if (mglevels[i]->smoothd != mglevels[i]->smoothu) PetscCall(KSPReset(mglevels[i]->smoothd)); 2289566063dSJacob Faibussowitsch PetscCall(KSPReset(mglevels[i]->smoothu)); 2299566063dSJacob Faibussowitsch if (mglevels[i]->cr) PetscCall(KSPReset(mglevels[i]->cr)); 2303aeaf226SBarry Smith } 231f3b08a26SMatthew G. Knepley mg->Nc = 0; 2323aeaf226SBarry Smith } 2333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2343aeaf226SBarry Smith } 2353aeaf226SBarry Smith 23641b6fd38SMatthew G. Knepley /* Implementing CR 23741b6fd38SMatthew G. Knepley 23841b6fd38SMatthew 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 23941b6fd38SMatthew G. Knepley 24041b6fd38SMatthew G. Knepley Inj^T (Inj Inj^T)^{-1} Inj 24141b6fd38SMatthew G. Knepley 24241b6fd38SMatthew G. Knepley and if Inj is a VecScatter, as it is now in PETSc, we have 24341b6fd38SMatthew G. Knepley 24441b6fd38SMatthew G. Knepley Inj^T Inj 24541b6fd38SMatthew G. Knepley 24641b6fd38SMatthew G. Knepley and 24741b6fd38SMatthew G. Knepley 24841b6fd38SMatthew G. Knepley S = I - Inj^T Inj 24941b6fd38SMatthew G. Knepley 25041b6fd38SMatthew G. Knepley since 25141b6fd38SMatthew G. Knepley 25241b6fd38SMatthew G. Knepley Inj S = Inj - (Inj Inj^T) Inj = 0. 25341b6fd38SMatthew G. Knepley 25441b6fd38SMatthew G. Knepley Brannick suggests 25541b6fd38SMatthew G. Knepley 25641b6fd38SMatthew G. Knepley A \to S^T A S \qquad\mathrm{and}\qquad M \to S^T M S 25741b6fd38SMatthew G. Knepley 25841b6fd38SMatthew 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 25941b6fd38SMatthew G. Knepley 26041b6fd38SMatthew G. Knepley M^{-1} A \to S M^{-1} A S 26141b6fd38SMatthew G. Knepley 26241b6fd38SMatthew G. Knepley In fact, since it is somewhat hard in PETSc to do the symmetric application, we will just apply S on the left. 26341b6fd38SMatthew G. Knepley 26441b6fd38SMatthew G. Knepley Check: || Inj P - I ||_F < tol 26541b6fd38SMatthew G. Knepley Check: In general, Inj Inj^T = I 26641b6fd38SMatthew G. Knepley */ 26741b6fd38SMatthew G. Knepley 26841b6fd38SMatthew G. Knepley typedef struct { 26941b6fd38SMatthew G. Knepley PC mg; /* The PCMG object */ 27041b6fd38SMatthew G. Knepley PetscInt l; /* The multigrid level for this solver */ 27141b6fd38SMatthew G. Knepley Mat Inj; /* The injection matrix */ 27241b6fd38SMatthew G. Knepley Mat S; /* I - Inj^T Inj */ 27341b6fd38SMatthew G. Knepley } CRContext; 27441b6fd38SMatthew G. Knepley 275d71ae5a4SJacob Faibussowitsch static PetscErrorCode CRSetup_Private(PC pc) 276d71ae5a4SJacob Faibussowitsch { 27741b6fd38SMatthew G. Knepley CRContext *ctx; 27841b6fd38SMatthew G. Knepley Mat It; 27941b6fd38SMatthew G. Knepley 28041b6fd38SMatthew G. Knepley PetscFunctionBeginUser; 2819566063dSJacob Faibussowitsch PetscCall(PCShellGetContext(pc, &ctx)); 2829566063dSJacob Faibussowitsch PetscCall(PCMGGetInjection(ctx->mg, ctx->l, &It)); 28328b400f6SJacob Faibussowitsch PetscCheck(It, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "CR requires that injection be defined for this PCMG"); 2849566063dSJacob Faibussowitsch PetscCall(MatCreateTranspose(It, &ctx->Inj)); 2859566063dSJacob Faibussowitsch PetscCall(MatCreateNormal(ctx->Inj, &ctx->S)); 2869566063dSJacob Faibussowitsch PetscCall(MatScale(ctx->S, -1.0)); 2879566063dSJacob Faibussowitsch PetscCall(MatShift(ctx->S, 1.0)); 2883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 28941b6fd38SMatthew G. Knepley } 29041b6fd38SMatthew G. Knepley 291d71ae5a4SJacob Faibussowitsch static PetscErrorCode CRApply_Private(PC pc, Vec x, Vec y) 292d71ae5a4SJacob Faibussowitsch { 29341b6fd38SMatthew G. Knepley CRContext *ctx; 29441b6fd38SMatthew G. Knepley 29541b6fd38SMatthew G. Knepley PetscFunctionBeginUser; 2969566063dSJacob Faibussowitsch PetscCall(PCShellGetContext(pc, &ctx)); 2979566063dSJacob Faibussowitsch PetscCall(MatMult(ctx->S, x, y)); 2983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 29941b6fd38SMatthew G. Knepley } 30041b6fd38SMatthew G. Knepley 301d71ae5a4SJacob Faibussowitsch static PetscErrorCode CRDestroy_Private(PC pc) 302d71ae5a4SJacob Faibussowitsch { 30341b6fd38SMatthew G. Knepley CRContext *ctx; 30441b6fd38SMatthew G. Knepley 30541b6fd38SMatthew G. Knepley PetscFunctionBeginUser; 3069566063dSJacob Faibussowitsch PetscCall(PCShellGetContext(pc, &ctx)); 3079566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx->Inj)); 3089566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx->S)); 3099566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 3109566063dSJacob Faibussowitsch PetscCall(PCShellSetContext(pc, NULL)); 3113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 31241b6fd38SMatthew G. Knepley } 31341b6fd38SMatthew G. Knepley 314d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateCR_Private(PC pc, PetscInt l, PC *cr) 315d71ae5a4SJacob Faibussowitsch { 31641b6fd38SMatthew G. Knepley CRContext *ctx; 31741b6fd38SMatthew G. Knepley 31841b6fd38SMatthew G. Knepley PetscFunctionBeginUser; 3199566063dSJacob Faibussowitsch PetscCall(PCCreate(PetscObjectComm((PetscObject)pc), cr)); 3209566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*cr, "S (complementary projector to injection)")); 3219566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(1, &ctx)); 32241b6fd38SMatthew G. Knepley ctx->mg = pc; 32341b6fd38SMatthew G. Knepley ctx->l = l; 3249566063dSJacob Faibussowitsch PetscCall(PCSetType(*cr, PCSHELL)); 3259566063dSJacob Faibussowitsch PetscCall(PCShellSetContext(*cr, ctx)); 3269566063dSJacob Faibussowitsch PetscCall(PCShellSetApply(*cr, CRApply_Private)); 3279566063dSJacob Faibussowitsch PetscCall(PCShellSetSetUp(*cr, CRSetup_Private)); 3289566063dSJacob Faibussowitsch PetscCall(PCShellSetDestroy(*cr, CRDestroy_Private)); 3293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 33041b6fd38SMatthew G. Knepley } 33141b6fd38SMatthew G. Knepley 332d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetLevels_MG(PC pc, PetscInt levels, MPI_Comm *comms) 333d71ae5a4SJacob Faibussowitsch { 334f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 335ce94432eSBarry Smith MPI_Comm comm; 3363aeaf226SBarry Smith PC_MG_Levels **mglevels = mg->levels; 33710eca3edSLisandro Dalcin PCMGType mgtype = mg->am; 33810167fecSBarry Smith PetscInt mgctype = (PetscInt)PC_MG_CYCLE_V; 339f3fbd535SBarry Smith PetscInt i; 340f3fbd535SBarry Smith PetscMPIInt size; 341f3fbd535SBarry Smith const char *prefix; 342f3fbd535SBarry Smith PC ipc; 3433aeaf226SBarry Smith PetscInt n; 3444b9ad928SBarry Smith 3454b9ad928SBarry Smith PetscFunctionBegin; 3460700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 347548767e0SJed Brown PetscValidLogicalCollectiveInt(pc, levels, 2); 3483ba16761SJacob Faibussowitsch if (mg->nlevels == levels) PetscFunctionReturn(PETSC_SUCCESS); 3499566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)pc, &comm)); 3503aeaf226SBarry Smith if (mglevels) { 35110eca3edSLisandro Dalcin mgctype = mglevels[0]->cycles; 3523aeaf226SBarry Smith /* changing the number of levels so free up the previous stuff */ 3539566063dSJacob Faibussowitsch PetscCall(PCReset_MG(pc)); 3543aeaf226SBarry Smith n = mglevels[0]->levels; 3553aeaf226SBarry Smith for (i = 0; i < n; i++) { 35648a46eb9SPierre Jolivet if (mglevels[i]->smoothd != mglevels[i]->smoothu) PetscCall(KSPDestroy(&mglevels[i]->smoothd)); 3579566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&mglevels[i]->smoothu)); 3589566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&mglevels[i]->cr)); 3599566063dSJacob Faibussowitsch PetscCall(PetscFree(mglevels[i])); 3603aeaf226SBarry Smith } 3619566063dSJacob Faibussowitsch PetscCall(PetscFree(mg->levels)); 3623aeaf226SBarry Smith } 363f3fbd535SBarry Smith 364f3fbd535SBarry Smith mg->nlevels = levels; 365f3fbd535SBarry Smith 3669566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(levels, &mglevels)); 367f3fbd535SBarry Smith 3689566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &prefix)); 369f3fbd535SBarry Smith 370a9db87a2SMatthew G Knepley mg->stageApply = 0; 371f3fbd535SBarry Smith for (i = 0; i < levels; i++) { 3724dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&mglevels[i])); 3732fa5cd67SKarl Rupp 374f3fbd535SBarry Smith mglevels[i]->level = i; 375f3fbd535SBarry Smith mglevels[i]->levels = levels; 37610eca3edSLisandro Dalcin mglevels[i]->cycles = mgctype; 377336babb1SJed Brown mg->default_smoothu = 2; 378336babb1SJed Brown mg->default_smoothd = 2; 37963e6d426SJed Brown mglevels[i]->eventsmoothsetup = 0; 38063e6d426SJed Brown mglevels[i]->eventsmoothsolve = 0; 38163e6d426SJed Brown mglevels[i]->eventresidual = 0; 38263e6d426SJed Brown mglevels[i]->eventinterprestrict = 0; 383f3fbd535SBarry Smith 384f3fbd535SBarry Smith if (comms) comm = comms[i]; 385d5a8f7e6SBarry Smith if (comm != MPI_COMM_NULL) { 3869566063dSJacob Faibussowitsch PetscCall(KSPCreate(comm, &mglevels[i]->smoothd)); 3873821be0aSBarry Smith PetscCall(KSPSetNestLevel(mglevels[i]->smoothd, pc->kspnestlevel)); 3889566063dSJacob Faibussowitsch PetscCall(KSPSetErrorIfNotConverged(mglevels[i]->smoothd, pc->erroriffailure)); 3899566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd, (PetscObject)pc, levels - i)); 3909566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(mglevels[i]->smoothd, prefix)); 3919566063dSJacob Faibussowitsch PetscCall(PetscObjectComposedDataSetInt((PetscObject)mglevels[i]->smoothd, PetscMGLevelId, mglevels[i]->level)); 3920ee9a668SBarry Smith if (i || levels == 1) { 3930ee9a668SBarry Smith char tprefix[128]; 3940ee9a668SBarry Smith 3959566063dSJacob Faibussowitsch PetscCall(KSPSetType(mglevels[i]->smoothd, KSPCHEBYSHEV)); 3969566063dSJacob Faibussowitsch PetscCall(KSPSetConvergenceTest(mglevels[i]->smoothd, KSPConvergedSkip, NULL, NULL)); 3979566063dSJacob Faibussowitsch PetscCall(KSPSetNormType(mglevels[i]->smoothd, KSP_NORM_NONE)); 3989566063dSJacob Faibussowitsch PetscCall(KSPGetPC(mglevels[i]->smoothd, &ipc)); 3999566063dSJacob Faibussowitsch PetscCall(PCSetType(ipc, PCSOR)); 4009566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(mglevels[i]->smoothd, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, mg->default_smoothd)); 401f3fbd535SBarry Smith 4029566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(tprefix, 128, "mg_levels_%d_", (int)i)); 4039566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(mglevels[i]->smoothd, tprefix)); 4040ee9a668SBarry Smith } else { 4059566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(mglevels[0]->smoothd, "mg_coarse_")); 406f3fbd535SBarry Smith 4079dbfc187SHong Zhang /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */ 4089566063dSJacob Faibussowitsch PetscCall(KSPSetType(mglevels[0]->smoothd, KSPPREONLY)); 4099566063dSJacob Faibussowitsch PetscCall(KSPGetPC(mglevels[0]->smoothd, &ipc)); 4109566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 411f3fbd535SBarry Smith if (size > 1) { 4129566063dSJacob Faibussowitsch PetscCall(PCSetType(ipc, PCREDUNDANT)); 413f3fbd535SBarry Smith } else { 4149566063dSJacob Faibussowitsch PetscCall(PCSetType(ipc, PCLU)); 415f3fbd535SBarry Smith } 4169566063dSJacob Faibussowitsch PetscCall(PCFactorSetShiftType(ipc, MAT_SHIFT_INBLOCKS)); 417f3fbd535SBarry Smith } 418d5a8f7e6SBarry Smith } 419f3fbd535SBarry Smith mglevels[i]->smoothu = mglevels[i]->smoothd; 42031567311SBarry Smith mg->rtol = 0.0; 42131567311SBarry Smith mg->abstol = 0.0; 42231567311SBarry Smith mg->dtol = 0.0; 42331567311SBarry Smith mg->ttol = 0.0; 42431567311SBarry Smith mg->cyclesperpcapply = 1; 425f3fbd535SBarry Smith } 426f3fbd535SBarry Smith mg->levels = mglevels; 4279566063dSJacob Faibussowitsch PetscCall(PCMGSetType(pc, mgtype)); 4283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4294b9ad928SBarry Smith } 4304b9ad928SBarry Smith 43197d33e41SMatthew G. Knepley /*@C 432f1580f4eSBarry Smith PCMGSetLevels - Sets the number of levels to use with `PCMG`. 433f1580f4eSBarry Smith Must be called before any other `PCMG` routine. 43497d33e41SMatthew G. Knepley 435c3339decSBarry Smith Logically Collective 43697d33e41SMatthew G. Knepley 43797d33e41SMatthew G. Knepley Input Parameters: 43897d33e41SMatthew G. Knepley + pc - the preconditioner context 43997d33e41SMatthew G. Knepley . levels - the number of levels 44097d33e41SMatthew G. Knepley - comms - optional communicators for each level; this is to allow solving the coarser problems 441d5a8f7e6SBarry Smith on smaller sets of processes. For processes that are not included in the computation 44220f4b53cSBarry Smith you must pass `MPI_COMM_NULL`. Use comms = `NULL` to specify that all processes 44305552d4bSJunchao Zhang should participate in each level of problem. 44497d33e41SMatthew G. Knepley 44597d33e41SMatthew G. Knepley Level: intermediate 44697d33e41SMatthew G. Knepley 44797d33e41SMatthew G. Knepley Notes: 44820f4b53cSBarry Smith If the number of levels is one then the multigrid uses the `-mg_levels` prefix 44920f4b53cSBarry Smith for setting the level options rather than the `-mg_coarse` prefix. 45097d33e41SMatthew G. Knepley 451d5a8f7e6SBarry Smith You can free the information in comms after this routine is called. 452d5a8f7e6SBarry Smith 453f1580f4eSBarry Smith The array of MPI communicators must contain `MPI_COMM_NULL` for those ranks that at each level 454d5a8f7e6SBarry Smith are not participating in the coarser solve. For example, with 2 levels and 1 and 2 ranks on 455d5a8f7e6SBarry Smith the two levels, rank 0 in the original communicator will pass in an array of 2 communicators 456d5a8f7e6SBarry Smith of size 2 and 1, while rank 1 in the original communicator will pass in array of 2 communicators 457f1580f4eSBarry Smith the first of size 2 and the second of value `MPI_COMM_NULL` since the rank 1 does not participate 458d5a8f7e6SBarry Smith in the coarse grid solve. 459d5a8f7e6SBarry Smith 460f1580f4eSBarry Smith Since each coarser level may have a new `MPI_Comm` with fewer ranks than the previous, one 461d5a8f7e6SBarry Smith must take special care in providing the restriction and interpolation operation. We recommend 462d5a8f7e6SBarry Smith providing these as two step operations; first perform a standard restriction or interpolation on 463d5a8f7e6SBarry Smith the full number of ranks for that level and then use an MPI call to copy the resulting vector 46405552d4bSJunchao Zhang array entries (after calls to VecGetArray()) to the smaller or larger number of ranks, note in both 465d5a8f7e6SBarry Smith cases the MPI calls must be made on the larger of the two communicators. Traditional MPI send and 46620f4b53cSBarry Smith receives or `MPI_AlltoAllv()` could be used to do the reshuffling of the vector entries. 467d5a8f7e6SBarry Smith 468feefa0e1SJacob Faibussowitsch Fortran Notes: 46920f4b53cSBarry Smith Use comms = `PETSC_NULL_MPI_COMM` as the equivalent of `NULL` in the C interface. Note `PETSC_NULL_MPI_COMM` 470f1580f4eSBarry Smith is not `MPI_COMM_NULL`. It is more like `PETSC_NULL_INTEGER`, `PETSC_NULL_REAL` etc. 471d5a8f7e6SBarry Smith 472*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMGSetType()`, `PCMGGetLevels()` 47397d33e41SMatthew G. Knepley @*/ 474d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetLevels(PC pc, PetscInt levels, MPI_Comm *comms) 475d71ae5a4SJacob Faibussowitsch { 47697d33e41SMatthew G. Knepley PetscFunctionBegin; 47797d33e41SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 4784f572ea9SToby Isaac if (comms) PetscAssertPointer(comms, 3); 479cac4c232SBarry Smith PetscTryMethod(pc, "PCMGSetLevels_C", (PC, PetscInt, MPI_Comm *), (pc, levels, comms)); 4803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 48197d33e41SMatthew G. Knepley } 48297d33e41SMatthew G. Knepley 483d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDestroy_MG(PC pc) 484d71ae5a4SJacob Faibussowitsch { 485a06653b4SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 486a06653b4SBarry Smith PC_MG_Levels **mglevels = mg->levels; 487a06653b4SBarry Smith PetscInt i, n; 488c07bf074SBarry Smith 489c07bf074SBarry Smith PetscFunctionBegin; 4909566063dSJacob Faibussowitsch PetscCall(PCReset_MG(pc)); 491a06653b4SBarry Smith if (mglevels) { 492a06653b4SBarry Smith n = mglevels[0]->levels; 493a06653b4SBarry Smith for (i = 0; i < n; i++) { 49448a46eb9SPierre Jolivet if (mglevels[i]->smoothd != mglevels[i]->smoothu) PetscCall(KSPDestroy(&mglevels[i]->smoothd)); 4959566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&mglevels[i]->smoothu)); 4969566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&mglevels[i]->cr)); 4979566063dSJacob Faibussowitsch PetscCall(PetscFree(mglevels[i])); 498a06653b4SBarry Smith } 4999566063dSJacob Faibussowitsch PetscCall(PetscFree(mg->levels)); 500a06653b4SBarry Smith } 5019566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 5029566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGetInterpolations_C", NULL)); 5039566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGetCoarseOperators_C", NULL)); 5042b3cbbdaSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetGalerkin_C", NULL)); 5052b3cbbdaSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetLevels_C", NULL)); 5062b3cbbdaSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetLevels_C", NULL)); 5072b3cbbdaSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGetInterpolations_C", NULL)); 5082b3cbbdaSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGetCoarseOperators_C", NULL)); 5092b3cbbdaSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetAdaptInterpolation_C", NULL)); 5102b3cbbdaSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetAdaptInterpolation_C", NULL)); 5112b3cbbdaSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetAdaptCR_C", NULL)); 5122b3cbbdaSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetAdaptCR_C", NULL)); 5132b3cbbdaSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetAdaptCoarseSpaceType_C", NULL)); 5142b3cbbdaSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetAdaptCoarseSpaceType_C", NULL)); 5153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 516f3fbd535SBarry Smith } 517f3fbd535SBarry Smith 518f3fbd535SBarry Smith /* 519f3fbd535SBarry Smith PCApply_MG - Runs either an additive, multiplicative, Kaskadic 520f3fbd535SBarry Smith or full cycle of multigrid. 521f3fbd535SBarry Smith 522f3fbd535SBarry Smith Note: 523f3fbd535SBarry Smith A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle(). 524f3fbd535SBarry Smith */ 525d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_MG_Internal(PC pc, Vec b, Vec x, Mat B, Mat X, PetscBool transpose) 526d71ae5a4SJacob Faibussowitsch { 527f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 528f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 529391689abSStefano Zampini PC tpc; 530f3fbd535SBarry Smith PetscInt levels = mglevels[0]->levels, i; 53130b0564aSStefano Zampini PetscBool changeu, changed, matapp; 532f3fbd535SBarry Smith 533f3fbd535SBarry Smith PetscFunctionBegin; 53430b0564aSStefano Zampini matapp = (PetscBool)(B && X); 5359566063dSJacob Faibussowitsch if (mg->stageApply) PetscCall(PetscLogStagePush(mg->stageApply)); 536e1d8e5deSBarry Smith /* When the DM is supplying the matrix then it will not exist until here */ 537298cc208SBarry Smith for (i = 0; i < levels; i++) { 538e1d8e5deSBarry Smith if (!mglevels[i]->A) { 5399566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[i]->smoothu, &mglevels[i]->A, NULL)); 5409566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mglevels[i]->A)); 541e1d8e5deSBarry Smith } 542e1d8e5deSBarry Smith } 543e1d8e5deSBarry Smith 5449566063dSJacob Faibussowitsch PetscCall(KSPGetPC(mglevels[levels - 1]->smoothd, &tpc)); 5459566063dSJacob Faibussowitsch PetscCall(PCPreSolveChangeRHS(tpc, &changed)); 5469566063dSJacob Faibussowitsch PetscCall(KSPGetPC(mglevels[levels - 1]->smoothu, &tpc)); 5479566063dSJacob Faibussowitsch PetscCall(PCPreSolveChangeRHS(tpc, &changeu)); 548391689abSStefano Zampini if (!changeu && !changed) { 54930b0564aSStefano Zampini if (matapp) { 5509566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[levels - 1]->B)); 55130b0564aSStefano Zampini mglevels[levels - 1]->B = B; 55230b0564aSStefano Zampini } else { 5539566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[levels - 1]->b)); 554f3fbd535SBarry Smith mglevels[levels - 1]->b = b; 55530b0564aSStefano Zampini } 556391689abSStefano Zampini } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */ 55730b0564aSStefano Zampini if (matapp) { 55830b0564aSStefano Zampini if (mglevels[levels - 1]->B) { 55930b0564aSStefano Zampini PetscInt N1, N2; 56030b0564aSStefano Zampini PetscBool flg; 56130b0564aSStefano Zampini 5629566063dSJacob Faibussowitsch PetscCall(MatGetSize(mglevels[levels - 1]->B, NULL, &N1)); 5639566063dSJacob Faibussowitsch PetscCall(MatGetSize(B, NULL, &N2)); 5649566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)mglevels[levels - 1]->B, ((PetscObject)B)->type_name, &flg)); 56548a46eb9SPierre Jolivet if (N1 != N2 || !flg) PetscCall(MatDestroy(&mglevels[levels - 1]->B)); 56630b0564aSStefano Zampini } 56730b0564aSStefano Zampini if (!mglevels[levels - 1]->B) { 5689566063dSJacob Faibussowitsch PetscCall(MatDuplicate(B, MAT_COPY_VALUES, &mglevels[levels - 1]->B)); 56930b0564aSStefano Zampini } else { 5709566063dSJacob Faibussowitsch PetscCall(MatCopy(B, mglevels[levels - 1]->B, SAME_NONZERO_PATTERN)); 57130b0564aSStefano Zampini } 57230b0564aSStefano Zampini } else { 573391689abSStefano Zampini if (!mglevels[levels - 1]->b) { 574391689abSStefano Zampini Vec *vec; 575391689abSStefano Zampini 5769566063dSJacob Faibussowitsch PetscCall(KSPCreateVecs(mglevels[levels - 1]->smoothd, 1, &vec, 0, NULL)); 577391689abSStefano Zampini mglevels[levels - 1]->b = *vec; 5789566063dSJacob Faibussowitsch PetscCall(PetscFree(vec)); 579391689abSStefano Zampini } 5809566063dSJacob Faibussowitsch PetscCall(VecCopy(b, mglevels[levels - 1]->b)); 581391689abSStefano Zampini } 58230b0564aSStefano Zampini } 5839371c9d4SSatish Balay if (matapp) { 5849371c9d4SSatish Balay mglevels[levels - 1]->X = X; 5859371c9d4SSatish Balay } else { 5869371c9d4SSatish Balay mglevels[levels - 1]->x = x; 5879371c9d4SSatish Balay } 58830b0564aSStefano Zampini 58930b0564aSStefano Zampini /* If coarser Xs are present, it means we have already block applied the PC at least once 59030b0564aSStefano Zampini Reset operators if sizes/type do no match */ 59130b0564aSStefano Zampini if (matapp && levels > 1 && mglevels[levels - 2]->X) { 59230b0564aSStefano Zampini PetscInt Xc, Bc; 59330b0564aSStefano Zampini PetscBool flg; 59430b0564aSStefano Zampini 5959566063dSJacob Faibussowitsch PetscCall(MatGetSize(mglevels[levels - 2]->X, NULL, &Xc)); 5969566063dSJacob Faibussowitsch PetscCall(MatGetSize(mglevels[levels - 1]->B, NULL, &Bc)); 5979566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)mglevels[levels - 2]->X, ((PetscObject)mglevels[levels - 1]->X)->type_name, &flg)); 59830b0564aSStefano Zampini if (Xc != Bc || !flg) { 5999566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[levels - 1]->R)); 60030b0564aSStefano Zampini for (i = 0; i < levels - 1; i++) { 6019566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[i]->R)); 6029566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[i]->B)); 6039566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[i]->X)); 60430b0564aSStefano Zampini } 60530b0564aSStefano Zampini } 60630b0564aSStefano Zampini } 607391689abSStefano Zampini 60831567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 6099566063dSJacob Faibussowitsch if (matapp) PetscCall(MatZeroEntries(X)); 6109566063dSJacob Faibussowitsch else PetscCall(VecZeroEntries(x)); 61148a46eb9SPierre Jolivet for (i = 0; i < mg->cyclesperpcapply; i++) PetscCall(PCMGMCycle_Private(pc, mglevels + levels - 1, transpose, matapp, NULL)); 6122fa5cd67SKarl Rupp } else if (mg->am == PC_MG_ADDITIVE) { 6139566063dSJacob Faibussowitsch PetscCall(PCMGACycle_Private(pc, mglevels, transpose, matapp)); 6142fa5cd67SKarl Rupp } else if (mg->am == PC_MG_KASKADE) { 6159566063dSJacob Faibussowitsch PetscCall(PCMGKCycle_Private(pc, mglevels, transpose, matapp)); 6162fa5cd67SKarl Rupp } else { 6179566063dSJacob Faibussowitsch PetscCall(PCMGFCycle_Private(pc, mglevels, transpose, matapp)); 618f3fbd535SBarry Smith } 6199566063dSJacob Faibussowitsch if (mg->stageApply) PetscCall(PetscLogStagePop()); 62030b0564aSStefano Zampini if (!changeu && !changed) { 6219371c9d4SSatish Balay if (matapp) { 6229371c9d4SSatish Balay mglevels[levels - 1]->B = NULL; 6239371c9d4SSatish Balay } else { 6249371c9d4SSatish Balay mglevels[levels - 1]->b = NULL; 6259371c9d4SSatish Balay } 62630b0564aSStefano Zampini } 6273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 628f3fbd535SBarry Smith } 629f3fbd535SBarry Smith 630d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_MG(PC pc, Vec b, Vec x) 631d71ae5a4SJacob Faibussowitsch { 632fcb023d4SJed Brown PetscFunctionBegin; 6339566063dSJacob Faibussowitsch PetscCall(PCApply_MG_Internal(pc, b, x, NULL, NULL, PETSC_FALSE)); 6343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 635fcb023d4SJed Brown } 636fcb023d4SJed Brown 637d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_MG(PC pc, Vec b, Vec x) 638d71ae5a4SJacob Faibussowitsch { 639fcb023d4SJed Brown PetscFunctionBegin; 6409566063dSJacob Faibussowitsch PetscCall(PCApply_MG_Internal(pc, b, x, NULL, NULL, PETSC_TRUE)); 6413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 64230b0564aSStefano Zampini } 64330b0564aSStefano Zampini 644d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_MG(PC pc, Mat b, Mat x) 645d71ae5a4SJacob Faibussowitsch { 64630b0564aSStefano Zampini PetscFunctionBegin; 6479566063dSJacob Faibussowitsch PetscCall(PCApply_MG_Internal(pc, NULL, NULL, b, x, PETSC_FALSE)); 6483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 649fcb023d4SJed Brown } 650f3fbd535SBarry Smith 651d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetFromOptions_MG(PC pc, PetscOptionItems *PetscOptionsObject) 652d71ae5a4SJacob Faibussowitsch { 653710315b6SLawrence Mitchell PetscInt levels, cycles; 654f3b08a26SMatthew G. Knepley PetscBool flg, flg2; 655f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 6563d3eaba7SBarry Smith PC_MG_Levels **mglevels; 657f3fbd535SBarry Smith PCMGType mgtype; 658f3fbd535SBarry Smith PCMGCycleType mgctype; 6592134b1e4SBarry Smith PCMGGalerkinType gtype; 6602b3cbbdaSStefano Zampini PCMGCoarseSpaceType coarseSpaceType; 661f3fbd535SBarry Smith 662f3fbd535SBarry Smith PetscFunctionBegin; 6631d743356SStefano Zampini levels = PetscMax(mg->nlevels, 1); 664d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Multigrid options"); 6659566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_mg_levels", "Number of Levels", "PCMGSetLevels", levels, &levels, &flg)); 6661a97d7b6SStefano Zampini if (!flg && !mg->levels && pc->dm) { 6679566063dSJacob Faibussowitsch PetscCall(DMGetRefineLevel(pc->dm, &levels)); 668eb3f98d2SBarry Smith levels++; 6693aeaf226SBarry Smith mg->usedmfornumberoflevels = PETSC_TRUE; 670eb3f98d2SBarry Smith } 6719566063dSJacob Faibussowitsch PetscCall(PCMGSetLevels(pc, levels, NULL)); 6723aeaf226SBarry Smith mglevels = mg->levels; 6733aeaf226SBarry Smith 674f3fbd535SBarry Smith mgctype = (PCMGCycleType)mglevels[0]->cycles; 6759566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_mg_cycle_type", "V cycle or for W-cycle", "PCMGSetCycleType", PCMGCycleTypes, (PetscEnum)mgctype, (PetscEnum *)&mgctype, &flg)); 6761baa6e33SBarry Smith if (flg) PetscCall(PCMGSetCycleType(pc, mgctype)); 6772134b1e4SBarry Smith gtype = mg->galerkin; 6789566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_mg_galerkin", "Use Galerkin process to compute coarser operators", "PCMGSetGalerkin", PCMGGalerkinTypes, (PetscEnum)gtype, (PetscEnum *)>ype, &flg)); 6791baa6e33SBarry Smith if (flg) PetscCall(PCMGSetGalerkin(pc, gtype)); 6802b3cbbdaSStefano Zampini coarseSpaceType = mg->coarseSpaceType; 6812b3cbbdaSStefano 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)); 6822b3cbbdaSStefano Zampini if (flg) PetscCall(PCMGSetAdaptCoarseSpaceType(pc, coarseSpaceType)); 6839566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_mg_adapt_interp_n", "Size of the coarse space for adaptive interpolation", "PCMGSetCoarseSpace", mg->Nc, &mg->Nc, &flg)); 6849566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_mg_mesp_monitor", "Monitor the multilevel eigensolver", "PCMGSetAdaptInterpolation", PETSC_FALSE, &mg->mespMonitor, &flg)); 68541b6fd38SMatthew G. Knepley flg2 = PETSC_FALSE; 6869566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_mg_adapt_cr", "Monitor coarse space quality using Compatible Relaxation (CR)", "PCMGSetAdaptCR", PETSC_FALSE, &flg2, &flg)); 6879566063dSJacob Faibussowitsch if (flg) PetscCall(PCMGSetAdaptCR(pc, flg2)); 688f56b1265SBarry Smith flg = PETSC_FALSE; 6899566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_mg_distinct_smoothup", "Create separate smoothup KSP and append the prefix _up", "PCMGSetDistinctSmoothUp", PETSC_FALSE, &flg, NULL)); 6901baa6e33SBarry Smith if (flg) PetscCall(PCMGSetDistinctSmoothUp(pc)); 69131567311SBarry Smith mgtype = mg->am; 6929566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_mg_type", "Multigrid type", "PCMGSetType", PCMGTypes, (PetscEnum)mgtype, (PetscEnum *)&mgtype, &flg)); 6931baa6e33SBarry Smith if (flg) PetscCall(PCMGSetType(pc, mgtype)); 69431567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 6959566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_mg_multiplicative_cycles", "Number of cycles for each preconditioner step", "PCMGMultiplicativeSetCycles", mg->cyclesperpcapply, &cycles, &flg)); 6961baa6e33SBarry Smith if (flg) PetscCall(PCMGMultiplicativeSetCycles(pc, cycles)); 697f3fbd535SBarry Smith } 698f3fbd535SBarry Smith flg = PETSC_FALSE; 6999566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_mg_log", "Log times for each multigrid level", "None", flg, &flg, NULL)); 700f3fbd535SBarry Smith if (flg) { 701f3fbd535SBarry Smith PetscInt i; 702f3fbd535SBarry Smith char eventname[128]; 7031a97d7b6SStefano Zampini 704f3fbd535SBarry Smith levels = mglevels[0]->levels; 705f3fbd535SBarry Smith for (i = 0; i < levels; i++) { 706a364092eSJacob Faibussowitsch PetscCall(PetscSNPrintf(eventname, PETSC_STATIC_ARRAY_LENGTH(eventname), "MGSetup Level %d", (int)i)); 7079566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(eventname, ((PetscObject)pc)->classid, &mglevels[i]->eventsmoothsetup)); 708a364092eSJacob Faibussowitsch PetscCall(PetscSNPrintf(eventname, PETSC_STATIC_ARRAY_LENGTH(eventname), "MGSmooth Level %d", (int)i)); 7099566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(eventname, ((PetscObject)pc)->classid, &mglevels[i]->eventsmoothsolve)); 710f3fbd535SBarry Smith if (i) { 711a364092eSJacob Faibussowitsch PetscCall(PetscSNPrintf(eventname, PETSC_STATIC_ARRAY_LENGTH(eventname), "MGResid Level %d", (int)i)); 7129566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(eventname, ((PetscObject)pc)->classid, &mglevels[i]->eventresidual)); 713a364092eSJacob Faibussowitsch PetscCall(PetscSNPrintf(eventname, PETSC_STATIC_ARRAY_LENGTH(eventname), "MGInterp Level %d", (int)i)); 7149566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(eventname, ((PetscObject)pc)->classid, &mglevels[i]->eventinterprestrict)); 715f3fbd535SBarry Smith } 716f3fbd535SBarry Smith } 717eec35531SMatthew G Knepley 7182611ad71SToby Isaac if (PetscDefined(USE_LOG)) { 7192611ad71SToby Isaac const char sname[] = "MG Apply"; 720eec35531SMatthew G Knepley 7212611ad71SToby Isaac PetscCall(PetscLogStageGetId(sname, &mg->stageApply)); 7222611ad71SToby Isaac if (mg->stageApply < 0) PetscCall(PetscLogStageRegister(sname, &mg->stageApply)); 723eec35531SMatthew G Knepley } 724f3fbd535SBarry Smith } 725d0609cedSBarry Smith PetscOptionsHeadEnd(); 726f3b08a26SMatthew G. Knepley /* Check option consistency */ 7279566063dSJacob Faibussowitsch PetscCall(PCMGGetGalerkin(pc, >ype)); 7289566063dSJacob Faibussowitsch PetscCall(PCMGGetAdaptInterpolation(pc, &flg)); 72908401ef6SPierre Jolivet PetscCheck(!flg || !(gtype >= PC_MG_GALERKIN_NONE), PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Must use Galerkin coarse operators when adapting the interpolator"); 7303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 731f3fbd535SBarry Smith } 732f3fbd535SBarry Smith 7330a545947SLisandro Dalcin const char *const PCMGTypes[] = {"MULTIPLICATIVE", "ADDITIVE", "FULL", "KASKADE", "PCMGType", "PC_MG", NULL}; 7340a545947SLisandro Dalcin const char *const PCMGCycleTypes[] = {"invalid", "v", "w", "PCMGCycleType", "PC_MG_CYCLE", NULL}; 7350a545947SLisandro Dalcin const char *const PCMGGalerkinTypes[] = {"both", "pmat", "mat", "none", "external", "PCMGGalerkinType", "PC_MG_GALERKIN", NULL}; 7362b3cbbdaSStefano Zampini const char *const PCMGCoarseSpaceTypes[] = {"none", "polynomial", "harmonic", "eigenvector", "generalized_eigenvector", "gdsw", "PCMGCoarseSpaceType", "PCMG_ADAPT_NONE", NULL}; 737f3fbd535SBarry Smith 7389804daf3SBarry Smith #include <petscdraw.h> 739d71ae5a4SJacob Faibussowitsch PetscErrorCode PCView_MG(PC pc, PetscViewer viewer) 740d71ae5a4SJacob Faibussowitsch { 741f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 742f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 743e3deeeafSJed Brown PetscInt levels = mglevels ? mglevels[0]->levels : 0, i; 744219b1045SBarry Smith PetscBool iascii, isbinary, isdraw; 745f3fbd535SBarry Smith 746f3fbd535SBarry Smith PetscFunctionBegin; 7479566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 7489566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary)); 7499566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw)); 750f3fbd535SBarry Smith if (iascii) { 751e3deeeafSJed Brown const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown"; 75263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " type is %s, levels=%" PetscInt_FMT " cycles=%s\n", PCMGTypes[mg->am], levels, cyclename)); 75348a46eb9SPierre Jolivet if (mg->am == PC_MG_MULTIPLICATIVE) PetscCall(PetscViewerASCIIPrintf(viewer, " Cycles per PCApply=%" PetscInt_FMT "\n", mg->cyclesperpcapply)); 7542134b1e4SBarry Smith if (mg->galerkin == PC_MG_GALERKIN_BOTH) { 7559566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Using Galerkin computed coarse grid matrices\n")); 7562134b1e4SBarry Smith } else if (mg->galerkin == PC_MG_GALERKIN_PMAT) { 7579566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Using Galerkin computed coarse grid matrices for pmat\n")); 7582134b1e4SBarry Smith } else if (mg->galerkin == PC_MG_GALERKIN_MAT) { 7599566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Using Galerkin computed coarse grid matrices for mat\n")); 7602134b1e4SBarry Smith } else if (mg->galerkin == PC_MG_GALERKIN_EXTERNAL) { 7619566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Using externally compute Galerkin coarse grid matrices\n")); 7624f66f45eSBarry Smith } else { 7639566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Not using Galerkin computed coarse grid matrices\n")); 764f3fbd535SBarry Smith } 7651baa6e33SBarry Smith if (mg->view) PetscCall((*mg->view)(pc, viewer)); 766f3fbd535SBarry Smith for (i = 0; i < levels; i++) { 76763a3b9bcSJacob Faibussowitsch if (i) { 76863a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Down solver (pre-smoother) on level %" PetscInt_FMT " -------------------------------\n", i)); 769f3fbd535SBarry Smith } else { 77063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Coarse grid solver -- level %" PetscInt_FMT " -------------------------------\n", i)); 771f3fbd535SBarry Smith } 7729566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 7739566063dSJacob Faibussowitsch PetscCall(KSPView(mglevels[i]->smoothd, viewer)); 7749566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 775f3fbd535SBarry Smith if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) { 7769566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Up solver (post-smoother) same as down solver (pre-smoother)\n")); 777f3fbd535SBarry Smith } else if (i) { 77863a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Up solver (post-smoother) on level %" PetscInt_FMT " -------------------------------\n", i)); 7799566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 7809566063dSJacob Faibussowitsch PetscCall(KSPView(mglevels[i]->smoothu, viewer)); 7819566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 782f3fbd535SBarry Smith } 78341b6fd38SMatthew G. Knepley if (i && mglevels[i]->cr) { 78463a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "CR solver on level %" PetscInt_FMT " -------------------------------\n", i)); 7859566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 7869566063dSJacob Faibussowitsch PetscCall(KSPView(mglevels[i]->cr, viewer)); 7879566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 78841b6fd38SMatthew G. Knepley } 789f3fbd535SBarry Smith } 7905b0b0462SBarry Smith } else if (isbinary) { 7915b0b0462SBarry Smith for (i = levels - 1; i >= 0; i--) { 7929566063dSJacob Faibussowitsch PetscCall(KSPView(mglevels[i]->smoothd, viewer)); 79348a46eb9SPierre Jolivet if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) PetscCall(KSPView(mglevels[i]->smoothu, viewer)); 7945b0b0462SBarry Smith } 795219b1045SBarry Smith } else if (isdraw) { 796219b1045SBarry Smith PetscDraw draw; 797b4375e8dSPeter Brune PetscReal x, w, y, bottom, th; 7989566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw)); 7999566063dSJacob Faibussowitsch PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y)); 8009566063dSJacob Faibussowitsch PetscCall(PetscDrawStringGetSize(draw, NULL, &th)); 801219b1045SBarry Smith bottom = y - th; 802219b1045SBarry Smith for (i = levels - 1; i >= 0; i--) { 803b4375e8dSPeter Brune if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) { 8049566063dSJacob Faibussowitsch PetscCall(PetscDrawPushCurrentPoint(draw, x, bottom)); 8059566063dSJacob Faibussowitsch PetscCall(KSPView(mglevels[i]->smoothd, viewer)); 8069566063dSJacob Faibussowitsch PetscCall(PetscDrawPopCurrentPoint(draw)); 807b4375e8dSPeter Brune } else { 808b4375e8dSPeter Brune w = 0.5 * PetscMin(1.0 - x, x); 8099566063dSJacob Faibussowitsch PetscCall(PetscDrawPushCurrentPoint(draw, x + w, bottom)); 8109566063dSJacob Faibussowitsch PetscCall(KSPView(mglevels[i]->smoothd, viewer)); 8119566063dSJacob Faibussowitsch PetscCall(PetscDrawPopCurrentPoint(draw)); 8129566063dSJacob Faibussowitsch PetscCall(PetscDrawPushCurrentPoint(draw, x - w, bottom)); 8139566063dSJacob Faibussowitsch PetscCall(KSPView(mglevels[i]->smoothu, viewer)); 8149566063dSJacob Faibussowitsch PetscCall(PetscDrawPopCurrentPoint(draw)); 815b4375e8dSPeter Brune } 8169566063dSJacob Faibussowitsch PetscCall(PetscDrawGetBoundingBox(draw, NULL, &bottom, NULL, NULL)); 8171cd381d6SBarry Smith bottom -= th; 818219b1045SBarry Smith } 8195b0b0462SBarry Smith } 8203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 821f3fbd535SBarry Smith } 822f3fbd535SBarry Smith 823af0996ceSBarry Smith #include <petsc/private/kspimpl.h> 824cab2e9ccSBarry Smith 825f3fbd535SBarry Smith /* 826f3fbd535SBarry Smith Calls setup for the KSP on each level 827f3fbd535SBarry Smith */ 828d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetUp_MG(PC pc) 829d71ae5a4SJacob Faibussowitsch { 830f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 831f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 8321a97d7b6SStefano Zampini PetscInt i, n; 83398e04984SBarry Smith PC cpc; 8343db492dfSBarry Smith PetscBool dump = PETSC_FALSE, opsset, use_amat, missinginterpolate = PETSC_FALSE; 835f3fbd535SBarry Smith Mat dA, dB; 836f3fbd535SBarry Smith Vec tvec; 837218a07d4SBarry Smith DM *dms; 8380a545947SLisandro Dalcin PetscViewer viewer = NULL; 83941b6fd38SMatthew G. Knepley PetscBool dAeqdB = PETSC_FALSE, needRestricts = PETSC_FALSE, doCR = PETSC_FALSE; 8402b3cbbdaSStefano Zampini PetscBool adaptInterpolation = mg->adaptInterpolation; 841f3fbd535SBarry Smith 842f3fbd535SBarry Smith PetscFunctionBegin; 84328b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels with PCMGSetLevels() before setting up"); 8441a97d7b6SStefano Zampini n = mglevels[0]->levels; 84501bc414fSDmitry Karpeev /* FIX: Move this to PCSetFromOptions_MG? */ 8463aeaf226SBarry Smith if (mg->usedmfornumberoflevels) { 8473aeaf226SBarry Smith PetscInt levels; 8489566063dSJacob Faibussowitsch PetscCall(DMGetRefineLevel(pc->dm, &levels)); 8493aeaf226SBarry Smith levels++; 8503aeaf226SBarry Smith if (levels > n) { /* the problem is now being solved on a finer grid */ 8519566063dSJacob Faibussowitsch PetscCall(PCMGSetLevels(pc, levels, NULL)); 8523aeaf226SBarry Smith n = levels; 8539566063dSJacob Faibussowitsch PetscCall(PCSetFromOptions(pc)); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */ 8543aeaf226SBarry Smith mglevels = mg->levels; 8553aeaf226SBarry Smith } 8563aeaf226SBarry Smith } 8579566063dSJacob Faibussowitsch PetscCall(KSPGetPC(mglevels[0]->smoothd, &cpc)); 8583aeaf226SBarry Smith 859f3fbd535SBarry Smith /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */ 860f3fbd535SBarry Smith /* so use those from global PC */ 861f3fbd535SBarry Smith /* Is this what we always want? What if user wants to keep old one? */ 8629566063dSJacob Faibussowitsch PetscCall(KSPGetOperatorsSet(mglevels[n - 1]->smoothd, NULL, &opsset)); 86398e04984SBarry Smith if (opsset) { 86498e04984SBarry Smith Mat mmat; 8659566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[n - 1]->smoothd, NULL, &mmat)); 86698e04984SBarry Smith if (mmat == pc->pmat) opsset = PETSC_FALSE; 86798e04984SBarry Smith } 8684a5f07a5SMark F. Adams 86941b6fd38SMatthew G. Knepley /* Create CR solvers */ 8709566063dSJacob Faibussowitsch PetscCall(PCMGGetAdaptCR(pc, &doCR)); 87141b6fd38SMatthew G. Knepley if (doCR) { 87241b6fd38SMatthew G. Knepley const char *prefix; 87341b6fd38SMatthew G. Knepley 8749566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &prefix)); 87541b6fd38SMatthew G. Knepley for (i = 1; i < n; ++i) { 87641b6fd38SMatthew G. Knepley PC ipc, cr; 87741b6fd38SMatthew G. Knepley char crprefix[128]; 87841b6fd38SMatthew G. Knepley 8799566063dSJacob Faibussowitsch PetscCall(KSPCreate(PetscObjectComm((PetscObject)pc), &mglevels[i]->cr)); 8803821be0aSBarry Smith PetscCall(KSPSetNestLevel(mglevels[i]->cr, pc->kspnestlevel)); 8819566063dSJacob Faibussowitsch PetscCall(KSPSetErrorIfNotConverged(mglevels[i]->cr, PETSC_FALSE)); 8829566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->cr, (PetscObject)pc, n - i)); 8839566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(mglevels[i]->cr, prefix)); 8849566063dSJacob Faibussowitsch PetscCall(PetscObjectComposedDataSetInt((PetscObject)mglevels[i]->cr, PetscMGLevelId, mglevels[i]->level)); 8859566063dSJacob Faibussowitsch PetscCall(KSPSetType(mglevels[i]->cr, KSPCHEBYSHEV)); 8869566063dSJacob Faibussowitsch PetscCall(KSPSetConvergenceTest(mglevels[i]->cr, KSPConvergedSkip, NULL, NULL)); 8879566063dSJacob Faibussowitsch PetscCall(KSPSetNormType(mglevels[i]->cr, KSP_NORM_PRECONDITIONED)); 8889566063dSJacob Faibussowitsch PetscCall(KSPGetPC(mglevels[i]->cr, &ipc)); 88941b6fd38SMatthew G. Knepley 8909566063dSJacob Faibussowitsch PetscCall(PCSetType(ipc, PCCOMPOSITE)); 8919566063dSJacob Faibussowitsch PetscCall(PCCompositeSetType(ipc, PC_COMPOSITE_MULTIPLICATIVE)); 8929566063dSJacob Faibussowitsch PetscCall(PCCompositeAddPCType(ipc, PCSOR)); 8939566063dSJacob Faibussowitsch PetscCall(CreateCR_Private(pc, i, &cr)); 8949566063dSJacob Faibussowitsch PetscCall(PCCompositeAddPC(ipc, cr)); 8959566063dSJacob Faibussowitsch PetscCall(PCDestroy(&cr)); 89641b6fd38SMatthew G. Knepley 8979566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(mglevels[i]->cr, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, mg->default_smoothd)); 8989566063dSJacob Faibussowitsch PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->cr, PETSC_TRUE)); 8999566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(crprefix, 128, "mg_levels_%d_cr_", (int)i)); 9009566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(mglevels[i]->cr, crprefix)); 90141b6fd38SMatthew G. Knepley } 90241b6fd38SMatthew G. Knepley } 90341b6fd38SMatthew G. Knepley 9044a5f07a5SMark F. Adams if (!opsset) { 9059566063dSJacob Faibussowitsch PetscCall(PCGetUseAmat(pc, &use_amat)); 9064a5f07a5SMark F. Adams if (use_amat) { 9079566063dSJacob 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")); 9089566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(mglevels[n - 1]->smoothd, pc->mat, pc->pmat)); 90969858f1bSStefano Zampini } else { 9109566063dSJacob 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")); 9119566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(mglevels[n - 1]->smoothd, pc->pmat, pc->pmat)); 9124a5f07a5SMark F. Adams } 9134a5f07a5SMark F. Adams } 914f3fbd535SBarry Smith 9159c7ffb39SBarry Smith for (i = n - 1; i > 0; i--) { 9169c7ffb39SBarry Smith if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) { 9179c7ffb39SBarry Smith missinginterpolate = PETSC_TRUE; 9182b3cbbdaSStefano Zampini break; 9199c7ffb39SBarry Smith } 9209c7ffb39SBarry Smith } 9212134b1e4SBarry Smith 9229566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[n - 1]->smoothd, &dA, &dB)); 9232134b1e4SBarry Smith if (dA == dB) dAeqdB = PETSC_TRUE; 9242b3cbbdaSStefano Zampini if (mg->galerkin == PC_MG_GALERKIN_NONE || ((mg->galerkin == PC_MG_GALERKIN_PMAT || mg->galerkin == PC_MG_GALERKIN_MAT) && !dAeqdB)) { 9252134b1e4SBarry Smith needRestricts = PETSC_TRUE; /* user must compute either mat, pmat, or both so must restrict x to coarser levels */ 9262134b1e4SBarry Smith } 9272134b1e4SBarry Smith 9282b3cbbdaSStefano Zampini if (pc->dm && !pc->setupcalled) { 9292b3cbbdaSStefano Zampini /* finest smoother also gets DM but it is not active, independent of whether galerkin==PC_MG_GALERKIN_EXTERNAL */ 9302b3cbbdaSStefano Zampini PetscCall(KSPSetDM(mglevels[n - 1]->smoothd, pc->dm)); 9312b3cbbdaSStefano Zampini PetscCall(KSPSetDMActive(mglevels[n - 1]->smoothd, PETSC_FALSE)); 9322b3cbbdaSStefano Zampini if (mglevels[n - 1]->smoothd != mglevels[n - 1]->smoothu) { 9332b3cbbdaSStefano Zampini PetscCall(KSPSetDM(mglevels[n - 1]->smoothu, pc->dm)); 9342b3cbbdaSStefano Zampini PetscCall(KSPSetDMActive(mglevels[n - 1]->smoothu, PETSC_FALSE)); 9352b3cbbdaSStefano Zampini } 9362b3cbbdaSStefano Zampini if (mglevels[n - 1]->cr) { 9372b3cbbdaSStefano Zampini PetscCall(KSPSetDM(mglevels[n - 1]->cr, pc->dm)); 9382b3cbbdaSStefano Zampini PetscCall(KSPSetDMActive(mglevels[n - 1]->cr, PETSC_FALSE)); 9392b3cbbdaSStefano Zampini } 9402b3cbbdaSStefano Zampini } 9412b3cbbdaSStefano Zampini 9429c7ffb39SBarry Smith /* 9439c7ffb39SBarry Smith Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS) 9442b3cbbdaSStefano Zampini Skipping for externally managed hierarchy (such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs? 9459c7ffb39SBarry Smith */ 94632fe681dSStefano Zampini if (missinginterpolate && mg->galerkin != PC_MG_GALERKIN_EXTERNAL && !pc->setupcalled) { 94732fe681dSStefano Zampini /* first see if we can compute a coarse space */ 94832fe681dSStefano Zampini if (mg->coarseSpaceType == PCMG_ADAPT_GDSW) { 94932fe681dSStefano Zampini for (i = n - 2; i > -1; i--) { 95032fe681dSStefano Zampini if (!mglevels[i + 1]->restrct && !mglevels[i + 1]->interpolate) { 95132fe681dSStefano Zampini PetscCall(PCMGComputeCoarseSpace_Internal(pc, i + 1, mg->coarseSpaceType, mg->Nc, NULL, &mglevels[i + 1]->coarseSpace)); 95232fe681dSStefano Zampini PetscCall(PCMGSetInterpolation(pc, i + 1, mglevels[i + 1]->coarseSpace)); 95332fe681dSStefano Zampini } 95432fe681dSStefano Zampini } 95532fe681dSStefano Zampini } else { /* construct the interpolation from the DMs */ 9562e499ae9SBarry Smith Mat p; 95773250ac0SBarry Smith Vec rscale; 9589566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &dms)); 959218a07d4SBarry Smith dms[n - 1] = pc->dm; 960ef1267afSMatthew G. Knepley /* Separately create them so we do not get DMKSP interference between levels */ 9619566063dSJacob Faibussowitsch for (i = n - 2; i > -1; i--) PetscCall(DMCoarsen(dms[i + 1], MPI_COMM_NULL, &dms[i])); 962218a07d4SBarry Smith for (i = n - 2; i > -1; i--) { 963942e3340SBarry Smith DMKSP kdm; 964eab52d2bSLawrence Mitchell PetscBool dmhasrestrict, dmhasinject; 9659566063dSJacob Faibussowitsch PetscCall(KSPSetDM(mglevels[i]->smoothd, dms[i])); 9669566063dSJacob Faibussowitsch if (!needRestricts) PetscCall(KSPSetDMActive(mglevels[i]->smoothd, PETSC_FALSE)); 967c27ee7a3SPatrick Farrell if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 9689566063dSJacob Faibussowitsch PetscCall(KSPSetDM(mglevels[i]->smoothu, dms[i])); 9699566063dSJacob Faibussowitsch if (!needRestricts) PetscCall(KSPSetDMActive(mglevels[i]->smoothu, PETSC_FALSE)); 970c27ee7a3SPatrick Farrell } 97141b6fd38SMatthew G. Knepley if (mglevels[i]->cr) { 9729566063dSJacob Faibussowitsch PetscCall(KSPSetDM(mglevels[i]->cr, dms[i])); 9739566063dSJacob Faibussowitsch if (!needRestricts) PetscCall(KSPSetDMActive(mglevels[i]->cr, PETSC_FALSE)); 97441b6fd38SMatthew G. Knepley } 9759566063dSJacob Faibussowitsch PetscCall(DMGetDMKSPWrite(dms[i], &kdm)); 976d1a3e677SJed Brown /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take 977d1a3e677SJed Brown * a bitwise OR of computing the matrix, RHS, and initial iterate. */ 9780298fd71SBarry Smith kdm->ops->computerhs = NULL; 9790298fd71SBarry Smith kdm->rhsctx = NULL; 98024c3aa18SBarry Smith if (!mglevels[i + 1]->interpolate) { 9819566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(dms[i], dms[i + 1], &p, &rscale)); 9829566063dSJacob Faibussowitsch PetscCall(PCMGSetInterpolation(pc, i + 1, p)); 9839566063dSJacob Faibussowitsch if (rscale) PetscCall(PCMGSetRScale(pc, i + 1, rscale)); 9849566063dSJacob Faibussowitsch PetscCall(VecDestroy(&rscale)); 9859566063dSJacob Faibussowitsch PetscCall(MatDestroy(&p)); 986218a07d4SBarry Smith } 9879566063dSJacob Faibussowitsch PetscCall(DMHasCreateRestriction(dms[i], &dmhasrestrict)); 9883ad4599aSBarry Smith if (dmhasrestrict && !mglevels[i + 1]->restrct) { 9899566063dSJacob Faibussowitsch PetscCall(DMCreateRestriction(dms[i], dms[i + 1], &p)); 9909566063dSJacob Faibussowitsch PetscCall(PCMGSetRestriction(pc, i + 1, p)); 9919566063dSJacob Faibussowitsch PetscCall(MatDestroy(&p)); 9923ad4599aSBarry Smith } 9939566063dSJacob Faibussowitsch PetscCall(DMHasCreateInjection(dms[i], &dmhasinject)); 994eab52d2bSLawrence Mitchell if (dmhasinject && !mglevels[i + 1]->inject) { 9959566063dSJacob Faibussowitsch PetscCall(DMCreateInjection(dms[i], dms[i + 1], &p)); 9969566063dSJacob Faibussowitsch PetscCall(PCMGSetInjection(pc, i + 1, p)); 9979566063dSJacob Faibussowitsch PetscCall(MatDestroy(&p)); 998eab52d2bSLawrence Mitchell } 99924c3aa18SBarry Smith } 10002d2b81a6SBarry Smith 10019566063dSJacob Faibussowitsch for (i = n - 2; i > -1; i--) PetscCall(DMDestroy(&dms[i])); 10029566063dSJacob Faibussowitsch PetscCall(PetscFree(dms)); 1003ce4cda84SJed Brown } 100432fe681dSStefano Zampini } 1005cab2e9ccSBarry Smith 10062134b1e4SBarry Smith if (mg->galerkin < PC_MG_GALERKIN_NONE) { 10072134b1e4SBarry Smith Mat A, B; 10082134b1e4SBarry Smith PetscBool doA = PETSC_FALSE, doB = PETSC_FALSE; 10092134b1e4SBarry Smith MatReuse reuse = MAT_INITIAL_MATRIX; 10102134b1e4SBarry Smith 10112b3cbbdaSStefano Zampini if (mg->galerkin == PC_MG_GALERKIN_PMAT || mg->galerkin == PC_MG_GALERKIN_BOTH) doB = PETSC_TRUE; 10122b3cbbdaSStefano Zampini if (mg->galerkin == PC_MG_GALERKIN_MAT || (mg->galerkin == PC_MG_GALERKIN_BOTH && dA != dB)) doA = PETSC_TRUE; 10132134b1e4SBarry Smith if (pc->setupcalled) reuse = MAT_REUSE_MATRIX; 1014f3fbd535SBarry Smith for (i = n - 2; i > -1; i--) { 10157827d75bSBarry 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"); 101648a46eb9SPierre Jolivet if (!mglevels[i + 1]->interpolate) PetscCall(PCMGSetInterpolation(pc, i + 1, mglevels[i + 1]->restrct)); 101748a46eb9SPierre Jolivet if (!mglevels[i + 1]->restrct) PetscCall(PCMGSetRestriction(pc, i + 1, mglevels[i + 1]->interpolate)); 101848a46eb9SPierre Jolivet if (reuse == MAT_REUSE_MATRIX) PetscCall(KSPGetOperators(mglevels[i]->smoothd, &A, &B)); 101948a46eb9SPierre Jolivet if (doA) PetscCall(MatGalerkin(mglevels[i + 1]->restrct, dA, mglevels[i + 1]->interpolate, reuse, 1.0, &A)); 102048a46eb9SPierre Jolivet if (doB) PetscCall(MatGalerkin(mglevels[i + 1]->restrct, dB, mglevels[i + 1]->interpolate, reuse, 1.0, &B)); 10212134b1e4SBarry Smith /* the management of the PetscObjectReference() and PetscObjecDereference() below is rather delicate */ 10222134b1e4SBarry Smith if (!doA && dAeqdB) { 10239566063dSJacob Faibussowitsch if (reuse == MAT_INITIAL_MATRIX) PetscCall(PetscObjectReference((PetscObject)B)); 10242134b1e4SBarry Smith A = B; 10252134b1e4SBarry Smith } else if (!doA && reuse == MAT_INITIAL_MATRIX) { 10269566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[i]->smoothd, &A, NULL)); 10279566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)A)); 1028b9367914SBarry Smith } 10292134b1e4SBarry Smith if (!doB && dAeqdB) { 10309566063dSJacob Faibussowitsch if (reuse == MAT_INITIAL_MATRIX) PetscCall(PetscObjectReference((PetscObject)A)); 10312134b1e4SBarry Smith B = A; 10322134b1e4SBarry Smith } else if (!doB && reuse == MAT_INITIAL_MATRIX) { 10339566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[i]->smoothd, NULL, &B)); 10349566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)B)); 10357d537d92SJed Brown } 10362134b1e4SBarry Smith if (reuse == MAT_INITIAL_MATRIX) { 10379566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(mglevels[i]->smoothd, A, B)); 10389566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)A)); 10399566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)B)); 10402134b1e4SBarry Smith } 10412134b1e4SBarry Smith dA = A; 1042f3fbd535SBarry Smith dB = B; 1043f3fbd535SBarry Smith } 1044f3fbd535SBarry Smith } 1045f3b08a26SMatthew G. Knepley 1046f3b08a26SMatthew G. Knepley /* Adapt interpolation matrices */ 10472b3cbbdaSStefano Zampini if (adaptInterpolation) { 1048f3b08a26SMatthew G. Knepley for (i = 0; i < n; ++i) { 104948a46eb9SPierre Jolivet if (!mglevels[i]->coarseSpace) PetscCall(PCMGComputeCoarseSpace_Internal(pc, i, mg->coarseSpaceType, mg->Nc, !i ? NULL : mglevels[i - 1]->coarseSpace, &mglevels[i]->coarseSpace)); 10502b3cbbdaSStefano Zampini if (i) PetscCall(PCMGAdaptInterpolator_Internal(pc, i, mglevels[i - 1]->smoothu, mglevels[i]->smoothu, mglevels[i - 1]->coarseSpace, mglevels[i]->coarseSpace)); 1051f3b08a26SMatthew G. Knepley } 105248a46eb9SPierre Jolivet for (i = n - 2; i > -1; --i) PetscCall(PCMGRecomputeLevelOperators_Internal(pc, i)); 1053f3b08a26SMatthew G. Knepley } 1054f3b08a26SMatthew G. Knepley 10552134b1e4SBarry Smith if (needRestricts && pc->dm) { 1056caa4e7f2SJed Brown for (i = n - 2; i >= 0; i--) { 1057ccceb50cSJed Brown DM dmfine, dmcoarse; 1058caa4e7f2SJed Brown Mat Restrict, Inject; 1059caa4e7f2SJed Brown Vec rscale; 10609566063dSJacob Faibussowitsch PetscCall(KSPGetDM(mglevels[i + 1]->smoothd, &dmfine)); 10619566063dSJacob Faibussowitsch PetscCall(KSPGetDM(mglevels[i]->smoothd, &dmcoarse)); 10629566063dSJacob Faibussowitsch PetscCall(PCMGGetRestriction(pc, i + 1, &Restrict)); 10639566063dSJacob Faibussowitsch PetscCall(PCMGGetRScale(pc, i + 1, &rscale)); 10649566063dSJacob Faibussowitsch PetscCall(PCMGGetInjection(pc, i + 1, &Inject)); 10659566063dSJacob Faibussowitsch PetscCall(DMRestrict(dmfine, Restrict, rscale, Inject, dmcoarse)); 1066caa4e7f2SJed Brown } 1067caa4e7f2SJed Brown } 1068f3fbd535SBarry Smith 1069f3fbd535SBarry Smith if (!pc->setupcalled) { 107048a46eb9SPierre Jolivet for (i = 0; i < n; i++) PetscCall(KSPSetFromOptions(mglevels[i]->smoothd)); 1071f3fbd535SBarry Smith for (i = 1; i < n; i++) { 107248a46eb9SPierre Jolivet if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) PetscCall(KSPSetFromOptions(mglevels[i]->smoothu)); 107348a46eb9SPierre Jolivet if (mglevels[i]->cr) PetscCall(KSPSetFromOptions(mglevels[i]->cr)); 1074f3fbd535SBarry Smith } 10753ad4599aSBarry Smith /* insure that if either interpolation or restriction is set the other other one is set */ 1076f3fbd535SBarry Smith for (i = 1; i < n; i++) { 10779566063dSJacob Faibussowitsch PetscCall(PCMGGetInterpolation(pc, i, NULL)); 10789566063dSJacob Faibussowitsch PetscCall(PCMGGetRestriction(pc, i, NULL)); 1079f3fbd535SBarry Smith } 1080f3fbd535SBarry Smith for (i = 0; i < n - 1; i++) { 1081f3fbd535SBarry Smith if (!mglevels[i]->b) { 1082f3fbd535SBarry Smith Vec *vec; 10839566063dSJacob Faibussowitsch PetscCall(KSPCreateVecs(mglevels[i]->smoothd, 1, &vec, 0, NULL)); 10849566063dSJacob Faibussowitsch PetscCall(PCMGSetRhs(pc, i, *vec)); 10859566063dSJacob Faibussowitsch PetscCall(VecDestroy(vec)); 10869566063dSJacob Faibussowitsch PetscCall(PetscFree(vec)); 1087f3fbd535SBarry Smith } 1088f3fbd535SBarry Smith if (!mglevels[i]->r && i) { 10899566063dSJacob Faibussowitsch PetscCall(VecDuplicate(mglevels[i]->b, &tvec)); 10909566063dSJacob Faibussowitsch PetscCall(PCMGSetR(pc, i, tvec)); 10919566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tvec)); 1092f3fbd535SBarry Smith } 1093f3fbd535SBarry Smith if (!mglevels[i]->x) { 10949566063dSJacob Faibussowitsch PetscCall(VecDuplicate(mglevels[i]->b, &tvec)); 10959566063dSJacob Faibussowitsch PetscCall(PCMGSetX(pc, i, tvec)); 10969566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tvec)); 1097f3fbd535SBarry Smith } 109841b6fd38SMatthew G. Knepley if (doCR) { 10999566063dSJacob Faibussowitsch PetscCall(VecDuplicate(mglevels[i]->b, &mglevels[i]->crx)); 11009566063dSJacob Faibussowitsch PetscCall(VecDuplicate(mglevels[i]->b, &mglevels[i]->crb)); 11019566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(mglevels[i]->crb)); 110241b6fd38SMatthew G. Knepley } 1103f3fbd535SBarry Smith } 1104f3fbd535SBarry Smith if (n != 1 && !mglevels[n - 1]->r) { 1105f3fbd535SBarry Smith /* PCMGSetR() on the finest level if user did not supply it */ 1106f3fbd535SBarry Smith Vec *vec; 11079566063dSJacob Faibussowitsch PetscCall(KSPCreateVecs(mglevels[n - 1]->smoothd, 1, &vec, 0, NULL)); 11089566063dSJacob Faibussowitsch PetscCall(PCMGSetR(pc, n - 1, *vec)); 11099566063dSJacob Faibussowitsch PetscCall(VecDestroy(vec)); 11109566063dSJacob Faibussowitsch PetscCall(PetscFree(vec)); 1111f3fbd535SBarry Smith } 111241b6fd38SMatthew G. Knepley if (doCR) { 11139566063dSJacob Faibussowitsch PetscCall(VecDuplicate(mglevels[n - 1]->r, &mglevels[n - 1]->crx)); 11149566063dSJacob Faibussowitsch PetscCall(VecDuplicate(mglevels[n - 1]->r, &mglevels[n - 1]->crb)); 11159566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(mglevels[n - 1]->crb)); 111641b6fd38SMatthew G. Knepley } 1117f3fbd535SBarry Smith } 1118f3fbd535SBarry Smith 111998e04984SBarry Smith if (pc->dm) { 112098e04984SBarry Smith /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */ 112198e04984SBarry Smith for (i = 0; i < n - 1; i++) { 112298e04984SBarry Smith if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX; 112398e04984SBarry Smith } 112498e04984SBarry Smith } 112508549ed4SJed Brown // We got here (PCSetUp_MG) because the matrix has changed, which means the smoother needs to be set up again (e.g., 112608549ed4SJed Brown // new diagonal for Jacobi). Setting it here allows it to be logged under PCSetUp rather than deep inside a PCApply. 112708549ed4SJed Brown if (mglevels[n - 1]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[n - 1]->smoothd->setupstage = KSP_SETUP_NEWMATRIX; 1128f3fbd535SBarry Smith 1129f3fbd535SBarry Smith for (i = 1; i < n; i++) { 11302a21e185SBarry Smith if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1) { 1131f3fbd535SBarry Smith /* if doing only down then initial guess is zero */ 11329566063dSJacob Faibussowitsch PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->smoothd, PETSC_TRUE)); 1133f3fbd535SBarry Smith } 11349566063dSJacob Faibussowitsch if (mglevels[i]->cr) PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->cr, PETSC_TRUE)); 11359566063dSJacob Faibussowitsch if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventBegin(mglevels[i]->eventsmoothsetup, 0, 0, 0, 0)); 11369566063dSJacob Faibussowitsch PetscCall(KSPSetUp(mglevels[i]->smoothd)); 1137d4645a75SStefano Zampini if (mglevels[i]->smoothd->reason) pc->failedreason = PC_SUBPC_ERROR; 11389566063dSJacob Faibussowitsch if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventEnd(mglevels[i]->eventsmoothsetup, 0, 0, 0, 0)); 1139d42688cbSBarry Smith if (!mglevels[i]->residual) { 1140d42688cbSBarry Smith Mat mat; 11419566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[i]->smoothd, &mat, NULL)); 11429566063dSJacob Faibussowitsch PetscCall(PCMGSetResidual(pc, i, PCMGResidualDefault, mat)); 1143d42688cbSBarry Smith } 1144fcb023d4SJed Brown if (!mglevels[i]->residualtranspose) { 1145fcb023d4SJed Brown Mat mat; 11469566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[i]->smoothd, &mat, NULL)); 11479566063dSJacob Faibussowitsch PetscCall(PCMGSetResidualTranspose(pc, i, PCMGResidualTransposeDefault, mat)); 1148fcb023d4SJed Brown } 1149f3fbd535SBarry Smith } 1150f3fbd535SBarry Smith for (i = 1; i < n; i++) { 1151f3fbd535SBarry Smith if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) { 1152f3fbd535SBarry Smith Mat downmat, downpmat; 1153f3fbd535SBarry Smith 1154f3fbd535SBarry Smith /* check if operators have been set for up, if not use down operators to set them */ 11559566063dSJacob Faibussowitsch PetscCall(KSPGetOperatorsSet(mglevels[i]->smoothu, &opsset, NULL)); 1156f3fbd535SBarry Smith if (!opsset) { 11579566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[i]->smoothd, &downmat, &downpmat)); 11589566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(mglevels[i]->smoothu, downmat, downpmat)); 1159f3fbd535SBarry Smith } 1160f3fbd535SBarry Smith 11619566063dSJacob Faibussowitsch PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->smoothu, PETSC_TRUE)); 11629566063dSJacob Faibussowitsch if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventBegin(mglevels[i]->eventsmoothsetup, 0, 0, 0, 0)); 11639566063dSJacob Faibussowitsch PetscCall(KSPSetUp(mglevels[i]->smoothu)); 1164d4645a75SStefano Zampini if (mglevels[i]->smoothu->reason) pc->failedreason = PC_SUBPC_ERROR; 11659566063dSJacob Faibussowitsch if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventEnd(mglevels[i]->eventsmoothsetup, 0, 0, 0, 0)); 1166f3fbd535SBarry Smith } 116741b6fd38SMatthew G. Knepley if (mglevels[i]->cr) { 116841b6fd38SMatthew G. Knepley Mat downmat, downpmat; 116941b6fd38SMatthew G. Knepley 117041b6fd38SMatthew G. Knepley /* check if operators have been set for up, if not use down operators to set them */ 11719566063dSJacob Faibussowitsch PetscCall(KSPGetOperatorsSet(mglevels[i]->cr, &opsset, NULL)); 117241b6fd38SMatthew G. Knepley if (!opsset) { 11739566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[i]->smoothd, &downmat, &downpmat)); 11749566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(mglevels[i]->cr, downmat, downpmat)); 117541b6fd38SMatthew G. Knepley } 117641b6fd38SMatthew G. Knepley 11779566063dSJacob Faibussowitsch PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->cr, PETSC_TRUE)); 11789566063dSJacob Faibussowitsch if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventBegin(mglevels[i]->eventsmoothsetup, 0, 0, 0, 0)); 11799566063dSJacob Faibussowitsch PetscCall(KSPSetUp(mglevels[i]->cr)); 1180d4645a75SStefano Zampini if (mglevels[i]->cr->reason) pc->failedreason = PC_SUBPC_ERROR; 11819566063dSJacob Faibussowitsch if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventEnd(mglevels[i]->eventsmoothsetup, 0, 0, 0, 0)); 118241b6fd38SMatthew G. Knepley } 1183f3fbd535SBarry Smith } 1184f3fbd535SBarry Smith 11859566063dSJacob Faibussowitsch if (mglevels[0]->eventsmoothsetup) PetscCall(PetscLogEventBegin(mglevels[0]->eventsmoothsetup, 0, 0, 0, 0)); 11869566063dSJacob Faibussowitsch PetscCall(KSPSetUp(mglevels[0]->smoothd)); 1187d4645a75SStefano Zampini if (mglevels[0]->smoothd->reason) pc->failedreason = PC_SUBPC_ERROR; 11889566063dSJacob Faibussowitsch if (mglevels[0]->eventsmoothsetup) PetscCall(PetscLogEventEnd(mglevels[0]->eventsmoothsetup, 0, 0, 0, 0)); 1189f3fbd535SBarry Smith 1190f3fbd535SBarry Smith /* 1191f3fbd535SBarry Smith Dump the interpolation/restriction matrices plus the 1192e3c5b3baSBarry Smith Jacobian/stiffness on each level. This allows MATLAB users to 1193f3fbd535SBarry Smith easily check if the Galerkin condition A_c = R A_f R^T is satisfied. 1194f3fbd535SBarry Smith 1195f3fbd535SBarry Smith Only support one or the other at the same time. 1196f3fbd535SBarry Smith */ 1197f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER) 11989566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_mg_dump_matlab", &dump, NULL)); 1199ce94432eSBarry Smith if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc)); 1200f3fbd535SBarry Smith dump = PETSC_FALSE; 1201f3fbd535SBarry Smith #endif 12029566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_mg_dump_binary", &dump, NULL)); 1203ce94432eSBarry Smith if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc)); 1204f3fbd535SBarry Smith 1205f3fbd535SBarry Smith if (viewer) { 120648a46eb9SPierre Jolivet for (i = 1; i < n; i++) PetscCall(MatView(mglevels[i]->restrct, viewer)); 1207f3fbd535SBarry Smith for (i = 0; i < n; i++) { 12089566063dSJacob Faibussowitsch PetscCall(KSPGetPC(mglevels[i]->smoothd, &pc)); 12099566063dSJacob Faibussowitsch PetscCall(MatView(pc->mat, viewer)); 1210f3fbd535SBarry Smith } 1211f3fbd535SBarry Smith } 12123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1213f3fbd535SBarry Smith } 1214f3fbd535SBarry Smith 1215d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetLevels_MG(PC pc, PetscInt *levels) 1216d71ae5a4SJacob Faibussowitsch { 121797d33e41SMatthew G. Knepley PC_MG *mg = (PC_MG *)pc->data; 121897d33e41SMatthew G. Knepley 121997d33e41SMatthew G. Knepley PetscFunctionBegin; 122097d33e41SMatthew G. Knepley *levels = mg->nlevels; 12213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 122297d33e41SMatthew G. Knepley } 122397d33e41SMatthew G. Knepley 12244b9ad928SBarry Smith /*@ 1225f1580f4eSBarry Smith PCMGGetLevels - Gets the number of levels to use with `PCMG`. 12264b9ad928SBarry Smith 12274b9ad928SBarry Smith Not Collective 12284b9ad928SBarry Smith 12294b9ad928SBarry Smith Input Parameter: 12304b9ad928SBarry Smith . pc - the preconditioner context 12314b9ad928SBarry Smith 1232feefa0e1SJacob Faibussowitsch Output Parameter: 12334b9ad928SBarry Smith . levels - the number of levels 12344b9ad928SBarry Smith 12354b9ad928SBarry Smith Level: advanced 12364b9ad928SBarry Smith 1237*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGSetLevels()` 12384b9ad928SBarry Smith @*/ 1239d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetLevels(PC pc, PetscInt *levels) 1240d71ae5a4SJacob Faibussowitsch { 12414b9ad928SBarry Smith PetscFunctionBegin; 12420700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 12434f572ea9SToby Isaac PetscAssertPointer(levels, 2); 124497d33e41SMatthew G. Knepley *levels = 0; 1245cac4c232SBarry Smith PetscTryMethod(pc, "PCMGGetLevels_C", (PC, PetscInt *), (pc, levels)); 12463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12474b9ad928SBarry Smith } 12484b9ad928SBarry Smith 12494b9ad928SBarry Smith /*@ 1250f1580f4eSBarry Smith PCMGGetGridComplexity - compute operator and grid complexity of the `PCMG` hierarchy 1251e7d4b4cbSMark Adams 1252e7d4b4cbSMark Adams Input Parameter: 1253e7d4b4cbSMark Adams . pc - the preconditioner context 1254e7d4b4cbSMark Adams 125590db8557SMark Adams Output Parameters: 125690db8557SMark Adams + gc - grid complexity = sum_i(n_i) / n_0 125790db8557SMark Adams - oc - operator complexity = sum_i(nnz_i) / nnz_0 1258e7d4b4cbSMark Adams 1259e7d4b4cbSMark Adams Level: advanced 1260e7d4b4cbSMark Adams 1261f1580f4eSBarry Smith Note: 1262f1580f4eSBarry Smith This is often call the operator complexity in multigrid literature 1263f1580f4eSBarry Smith 1264*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGGetLevels()`, `PCMGSetLevels()` 1265e7d4b4cbSMark Adams @*/ 1266d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetGridComplexity(PC pc, PetscReal *gc, PetscReal *oc) 1267d71ae5a4SJacob Faibussowitsch { 1268e7d4b4cbSMark Adams PC_MG *mg = (PC_MG *)pc->data; 1269e7d4b4cbSMark Adams PC_MG_Levels **mglevels = mg->levels; 1270e7d4b4cbSMark Adams PetscInt lev, N; 1271e7d4b4cbSMark Adams PetscLogDouble nnz0 = 0, sgc = 0, soc = 0, n0 = 0; 1272e7d4b4cbSMark Adams MatInfo info; 1273e7d4b4cbSMark Adams 1274e7d4b4cbSMark Adams PetscFunctionBegin; 127590db8557SMark Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 12764f572ea9SToby Isaac if (gc) PetscAssertPointer(gc, 2); 12774f572ea9SToby Isaac if (oc) PetscAssertPointer(oc, 3); 1278e7d4b4cbSMark Adams if (!pc->setupcalled) { 127990db8557SMark Adams if (gc) *gc = 0; 128090db8557SMark Adams if (oc) *oc = 0; 12813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1282e7d4b4cbSMark Adams } 1283e7d4b4cbSMark Adams PetscCheck(mg->nlevels > 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MG has no levels"); 1284e7d4b4cbSMark Adams for (lev = 0; lev < mg->nlevels; lev++) { 1285e7d4b4cbSMark Adams Mat dB; 12869566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[lev]->smoothd, NULL, &dB)); 12879566063dSJacob Faibussowitsch PetscCall(MatGetInfo(dB, MAT_GLOBAL_SUM, &info)); /* global reduction */ 12889566063dSJacob Faibussowitsch PetscCall(MatGetSize(dB, &N, NULL)); 1289e7d4b4cbSMark Adams sgc += N; 1290e7d4b4cbSMark Adams soc += info.nz_used; 12919371c9d4SSatish Balay if (lev == mg->nlevels - 1) { 12929371c9d4SSatish Balay nnz0 = info.nz_used; 12939371c9d4SSatish Balay n0 = N; 12949371c9d4SSatish Balay } 1295e7d4b4cbSMark Adams } 12960fdf79fbSJacob Faibussowitsch PetscCheck(n0 > 0 && gc, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number for grid points on finest level is not available"); 12970fdf79fbSJacob Faibussowitsch *gc = (PetscReal)(sgc / n0); 129890db8557SMark Adams if (nnz0 > 0 && oc) *oc = (PetscReal)(soc / nnz0); 12993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1300e7d4b4cbSMark Adams } 1301e7d4b4cbSMark Adams 1302e7d4b4cbSMark Adams /*@ 130304c3f3b8SBarry Smith PCMGSetType - Determines the form of multigrid to use, either 13044b9ad928SBarry Smith multiplicative, additive, full, or the Kaskade algorithm. 13054b9ad928SBarry Smith 1306c3339decSBarry Smith Logically Collective 13074b9ad928SBarry Smith 13084b9ad928SBarry Smith Input Parameters: 13094b9ad928SBarry Smith + pc - the preconditioner context 1310f1580f4eSBarry Smith - form - multigrid form, one of `PC_MG_MULTIPLICATIVE`, `PC_MG_ADDITIVE`, `PC_MG_FULL`, `PC_MG_KASKADE` 13114b9ad928SBarry Smith 13124b9ad928SBarry Smith Options Database Key: 131320f4b53cSBarry Smith . -pc_mg_type <form> - Sets <form>, one of multiplicative, additive, full, kaskade 13144b9ad928SBarry Smith 13154b9ad928SBarry Smith Level: advanced 13164b9ad928SBarry Smith 1317*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMGType`, `PCMG`, `PCMGGetLevels()`, `PCMGSetLevels()`, `PCMGGetType()`, `PCMGCycleType` 13184b9ad928SBarry Smith @*/ 1319d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetType(PC pc, PCMGType form) 1320d71ae5a4SJacob Faibussowitsch { 1321f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 13224b9ad928SBarry Smith 13234b9ad928SBarry Smith PetscFunctionBegin; 13240700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1325c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc, form, 2); 132631567311SBarry Smith mg->am = form; 13279dcbbd2bSBarry Smith if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG; 1328c60c7ad4SBarry Smith else pc->ops->applyrichardson = NULL; 13293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1330c60c7ad4SBarry Smith } 1331c60c7ad4SBarry Smith 1332c60c7ad4SBarry Smith /*@ 1333f1580f4eSBarry Smith PCMGGetType - Finds the form of multigrid the `PCMG` is using multiplicative, additive, full, or the Kaskade algorithm. 1334c60c7ad4SBarry Smith 1335c3339decSBarry Smith Logically Collective 1336c60c7ad4SBarry Smith 1337c60c7ad4SBarry Smith Input Parameter: 1338c60c7ad4SBarry Smith . pc - the preconditioner context 1339c60c7ad4SBarry Smith 1340c60c7ad4SBarry Smith Output Parameter: 1341f1580f4eSBarry Smith . type - one of `PC_MG_MULTIPLICATIVE`, `PC_MG_ADDITIVE`, `PC_MG_FULL`, `PC_MG_KASKADE`, `PCMGCycleType` 1342c60c7ad4SBarry Smith 1343c60c7ad4SBarry Smith Level: advanced 1344c60c7ad4SBarry Smith 1345*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMGType`, `PCMG`, `PCMGGetLevels()`, `PCMGSetLevels()`, `PCMGSetType()` 1346c60c7ad4SBarry Smith @*/ 1347d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetType(PC pc, PCMGType *type) 1348d71ae5a4SJacob Faibussowitsch { 1349c60c7ad4SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 1350c60c7ad4SBarry Smith 1351c60c7ad4SBarry Smith PetscFunctionBegin; 1352c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1353c60c7ad4SBarry Smith *type = mg->am; 13543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13554b9ad928SBarry Smith } 13564b9ad928SBarry Smith 13574b9ad928SBarry Smith /*@ 1358f1580f4eSBarry Smith PCMGSetCycleType - Sets the type cycles to use. Use `PCMGSetCycleTypeOnLevel()` for more 13594b9ad928SBarry Smith complicated cycling. 13604b9ad928SBarry Smith 1361c3339decSBarry Smith Logically Collective 13624b9ad928SBarry Smith 13634b9ad928SBarry Smith Input Parameters: 1364c2be2410SBarry Smith + pc - the multigrid context 1365f1580f4eSBarry Smith - n - either `PC_MG_CYCLE_V` or `PC_MG_CYCLE_W` 13664b9ad928SBarry Smith 13674b9ad928SBarry Smith Options Database Key: 1368c1cbb1deSBarry Smith . -pc_mg_cycle_type <v,w> - provide the cycle desired 13694b9ad928SBarry Smith 13704b9ad928SBarry Smith Level: advanced 13714b9ad928SBarry Smith 1372*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGSetCycleTypeOnLevel()`, `PCMGType`, `PCMGCycleType` 13734b9ad928SBarry Smith @*/ 1374d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetCycleType(PC pc, PCMGCycleType n) 1375d71ae5a4SJacob Faibussowitsch { 1376f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 1377f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 137879416396SBarry Smith PetscInt i, levels; 13794b9ad928SBarry Smith 13804b9ad928SBarry Smith PetscFunctionBegin; 13810700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 138269fbec6eSBarry Smith PetscValidLogicalCollectiveEnum(pc, n, 2); 138328b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ORDER, "Must set MG levels with PCMGSetLevels() before calling"); 1384f3fbd535SBarry Smith levels = mglevels[0]->levels; 13852fa5cd67SKarl Rupp for (i = 0; i < levels; i++) mglevels[i]->cycles = n; 13863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13874b9ad928SBarry Smith } 13884b9ad928SBarry Smith 13898cc2d5dfSBarry Smith /*@ 13908cc2d5dfSBarry Smith PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step 1391f1580f4eSBarry Smith of multigrid when `PCMGType` is `PC_MG_MULTIPLICATIVE` 13928cc2d5dfSBarry Smith 1393c3339decSBarry Smith Logically Collective 13948cc2d5dfSBarry Smith 13958cc2d5dfSBarry Smith Input Parameters: 13968cc2d5dfSBarry Smith + pc - the multigrid context 13978cc2d5dfSBarry Smith - n - number of cycles (default is 1) 13988cc2d5dfSBarry Smith 13998cc2d5dfSBarry Smith Options Database Key: 1400147403d9SBarry Smith . -pc_mg_multiplicative_cycles n - set the number of cycles 14018cc2d5dfSBarry Smith 14028cc2d5dfSBarry Smith Level: advanced 14038cc2d5dfSBarry Smith 1404f1580f4eSBarry Smith Note: 1405f1580f4eSBarry Smith This is not associated with setting a v or w cycle, that is set with `PCMGSetCycleType()` 14068cc2d5dfSBarry Smith 1407*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMGSetCycleTypeOnLevel()`, `PCMGSetCycleType()`, `PCMGCycleType`, `PCMGType` 14088cc2d5dfSBarry Smith @*/ 1409d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGMultiplicativeSetCycles(PC pc, PetscInt n) 1410d71ae5a4SJacob Faibussowitsch { 1411f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 14128cc2d5dfSBarry Smith 14138cc2d5dfSBarry Smith PetscFunctionBegin; 14140700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1415c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc, n, 2); 14162a21e185SBarry Smith mg->cyclesperpcapply = n; 14173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14188cc2d5dfSBarry Smith } 14198cc2d5dfSBarry Smith 142066976f2fSJacob Faibussowitsch static PetscErrorCode PCMGSetGalerkin_MG(PC pc, PCMGGalerkinType use) 1421d71ae5a4SJacob Faibussowitsch { 1422fb15c04eSBarry Smith PC_MG *mg = (PC_MG *)pc->data; 1423fb15c04eSBarry Smith 1424fb15c04eSBarry Smith PetscFunctionBegin; 14252134b1e4SBarry Smith mg->galerkin = use; 14263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1427fb15c04eSBarry Smith } 1428fb15c04eSBarry Smith 1429c2be2410SBarry Smith /*@ 143097177400SBarry Smith PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the 143182b87d87SMatthew G. Knepley finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i 1432c2be2410SBarry Smith 1433c3339decSBarry Smith Logically Collective 1434c2be2410SBarry Smith 1435c2be2410SBarry Smith Input Parameters: 1436c91913d3SJed Brown + pc - the multigrid context 1437f1580f4eSBarry Smith - use - one of `PC_MG_GALERKIN_BOTH`, `PC_MG_GALERKIN_PMAT`, `PC_MG_GALERKIN_MAT`, or `PC_MG_GALERKIN_NONE` 1438c2be2410SBarry Smith 1439c2be2410SBarry Smith Options Database Key: 1440147403d9SBarry Smith . -pc_mg_galerkin <both,pmat,mat,none> - set the matrices to form via the Galerkin process 1441c2be2410SBarry Smith 1442c2be2410SBarry Smith Level: intermediate 1443c2be2410SBarry Smith 1444f1580f4eSBarry Smith Note: 1445f1580f4eSBarry Smith Some codes that use `PCMG` such as `PCGAMG` use Galerkin internally while constructing the hierarchy and thus do not 1446f1580f4eSBarry Smith use the `PCMG` construction of the coarser grids. 14471c1aac46SBarry Smith 1448*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGGetGalerkin()`, `PCMGGalerkinType` 1449c2be2410SBarry Smith @*/ 1450d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetGalerkin(PC pc, PCMGGalerkinType use) 1451d71ae5a4SJacob Faibussowitsch { 1452c2be2410SBarry Smith PetscFunctionBegin; 14530700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1454cac4c232SBarry Smith PetscTryMethod(pc, "PCMGSetGalerkin_C", (PC, PCMGGalerkinType), (pc, use)); 14553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1456c2be2410SBarry Smith } 1457c2be2410SBarry Smith 14583fc8bf9cSBarry Smith /*@ 1459f1580f4eSBarry Smith PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. A_i-1 = r_i * A_i * p_i 14603fc8bf9cSBarry Smith 14613fc8bf9cSBarry Smith Not Collective 14623fc8bf9cSBarry Smith 14633fc8bf9cSBarry Smith Input Parameter: 14643fc8bf9cSBarry Smith . pc - the multigrid context 14653fc8bf9cSBarry Smith 14663fc8bf9cSBarry Smith Output Parameter: 1467f1580f4eSBarry 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` 14683fc8bf9cSBarry Smith 14693fc8bf9cSBarry Smith Level: intermediate 14703fc8bf9cSBarry Smith 1471*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGSetGalerkin()`, `PCMGGalerkinType` 14723fc8bf9cSBarry Smith @*/ 1473d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetGalerkin(PC pc, PCMGGalerkinType *galerkin) 1474d71ae5a4SJacob Faibussowitsch { 1475f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 14763fc8bf9cSBarry Smith 14773fc8bf9cSBarry Smith PetscFunctionBegin; 14780700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 14792134b1e4SBarry Smith *galerkin = mg->galerkin; 14803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14813fc8bf9cSBarry Smith } 14823fc8bf9cSBarry Smith 148366976f2fSJacob Faibussowitsch static PetscErrorCode PCMGSetAdaptInterpolation_MG(PC pc, PetscBool adapt) 1484d71ae5a4SJacob Faibussowitsch { 1485f3b08a26SMatthew G. Knepley PC_MG *mg = (PC_MG *)pc->data; 1486f3b08a26SMatthew G. Knepley 1487f3b08a26SMatthew G. Knepley PetscFunctionBegin; 1488f3b08a26SMatthew G. Knepley mg->adaptInterpolation = adapt; 14893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1490f3b08a26SMatthew G. Knepley } 1491f3b08a26SMatthew G. Knepley 149266976f2fSJacob Faibussowitsch static PetscErrorCode PCMGGetAdaptInterpolation_MG(PC pc, PetscBool *adapt) 1493d71ae5a4SJacob Faibussowitsch { 1494f3b08a26SMatthew G. Knepley PC_MG *mg = (PC_MG *)pc->data; 1495f3b08a26SMatthew G. Knepley 1496f3b08a26SMatthew G. Knepley PetscFunctionBegin; 1497f3b08a26SMatthew G. Knepley *adapt = mg->adaptInterpolation; 14983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1499f3b08a26SMatthew G. Knepley } 1500f3b08a26SMatthew G. Knepley 150166976f2fSJacob Faibussowitsch static PetscErrorCode PCMGSetAdaptCoarseSpaceType_MG(PC pc, PCMGCoarseSpaceType ctype) 1502d71ae5a4SJacob Faibussowitsch { 15032b3cbbdaSStefano Zampini PC_MG *mg = (PC_MG *)pc->data; 15042b3cbbdaSStefano Zampini 15052b3cbbdaSStefano Zampini PetscFunctionBegin; 15062b3cbbdaSStefano Zampini mg->adaptInterpolation = ctype != PCMG_ADAPT_NONE ? PETSC_TRUE : PETSC_FALSE; 15072b3cbbdaSStefano Zampini mg->coarseSpaceType = ctype; 15083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15092b3cbbdaSStefano Zampini } 15102b3cbbdaSStefano Zampini 151166976f2fSJacob Faibussowitsch static PetscErrorCode PCMGGetAdaptCoarseSpaceType_MG(PC pc, PCMGCoarseSpaceType *ctype) 1512d71ae5a4SJacob Faibussowitsch { 15132b3cbbdaSStefano Zampini PC_MG *mg = (PC_MG *)pc->data; 15142b3cbbdaSStefano Zampini 15152b3cbbdaSStefano Zampini PetscFunctionBegin; 15162b3cbbdaSStefano Zampini *ctype = mg->coarseSpaceType; 15173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15182b3cbbdaSStefano Zampini } 15192b3cbbdaSStefano Zampini 152066976f2fSJacob Faibussowitsch static PetscErrorCode PCMGSetAdaptCR_MG(PC pc, PetscBool cr) 1521d71ae5a4SJacob Faibussowitsch { 152241b6fd38SMatthew G. Knepley PC_MG *mg = (PC_MG *)pc->data; 152341b6fd38SMatthew G. Knepley 152441b6fd38SMatthew G. Knepley PetscFunctionBegin; 152541b6fd38SMatthew G. Knepley mg->compatibleRelaxation = cr; 15263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 152741b6fd38SMatthew G. Knepley } 152841b6fd38SMatthew G. Knepley 152966976f2fSJacob Faibussowitsch static PetscErrorCode PCMGGetAdaptCR_MG(PC pc, PetscBool *cr) 1530d71ae5a4SJacob Faibussowitsch { 153141b6fd38SMatthew G. Knepley PC_MG *mg = (PC_MG *)pc->data; 153241b6fd38SMatthew G. Knepley 153341b6fd38SMatthew G. Knepley PetscFunctionBegin; 153441b6fd38SMatthew G. Knepley *cr = mg->compatibleRelaxation; 15353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 153641b6fd38SMatthew G. Knepley } 153741b6fd38SMatthew G. Knepley 15382b3cbbdaSStefano Zampini /*@C 15392b3cbbdaSStefano Zampini PCMGSetAdaptCoarseSpaceType - Set the type of adaptive coarse space. 15402b3cbbdaSStefano Zampini 15412b3cbbdaSStefano 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. 15422b3cbbdaSStefano Zampini 1543c3339decSBarry Smith Logically Collective 15442b3cbbdaSStefano Zampini 15452b3cbbdaSStefano Zampini Input Parameters: 15462b3cbbdaSStefano Zampini + pc - the multigrid context 15472b3cbbdaSStefano Zampini - ctype - the type of coarse space 15482b3cbbdaSStefano Zampini 15492b3cbbdaSStefano Zampini Options Database Keys: 15502b3cbbdaSStefano Zampini + -pc_mg_adapt_interp_n <int> - The number of modes to use 15512b3cbbdaSStefano Zampini - -pc_mg_adapt_interp_coarse_space <type> - The type of coarse space: none, polynomial, harmonic, eigenvector, generalized_eigenvector, gdsw 15522b3cbbdaSStefano Zampini 15532b3cbbdaSStefano Zampini Level: intermediate 15542b3cbbdaSStefano Zampini 1555*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGCoarseSpaceType`, `PCMGGetAdaptCoarseSpaceType()`, `PCMGSetGalerkin()`, `PCMGSetAdaptInterpolation()` 15562b3cbbdaSStefano Zampini @*/ 1557d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetAdaptCoarseSpaceType(PC pc, PCMGCoarseSpaceType ctype) 1558d71ae5a4SJacob Faibussowitsch { 15592b3cbbdaSStefano Zampini PetscFunctionBegin; 15602b3cbbdaSStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 15612b3cbbdaSStefano Zampini PetscValidLogicalCollectiveEnum(pc, ctype, 2); 15622b3cbbdaSStefano Zampini PetscTryMethod(pc, "PCMGSetAdaptCoarseSpaceType_C", (PC, PCMGCoarseSpaceType), (pc, ctype)); 15633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15642b3cbbdaSStefano Zampini } 15652b3cbbdaSStefano Zampini 15662b3cbbdaSStefano Zampini /*@C 15672b3cbbdaSStefano Zampini PCMGGetAdaptCoarseSpaceType - Get the type of adaptive coarse space. 15682b3cbbdaSStefano Zampini 15692b3cbbdaSStefano Zampini Not Collective 15702b3cbbdaSStefano Zampini 15712b3cbbdaSStefano Zampini Input Parameter: 15722b3cbbdaSStefano Zampini . pc - the multigrid context 15732b3cbbdaSStefano Zampini 15742b3cbbdaSStefano Zampini Output Parameter: 15752b3cbbdaSStefano Zampini . ctype - the type of coarse space 15762b3cbbdaSStefano Zampini 15772b3cbbdaSStefano Zampini Level: intermediate 15782b3cbbdaSStefano Zampini 1579*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGCoarseSpaceType`, `PCMGSetAdaptCoarseSpaceType()`, `PCMGSetGalerkin()`, `PCMGSetAdaptInterpolation()` 15802b3cbbdaSStefano Zampini @*/ 1581d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetAdaptCoarseSpaceType(PC pc, PCMGCoarseSpaceType *ctype) 1582d71ae5a4SJacob Faibussowitsch { 15832b3cbbdaSStefano Zampini PetscFunctionBegin; 15842b3cbbdaSStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 15854f572ea9SToby Isaac PetscAssertPointer(ctype, 2); 15862b3cbbdaSStefano Zampini PetscUseMethod(pc, "PCMGGetAdaptCoarseSpaceType_C", (PC, PCMGCoarseSpaceType *), (pc, ctype)); 15873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15882b3cbbdaSStefano Zampini } 15892b3cbbdaSStefano Zampini 15902b3cbbdaSStefano Zampini /* MATT: REMOVE? */ 1591f3b08a26SMatthew G. Knepley /*@ 1592f3b08a26SMatthew 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. 1593f3b08a26SMatthew G. Knepley 1594c3339decSBarry Smith Logically Collective 1595f3b08a26SMatthew G. Knepley 1596f3b08a26SMatthew G. Knepley Input Parameters: 1597f3b08a26SMatthew G. Knepley + pc - the multigrid context 1598f3b08a26SMatthew G. Knepley - adapt - flag for adaptation of the interpolator 1599f3b08a26SMatthew G. Knepley 1600f3b08a26SMatthew G. Knepley Options Database Keys: 1601f3b08a26SMatthew G. Knepley + -pc_mg_adapt_interp - Turn on adaptation 1602f3b08a26SMatthew G. Knepley . -pc_mg_adapt_interp_n <int> - The number of modes to use, should be divisible by dimension 1603f3b08a26SMatthew G. Knepley - -pc_mg_adapt_interp_coarse_space <type> - The type of coarse space: polynomial, harmonic, eigenvector, generalized_eigenvector 1604f3b08a26SMatthew G. Knepley 1605f3b08a26SMatthew G. Knepley Level: intermediate 1606f3b08a26SMatthew G. Knepley 1607*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGGetAdaptInterpolation()`, `PCMGSetGalerkin()`, `PCMGGetAdaptCoarseSpaceType()`, `PCMGSetAdaptCoarseSpaceType()` 1608f3b08a26SMatthew G. Knepley @*/ 1609d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetAdaptInterpolation(PC pc, PetscBool adapt) 1610d71ae5a4SJacob Faibussowitsch { 1611f3b08a26SMatthew G. Knepley PetscFunctionBegin; 1612f3b08a26SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1613cac4c232SBarry Smith PetscTryMethod(pc, "PCMGSetAdaptInterpolation_C", (PC, PetscBool), (pc, adapt)); 16143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1615f3b08a26SMatthew G. Knepley } 1616f3b08a26SMatthew G. Knepley 1617f3b08a26SMatthew G. Knepley /*@ 1618f1580f4eSBarry Smith PCMGGetAdaptInterpolation - Get the flag to adapt the interpolator based upon a vector space which should be accurately captured by the next coarser mesh, 1619f1580f4eSBarry Smith and thus accurately interpolated. 1620f3b08a26SMatthew G. Knepley 162120f4b53cSBarry Smith Not Collective 1622f3b08a26SMatthew G. Knepley 1623f3b08a26SMatthew G. Knepley Input Parameter: 1624f3b08a26SMatthew G. Knepley . pc - the multigrid context 1625f3b08a26SMatthew G. Knepley 1626f3b08a26SMatthew G. Knepley Output Parameter: 1627f3b08a26SMatthew G. Knepley . adapt - flag for adaptation of the interpolator 1628f3b08a26SMatthew G. Knepley 1629f3b08a26SMatthew G. Knepley Level: intermediate 1630f3b08a26SMatthew G. Knepley 1631*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGSetAdaptInterpolation()`, `PCMGSetGalerkin()`, `PCMGGetAdaptCoarseSpaceType()`, `PCMGSetAdaptCoarseSpaceType()` 1632f3b08a26SMatthew G. Knepley @*/ 1633d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetAdaptInterpolation(PC pc, PetscBool *adapt) 1634d71ae5a4SJacob Faibussowitsch { 1635f3b08a26SMatthew G. Knepley PetscFunctionBegin; 1636f3b08a26SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 16374f572ea9SToby Isaac PetscAssertPointer(adapt, 2); 1638cac4c232SBarry Smith PetscUseMethod(pc, "PCMGGetAdaptInterpolation_C", (PC, PetscBool *), (pc, adapt)); 16393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1640f3b08a26SMatthew G. Knepley } 1641f3b08a26SMatthew G. Knepley 16424b9ad928SBarry Smith /*@ 164341b6fd38SMatthew G. Knepley PCMGSetAdaptCR - Monitor the coarse space quality using an auxiliary solve with compatible relaxation. 164441b6fd38SMatthew G. Knepley 1645c3339decSBarry Smith Logically Collective 164641b6fd38SMatthew G. Knepley 164741b6fd38SMatthew G. Knepley Input Parameters: 164841b6fd38SMatthew G. Knepley + pc - the multigrid context 164941b6fd38SMatthew G. Knepley - cr - flag for compatible relaxation 165041b6fd38SMatthew G. Knepley 1651f1580f4eSBarry Smith Options Database Key: 165241b6fd38SMatthew G. Knepley . -pc_mg_adapt_cr - Turn on compatible relaxation 165341b6fd38SMatthew G. Knepley 165441b6fd38SMatthew G. Knepley Level: intermediate 165541b6fd38SMatthew G. Knepley 1656*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGGetAdaptCR()`, `PCMGSetAdaptInterpolation()`, `PCMGSetGalerkin()`, `PCMGGetAdaptCoarseSpaceType()`, `PCMGSetAdaptCoarseSpaceType()` 165741b6fd38SMatthew G. Knepley @*/ 1658d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetAdaptCR(PC pc, PetscBool cr) 1659d71ae5a4SJacob Faibussowitsch { 166041b6fd38SMatthew G. Knepley PetscFunctionBegin; 166141b6fd38SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1662cac4c232SBarry Smith PetscTryMethod(pc, "PCMGSetAdaptCR_C", (PC, PetscBool), (pc, cr)); 16633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 166441b6fd38SMatthew G. Knepley } 166541b6fd38SMatthew G. Knepley 166641b6fd38SMatthew G. Knepley /*@ 166741b6fd38SMatthew G. Knepley PCMGGetAdaptCR - Get the flag to monitor coarse space quality using an auxiliary solve with compatible relaxation. 166841b6fd38SMatthew G. Knepley 166920f4b53cSBarry Smith Not Collective 167041b6fd38SMatthew G. Knepley 167141b6fd38SMatthew G. Knepley Input Parameter: 167241b6fd38SMatthew G. Knepley . pc - the multigrid context 167341b6fd38SMatthew G. Knepley 167441b6fd38SMatthew G. Knepley Output Parameter: 167541b6fd38SMatthew G. Knepley . cr - flag for compatible relaxaion 167641b6fd38SMatthew G. Knepley 167741b6fd38SMatthew G. Knepley Level: intermediate 167841b6fd38SMatthew G. Knepley 1679*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMGSetAdaptCR()`, `PCMGGetAdaptInterpolation()`, `PCMGSetGalerkin()`, `PCMGGetAdaptCoarseSpaceType()`, `PCMGSetAdaptCoarseSpaceType()` 168041b6fd38SMatthew G. Knepley @*/ 1681d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetAdaptCR(PC pc, PetscBool *cr) 1682d71ae5a4SJacob Faibussowitsch { 168341b6fd38SMatthew G. Knepley PetscFunctionBegin; 168441b6fd38SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 16854f572ea9SToby Isaac PetscAssertPointer(cr, 2); 1686cac4c232SBarry Smith PetscUseMethod(pc, "PCMGGetAdaptCR_C", (PC, PetscBool *), (pc, cr)); 16873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 168841b6fd38SMatthew G. Knepley } 168941b6fd38SMatthew G. Knepley 169041b6fd38SMatthew G. Knepley /*@ 169106792cafSBarry Smith PCMGSetNumberSmooth - Sets the number of pre and post-smoothing steps to use 1692f1580f4eSBarry Smith on all levels. Use `PCMGDistinctSmoothUp()` to create separate up and down smoothers if you want different numbers of 1693710315b6SLawrence Mitchell pre- and post-smoothing steps. 169406792cafSBarry Smith 1695c3339decSBarry Smith Logically Collective 169606792cafSBarry Smith 169706792cafSBarry Smith Input Parameters: 1698feefa0e1SJacob Faibussowitsch + pc - the multigrid context 169906792cafSBarry Smith - n - the number of smoothing steps 170006792cafSBarry Smith 170106792cafSBarry Smith Options Database Key: 1702a2b725a8SWilliam Gropp . -mg_levels_ksp_max_it <n> - Sets number of pre and post-smoothing steps 170306792cafSBarry Smith 170406792cafSBarry Smith Level: advanced 170506792cafSBarry Smith 1706f1580f4eSBarry Smith Note: 1707f1580f4eSBarry 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. 170806792cafSBarry Smith 1709*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGSetDistinctSmoothUp()` 171006792cafSBarry Smith @*/ 1711d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetNumberSmooth(PC pc, PetscInt n) 1712d71ae5a4SJacob Faibussowitsch { 171306792cafSBarry Smith PC_MG *mg = (PC_MG *)pc->data; 171406792cafSBarry Smith PC_MG_Levels **mglevels = mg->levels; 171506792cafSBarry Smith PetscInt i, levels; 171606792cafSBarry Smith 171706792cafSBarry Smith PetscFunctionBegin; 171806792cafSBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 171906792cafSBarry Smith PetscValidLogicalCollectiveInt(pc, n, 2); 172028b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ORDER, "Must set MG levels with PCMGSetLevels() before calling"); 172106792cafSBarry Smith levels = mglevels[0]->levels; 172206792cafSBarry Smith 172306792cafSBarry Smith for (i = 1; i < levels; i++) { 17249566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(mglevels[i]->smoothu, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, n)); 17259566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(mglevels[i]->smoothd, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, n)); 172606792cafSBarry Smith mg->default_smoothu = n; 172706792cafSBarry Smith mg->default_smoothd = n; 172806792cafSBarry Smith } 17293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 173006792cafSBarry Smith } 173106792cafSBarry Smith 1732f442ab6aSBarry Smith /*@ 1733f1580f4eSBarry Smith PCMGSetDistinctSmoothUp - sets the up (post) smoother to be a separate `KSP` from the down (pre) smoother on all levels 1734710315b6SLawrence Mitchell and adds the suffix _up to the options name 1735f442ab6aSBarry Smith 1736c3339decSBarry Smith Logically Collective 1737f442ab6aSBarry Smith 1738f1580f4eSBarry Smith Input Parameter: 1739f442ab6aSBarry Smith . pc - the preconditioner context 1740f442ab6aSBarry Smith 1741f442ab6aSBarry Smith Options Database Key: 1742147403d9SBarry Smith . -pc_mg_distinct_smoothup <bool> - use distinct smoothing objects 1743f442ab6aSBarry Smith 1744f442ab6aSBarry Smith Level: advanced 1745f442ab6aSBarry Smith 1746f1580f4eSBarry Smith Note: 1747f1580f4eSBarry 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. 1748f442ab6aSBarry Smith 1749*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGSetNumberSmooth()` 1750f442ab6aSBarry Smith @*/ 1751d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetDistinctSmoothUp(PC pc) 1752d71ae5a4SJacob Faibussowitsch { 1753f442ab6aSBarry Smith PC_MG *mg = (PC_MG *)pc->data; 1754f442ab6aSBarry Smith PC_MG_Levels **mglevels = mg->levels; 1755f442ab6aSBarry Smith PetscInt i, levels; 1756f442ab6aSBarry Smith KSP subksp; 1757f442ab6aSBarry Smith 1758f442ab6aSBarry Smith PetscFunctionBegin; 1759f442ab6aSBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 176028b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ORDER, "Must set MG levels with PCMGSetLevels() before calling"); 1761f442ab6aSBarry Smith levels = mglevels[0]->levels; 1762f442ab6aSBarry Smith 1763f442ab6aSBarry Smith for (i = 1; i < levels; i++) { 1764710315b6SLawrence Mitchell const char *prefix = NULL; 1765f442ab6aSBarry Smith /* make sure smoother up and down are different */ 17669566063dSJacob Faibussowitsch PetscCall(PCMGGetSmootherUp(pc, i, &subksp)); 17679566063dSJacob Faibussowitsch PetscCall(KSPGetOptionsPrefix(mglevels[i]->smoothd, &prefix)); 17689566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(subksp, prefix)); 17699566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(subksp, "up_")); 1770f442ab6aSBarry Smith } 17713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1772f442ab6aSBarry Smith } 1773f442ab6aSBarry Smith 177407a4832bSFande Kong /* No new matrices are created, and the coarse operator matrices are the references to the original ones */ 177566976f2fSJacob Faibussowitsch static PetscErrorCode PCGetInterpolations_MG(PC pc, PetscInt *num_levels, Mat *interpolations[]) 1776d71ae5a4SJacob Faibussowitsch { 177707a4832bSFande Kong PC_MG *mg = (PC_MG *)pc->data; 177807a4832bSFande Kong PC_MG_Levels **mglevels = mg->levels; 177907a4832bSFande Kong Mat *mat; 178007a4832bSFande Kong PetscInt l; 178107a4832bSFande Kong 178207a4832bSFande Kong PetscFunctionBegin; 178328b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling"); 17849566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(mg->nlevels, &mat)); 178507a4832bSFande Kong for (l = 1; l < mg->nlevels; l++) { 178607a4832bSFande Kong mat[l - 1] = mglevels[l]->interpolate; 17879566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat[l - 1])); 178807a4832bSFande Kong } 178907a4832bSFande Kong *num_levels = mg->nlevels; 179007a4832bSFande Kong *interpolations = mat; 17913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 179207a4832bSFande Kong } 179307a4832bSFande Kong 179407a4832bSFande Kong /* No new matrices are created, and the coarse operator matrices are the references to the original ones */ 179566976f2fSJacob Faibussowitsch static PetscErrorCode PCGetCoarseOperators_MG(PC pc, PetscInt *num_levels, Mat *coarseOperators[]) 1796d71ae5a4SJacob Faibussowitsch { 179707a4832bSFande Kong PC_MG *mg = (PC_MG *)pc->data; 179807a4832bSFande Kong PC_MG_Levels **mglevels = mg->levels; 179907a4832bSFande Kong PetscInt l; 180007a4832bSFande Kong Mat *mat; 180107a4832bSFande Kong 180207a4832bSFande Kong PetscFunctionBegin; 180328b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling"); 18049566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(mg->nlevels, &mat)); 180507a4832bSFande Kong for (l = 0; l < mg->nlevels - 1; l++) { 18069566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(mglevels[l]->smoothd, NULL, &(mat[l]))); 18079566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat[l])); 180807a4832bSFande Kong } 180907a4832bSFande Kong *num_levels = mg->nlevels; 181007a4832bSFande Kong *coarseOperators = mat; 18113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 181207a4832bSFande Kong } 181307a4832bSFande Kong 1814f3b08a26SMatthew G. Knepley /*@C 1815f1580f4eSBarry Smith PCMGRegisterCoarseSpaceConstructor - Adds a method to the `PCMG` package for coarse space construction. 1816f3b08a26SMatthew G. Knepley 181720f4b53cSBarry Smith Not Collective 1818f3b08a26SMatthew G. Knepley 1819f3b08a26SMatthew G. Knepley Input Parameters: 1820f3b08a26SMatthew G. Knepley + name - name of the constructor 1821f3b08a26SMatthew G. Knepley - function - constructor routine 1822f3b08a26SMatthew G. Knepley 182320f4b53cSBarry Smith Calling sequence of `function`: 182420f4b53cSBarry Smith + pc - The `PC` object 1825f1580f4eSBarry Smith . l - The multigrid level, 0 is the coarse level 182620f4b53cSBarry Smith . dm - The `DM` for this level 1827f1580f4eSBarry Smith . smooth - The level smoother 1828f1580f4eSBarry Smith . Nc - The size of the coarse space 1829f1580f4eSBarry Smith . initGuess - Basis for an initial guess for the space 1830f1580f4eSBarry Smith - coarseSp - A basis for the computed coarse space 1831f3b08a26SMatthew G. Knepley 1832f3b08a26SMatthew G. Knepley Level: advanced 1833f3b08a26SMatthew G. Knepley 1834feefa0e1SJacob Faibussowitsch Developer Notes: 1835f1580f4eSBarry Smith How come this is not used by `PCGAMG`? 1836f1580f4eSBarry Smith 1837*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGGetCoarseSpaceConstructor()`, `PCRegister()` 1838f3b08a26SMatthew G. Knepley @*/ 183904c3f3b8SBarry Smith PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char name[], PetscErrorCode (*function)(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, Mat initGuess, Mat *coarseSp)) 1840d71ae5a4SJacob Faibussowitsch { 1841f3b08a26SMatthew G. Knepley PetscFunctionBegin; 18429566063dSJacob Faibussowitsch PetscCall(PCInitializePackage()); 18439566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&PCMGCoarseList, name, function)); 18443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1845f3b08a26SMatthew G. Knepley } 1846f3b08a26SMatthew G. Knepley 1847f3b08a26SMatthew G. Knepley /*@C 1848f3b08a26SMatthew G. Knepley PCMGGetCoarseSpaceConstructor - Returns the given coarse space construction method. 1849f3b08a26SMatthew G. Knepley 185020f4b53cSBarry Smith Not Collective 1851f3b08a26SMatthew G. Knepley 1852f3b08a26SMatthew G. Knepley Input Parameter: 1853f3b08a26SMatthew G. Knepley . name - name of the constructor 1854f3b08a26SMatthew G. Knepley 1855f3b08a26SMatthew G. Knepley Output Parameter: 1856f3b08a26SMatthew G. Knepley . function - constructor routine 1857f3b08a26SMatthew G. Knepley 1858f3b08a26SMatthew G. Knepley Level: advanced 1859f3b08a26SMatthew G. Knepley 1860*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGRegisterCoarseSpaceConstructor()`, `PCRegister()` 1861f3b08a26SMatthew G. Knepley @*/ 1862d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetCoarseSpaceConstructor(const char name[], PetscErrorCode (**function)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat *)) 1863d71ae5a4SJacob Faibussowitsch { 1864f3b08a26SMatthew G. Knepley PetscFunctionBegin; 18659566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(PCMGCoarseList, name, function)); 18663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1867f3b08a26SMatthew G. Knepley } 1868f3b08a26SMatthew G. Knepley 18693b09bd56SBarry Smith /*MC 1870ccb205f8SBarry Smith PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional 18713b09bd56SBarry Smith information about the coarser grid matrices and restriction/interpolation operators. 18723b09bd56SBarry Smith 18733b09bd56SBarry Smith Options Database Keys: 18743b09bd56SBarry Smith + -pc_mg_levels <nlevels> - number of levels including finest 1875391689abSStefano Zampini . -pc_mg_cycle_type <v,w> - provide the cycle desired 18768c1c2452SJed Brown . -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default 18773b09bd56SBarry Smith . -pc_mg_log - log information about time spent on each level of the solver 1878710315b6SLawrence Mitchell . -pc_mg_distinct_smoothup - configure up (after interpolation) and down (before restriction) smoothers separately (with different options prefixes) 18792134b1e4SBarry Smith . -pc_mg_galerkin <both,pmat,mat,none> - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R' 18808cf6037eSBarry Smith . -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1) 18818cf6037eSBarry Smith . -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 1882e3c5b3baSBarry Smith to the Socket viewer for reading from MATLAB. 18838cf6037eSBarry Smith - -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices 18848cf6037eSBarry Smith to the binary output file called binaryoutput 18853b09bd56SBarry Smith 188620f4b53cSBarry Smith Level: intermediate 188720f4b53cSBarry Smith 188895452b02SPatrick Sanan Notes: 188904c3f3b8SBarry Smith The Krylov solver (if any) and preconditioner (smoother) and their parameters are controlled from the options database with the standard 189004c3f3b8SBarry Smith options database keywords prefixed with `-mg_levels_` to affect all the levels but the coarsest, which is controlled with `-mg_coarse_`. 189104c3f3b8SBarry Smith One can set different preconditioners etc on specific levels with the prefix `-mg_levels_n_` where `n` is the level number (zero being 189204c3f3b8SBarry Smith the coarse level. For example 189304c3f3b8SBarry Smith .vb 189404c3f3b8SBarry Smith -mg_levels_ksp_type gmres -mg_levels_pc_type bjacobi -mg_coarse_pc_type svd -mg_levels_2_pc_type sor 189504c3f3b8SBarry Smith .ve 189604c3f3b8SBarry Smith These options also work for controlling the smoothers etc inside `PCGAMG` 189704c3f3b8SBarry Smith 1898f1580f4eSBarry 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 18993b09bd56SBarry Smith 19008cf6037eSBarry Smith When run with a single level the smoother options are used on that level NOT the coarse grid solver options 19018cf6037eSBarry Smith 1902f1580f4eSBarry Smith When run with `KSPRICHARDSON` the convergence test changes slightly if monitor is turned on. The iteration count may change slightly. This 190323067569SBarry Smith is because without monitoring the residual norm is computed WITHIN each multigrid cycle on the finest level after the pre-smoothing 190423067569SBarry Smith (because the residual has just been computed for the multigrid algorithm and is hence available for free) while with monitoring the 190523067569SBarry Smith residual is computed at the end of each cycle. 190623067569SBarry Smith 190704c3f3b8SBarry Smith .seealso: [](sec_mg), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCMGType`, `PCEXOTIC`, `PCGAMG`, `PCML`, `PCHYPRE` 1908db781477SPatrick Sanan `PCMGSetLevels()`, `PCMGGetLevels()`, `PCMGSetType()`, `PCMGSetCycleType()`, 1909db781477SPatrick Sanan `PCMGSetDistinctSmoothUp()`, `PCMGGetCoarseSolve()`, `PCMGSetResidual()`, `PCMGSetInterpolation()`, 1910db781477SPatrick Sanan `PCMGSetRestriction()`, `PCMGGetSmoother()`, `PCMGGetSmootherUp()`, `PCMGGetSmootherDown()`, 1911f1580f4eSBarry Smith `PCMGSetCycleTypeOnLevel()`, `PCMGSetRhs()`, `PCMGSetX()`, `PCMGSetR()`, 1912f1580f4eSBarry Smith `PCMGSetAdaptCR()`, `PCMGGetAdaptInterpolation()`, `PCMGSetGalerkin()`, `PCMGGetAdaptCoarseSpaceType()`, `PCMGSetAdaptCoarseSpaceType()` 19133b09bd56SBarry Smith M*/ 19143b09bd56SBarry Smith 1915d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc) 1916d71ae5a4SJacob Faibussowitsch { 1917f3fbd535SBarry Smith PC_MG *mg; 1918f3fbd535SBarry Smith 19194b9ad928SBarry Smith PetscFunctionBegin; 19204dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&mg)); 19213ec1f749SStefano Zampini pc->data = mg; 1922f3fbd535SBarry Smith mg->nlevels = -1; 192310eca3edSLisandro Dalcin mg->am = PC_MG_MULTIPLICATIVE; 19242134b1e4SBarry Smith mg->galerkin = PC_MG_GALERKIN_NONE; 1925f3b08a26SMatthew G. Knepley mg->adaptInterpolation = PETSC_FALSE; 1926f3b08a26SMatthew G. Knepley mg->Nc = -1; 1927f3b08a26SMatthew G. Knepley mg->eigenvalue = -1; 1928f3fbd535SBarry Smith 192937a44384SMark Adams pc->useAmat = PETSC_TRUE; 193037a44384SMark Adams 19314b9ad928SBarry Smith pc->ops->apply = PCApply_MG; 1932fcb023d4SJed Brown pc->ops->applytranspose = PCApplyTranspose_MG; 193330b0564aSStefano Zampini pc->ops->matapply = PCMatApply_MG; 19344b9ad928SBarry Smith pc->ops->setup = PCSetUp_MG; 1935a06653b4SBarry Smith pc->ops->reset = PCReset_MG; 19364b9ad928SBarry Smith pc->ops->destroy = PCDestroy_MG; 19374b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_MG; 19384b9ad928SBarry Smith pc->ops->view = PCView_MG; 1939fb15c04eSBarry Smith 19409566063dSJacob Faibussowitsch PetscCall(PetscObjectComposedDataRegister(&mg->eigenvalue)); 19419566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetGalerkin_C", PCMGSetGalerkin_MG)); 19429566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetLevels_C", PCMGGetLevels_MG)); 19439566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetLevels_C", PCMGSetLevels_MG)); 19449566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGetInterpolations_C", PCGetInterpolations_MG)); 19459566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGetCoarseOperators_C", PCGetCoarseOperators_MG)); 19469566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetAdaptInterpolation_C", PCMGSetAdaptInterpolation_MG)); 19479566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetAdaptInterpolation_C", PCMGGetAdaptInterpolation_MG)); 19489566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetAdaptCR_C", PCMGSetAdaptCR_MG)); 19499566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetAdaptCR_C", PCMGGetAdaptCR_MG)); 19502b3cbbdaSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetAdaptCoarseSpaceType_C", PCMGSetAdaptCoarseSpaceType_MG)); 19512b3cbbdaSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetAdaptCoarseSpaceType_C", PCMGGetAdaptCoarseSpaceType_MG)); 19523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 19534b9ad928SBarry Smith } 1954